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