Salome HOME
c6fb501f9575e722992e13a3ea74bbe7de773972
[tools/medcoupling.git] / src / MEDCoupling_Swig / UsersGuideExamplesTest.py
1 #  -*- coding: iso-8859-1 -*-
2 # Copyright (C) 2007-2016  CEA/DEN, EDF R&D
3 #
4 # This library is free software; you can redistribute it and/or
5 # modify it under the terms of the GNU Lesser General Public
6 # License as published by the Free Software Foundation; either
7 # version 2.1 of the License, or (at your option) any later version.
8 #
9 # This library is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12 # Lesser General Public License for more details.
13 #
14 # You should have received a copy of the GNU Lesser General Public
15 # License along with this library; if not, write to the Free Software
16 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
17 #
18 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 #
20
21 from MEDCoupling import *
22 from math import pi, sqrt
23
24 # ! [PySnippetUMeshStdBuild1_1]
25 coords=[-0.3,-0.3,0.,   0.2,-0.3,0.,   0.7,-0.3,0.,   -0.3,0.2,0.,   0.2,0.2,0.,
26         0.7,0.2,0.,    -0.3,0.7,0.,    0.2,0.7,0.,     0.7,0.7,0. ]
27 nodalConnPerCell=[0,3,4,1, 1,4,2, 4,5,2, 6,7,4,3, 7,8,5,4]
28 # ! [PySnippetUMeshStdBuild1_1]
29 # ! [PySnippetUMeshStdBuild1_2]
30 mesh=MEDCouplingUMesh("My2DMesh",2)
31 # ! [PySnippetUMeshStdBuild1_2]
32 # ! [PySnippetUMeshStdBuild1_3]
33 mesh.allocateCells(5)#You can put more than 5 if you want but not less.
34 # adding cells
35 mesh.insertNextCell(NORM_QUAD4,nodalConnPerCell[:4])
36 mesh.insertNextCell(NORM_TRI3,nodalConnPerCell[4:7])
37 mesh.insertNextCell(NORM_TRI3,nodalConnPerCell[7:10])
38 mesh.insertNextCell(NORM_QUAD4,nodalConnPerCell[10:14])
39 mesh.insertNextCell(NORM_QUAD4,nodalConnPerCell[14:])
40 # compacting
41 mesh.finishInsertingCells()
42 # ! [PySnippetUMeshStdBuild1_3]
43 # ! [PySnippetUMeshStdBuild1_4]
44 coordsArr=DataArrayDouble(coords,9,3)#here coordsArr are declared to have 3 components, mesh will deduce that its spaceDim==3.
45 mesh.setCoords(coordsArr)#coordsArr contains 9 tuples, that is to say mesh contains 9 nodes.
46 # ! [PySnippetUMeshStdBuild1_4]
47 # ! [PySnippetUMeshStdBuild1_5]
48 mesh.checkConsistencyLight()
49 # ! [PySnippetUMeshStdBuild1_5]
50
51 # ! [PySnippetCMeshStdBuild1_1]
52 XCoords=[-0.3,0.,0.1,0.3,0.45,0.47,0.49,1.,1.22] # 9 values along X
53 YCoords=[0.,0.1,0.37,0.45,0.47,0.49,1.007] # 7 values along Y
54 arrX=DataArrayDouble(XCoords)
55 arrX.setInfoOnComponent(0,"X [m]")
56 arrY=DataArrayDouble(YCoords)
57 arrY.setInfoOnComponent(0,"Y [m]")
58 # ! [PySnippetCMeshStdBuild1_1]
59 # ! [PySnippetCMeshStdBuild1_2]
60 mesh=MEDCouplingCMesh("My2D_CMesh")
61 mesh.setCoords(arrX,arrY)
62 # ! [PySnippetCMeshStdBuild1_2]
63
64 nodalConnPerCell=range(4*4)
65 # ! [GU_MEDCoupling1SGTUMesh_0]
66 mesh=MEDCoupling1SGTUMesh("myQuadMesh",NORM_QUAD4)
67 mesh.allocateCells(3)
68 mesh.insertNextCell(nodalConnPerCell[:4])
69 mesh.insertNextCell(nodalConnPerCell[4:8])
70 mesh.insertNextCell(nodalConnPerCell[8:12])
71 # ! [GU_MEDCoupling1SGTUMesh_0]
72
73 # ! [GU_MEDCoupling1SGTUMesh_1]
74 polymesh=MEDCoupling1DGTUMesh("myPolyhedra",NORM_POLYHED)
75 polymesh.allocateCells(1)
76 polymesh.insertNextCell([0,1,2,3,-1,7,6,5,4,-1,0,4,5,1,-1,1,5,6,2,-1,3,2,6,7,-1,0,3,7,4])
77 # ! [GU_MEDCoupling1SGTUMesh_1]
78
79 #! [UG_DataArrayDouble_0]
80 d=DataArrayDouble([1,2,3,4,5,6],3,2)
81 #! [UG_DataArrayDouble_0]
82 #! [UG_DataArrayDouble_1]
83 d=DataArrayDouble([(1,2),(3,4),(5,6)])
84 #! [UG_DataArrayDouble_1]
85 #! [UG_DataArrayDouble_2]
86 d=DataArrayDouble([(1,2,3),(4,5,6)])
87 #! [UG_DataArrayDouble_2]
88 #! [UG_DataArrayDouble_3]
89 d.rearrange(2)
90 #! [UG_DataArrayDouble_3]
91 #! [UG_DataArrayDouble_4]
92 i=DataArrayInt([(1,2,3),(4,5,6)])
93 f=DataArrayFloat([(1,2,3),(4,5,6)])
94 #! [UG_DataArrayDouble_4]
95
96 #! [UG_MEDCouplingCurveLinearMesh_0]
97 m=MEDCouplingCurveLinearMesh("myCurveLinearMesh")
98 m.setNodeGridStructure([2,3])
99 #! [UG_MEDCouplingCurveLinearMesh_0]
100 #! [UG_MEDCouplingCurveLinearMesh_1]
101 coords=DataArrayDouble([0.,0., 2.,0., 0.,1., 1.9,1.1, 0.3,1.9, 2.2,2.1],6,2)
102 coords.setInfoOnComponents(["X [m]","Y [m]"])
103 m.setCoords(coords)
104 #! [UG_MEDCouplingCurveLinearMesh_1]
105 #! [UG_MEDCouplingCurveLinearMesh_2]
106 m.checkConsistencyLight()
107 #! [UG_MEDCouplingCurveLinearMesh_2]
108
109 XCoords=[-0.3,0.,0.1,0.3,0.45,0.47,0.49,1.,1.22];  arrX=DataArrayDouble(XCoords)
110 YCoords=[0.,0.1,0.37,0.45,0.47,0.49,1.007];  arrY=DataArrayDouble(YCoords)
111 mesh=MEDCouplingCMesh("My2D_CMesh")
112 mesh.setCoords(arrX,arrY)
113 #! [UG_MEDCouplingFieldDouble_0]
114 fieldOnCells=MEDCouplingFieldDouble(ON_CELLS,ONE_TIME)
115 fieldOnCells.setName("MyTensorFieldOnCellOneTime")
116 fieldOnCells.setMesh(mesh)
117 #! [UG_MEDCouplingFieldDouble_0]
118 #! [UG_MEDCouplingFieldDouble_1]
119 fieldOnCells.setTimeUnit("ms") # Time unit is ms.
120 fieldOnCells.setTime(4.22,2,-1) # Time attached is 4.22 ms, iteration id is 2 and order id (or sub iteration id) is -1
121 #! [UG_MEDCouplingFieldDouble_1]
122 #! [UG_MEDCouplingFieldDouble_2]
123 array=DataArrayDouble()
124 array.alloc(fieldOnCells.getMesh().getNumberOfCells(),2) # Implicitly fieldOnCells will be a 2 components field.
125 array.fillWithValue(7.)
126 fieldOnCells.setArray(array)
127 fieldOnCells.checkConsistencyLight()
128 # fieldOnCells is now usable
129 # ...
130 #! [UG_MEDCouplingFieldDouble_2]
131 #
132 #! [UG_MEDCouplingFieldDouble_3]
133 fieldOnNodes=MEDCouplingFieldDouble(ON_NODES,ONE_TIME)
134 fieldOnNodes.setName("MyScalarFieldOnNodeOneTime")
135 fieldOnNodes.setMesh(mesh)
136 fieldOnNodes.setTimeUnit("ms") # Time unit is ms.
137 fieldOnNodes.setTime(4.22,2,-1) # Time attached is 4.22 ms, iteration id is 2 and order id (or sub iteration id) is -1
138 array=DataArrayDouble()
139 array.alloc(fieldOnNodes.getMesh().getNumberOfNodes(),1) # Implicitly fieldOnNodes will be a 1 component field.
140 array.fillWithValue(7.)
141 fieldOnNodes.setArray(array)
142 fieldOnNodes.checkConsistencyLight()
143 # fieldOnNodes is now usable
144 # ...
145 #! [UG_MEDCouplingFieldDouble_3]
146 field=fieldOnCells
147 #! [UG_MEDCouplingFieldDouble_4]
148 print(mesh.getHeapMemorySizeStr())
149 print(field.getHeapMemorySizeStr())
150 print(array.getHeapMemorySizeStr())
151 #! [UG_MEDCouplingFieldDouble_4]
152
153 f=fieldOnCells
154 nbComp=3
155 val=5.
156 #! [UG_MEDCouplingFieldDouble_5]
157 f.applyFunc(nbComp,val)
158 #! [UG_MEDCouplingFieldDouble_5]
159 #! [UG_MEDCouplingFieldDouble_6]
160 f.applyFunc(1,"sqrt(X*X+Y*Y+Z*Z)")
161 #! [UG_MEDCouplingFieldDouble_6]
162 f.applyFunc(nbComp,val)
163 #! [UG_MEDCouplingFieldDouble_7]
164 f.applyFunc(4,"IVec*y+JVec*x+KVec*z+LVec*sqrt(x*x+y*y+z*z)")
165 #! [UG_MEDCouplingFieldDouble_7]
166
167 field=fieldOnCells
168 #! [UG_MEDCouplingFieldDouble_8]
169 points=DataArrayDouble([(0.,0.),(1,1)])
170 values=field.getValueOnMulti(points)
171 #! [UG_MEDCouplingFieldDouble_8]
172 #! [UG_MEDCouplingFieldDouble_9]
173 field.integral(True)
174 field.integral(0,True)
175 #! [UG_MEDCouplingFieldDouble_9]
176 field.applyFunc(6,1.)
177 #! [UG_MEDCouplingFieldDouble_10]
178 assert(field.getNumberOfComponents()==6)
179 diviatorfield=field.deviator()
180 #! [UG_MEDCouplingFieldDouble_10]
181
182
183 from MEDCouplingDataForTest import MEDCouplingDataForTest
184 mesh=MEDCouplingDataForTest.build2DTargetMesh_1();
185 #! [UG_MEDCouplingGaussPointField_0]
186 fieldGauss=MEDCouplingFieldDouble.New(ON_GAUSS_PT,ONE_TIME);
187 fieldGauss.setMesh(mesh);
188 fieldGauss.setName("MyFirstFieldOnGaussPoint");
189 fieldGauss.setTimeUnit("ms") # Time unit is ms.
190 fieldGauss.setTime(4.22,2,-1) # Time attached is 4.22 ms, iteration id is 2 and order id (or sub iteration id) is -1
191 #! [UG_MEDCouplingGaussPointField_0]
192 #! [UG_MEDCouplingGaussPointField_1]
193 tria3CooRef=[ 0.0, 0.0, 1.0 , 0.0, 0.0, 1.0 ]
194 tria3CooGauss=[ 0.1, 0.8, 0.2, 0.7 ]
195 wg3=[0.3,0.3];
196 fieldGauss.setGaussLocalizationOnType(NORM_TRI3,tria3CooRef,tria3CooGauss,wg3);
197 #
198 quad4CooRef=[-1.0, 1.0, -1.0, -1.0, 1.0, -1.0, 1.0, 1.0]
199 quad4CooGauss=[ 0.3, 0.2, 0.2, 0.1, 0.2, 0.4, 0.15, 0.27 ]
200 wg4=[0.3,0.3,0.3,0.3];
201 fieldGauss.setGaussLocalizationOnType(NORM_QUAD4,quad4CooRef,quad4CooGauss,wg4);
202 #! [UG_MEDCouplingGaussPointField_1]
203 #! [UG_MEDCouplingGaussPointField_2]
204 nbTuples=mesh.getNumberOfCellsWithType(NORM_TRI3)*len(wg3)+mesh.getNumberOfCellsWithType(NORM_QUAD4)*len(wg4)
205 array=DataArrayDouble.New();
206 values=[float(i+1) for i in range(nbTuples)]
207 array.setValues(values,nbTuples,1);
208 fieldGauss.setArray(array);
209 fieldGauss.checkConsistencyLight();
210 #! [UG_MEDCouplingGaussPointField_2]
211
212 from MEDCouplingDataForTest import MEDCouplingDataForTest
213 m=MEDCouplingDataForTest.build2DTargetMesh_1();
214 Ids=range(1,3)
215 #! [UG_ExtractForMeshes_0]
216 part=m[Ids]
217 #! [UG_ExtractForMeshes_0]
218 #! [UG_ExtractForMeshes_1]
219 subNodeIds=part.computeFetchedNodeIds()
220 #! [UG_ExtractForMeshes_1]
221 #! [UG_ExtractForMeshes_2]
222 m.getCoords()[subNodeIds]
223 #! [UG_ExtractForMeshes_2]
224 #! [UG_ExtractForMeshes_3]
225 part.zipCoords()
226 #! [UG_ExtractForMeshes_3]
227 #! [UG_ExtractForMeshes_4]
228 o2n=part.zipCoordsTraducer()
229 #! [UG_ExtractForMeshes_4]
230
231 m2=MEDCouplingDataForTest.build3DExtrudedUMesh_1()[0]
232 #! [UG_ExtractForMeshes_5]
233 bn = m2.findBoundaryNodes()
234 #! [UG_ExtractForMeshes_5]
235 #! [UG_ExtractForMeshes_6]
236 bc = m2.getCellIdsLyingOnNodes(bn,False)
237 #! [UG_ExtractForMeshes_6]
238
239 #! [UG_ExtractForMeshes_7]
240 m2.translate([1.,2.,3.])
241 #! [UG_ExtractForMeshes_7]
242 #! [UG_ExtractForMeshes_8]
243 m2.getCoords()[:]+=DataArrayDouble([1.,2.,3.],1,3)
244 #! [UG_ExtractForMeshes_8]
245 import math
246 #! [UG_ExtractForMeshes_9]
247 m2.rotate([1,2,1],[0,1,0],math.pi/3)
248 #! [UG_ExtractForMeshes_9]
249 #! [UG_ExtractForMeshes_10]
250 MEDCouplingPointSet.Rotate3DAlg([1,2,1],[0,1,0],math.pi/3,m2.getCoords())
251 #! [UG_ExtractForMeshes_10]
252
253 #! [UG_ExtractForMeshes_11]
254 volPerCell=m2.getMeasureField(True)
255 #! [UG_ExtractForMeshes_11]
256 #! [UG_ExtractForMeshes_12]
257 volPerCell.getArray().accumulate()
258 #! [UG_ExtractForMeshes_12]
259 t1=-1
260 #! [UG_ExtractForMeshes_13]
261 part=volPerCell.getArray().findIdsGreaterOrEqualTo(t1)
262 #! [UG_ExtractForMeshes_13]
263 #! [UG_ExtractForMeshes_14]
264 m2[part]
265 #! [UG_ExtractForMeshes_14]
266 #! [UG_ExtractForMeshes_15]
267 centers=m2.computeCellCenterOfMass()
268 #! [UG_ExtractForMeshes_15]
269 #! [UG_ExtractForMeshes_16]
270 (centers*volPerCell.getArray()).accumulate()/DataArrayDouble(volPerCell.accumulate())
271 #! [UG_ExtractForMeshes_16]
272
273 #! [UG_ExtractForMeshes_17]
274 m2.scale( [1,2,4], 6. )
275 #! [UG_ExtractForMeshes_17]
276
277 #! [UG_ExtractForMeshes_18]
278 ortho_field=m.buildOrthogonalField()
279 #! [UG_ExtractForMeshes_18]
280
281 #! [UG_ExtractForMeshes_19]
282 ibc=m2.computeIsoBarycenterOfNodesPerCell()
283 #! [UG_ExtractForMeshes_19]
284
285 #! [UG_ExtractForMeshes_20]
286 # make a structured mesh 1x5
287 coords=DataArrayDouble(range(6))
288 cmesh=MEDCouplingCMesh("cmesh")
289 cmesh.setCoords(coords,coords[:2])
290
291 # make a mesh with two zones
292 zmesh=cmesh.buildUnstructured()[0,1,3,4]
293
294 # get cells ids of zones
295 zoneArrays=zmesh.partitionBySpreadZone()
296 print([ ids.getValues() for ids in zoneArrays])
297 #! [UG_ExtractForMeshes_20]
298
299 coordsArr=DataArrayDouble(range(6))
300 mesh2d=MEDCouplingCMesh("mesh2d")
301 mesh2d.setCoords(coordsArr,coordsArr[:2])
302 mesh2d=mesh2d.buildUnstructured()
303 mesh1d=MEDCouplingCMesh("mesh1d")
304 mesh1d.setCoords(coordsArr,coordsArr[:1])
305 mesh1d=mesh1d.buildUnstructured()
306 mesh1d.rotate( [2.3,0], math.radians( 25 ))
307 mesh1d.translate( [0.2,0.4] )
308 #! [UG_ExtractForMeshes_21]
309 m2d,m1d,a2d,a1d=MEDCouplingUMesh.Intersect2DMeshWith1DLine( mesh2d, mesh1d, 1e-12 )
310 #! [UG_ExtractForMeshes_21]
311 #print ("a2d",a2d.getValues())
312 #print (a1d.getValues())
313
314 mesh1=mesh2d.deepCopy()
315 mesh1.rotate( [2.3,0], math.radians( 25 ))
316 mesh2=mesh2d
317 #! [UG_ExtractForMeshes_22]
318 m,a1,a2=MEDCouplingUMesh.Intersect2DMeshes( mesh1, mesh2, 1e-12 )
319 #! [UG_ExtractForMeshes_22]
320
321
322 points=DataArrayDouble([(0,0,0)])
323 mesh=m2.computeSkin()
324 #! [UG_ExtractForMeshes_23]
325 dist,cells=mesh.distanceToPoints(points)
326 #! [UG_ExtractForMeshes_23]
327
328 arr=DataArrayDouble([1,2,3,4,5,6],3,2)
329 a,b=2,5
330 #! [UG_ExtractForArrays_0]
331 tupleIds = arr[:,0].findIdsInRange(a,b)
332 #! [UG_ExtractForArrays_0]
333 c,d=0,7
334 #! [UG_ExtractForArrays_1]
335 tupleIds1 = arr.magnitude().findIdsInRange(c,d)
336 #! [UG_ExtractForArrays_1]
337 #! [UG_ExtractForArrays_2]
338 tupleIds2 = DataArrayInt.buildSubstraction(tupleIds,tupleIds1)
339 #! [UG_ExtractForArrays_2]
340
341 valsArr1=DataArrayDouble(range(9*2),9,2)
342 field4 = MEDCouplingFieldDouble(ON_NODES)
343 field4.setArray(valsArr1)
344 mesh=MEDCouplingCMesh("My2D_CMesh")
345 coo=DataArrayDouble([0,1,2])
346 mesh.setCoords(coo,coo)
347 field4.setMesh(mesh)
348 ids4=[1,2]
349 #! [UG_ExtractForFields_0]
350 subField = field4[ids4]
351 #! [UG_ExtractForFields_0]
352
353 m4=MEDCouplingCMesh("box")
354 coo=DataArrayDouble(range(7))
355 m4.setCoords(coo[:5],coo[:5],coo)
356 m4=m4.buildUnstructured()
357 valsArr1=m4.computeCellCenterOfMass()
358 valsArr1.applyFunc(1,"sqrt(X*X+Y*Y+Z*Z)")
359 field5 = MEDCouplingFieldDouble(ON_CELLS)
360 field5.setArray(valsArr1)
361 field5.setMesh(m4)
362 #! [UG_ExtractForFields_1]
363 origin=[0,0,2]
364 normvec=[-1,-1,6]
365 slice5=field5.extractSlice3D(origin,normvec,1e-10)
366 #! [UG_ExtractForFields_1]
367
368 from MEDCouplingDataForTest import MEDCouplingDataForTest
369 m1=MEDCouplingDataForTest.build2DTargetMesh_1();
370 m2=m1
371 eps=1e-12
372 #! [UG_MeshComparison_0]
373 m1.isEqual(m2,eps)
374 #! [UG_MeshComparison_0]
375 #! [UG_MeshComparison_1]
376 m1.isEqualWithoutConsideringStr(m2,eps)
377 #! [UG_MeshComparison_1]
378
379 arr=DataArrayDouble(5) ; arr.iota()
380 m=MEDCouplingCMesh() ; m.setCoords(arr,arr)
381 m=m.buildUnstructured()
382 m1=m.deepCopy()
383 arr=DataArrayInt([0,3,1,4,2])
384 m2=m1[arr]
385 #! [UG_MeshComparison_2]
386 a,_=m1.checkGeoEquivalWith(m2,20,1e-12)
387 assert(m1[a].isEqualWithoutConsideringStr(m2,1e-12))
388 #! [UG_MeshComparison_2]
389 m3=m1
390 m3.translate([100.,0])
391 m1=MEDCouplingUMesh.MergeUMeshes([m3,m])
392 m2=MEDCouplingUMesh.MergeUMeshes([m,m3])
393 #! [UG_MeshComparison_3]
394 a,b=m1.checkGeoEquivalWith(m2,12,1e-12)
395 m2.renumberNodes(b,len(b))
396 assert(m1[a].isEqualWithoutConsideringStr(m2,1e-12))
397 #! [UG_MeshComparison_3]
398
399
400 coords = [0.,2.,4.]
401 coordsArr=DataArrayDouble(coords,3,1)
402 m4=MEDCouplingCMesh()
403 m4.setCoords(coordsArr,coordsArr,coordsArr)
404 m4=m4.buildUnstructured()
405 field1 = m4.fillFromAnalytic(ON_NODES,1,"x+y")
406 pts = [0.5,0.5,0.5,1.,1,1.]
407 #! [UG_CommonHandlingMesh_0]
408 assert(field1.getTypeOfField()==ON_NODES)
409 field1.getMesh().simplexize(PLANAR_FACE_5)
410 field1.getValueOnMulti(pts)
411 #! [UG_CommonHandlingMesh_0]
412
413 from MEDCouplingDataForTest import MEDCouplingDataForTest
414 m2=MEDCouplingDataForTest.build2DTargetMesh_1();
415 m2.changeSpaceDimension(3);
416 m1=MEDCouplingDataForTest.buildCU1DMesh_U();
417 m1.changeSpaceDimension(3);
418 center=[0.,0.,0.]
419 vector=[0.,1.,0.]
420 m1.rotate(center,vector,-pi/2.);
421 #! [UG_CommonHandlingMesh_1]
422 m3=m2.buildExtrudedMesh(m1,0);
423 #! [UG_CommonHandlingMesh_1]
424
425 #! [UG_CommonHandlingMesh_2]
426 m5=MEDCouplingUMesh.MergeUMeshes([m3,m4])
427 #! [UG_CommonHandlingMesh_2]
428
429 #! [UG_CommonHandlingMesh_3]
430 m1.convertLinearCellsToQuadratic(0)
431 #! [UG_CommonHandlingMesh_3]
432
433 #! [UG_CommonHandlingMesh_4]
434 skin = m1.computeSkin()
435 #! [UG_CommonHandlingMesh_4]
436 #! [UG_CommonHandlingMesh_5]
437 #bc = m1.getCellIdsLyingOnNodes(bn,False)
438 #! [UG_CommonHandlingMesh_5]
439
440 mesh3d=m3
441 #! [UG_CommonHandlingMesh_6]
442 mesh1d,d,di,r,ri=mesh3d.explodeIntoEdges()
443 #! [UG_CommonHandlingMesh_6]
444 #! [UG_CommonHandlingMesh_7]
445 mesh2d,d,di,r,ri=mesh3d.buildDescendingConnectivity()
446 #! [UG_CommonHandlingMesh_7]
447
448 mesh2d=MEDCouplingCMesh()
449 mesh2d.setCoords(coordsArr,coordsArr)
450 mesh2d=mesh2d.buildUnstructured()
451 eps=1e-12
452 #! [UG_CommonHandlingMesh_8]
453 changedCells=mesh2d.conformize2D(eps)
454 #! [UG_CommonHandlingMesh_8]
455
456 #! [UG_CommonHandlingMesh_9]
457 mesh2d.duplicateNodes([3,4])
458 #! [UG_CommonHandlingMesh_9]
459
460 m1d=MEDCouplingCMesh()
461 m1d.setCoords(coordsArr)
462 m1d=m1d.buildUnstructured()
463 #! [UG_CommonHandlingMesh_10]
464 m1d.renumberCells(m1d.orderConsecutiveCells1D().invertArrayN2O2O2N(m1d.getNumberOfCells()))
465 #! [UG_CommonHandlingMesh_10]
466
467 skin=m3.computeSkin()
468 #! [UG_CommonHandlingMesh_11]
469 vec=[0,0,-1]
470 skin.orientCorrectly2DCells(vec,False)
471 #! [UG_CommonHandlingMesh_11]
472
473 #! [UG_CommonHandlingMesh_12]
474 m3.orientCorrectlyPolyhedrons()
475 #! [UG_CommonHandlingMesh_12]
476
477 #! [UG_CommonHandlingMesh_13]
478 m3.renumberCells(m3.rearrange2ConsecutiveCellTypes())
479 m3.sortCellsInMEDFileFrmt()
480 #! [UG_CommonHandlingMesh_13]
481
482 m=MEDCouplingCMesh()
483 m.setCoords(coordsArr[:2],coordsArr[:2])
484 m=m.buildUnstructured()
485 #! [UG_CommonHandlingMesh_14]
486 m.renumberNodes([2,1,0,-1],3)
487 #! [UG_CommonHandlingMesh_14]
488
489 #! [UG_CommonHandlingMesh_15]
490 mtet,n2ocells,np=m3.tetrahedrize(PLANAR_FACE_5)
491 #! [UG_CommonHandlingMesh_15]
492 #! [UG_CommonHandlingMesh_16]
493 m.mergeNodes(1e-12)
494 #! [UG_CommonHandlingMesh_16]
495
496
497 #! [UG_Projection_0]
498 srcCoo=DataArrayDouble([(0,0),(1,0),(3,0),(0,1),(1,1),(3,1)])
499 src=MEDCouplingUMesh("src",2)
500 src.setCoords(srcCoo)
501 src.allocateCells()
502 src.insertNextCell(NORM_QUAD4,[0,3,4,1])
503 src.insertNextCell(NORM_QUAD4,[1,4,5,2])
504 #
505 trgCoo=DataArrayDouble([(0.5,0.5),(1.5,0.5),(1.5,1.5)])
506 trg=MEDCouplingUMesh("trg",2)
507 trg.setCoords(trgCoo)
508 trg.allocateCells()
509 trg.insertNextCell(NORM_TRI3,[0,2,1])
510 #! [UG_Projection_0]
511 from MEDCouplingRemapper import MEDCouplingRemapper
512 #! [UG_Projection_1]
513 rem=MEDCouplingRemapper()
514 rem.prepare(src,trg,"P0P0")
515 print(rem.getCrudeMatrix())
516 #! [UG_Projection_1]
517 #! [UG_Projection_2]
518 srcF=MEDCouplingFieldDouble(ON_CELLS)
519 srcF.setMesh(src)
520 srcF.setArray(DataArrayDouble([3,4]))
521 srcF.setNature(IntensiveMaximum)
522 #
523 trgF=rem.transferField(srcF,-1)
524 #! [UG_Projection_2]
525 #! [UG_Projection_3]
526 rem=MEDCouplingRemapper()
527 rem.prepare(src,trg,"P0P1")
528 print(rem.getCrudeMatrix())
529 #! [UG_Projection_3]
530 #! [UG_Projection_4]
531 rem=MEDCouplingRemapper()
532 rem.prepare(src,trg,"P1P0")
533 print(rem.getCrudeMatrix())
534 #! [UG_Projection_4]
535
536
537 #! [UG_Projection_5]
538 coo=DataArrayDouble([(0.,0.,0.), (1,0,0), (0,1,0)])
539 src=MEDCouplingUMesh("src",2)
540 src.setCoords(coo)
541 src.allocateCells(1)
542 src.insertNextCell(NORM_TRI3,[0,1,2])
543 tgt = src.deepCopy()
544 rem=MEDCouplingRemapper()
545 rem.prepare(src,tgt,"P0P0")
546 print(rem.getCrudeMatrix())
547 #! [UG_Projection_5]
548 #! [UG_Projection_6]
549 src.translate([0,0,1e-3])
550 rem.prepare(src,tgt,"P0P0")
551 print(rem.getCrudeMatrix())
552 #! [UG_Projection_6]
553 #! [UG_Projection_7]
554 rem.setBoundingBoxAdjustmentAbs( 1e-3 )
555 rem.prepare(src,tgt,"P0P0")
556 print(rem.getCrudeMatrix())
557 #! [UG_Projection_7]
558
559 import math
560 #! [UG_Projection_8]
561 src.rotate([0,0,0],[0,1,0],math.pi/4)
562 rem.prepare(src,tgt,"P0P0")
563 print(rem.getCrudeMatrix())
564 #! [UG_Projection_8]
565 #! [UG_Projection_9]
566 rem.setMaxDistance3DSurfIntersect( 0.1 )
567 rem.prepare(src,tgt,"P0P0")
568 print(rem.getCrudeMatrix())
569 #! [UG_Projection_9]
570
571 #! [UG_Projection_10]
572 rem.setMaxDistance3DSurfIntersect( -1 ) # switch it off
573 rem.setMinDotBtwPlane3DSurfIntersect( 0.8 )
574 rem.prepare(src,tgt,"P0P0")
575 print(rem.getCrudeMatrix())
576 #! [UG_Projection_10]
577
578 from MEDCouplingDataForTest import MEDCouplingDataForTest
579 m=MEDCouplingDataForTest.build2DTargetMesh_1();
580 #! [UG_Optimization_0]
581 from MEDRenumber import RenumberingFactory
582 ren=RenumberingFactory("BOOST")
583 a,b=m.computeNeighborsOfCells()
584 n2o,_=ren.renumber(a,b)
585 mrenum=m[n2o]
586 #! [UG_Optimization_0]
587
588 #! [UG_Optimization_1]
589 from MEDCoupling import MEDCouplingSkyLineArray
590 import MEDPartitioner
591 # prepare a MEDPartitioner
592 a,b=m.computeNeighborsOfCells()
593 sk=MEDCouplingSkyLineArray(b,a)
594 g=MEDPartitioner.MEDPartitioner.Graph(sk)
595 # compute partitioning into 4 parts
596 g.partGraph(4)
597 # get the 1st of parts of m
598 procIdOnCells=g.getPartition().getValuesArray()
599 p0=procIdOnCells.findIdsEqual(0)
600 part0=m[p0]
601 #! [UG_Optimization_1]
602 #! [UG_Optimization_2]
603 boundary_nodes_part0=part0.findBoundaryNodes()
604 boundary_cells_part0=p0[part0.getCellIdsLyingOnNodes(boundary_nodes_part0,False)]
605 # starting from knowledge of neighborhood it s possible to know the neighbors of boundary_cells_part0
606 neighbors_boundary_cells_part0=MEDCouplingUMesh.ExtractFromIndexedArrays(boundary_cells_part0,a,b)[0]
607 neighbors_boundary_cells_part0.sort()
608 neighbors_boundary_cells_part0=neighbors_boundary_cells_part0.buildUnique()
609 #
610 layer_of_part0=neighbors_boundary_cells_part0.buildSubstraction(p0)
611 #
612 whole_part_with_layer=DataArrayInt.Aggregate([p0,layer_of_part0])
613 whole_part_with_layer.sort()
614 part0_with_layer=m[whole_part_with_layer]
615 #! [UG_Optimization_2]
616