Salome HOME
042d45b946266dd73a1875fc6fe40f559b355207
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_StructuredCoincidentDEC.cxx
1 // Copyright (C) 2007-2020  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 "ICoCoMEDField.hxx"
33 #include "MEDLoader.hxx"
34 #include "TestInterpKernelUtils.hxx"
35
36 #include <string>
37
38 // use this define to enable lines, execution of which leads to Segmentation Fault
39 #define ENABLE_FAULTS
40
41 // use this define to enable CPPUNIT asserts and fails, showing bugs
42 #define ENABLE_FORCED_FAILURES
43
44 using namespace std;
45 using namespace MEDCoupling;
46
47 /*
48  * Check methods defined in StructuredCoincidentDEC.hxx
49  *
50  StructuredCoincidentDEC();
51  StructuredCoincidentDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
52  virtual ~StructuredCoincidentDEC();
53  void synchronize();
54  void recvData();
55  void sendData();
56 */
57
58 void ParaMEDMEMTest::testStructuredCoincidentDEC() {
59   string testname="ParaMEDMEM - testStructured CoincidentDEC";
60   //  MPI_Init(&argc, &argv); 
61   int size;
62   int rank;
63   MPI_Comm_size(MPI_COMM_WORLD, &size);
64   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
65   if (size<4) {
66     return;
67   }
68
69   MEDCoupling::CommInterface interface;
70
71   MEDCoupling::MPIProcessorGroup self_group (interface,rank,rank);
72   MEDCoupling::MPIProcessorGroup target_group(interface,3,size-1);
73   MEDCoupling::MPIProcessorGroup source_group (interface,0,2);
74
75   MEDCoupling::MEDCouplingUMesh* mesh = nullptr;
76   MEDCoupling::ParaMESH* paramesh = nullptr;
77   MEDCoupling::ParaFIELD* parafield = nullptr;
78
79   string filename_xml1 = INTERP_TEST::getResourceFile("square1_split");
80   string filename_2    = INTERP_TEST::getResourceFile("square1.med");
81   //string filename_seq_wr  = makeTmpFile("");
82   //string filename_seq_med = makeTmpFile("myWrField_seq_pointe221.med");
83
84   // To remove tmp files from disk
85   ParaMEDMEMTest_TmpFilesRemover aRemover;
86
87   //loading the geometry for the source group
88
89   MEDCoupling::StructuredCoincidentDEC dec(source_group, target_group);
90
91   MPI_Barrier(MPI_COMM_WORLD);
92   if (source_group.containsMyRank()) {
93     string master = filename_xml1;
94
95     ostringstream strstream;
96     strstream <<master<<rank+1<<".med";
97     ostringstream meshname;
98     meshname<< "Mesh_2_"<< rank+1;
99
100     mesh=ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
101     
102
103     paramesh=new ParaMESH (mesh,source_group,"source mesh");
104
105     MEDCoupling::ComponentTopology comptopo(6);
106     parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
107
108     int nb_local=mesh->getNumberOfCells();
109     const mcIdType* global_numbering = paramesh->getGlobalNumberingCell();
110     
111     double *value=parafield->getField()->getArray()->getPointer();
112     for(int ielem=0; ielem<nb_local;ielem++)
113       for (int icomp=0; icomp<6; icomp++)
114         value[ielem*6+icomp]=global_numbering[ielem]*6+icomp;
115
116     //ICoCo::Field* icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
117
118     dec.attachLocalField(parafield);
119     dec.synchronize();
120     dec.sendData();
121     //delete icocofield;
122   }
123
124   //loading the geometry for the target group
125   if (target_group.containsMyRank()) {
126
127     string meshname2("Mesh_2");
128     mesh = ReadUMeshFromFile(filename_2.c_str(),meshname2.c_str(),0);
129     
130     paramesh=new ParaMESH (mesh,self_group,"target mesh");
131     MEDCoupling::ComponentTopology comptopo(6, &target_group);
132
133     parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
134
135     int nb_local=mesh->getNumberOfCells();
136     double *value=parafield->getField()->getArray()->getPointer();
137     for (int ielem=0; ielem<nb_local; ielem++)
138       for (int icomp=0; icomp<comptopo.nbLocalComponents(); icomp++)
139         value[ielem*comptopo.nbLocalComponents()+icomp]=0.0;
140     //ICoCo::Field* icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
141
142     dec.attachLocalField(parafield);
143     dec.synchronize();
144     dec.recvData();
145
146     //checking validity of field
147     const double* recv_value = parafield->getField()->getArray()->getPointer();
148     for (int i=0; i< nb_local; i++) {
149       int first = comptopo.firstLocalComponent();
150       for (int icomp = 0; icomp < comptopo.nbLocalComponents(); icomp++)
151         CPPUNIT_ASSERT_DOUBLES_EQUAL(recv_value[i*comptopo.nbLocalComponents()+icomp],(double)(i*6+icomp+first),1e-12);
152     }
153     //delete icocofield;
154   }
155   delete parafield;
156   delete paramesh;
157   mesh->decrRef();
158
159   //  MPI_Barrier(MPI_COMM_WORLD);
160
161 }