Salome HOME
Copyright update 2022
[modules/paravis.git] / src / Plugins / ParaMEDCorba / plugin / ParaMEDMEM2VTK / VTKMEDCouplingCurveLinearMeshClient.cxx
1 // Copyright (C) 2010-2022  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "VTKMEDCouplingCurveLinearMeshClient.hxx"
21
22 #include "vtkPoints.h"
23 #include "vtkErrorCode.h"
24 #include "vtkDoubleArray.h"
25 #include "vtkStructuredGrid.h"
26
27 #include <vector>
28 #include <string>
29
30 void ParaMEDMEM2VTK::FillMEDCouplingCurveLinearMeshInstanceFrom(SALOME_MED::MEDCouplingCurveLinearMeshCorbaInterface_ptr meshPtr, vtkStructuredGrid *ret)
31 {
32   meshPtr->Register();
33   SALOME_TYPES::ListOfDouble *tinyD;
34   SALOME_TYPES::ListOfLong *tinyI;
35   SALOME_TYPES::ListOfString *tinyS;
36   meshPtr->getTinyInfo(tinyD,tinyI,tinyS);
37   int meshStr[3]={1,1,1};
38   int meshDim=(*tinyI)[2];
39   if(meshDim<0 || meshDim>3)
40     vtkErrorWithObjectMacro(ret,"Internal error in ParaMEDCorba plugin : distant curvilinear mesh has a mesh dimension not in [0,3] !");
41   for(int i=0;i<meshDim;i++)
42     meshStr[i]=(*tinyI)[3+i];
43   ret->SetDimensions(meshStr);
44   int nbOfNodes=(*tinyI)[3+meshDim];
45   int spaceDim=(*tinyI)[3+meshDim+1];
46   delete tinyD;
47   delete tinyI;
48   delete tinyS;
49   //
50   vtkDoubleArray *da=vtkDoubleArray::New();
51   da->SetNumberOfComponents(3);
52   da->SetNumberOfTuples(nbOfNodes);
53   double *pt=da->GetPointer(0);
54   SALOME_TYPES::ListOfDouble *bigD;
55   meshPtr->getSerialisationData(tinyI,bigD);
56   delete tinyI;
57   if((int)bigD->length()!=nbOfNodes*spaceDim)
58     vtkErrorWithObjectMacro(ret,"Internal error in ParaMEDCorba plugin : distant curvilinear mesh, mismatch between informations ! Internal error !");
59   for(int i=0;i<nbOfNodes;i++)
60     {
61       for(int j=0;j<spaceDim;j++,pt++)
62         *pt=(*bigD)[spaceDim*i+j];
63       for(int j=spaceDim;j<3;j++,pt++)
64         *pt=0.;
65     }
66   delete bigD;
67   vtkPoints *points=vtkPoints::New();
68   ret->SetPoints(points);
69   points->SetData(da);
70   points->Delete();
71   da->Delete();
72   meshPtr->UnRegister();
73 }