Salome HOME
Add missing folder from previous commit
[tools/solverlab.git] / CoreFlows / examples / Python / MPI4PY / testSendRecvFieldSameMesh.py
1 #!/usr/bin/env python3
2 # -*-coding:utf-8 -*
3
4 #===============================================================================================================================
5 # Name        : Tests of sending and receiving a 3D MEDCoupling field on cells (P0) lying on the same mesh between two processors
6 # Author      : MichaĆ«l Ndjinga
7 # Copyright   : CEA Saclay 2021
8 # Description : Use of the parallel Data Exchange Channel StructuredCoincidentDEC 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.StructuredCoincidentDEC(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,xarr)
36 mesh=cmesh.buildUnstructured()
37 mesh.setName("RegularSquare")
38
39 #Create a field by application of an analytic function 
40 if source_group.containsMyRank():
41         field=mesh.fillFromAnalytic(mc.ON_CELLS,1,"(x-5.)*(x-5.)+(y-5.)*(y-5.)+(z-5.)*(z-5.)")
42         field.setName("SourceField")
43         mc.WriteField("source_field.med", field, True)
44         print("Processor ", rank, " has created and saved the source field")
45 else:
46         field=mesh.fillFromAnalytic(mc.ON_CELLS,1,"0")
47         field.setName("TargetField")
48         print("Processor ", rank, " has created the target field")
49         
50 dec.attachLocalField(field)
51 dec.synchronize()
52
53 if source_group.containsMyRank():
54         dec.sendData()
55         print("Processor ", rank, " has sent the source field")
56 else:
57         dec.recvData()
58         print("Processor ", rank, " has received the source field on the target mesh")
59         exact_field=mesh.fillFromAnalytic(mc.ON_CELLS,1,"(x-5.)*(x-5.)+(y-5.)*(y-5.)+(z-5.)*(z-5.)")
60         exact_field.setName("ExactField")
61         error=(field-exact_field).normL2()[0]
62         print("Processor ", rank, " received source field differs from theoretical value by ", error, " (L2 norm on cells)" )
63         assert abs(error)<1.e-6
64         mc.WriteField("target_field.med", field, True)
65         mc.WriteField("exact_field.med", exact_field, True)