Salome HOME
Doc: indicating how to pass MPI_Comm from mpi4py
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_StructuredCoincidentDEC.cxx
1 // Copyright (C) 2007-2022  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
22
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "Topology.hxx"
27 #include "DEC.hxx"
28 #include "StructuredCoincidentDEC.hxx"
29 #include "ParaMESH.hxx"
30 #include "ParaFIELD.hxx"
31 #include "ComponentTopology.hxx"
32 #include "ICoCoMEDDoubleField.hxx"
33 #include "MEDLoader.hxx"
34 #include "MEDCouplingUMesh.hxx"
35 #include "TestInterpKernelUtils.hxx"
36
37 #include <string>
38
39 // use this define to enable lines, execution of which leads to Segmentation Fault
40 #define ENABLE_FAULTS
41
42 // use this define to enable CPPUNIT asserts and fails, showing bugs
43 #define ENABLE_FORCED_FAILURES
44
45 using namespace std;
46 using namespace MEDCoupling;
47
48 /*
49  * Check methods defined in StructuredCoincidentDEC.hxx
50  *
51  StructuredCoincidentDEC();
52  StructuredCoincidentDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
53  virtual ~StructuredCoincidentDEC();
54  void synchronize();
55  void recvData();
56  void sendData();
57 */
58
59 void ParaMEDMEMTest::testStructuredCoincidentDEC() {
60   string testname="ParaMEDMEM - testStructured CoincidentDEC";
61   //  MPI_Init(&argc, &argv); 
62   int size;
63   int rank;
64   MPI_Comm_size(MPI_COMM_WORLD, &size);
65   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
66   if (size<4) {
67     return;
68   }
69
70   MEDCoupling::CommInterface interface;
71
72   MEDCoupling::MPIProcessorGroup self_group (interface,rank,rank);
73   MEDCoupling::MPIProcessorGroup target_group(interface,3,size-1);
74   MEDCoupling::MPIProcessorGroup source_group (interface,0,2);
75
76   MEDCoupling::MEDCouplingUMesh* mesh = nullptr;
77   MEDCoupling::ParaMESH* paramesh = nullptr;
78   MEDCoupling::ParaFIELD* parafield = nullptr;
79
80   string filename_xml1 = INTERP_TEST::getResourceFile("square1_split");
81   string filename_2    = INTERP_TEST::getResourceFile("square1.med");
82   //string filename_seq_wr  = makeTmpFile("");
83   //string filename_seq_med = makeTmpFile("myWrField_seq_pointe221.med");
84
85   // To remove tmp files from disk
86   ParaMEDMEMTest_TmpFilesRemover aRemover;
87
88   //loading the geometry for the source group
89
90   MEDCoupling::StructuredCoincidentDEC dec(source_group, target_group);
91
92   MPI_Barrier(MPI_COMM_WORLD);
93   if (source_group.containsMyRank()) {
94     string master = filename_xml1;
95
96     ostringstream strstream;
97     strstream <<master<<rank+1<<".med";
98     ostringstream meshname;
99     meshname<< "Mesh_2_"<< rank+1;
100
101     mesh=ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
102     
103
104     paramesh=new ParaMESH (mesh,source_group,"source mesh");
105
106     MEDCoupling::ComponentTopology comptopo(6);
107     parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
108
109     int nb_local=mesh->getNumberOfCells();
110     const mcIdType* global_numbering = paramesh->getGlobalNumberingCell();
111     
112     double *value=parafield->getField()->getArray()->getPointer();
113     for(int ielem=0; ielem<nb_local;ielem++)
114       for (int icomp=0; icomp<6; icomp++)
115         value[ielem*6+icomp]=global_numbering[ielem]*6+icomp;
116
117     //ICoCo::Field* icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
118
119     dec.attachLocalField(parafield);
120     dec.synchronize();
121     dec.sendData();
122     //delete icocofield;
123   }
124
125   //loading the geometry for the target group
126   if (target_group.containsMyRank()) {
127
128     string meshname2("Mesh_2");
129     mesh = ReadUMeshFromFile(filename_2.c_str(),meshname2.c_str(),0);
130     
131     paramesh=new ParaMESH (mesh,self_group,"target mesh");
132     MEDCoupling::ComponentTopology comptopo(6, &target_group);
133
134     parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
135
136     int nb_local=mesh->getNumberOfCells();
137     double *value=parafield->getField()->getArray()->getPointer();
138     for (int ielem=0; ielem<nb_local; ielem++)
139       for (int icomp=0; icomp<comptopo.nbLocalComponents(); icomp++)
140         value[ielem*comptopo.nbLocalComponents()+icomp]=0.0;
141     //ICoCo::Field* icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
142
143     dec.attachLocalField(parafield);
144     dec.synchronize();
145     dec.recvData();
146
147     //checking validity of field
148     const double* recv_value = parafield->getField()->getArray()->getPointer();
149     for (int i=0; i< nb_local; i++) {
150       int first = comptopo.firstLocalComponent();
151       for (int icomp = 0; icomp < comptopo.nbLocalComponents(); icomp++)
152         CPPUNIT_ASSERT_DOUBLES_EQUAL(recv_value[i*comptopo.nbLocalComponents()+icomp],(double)(i*6+icomp+first),1e-12);
153     }
154     //delete icocofield;
155   }
156   delete parafield;
157   delete paramesh;
158   mesh->decrRef();
159
160   //  MPI_Barrier(MPI_COMM_WORLD);
161
162 }