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