Salome HOME
514251c47bfb6da06cf2663ee074c383ad38807a
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_FabienAPI.cxx
1 // Copyright (C) 2007-2019  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 "CommInterface.hxx"
22 #include "ProcessorGroup.hxx"
23 #include "MPIProcessorGroup.hxx"
24 #include "InterpKernelDEC.hxx"
25 #include "MEDCouplingUMesh.hxx"
26 #include "ParaMESH.hxx"
27 #include "ParaFIELD.hxx"
28 #include "ComponentTopology.hxx"
29
30 #include <set>
31
32 using namespace MEDCoupling;
33
34 void ParaMEDMEMTest::testFabienAPI1()
35 {
36   int size;
37   int rank;
38   MPI_Comm_size(MPI_COMM_WORLD,&size);
39   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
40   //
41   if(size!=3)
42     return ;
43   int procs_source_c[1]={0};
44   std::set<int> procs_source(procs_source_c,procs_source_c+1);
45   int procs_target_c[1]={1};
46   std::set<int> procs_target(procs_target_c,procs_target_c+1);
47   //
48   MEDCoupling::MEDCouplingUMesh *mesh=0;
49   MEDCoupling::ParaMESH *paramesh=0;
50   MEDCoupling::ParaFIELD *parafield=0;
51   //
52   MEDCoupling::CommInterface interface;
53   //
54   MPI_Barrier(MPI_COMM_WORLD);
55   double targetCoords[8]={ 0.,0., 1., 0., 0., 1., 1., 1. };
56   CommInterface comm;
57   //
58   MEDCoupling::InterpKernelDEC *dec=new MEDCoupling::InterpKernelDEC(procs_source,procs_target);
59   if(dec->isInSourceSide())
60     {    
61       mesh=MEDCouplingUMesh::New();
62       mesh->setMeshDimension(2);
63       DataArrayDouble *myCoords=DataArrayDouble::New();
64       myCoords->alloc(4,2);
65       std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
66       mesh->setCoords(myCoords);
67       myCoords->decrRef();
68       mcIdType targetConn[4]={0,2,3,1};
69       mesh->allocateCells(1);
70       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn);
71       mesh->finishInsertingCells();
72       MEDCoupling::ComponentTopology comptopo;
73       paramesh=new ParaMESH(mesh,*dec->getSourceGrp(),"source mesh");
74       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
75       parafield->getField()->setNature(IntensiveMaximum);
76       double *vals=parafield->getField()->getArray()->getPointer();
77       vals[0]=7.;
78     }
79   if(dec->isInTargetSide())
80     {
81       mesh=MEDCouplingUMesh::New();
82       mesh->setMeshDimension(2);
83       DataArrayDouble *myCoords=DataArrayDouble::New();
84       myCoords->alloc(4,2);
85       std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
86       mesh->setCoords(myCoords);
87       myCoords->decrRef();
88       mcIdType targetConn[6]={0,2,1,2,3,1};
89       mesh->allocateCells(2);
90       mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
91       mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
92       mesh->finishInsertingCells();
93       MEDCoupling::ComponentTopology comptopo;
94       paramesh=new ParaMESH(mesh,*dec->getTargetGrp(),"target mesh");
95       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
96       parafield->getField()->setNature(IntensiveMaximum);
97     }
98   dec->attachLocalField(parafield);
99   dec->synchronize();
100   dec->sendRecvData();
101   if(dec->isInTargetSide())
102     {
103       const double *valsToTest=parafield->getField()->getArray()->getConstPointer();
104       CPPUNIT_ASSERT_DOUBLES_EQUAL(valsToTest[0],7.,1e-14);
105       CPPUNIT_ASSERT_DOUBLES_EQUAL(valsToTest[1],7.,1e-14);
106     }
107   //
108   delete parafield;
109   delete paramesh;
110   if(mesh)
111     mesh->decrRef();
112   delete dec;
113   MPI_Barrier(MPI_COMM_WORLD);
114 }
115
116 /*!
117  * Idem testFabienAPI1 except that procs are shuffled. Test of the good management of group translation in newly created communicator.
118  */
119 void ParaMEDMEMTest::testFabienAPI2()
120 {
121   int size;
122   int rank;
123   MPI_Comm_size(MPI_COMM_WORLD,&size);
124   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
125   //
126   if(size!=3)
127     return ;
128   int procs_source_c[1]={2};//difference with testFabienAPI1
129   std::set<int> procs_source(procs_source_c,procs_source_c+1);
130   int procs_target_c[1]={1};
131   std::set<int> procs_target(procs_target_c,procs_target_c+1);
132   //
133   MEDCoupling::MEDCouplingUMesh *mesh=0;
134   MEDCoupling::ParaMESH *paramesh=0;
135   MEDCoupling::ParaFIELD *parafield=0;
136   //
137   MEDCoupling::CommInterface interface;
138   //
139   MPI_Barrier(MPI_COMM_WORLD);
140   double targetCoords[8]={ 0.,0., 1., 0., 0., 1., 1., 1. };
141   CommInterface comm;
142   //
143   MEDCoupling::InterpKernelDEC *dec=new MEDCoupling::InterpKernelDEC(procs_source,procs_target);
144   if(dec->isInSourceSide())
145     {    
146       mesh=MEDCouplingUMesh::New();
147       mesh->setMeshDimension(2);
148       DataArrayDouble *myCoords=DataArrayDouble::New();
149       myCoords->alloc(4,2);
150       std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
151       mesh->setCoords(myCoords);
152       myCoords->decrRef();
153       mcIdType targetConn[4]={0,2,3,1};
154       mesh->allocateCells(1);
155       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn);
156       mesh->finishInsertingCells();
157       MEDCoupling::ComponentTopology comptopo;
158       paramesh=new ParaMESH(mesh,*dec->getSourceGrp(),"source mesh");
159       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
160       parafield->getField()->setNature(IntensiveMaximum);
161       double *vals=parafield->getField()->getArray()->getPointer();
162       vals[0]=7.;
163     }
164   if(dec->isInTargetSide())
165     {
166       mesh=MEDCouplingUMesh::New();
167       mesh->setMeshDimension(2);
168       DataArrayDouble *myCoords=DataArrayDouble::New();
169       myCoords->alloc(4,2);
170       std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
171       mesh->setCoords(myCoords);
172       myCoords->decrRef();
173       mcIdType targetConn[6]={0,2,1,2,3,1};
174       mesh->allocateCells(2);
175       mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
176       mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
177       mesh->finishInsertingCells();
178       MEDCoupling::ComponentTopology comptopo;
179       paramesh=new ParaMESH(mesh,*dec->getTargetGrp(),"target mesh");
180       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
181       parafield->getField()->setNature(IntensiveMaximum);
182     }
183   dec->attachLocalField(parafield);
184   dec->synchronize();
185   dec->sendRecvData();
186   if(dec->isInTargetSide())
187     {
188       const double *valsToTest=parafield->getField()->getArray()->getConstPointer();
189       CPPUNIT_ASSERT_DOUBLES_EQUAL(valsToTest[0],7.,1e-14);
190       CPPUNIT_ASSERT_DOUBLES_EQUAL(valsToTest[1],7.,1e-14);
191     }
192   //
193   delete parafield;
194   delete paramesh;
195   if(mesh)
196     mesh->decrRef();
197   delete dec;
198   MPI_Barrier(MPI_COMM_WORLD);
199 }