]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/ParaMEDCorba/VTKMEDCouplingMeshClient.cxx
Salome HOME
Merge from V6_main_20120808 08Aug12
[modules/paravis.git] / src / Plugins / ParaMEDCorba / VTKMEDCouplingMeshClient.cxx
1 // Copyright (C) 2010-2012  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.
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 "VTKMEDCouplingMeshClient.hxx"
21 #include "VTKMEDCouplingUMeshClient.hxx"
22 #include "VTKMEDCouplingCMeshClient.hxx"
23
24 #include "vtkErrorCode.h"
25 #include "vtkUnstructuredGrid.h"
26 #include "vtkRectilinearGrid.h"
27
28 #include <vector>
29 #include <string>
30
31 //vtkErrorMacro("Cannot call Start() twice before calling Finish().");
32 void ParaMEDMEM2VTK::FillMEDCouplingMeshInstanceFrom(SALOME_MED::MEDCouplingMeshCorbaInterface_ptr meshPtr, vtkDataSet *ret)
33 {
34   SALOME_MED::MEDCouplingUMeshCorbaInterface_var umeshPtr=SALOME_MED::MEDCouplingUMeshCorbaInterface::_narrow(meshPtr);
35   if(!CORBA::is_nil(umeshPtr))
36     {
37       vtkUnstructuredGrid *ret1=vtkUnstructuredGrid::SafeDownCast(ret);
38       if(!ret1)
39         {
40           vtkErrorWithObjectMacro(ret,"Internal error in ParaMEDCorba plugin : mismatch between VTK type and CORBA type UMesh !");
41           return ;
42         }
43       bool dummy;//VTK bug
44       ParaMEDMEM2VTK::FillMEDCouplingUMeshInstanceFrom(umeshPtr,ret1,dummy);//VTK bug
45       return ;
46     }
47   SALOME_MED::MEDCouplingCMeshCorbaInterface_var cmeshPtr=SALOME_MED::MEDCouplingCMeshCorbaInterface::_narrow(meshPtr);
48   if(!CORBA::is_nil(cmeshPtr))
49     {
50       vtkRectilinearGrid *ret1=vtkRectilinearGrid::SafeDownCast(ret);
51       if(!ret1)
52         {
53           vtkErrorWithObjectMacro(ret,"Internal error in ParaMEDCorba plugin : mismatch between VTK type and CORBA type CMesh !");
54           return ;
55         }
56       ParaMEDMEM2VTK::FillMEDCouplingCMeshInstanceFrom(cmeshPtr,ret1);
57       return ;
58     }
59   vtkErrorWithObjectMacro(ret,"Error : CORBA mesh type ! Mesh type not managed !");
60 }
61
62 vtkDataSet *ParaMEDMEM2VTK::BuildFromMEDCouplingMeshInstance(SALOME_MED::MEDCouplingMeshCorbaInterface_ptr meshPtr, bool& isPolyh)
63 {
64   SALOME_MED::MEDCouplingUMeshCorbaInterface_var umeshPtr=SALOME_MED::MEDCouplingUMeshCorbaInterface::_narrow(meshPtr);
65   if(!CORBA::is_nil(umeshPtr))
66     {
67       vtkUnstructuredGrid *ret1=vtkUnstructuredGrid::New();
68       ParaMEDMEM2VTK::FillMEDCouplingUMeshInstanceFrom(umeshPtr,ret1,isPolyh);
69       return ret1;
70     }
71   SALOME_MED::MEDCouplingCMeshCorbaInterface_var cmeshPtr=SALOME_MED::MEDCouplingCMeshCorbaInterface::_narrow(meshPtr);
72   if(!CORBA::is_nil(cmeshPtr))
73     {
74       vtkRectilinearGrid *ret1=vtkRectilinearGrid::New();
75       ParaMEDMEM2VTK::FillMEDCouplingCMeshInstanceFrom(cmeshPtr,ret1);
76       return ret1;
77     }
78   vtkOutputWindowDisplayErrorText("Error : CORBA mesh type ! Mesh type not managed #2 !");
79   return 0;
80 }