]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/ParaMEDCorba/VTKMEDCouplingCMeshClient.cxx
Salome HOME
Add missing MEDFile dependency in DevelopedSurface filter CMakeLists.txt
[modules/paravis.git] / src / Plugins / ParaMEDCorba / VTKMEDCouplingCMeshClient.cxx
1 // Copyright (C) 2010-2016  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 "VTKMEDCouplingCMeshClient.hxx"
21
22 #include "vtkErrorCode.h"
23 #include "vtkDoubleArray.h"
24 #include "vtkRectilinearGrid.h"
25
26 #include <vector>
27 #include <string>
28
29 void ParaMEDMEM2VTK::FillMEDCouplingCMeshInstanceFrom(SALOME_MED::MEDCouplingCMeshCorbaInterface_ptr meshPtr, vtkRectilinearGrid *ret)
30 {
31   meshPtr->Register();
32   SALOME_TYPES::ListOfDouble *tinyD;
33   SALOME_TYPES::ListOfLong *tinyI;
34   SALOME_TYPES::ListOfString *tinyS;
35   meshPtr->getTinyInfo(tinyD,tinyI,tinyS);
36   int sizePerAxe[3]={1,1,1};
37   bool isOK[3]={false,false,false};
38   if((*tinyI)[0]>0)
39     { sizePerAxe[0]=(*tinyI)[0]; isOK[0]=true; }
40   if((*tinyI)[1]>0)
41     { sizePerAxe[1]=(*tinyI)[1]; isOK[1]=true; }
42   if((*tinyI)[2]>0)
43     { sizePerAxe[2]=(*tinyI)[2]; isOK[2]=true; }
44   ret->SetDimensions(sizePerAxe[0],sizePerAxe[1],sizePerAxe[2]);
45   delete tinyI;
46   delete tinyS;
47   delete tinyD;
48   SALOME_TYPES::ListOfDouble *bigD;
49   meshPtr->getSerialisationData(tinyI,bigD);
50   delete tinyI;
51   int offset=0;
52   vtkDoubleArray *da=0;
53   if(isOK[0])
54     {
55       da=vtkDoubleArray::New();
56       da->SetNumberOfTuples(sizePerAxe[0]);
57       da->SetNumberOfComponents(1);
58       double *pt=da->GetPointer(0);
59       for(int i=0;i<sizePerAxe[0];i++)
60         pt[i]=(*bigD)[i];
61       ret->SetXCoordinates(da);
62       da->Delete();
63       offset+=sizePerAxe[0];
64     }
65   if(isOK[1])
66     {
67       da=vtkDoubleArray::New();
68       da->SetNumberOfTuples(sizePerAxe[1]);
69       da->SetNumberOfComponents(1);
70       double *pt=da->GetPointer(0);
71       for(int i=0;i<sizePerAxe[1];i++)
72         pt[i]=(*bigD)[offset+i];
73       ret->SetYCoordinates(da);
74       da->Delete();
75       offset+=sizePerAxe[1];
76     }
77   else
78     {
79       da=vtkDoubleArray::New(); da->SetNumberOfTuples(1); da->SetNumberOfComponents(1);
80       double *pt=da->GetPointer(0); *pt=0.; ret->SetYCoordinates(da); da->Delete();
81     }
82   if(isOK[2])
83     {
84       da=vtkDoubleArray::New();
85       da->SetNumberOfTuples(sizePerAxe[2]);
86       da->SetNumberOfComponents(1);
87       double *pt=da->GetPointer(0);
88       for(int i=0;i<sizePerAxe[2];i++)
89         pt[i]=(*bigD)[offset+i];
90       ret->SetZCoordinates(da);
91       da->Delete();
92     }
93   else
94     {
95       da=vtkDoubleArray::New(); da->SetNumberOfTuples(1); da->SetNumberOfComponents(1);
96       double *pt=da->GetPointer(0); *pt=0.; ret->SetZCoordinates(da); da->Delete();
97     }
98   delete bigD;
99   meshPtr->UnRegister();
100 }