From: Anthony Geay Date: Wed, 13 Nov 2024 07:01:03 +0000 (+0100) Subject: [EDF31315] : Quick fix in python X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=b93e6ebe3fb4f504e999bfa9a89acfa1254e0074;p=tools%2Fmedcoupling.git [EDF31315] : Quick fix in python --- diff --git a/src/MEDCoupling_Swig/vtk2medcoupling.py b/src/MEDCoupling_Swig/vtk2medcoupling.py index 887db5729..1168102f4 100644 --- a/src/MEDCoupling_Swig/vtk2medcoupling.py +++ b/src/MEDCoupling_Swig/vtk2medcoupling.py @@ -25,6 +25,39 @@ from vtk.util import numpy_support import medcoupling as mc import numpy as np +def patchVTKArr( arrLi ): + i = 0 + nbFaces = arrLi[i] ; i += 1 + ret = [None for elt in range(nbFaces)] + for iFace in range(nbFaces): + nbPtsInFace = arrLi[i] ; i += 1 + if iFace > 0: + ret[iFace] = [-1] + else: + ret[iFace] = [] + ret[iFace] += arrLi[i:i+nbPtsInFace] ; i += nbPtsInFace + ret = [ mc.NORM_POLYHED] + sum( ret, [] ) + return mc.DataArrayInt( ret ) + +def patchForPolyedra(polyhedCellIds, ug, mesh): + """ + Method in charge to change the connectivity of polyedra contained in mesh using ug vtkUnstructuredGrid. + + :param in polyhedCellIds: mc.DataArrayInt of cells ids in mesh to be patched + :param in ug: vtkUnstructuredGrid containing polyhedra + :param in-out mesh: mc.MEDCouplingUMesh. 3D Mesh whose polyedra cells connectivity will be modified + """ + c, ci = mesh.getNodalConnectivity(), mesh.getNodalConnectivityIndex() + facesLoc = mc.DataArrayInt( numpy_support.vtk_to_numpy( ug.GetFaceLocations() ) ) + faces = mc.DataArrayInt( numpy_support.vtk_to_numpy( ug.GetFaces() ) ) + facesLoc = mc.DataArrayInt.Aggregate( [ facesLoc, mc.DataArrayInt([ len(faces) ]) ] ) + connForPoly = mc.DataArrayInt.Aggregate( [patchVTKArr( ( faces[facesLoc[i]:facesLoc[i+1]] ).getValues() ) for i in range(ug.GetNumberOfCells()) ] ) + tmpDa = mc.DataArrayInt( len(facesLoc) ) ; tmpDa.iota() + facesLoc -= tmpDa + meshPoly = mc.MEDCouplingUMesh( mesh.getName(), 3) ; meshPoly.setCoords( mesh.getCoords() ) ; meshPoly.setConnectivity(connForPoly,facesLoc,True) + mesh[polyhedCellIds] = meshPoly + pass + def mesh_convertor(fileName): #vtk.vtkDataSetReader() reader = vtk.vtkXMLUnstructuredGridReader() @@ -69,4 +102,9 @@ def mesh_convertor_mem(ug): m.setConnectivity(c,ci,True) m.checkConsistencyLight() # + if m.getMeshDimension() == 3: + polyhedCellIds = ct.findIdsEqual(mc.NORM_POLYHED) + if not polyhedCellIds.empty(): + patchForPolyedra(polyhedCellIds, ug, m) + # return m