4 #===============================================================================================================================
5 # Name : Tests of using a subcommnicator for sending and receiving a 3D MEDCoupling field on cells (P0) lying on the same mesh between two groups of 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 #================================================================================================================================
11 from mpi4py import MPI
12 import medcoupling as mc
15 size = MPI.COMM_WORLD.Get_size()
16 rank = MPI.COMM_WORLD.Get_rank()
19 raise ValueError("Processor ", rank, " : aborting.\n Simulation should done on four processors.\n", size, " processors given")
23 print("My rank is ", rank, " among ", size, "processors, my color is ", color)
25 sub_comm = MPI.COMM_WORLD.Split(color, rank); # two groups (0,2) and (1,3)
26 sub_rank = sub_comm.Get_rank();
27 sub_size = sub_comm.Get_size();
32 print("WORLD RANK/SIZE: ",rank,"/",size," \t subcommunicator RANK/SIZE: ",sub_rank,"/",sub_size,"\n")
34 interface = mc.CommInterface()
35 source_group = mc.MPIProcessorGroup(interface, procs_source,sub_comm)
36 target_group = mc.MPIProcessorGroup(interface, procs_target,sub_comm)
37 dec = mc.StructuredCoincidentDEC(source_group, target_group)
39 # Create a MEDCouplingUMesh from a 3D cartesian mesh
40 xarr=mc.DataArrayDouble.New(11,1)
42 cmesh=mc.MEDCouplingCMesh.New()
43 cmesh.setCoords(xarr,xarr,xarr)
44 mesh=cmesh.buildUnstructured()
45 mesh.setName("RegularSquare")
47 #Create a field by application of an analytic function
48 if source_group.containsMyRank():
49 field=mesh.fillFromAnalytic(mc.ON_CELLS,1,"(x-5.)*(x-5.)+(y-5.)*(y-5.)+(z-5.)*(z-5.)")
50 field.setName("SourceField")
51 mc.WriteField("source_field"+str(rank)+".med", field, True)
52 print("Processor ", rank, " has created and saved the source field")
54 field=mesh.fillFromAnalytic(mc.ON_CELLS,1,"0")
55 field.setName("TargetField")
56 print("Processor ", rank, " has created the target field")
58 dec.attachLocalField(field)
61 if source_group.containsMyRank():
63 print("Processor ", rank, " has sent the source field")
64 elif target_group.containsMyRank():
66 print("Processor ", rank, " has received the source field on the target mesh")
67 exact_field=mesh.fillFromAnalytic(mc.ON_CELLS,1,"(x-5.)*(x-5.)+(y-5.)*(y-5.)+(z-5.)*(z-5.)")
68 exact_field.setName("ExactField")
69 error=(field-exact_field).normMax()[0]/exact_field.normMax()[0]
70 print("Processor ", rank, " received source field that differs from theoretical value by ", error, " (maximum relative norm on cell values)" )
71 assert abs(error)<1.e-6
72 mc.WriteField("target_field"+str(rank)+".med", field, True)
73 mc.WriteField("exact_field"+str(rank)+".med", exact_field, True)
75 print("Processor ", rank, " did nothing" )