1 // Copyright (C) 2010-2016 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 #include "VTKMEDCouplingUMeshClient.hxx"
22 #include "vtkPoints.h"
23 #include "vtkCellArray.h"
24 #include "vtkDoubleArray.h"
25 #include "vtkSmartPointer.h"
26 #include "vtkUnstructuredGrid.h"
33 static const int ParaMEDMEM2VTKTypeTraducer[34]={1,3,21,5,9,7,22,34,23,28,-1,-1,-1,-1,10,14,13,-1,12,-1,24,-1,16,27,-1,26,-1,29,-1,-1,25,42,-1,4};
35 void ParaMEDMEM2VTK::FillMEDCouplingUMeshInstanceFrom(SALOME_MED::MEDCouplingUMeshCorbaInterface_ptr meshPtr, vtkUnstructuredGrid *ret, bool& isPolyh)
39 SALOME_TYPES::ListOfDouble *tinyD;
40 SALOME_TYPES::ListOfLong *tinyI;
41 SALOME_TYPES::ListOfString *tinyS;
42 meshPtr->getTinyInfo(tinyD,tinyI,tinyS);
44 int spaceDim=(*tinyI)[1];
45 int nbOfNodes=(*tinyI)[2];
46 int meshDim=(*tinyI)[5];
47 int nbOfCells=(*tinyI)[6];
48 int meshLength=(*tinyI)[7];
49 std::string name((*tinyS)[0]);
50 //std::vector<std::string> compoNames(spaceDim);
51 //for(int i=0;i<spaceDim;i++)
52 // compoNames[i]=(*tinyS)[i+4];
58 ret->Allocate(nbOfCells);
59 vtkPoints *points=vtkPoints::New();
60 vtkDoubleArray *da=vtkDoubleArray::New();
61 da->SetNumberOfComponents(3);
62 da->SetNumberOfTuples(nbOfNodes);
63 double *pts=da->GetPointer(0);
65 SALOME_TYPES::ListOfLong *a1Corba;
66 SALOME_TYPES::ListOfDouble *a2Corba;
67 meshPtr->getSerialisationData(a1Corba,a2Corba);
70 int myLgth=a2Corba->length();
71 for(int i=0;i<myLgth;i++)
76 int offset=3-spaceDim;
77 for(int i=0;i<nbOfNodes;i++)
79 for(int j=0;j<spaceDim;j++)
80 *pts++=(*a2Corba)[spaceDim*i+j];
81 std::fill(pts,pts+offset,0.);
86 vtkIdType *tmp=new vtkIdType[1000];
88 for(int i=0;i<nbOfCells;i++)
90 int pos=(*a1Corba)[i];
91 int pos1=(*a1Corba)[i+1];
92 int nbOfNodeInCurCell=pos1-pos-1;
93 int typeOfCell=(*a1Corba)[pos+nbOfCells+1];
94 for(int j=0;j<nbOfNodeInCurCell;j++)
95 tmp[j]=(*a1Corba)[pos+1+j+nbOfCells+1];
96 int vtkType=ParaMEDMEM2VTKTypeTraducer[typeOfCell];
98 ret->InsertNextCell(vtkType,nbOfNodeInCurCell,tmp);
102 std::set<vtkIdType> s(tmp,tmp+nbOfNodeInCurCell);
103 vtkSmartPointer<vtkCellArray> faces=vtkSmartPointer<vtkCellArray>::New();
104 int nbOfFaces=std::count(tmp,tmp+nbOfNodeInCurCell,-1)+1;
106 for(int i=0;i<nbOfFaces;i++)
108 vtkIdType *work2=std::find(work,tmp+nbOfNodeInCurCell,-1);
109 int nbOfNodesInFace=std::distance(work,work2);
110 faces->InsertNextCell(nbOfNodesInFace,work);
114 std::vector<vtkIdType> v(s.begin(),s.end());
115 ret->InsertNextCell(VTK_POLYHEDRON,v.size(),&v[0],nbOfFaces,faces->GetPointer());
123 ret->SetPoints(points);
128 meshPtr->UnRegister();