]> SALOME platform Git repositories - tools/medcoupling.git/blob - src/ParaMEDMEMTest/ParaMEDMEMTestMPI2_1.cxx
Salome HOME
Merge remote branch 'origin/abn/icoco_fix' into V7_5_BR
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTestMPI2_1.cxx
1 // Copyright (C) 2007-2014  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 <cppunit/extensions/HelperMacros.h>
21
22 #include "MPI2Connector.hxx"
23 #include "ParaMESH.hxx"
24 #include "ParaFIELD.hxx"
25 #include "MEDCouplingUMesh.hxx"
26 #include "MEDCouplingFieldDouble.hxx"
27 #include "InterpKernelDEC.hxx"
28 #include "MPIProcessorGroup.hxx"
29 #include "CommInterface.hxx"
30
31 #include <mpi.h>
32 #include <iostream>
33 #include <stdlib.h>
34
35 class MPI2ParaMEDMEMTest : public CppUnit::TestFixture
36 {
37   CPPUNIT_TEST_SUITE( MPI2ParaMEDMEMTest );
38   CPPUNIT_TEST( testBasicMPI2_1 );
39   CPPUNIT_TEST_SUITE_END();
40 public:
41   void testBasicMPI2_1();
42 };
43
44 using namespace ParaMEDMEM;
45
46 void MPI2ParaMEDMEMTest::testBasicMPI2_1()
47 {
48   int lsize, lrank, gsize, grank;
49   MPI_Comm gcom;
50   std::string service = "SERVICE";
51   std::ostringstream meshfilename, meshname;
52   ParaMEDMEM::ParaMESH *paramesh=0;
53   ParaMEDMEM::MEDCouplingUMesh *mesh;
54   ParaMEDMEM::ParaFIELD *parafield=0;
55   ParaMEDMEM::CommInterface *interface;
56   ParaMEDMEM::MPIProcessorGroup *source, *target;
57
58   MPI_Comm_size( MPI_COMM_WORLD, &lsize );
59   MPI_Comm_rank( MPI_COMM_WORLD, &lrank );
60   if(lsize!=2)
61     {
62       CPPUNIT_ASSERT(false);
63       return;
64     }
65
66   /* Connection to remote programm */
67   MPI2Connector *mpio = new MPI2Connector;
68   gcom = mpio->remoteMPI2Connect(service);
69   MPI_Comm_size( gcom, &gsize );
70   MPI_Comm_rank( gcom, &grank );
71   if(gsize!=5)
72     {
73       CPPUNIT_ASSERT(false);
74       return;
75     }
76   interface = new ParaMEDMEM::CommInterface;
77   source = new ParaMEDMEM::MPIProcessorGroup(*interface,0,lsize-1,gcom);
78   target = new ParaMEDMEM::MPIProcessorGroup(*interface,lsize,gsize-1,gcom);
79
80   const double sourceCoordsAll[2][8]={{0.4,0.5,0.4,1.5,1.6,1.5,1.6,0.5},
81                                       {0.3,-0.5,1.6,-0.5,1.6,-1.5,0.3,-1.5}};
82   
83   int conn4All[8]={0,1,2,3,4,5,6,7};
84   
85   std::ostringstream stream; stream << "sourcemesh2D proc " << grank;
86   mesh=MEDCouplingUMesh::New(stream.str().c_str(),2);
87   mesh->allocateCells(2);
88   mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All);
89   mesh->finishInsertingCells();
90   DataArrayDouble *myCoords=DataArrayDouble::New();
91   myCoords->alloc(4,2);
92   const double *sourceCoords=sourceCoordsAll[grank];
93   std::copy(sourceCoords,sourceCoords+8,myCoords->getPointer());
94   mesh->setCoords(myCoords);
95   myCoords->decrRef();
96   paramesh=new ParaMESH(mesh,*source,"source mesh");
97   ParaMEDMEM::ComponentTopology comptopo;
98   parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
99   double *value=parafield->getField()->getArray()->getPointer();
100   value[0]=34+13*((double)grank);
101
102   ParaMEDMEM::InterpKernelDEC dec(*source,*target);
103   parafield->getField()->setNature(ConservativeVolumic);
104
105
106   dec.setMethod("P0");
107   dec.attachLocalField(parafield);
108   dec.synchronize();
109   dec.setForcedRenormalization(false);
110   dec.sendData();
111   /* Deconnection of remote programm */
112   mpio->remoteMPI2Disconnect(service);
113   /* clean-up */
114   delete mpio;
115   delete parafield;
116   mesh->decrRef();
117   delete paramesh;
118   delete source;
119   delete target;
120   delete interface;
121 }
122
123 CPPUNIT_TEST_SUITE_REGISTRATION( MPI2ParaMEDMEMTest );
124
125 #include "MPIMainTest.hxx"