]> SALOME platform Git repositories - tools/solverlab.git/blob - CoreFlows/examples/Python/MPI4PY/testSendRecvFieldDifferentMeshes.py
Salome HOME
Add missing folder from previous commit
[tools/solverlab.git] / CoreFlows / examples / Python / MPI4PY / testSendRecvFieldDifferentMeshes.py
1 #!/usr/bin/env python3
2 # -*-coding:utf-8 -*
3
4 #===============================================================================================================================
5 # Name        : Tests of sending and receiving 2D MEDCoupling fields on nodes (P1) lying on different meshes between two processors
6 # Author      : Michaël Ndjinga
7 # Copyright   : CEA Saclay 2021
8 # Description : Use of the parallel Data Exchange Channel InterpKernelDEC of MEDCoupling
9 #================================================================================================================================
10
11 from mpi4py import MPI
12 import medcoupling as mc
13
14 comm = MPI.COMM_WORLD
15 size = comm.Get_size()
16 rank = comm.Get_rank()
17
18 if(size!=2):
19         raise ValueError("Processor ", rank, " : aborting.\n Simulation should done on two processors.\n", size, " processors given")
20         
21 print("My rank is ", rank, " among ", size, "processors")
22
23 procs_source = [0]
24 procs_target = [1]
25
26 interface = mc.CommInterface()
27 source_group = mc.MPIProcessorGroup(interface, procs_source)
28 target_group = mc.MPIProcessorGroup(interface, procs_target)
29 dec = mc.InterpKernelDEC(source_group, target_group)
30
31 # Create a MEDCouplingUMesh from a 3D cartesian mesh
32 xarr=mc.DataArrayDouble.New(11,1)
33 xarr.iota(0.)
34 cmesh=mc.MEDCouplingCMesh.New()
35 cmesh.setCoords(xarr,xarr)
36 mesh=cmesh.buildUnstructured()
37 mesh.setName("RegularSquare")
38 mesh.simplexize(rank)#The squares are cut in two right triangles according to one of the two possible diagonals
39
40 #Create a field by application of an analytic function 
41 if source_group.containsMyRank():
42         field=mesh.fillFromAnalytic(mc.ON_NODES,1,"(x-5.)*(x-5.)+(y-5.)*(y-5.)")
43         field.setName("SourceField")
44         #field.setNature(mc.ExtensiveConservation)
45         field.setNature(mc.IntensiveMaximum)
46         mc.WriteField("source_field.med", field, True)
47         print("Processor ", rank, " has created and saved the source field")
48 else:
49         field=mesh.fillFromAnalytic(mc.ON_NODES,1,"0")
50         field.setName("TargetField")
51         #field.setNature(mc.ExtensiveConservation)
52         field.setNature(mc.IntensiveMaximum)
53         print("Processor ", rank, " has created the target field")
54         
55 dec.setMethod("P1")
56 dec.attachLocalField(field)
57 dec.synchronize()
58
59 if source_group.containsMyRank():
60         dec.sendData()
61         print("Processor ", rank, " has sent the source field")
62 else:
63         dec.recvData()
64         print("Processor ", rank, " has received the source field on the target mesh")
65         exact_field=mesh.fillFromAnalytic(mc.ON_NODES,1,"(x-5.)*(x-5.)+(y-5.)*(y-5.)")
66         exact_field.setName("ExactField")
67         #Computing maximum error
68         coordsArr = mesh.getCoords()
69         values=field.getValueOnMulti(coordsArr)
70         exact_values=exact_field.getValueOnMulti(coordsArr)
71         error=(values-exact_values).normMax()/exact_values.normMax()
72         print("Processor ", rank, " received source field differs from theoretical value by less than", error, " (maximum relative norm on node values)" )
73         assert error<1.e-1
74         mc.WriteField("target_field.med", field, True)
75         mc.WriteField("exact_field.med", exact_field, True)