Salome HOME
activate paramedmem test with "salome test" mechanism
[tools/medcoupling.git] / src / MEDCoupling_Swig / vtk2medcoupling.py
1 #! /usr/bin/env python3
2 #  -*- coding: utf-8 -*-
3 # Copyright (C) 2020  CEA/DEN, EDF R&D
4 #
5 # This library is free software; you can redistribute it and/or
6 # modify it under the terms of the GNU Lesser General Public
7 # License as published by the Free Software Foundation; either
8 # version 2.1 of the License, or (at your option) any later version.
9 #
10 # This library is distributed in the hope that it will be useful,
11 # but WITHOUT ANY WARRANTY; without even the implied warranty of
12 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
13 # Lesser General Public License for more details.
14 #
15 # You should have received a copy of the GNU Lesser General Public
16 # License along with this library; if not, write to the Free Software
17 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
18 #
19 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 #
21 # Author : Anthony GEAY (EDF R&D)
22
23 import vtk
24 from vtk.util import numpy_support
25 import medcoupling as mc
26 import numpy as np
27
28 def mesh_convertor(fileName):
29     #vtk.vtkDataSetReader()
30     reader = vtk.vtkXMLUnstructuredGridReader()
31     reader.SetFileName(fileName)
32     reader.Update()
33     ug = reader.GetOutput()
34     return mesh_convertor_mem(ug),ug
35     
36 def mesh_convertor_mem(ug):
37     from distutils.version import LooseVersion
38     #
39     pts = numpy_support.vtk_to_numpy(ug.GetPoints().GetData())
40     #
41     cla = numpy_support.vtk_to_numpy(ug.GetCellLocationsArray())
42     ctvtk = numpy_support.vtk_to_numpy(ug.GetCellTypesArray())
43     conn = numpy_support.vtk_to_numpy(ug.GetCells().GetData())
44     #
45     ct=mc.DataArrayInt(np.array(ctvtk,dtype="int{}".format(mc.MEDCouplingSizeOfIDs())))[:]
46     c=mc.DataArrayInt(conn)[:]
47     ci=mc.DataArrayInt(cla)[:]
48     # for pv580
49     if LooseVersion(vtk.VTK_VERSION) >= LooseVersion("8.90.0"):
50         ci = ci.deltaShiftIndex()+1
51         ci.computeOffsetsFull()
52     #
53     vtk2med = mc.DataArrayInt(mc.vtk2med_cell_types())
54     #
55     ct.transformWithIndArr(vtk2med)
56     c[ci]=ct
57     ci = mc.DataArrayInt.Aggregate([ci,mc.DataArrayInt([len(c)])])
58     #
59     gtv = ct.getDifferentValues().getValues()
60     dimArray = mc.DataArrayInt(len(gtv)) ; dimArray[:] = -1
61     for i,gt in enumerate(gtv):
62         dimArray[i] = mc.MEDCouplingUMesh.GetDimensionOfGeometricType(gt)
63     dim = dimArray.getMaxValueInArray()
64     if not dimArray.isUniform(dim):
65         raise RuntimeError("The input vtkUnstructuredGrid instance is not a single level mesh ! need to split it !")
66     #
67     m=mc.MEDCouplingUMesh("mesh",dim)
68     m.setCoords(mc.DataArrayDouble(np.array(pts,dtype=np.float64)))
69     m.setConnectivity(c,ci,True)
70     m.checkConsistencyLight()
71     #
72     return m