Salome HOME
50% of work performed of porting tests.
[modules/med.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_Gauthier1.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 "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
22
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "DEC.hxx"
27 #include "InterpKernelDEC.hxx"
28 #include "ICoCoTrioField.hxx"
29 #include "MEDCouplingUMesh.hxx"
30 #include "MEDCouplingFieldDouble.hxx"
31 #include "ParaMESH.hxx"
32 #include "ParaFIELD.hxx"
33 #include "ComponentTopology.hxx"
34 #include "BlockTopology.hxx"
35
36 #include <set>
37 #include <time.h>
38 #include <iostream>
39 #include <assert.h>
40 #include <string>
41 #include <math.h>
42
43 using namespace std;
44 using namespace ParaMEDMEM;
45 using namespace ICoCo;
46
47 void afficheGauthier1(const ParaFIELD& field, const double *vals, int lgth)
48 {
49   const DataArrayDouble *valsOfField(field.getField()->getArray());
50   CPPUNIT_ASSERT_EQUAL(lgth,valsOfField->getNumberOfTuples());
51   for (int ele=0;ele<valsOfField->getNumberOfTuples();ele++)
52     CPPUNIT_ASSERT_DOUBLES_EQUAL(vals[ele],valsOfField->getIJ(ele,0),1e-12);
53 }
54
55 MEDCouplingUMesh *init_quadGauthier1(int is_master)
56 {
57   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> m(MEDCouplingUMesh::New("champ_quad",2));
58   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coo(DataArrayDouble::New());
59   if(is_master)
60     {
61       const double dataCoo[24]={0,0,0,1,0,0,0,0,1,1,0,1,0,1,0,1,1,0,0,1,1,1,1,1};
62       coo->alloc(8,3);
63       std::copy(dataCoo,dataCoo+24,coo->getPointer());
64       const int conn[8]={0,1,3,2,4,5,7,6};
65       m->allocateCells(2);
66       m->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
67       m->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+4);
68     }
69   else
70     {
71       coo->alloc(0,3);
72       m->allocateCells(0);
73     }
74   m->setCoords(coo);
75   return m.retn();
76 }
77
78 MEDCouplingUMesh *init_triangleGauthier1(int is_master)
79 {
80   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> m(MEDCouplingUMesh::New("champ_triangle",2));
81   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coo(DataArrayDouble::New());
82   if(is_master)
83     {
84       const double dataCoo[24]={0,0,0,1,0,0,0,0,1,1,0,1,0,1,0,1,1,0,0,1,1,1,1,1};
85       coo->alloc(8,3);
86       std::copy(dataCoo,dataCoo+24,coo->getPointer());
87       const int conn[12]={0,1,2,1,2,3,4,5,7,4,6,7};
88       m->allocateCells(2);
89       for(int i=0;i<4;i++)
90         m->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+3*i);
91     }
92   else
93     {
94       coo->alloc(0,3);
95       m->allocateCells(0);
96     }
97   m->setCoords(coo);
98   return m.retn();
99 }
100
101
102 void ParaMEDMEMTest::testGauthier1()
103 {
104   int num_cas=0;
105   int rank, size;
106   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
107   MPI_Comm_size(MPI_COMM_WORLD,&size);
108   
109   int is_master=0;
110
111   CommInterface comm;
112   set<int> emetteur_ids;
113   set<int> recepteur_ids;
114   emetteur_ids.insert(0);
115   if(size!=4)
116     return;
117   recepteur_ids.insert(1);
118   if (size >2) 
119     recepteur_ids.insert(2);
120   if (size >2) 
121     emetteur_ids.insert(3);
122   if ((rank==0)||(rank==1)) 
123     is_master=1;
124   
125   MPIProcessorGroup recepteur_group(comm,recepteur_ids);
126   MPIProcessorGroup emetteur_group(comm,emetteur_ids);
127
128   string cas;
129   if (recepteur_group.containsMyRank())
130     {
131       cas="recepteur";
132       //freopen("recpeteur.out","w",stdout);
133       //freopen("recepteur.err","w",stderr);
134     }
135   else
136     {
137       cas="emetteur";
138       // freopen("emetteur.out","w",stdout);
139       //freopen("emetteur.err","w",stderr);
140     }
141   double expected[8][4]={
142     {1.,1.,1.,1.},
143     {40., 40., 1., 1.},
144     {1.,1.,1e200,1e200},
145     {40.,1.,1e200,1e200},
146     {1.,1.,1.,1.},
147     {40.,1.,1.,1.},
148     {1.,1.,1e200,1e200},
149     {20.5,1.,1e200,1e200}
150   };
151   int expectedLgth[8]={4,4,2,2,4,4,2,2};
152   
153   for (int send=0;send<2;send++)
154     for (int rec=0;rec<2;rec++)
155       {
156         InterpKernelDEC dec_emetteur(emetteur_group, recepteur_group);
157         ParaMEDMEM::ParaFIELD *champ_emetteur(0),*champ_recepteur(0);
158         ParaMEDMEM::ParaMESH *paramesh(0);
159         MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh;
160         dec_emetteur.setOrientation(2);
161         if (send==0)
162           {
163             mesh=init_quadGauthier1(is_master);
164           }
165         else
166           {
167             mesh=init_triangleGauthier1(is_master);
168           }
169         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"emetteur mesh");
170         ParaMEDMEM::ComponentTopology comptopo;
171         champ_emetteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
172         champ_emetteur->getField()->setNature(ConservativeVolumic);
173         champ_emetteur->setOwnSupport(true);
174         if (rec==0)
175           {
176             mesh=init_triangleGauthier1(is_master);
177           }
178         else
179           {
180             mesh=init_quadGauthier1(is_master);
181           }
182         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"recepteur mesh");
183         champ_recepteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
184         champ_recepteur->getField()->setNature(ConservativeVolumic);
185         champ_recepteur->setOwnSupport(true);
186         if (cas=="emetteur") 
187           {
188             champ_emetteur->getField()->getArray()->fillWithValue(1.);
189           }
190   
191   
192         MPI_Barrier(MPI_COMM_WORLD);
193
194         //clock_t clock0= clock ();
195         int compti=0;
196
197         bool init=true; // first time step ??
198         bool stop=false;
199         //boucle sur les pas de quads
200         while (!stop) {
201   
202           compti++;
203           //clock_t clocki= clock ();
204           //cout << compti << " CLOCK " << (clocki-clock0)*1.e-6 << endl; 
205           for (int non_unif=0;non_unif<2;non_unif++)
206             {
207               if (cas=="emetteur") 
208                 {
209                   if (non_unif)
210                     if(rank!=3)
211                       champ_emetteur->getField()->getArray()->setIJ(0,0,40);
212                 }
213               //bool ok=false; // Is the time interval successfully solved ?
214     
215               // Loop on the time interval tries
216               if(1) {
217       
218
219                 if (cas=="emetteur")
220                   dec_emetteur.attachLocalField(champ_emetteur);
221                 else
222                   dec_emetteur.attachLocalField(champ_recepteur);
223
224
225                 if(init) dec_emetteur.synchronize();
226                 init=false;
227
228                 if (cas=="emetteur") {
229                   //    affiche(champ_emetteur);
230                   dec_emetteur.sendData();
231                 }
232                 else if (cas=="recepteur")
233                   {
234                     dec_emetteur.recvData();
235                     if (is_master)
236                       afficheGauthier1(*champ_recepteur,expected[num_cas],expectedLgth[num_cas]);
237                   }
238                 else
239                   throw 0;
240                 MPI_Barrier(MPI_COMM_WORLD);
241               }
242               stop=true;
243               num_cas++;
244             }
245         }
246         delete champ_emetteur;
247         delete champ_recepteur;
248       }
249 }
250
251 void ParaMEDMEMTest::testGauthier2()
252 {
253   const char save_vit_in_2[]="VITESSE_P1_OUT\n1\n2\n3\n63\n3\n80\n0\n 0 1 2\n 3 4 5\n 6 7 8\n 9 10 11\n 12 13 14\n 15 16 17\n 18 19 20\n 21 22 23\n 24 25 26\n 27 28 29\n 30 2 1\n 31 5 4\n 32 8 7\n 33 11 10\n 34 14 13\n 35 17 16\n 36 20 19\n 37 23 22\n 38 26 25\n 39 29 28\n 30 40 2\n 31 41 5\n 32 42 8\n 33 43 11\n 34 44 14\n 35 45 17\n 36 46 20\n 37 47 23\n 38 48 26\n 39 49 29\n 31 2 40\n 32 5 41\n 33 8 42\n 34 11 43\n 35 14 44\n 36 17 45\n 37 20 46\n 38 23 47\n 39 26 48\n 50 29 49\n 3 2 4\n 6 5 7\n 9 8 10\n 12 11 13\n 15 14 16\n 18 17 19\n 21 20 22\n 24 23 25\n 27 26 28\n 51 29 52\n 31 4 2\n 32 7 5\n 33 10 8\n 34 13 11\n 35 16 14\n 36 19 17\n 37 22 20\n 38 25 23\n 39 28 26\n 50 52 29\n 0 2 53\n 3 5 54\n 6 8 55\n 9 11 56\n 12 14 57\n 15 17 58\n 18 20 59\n 21 23 60\n 24 26 61\n 27 29 62\n 3 53 2\n 6 54 5\n 9 55 8\n 12 56 11\n 15 57 14\n 18 58 17\n 21 59 20\n 24 60 23\n 27 61 26\n 51 62 29\n 0 0 0\n 0.5 0 0\n 0.5 0.05 0\n 0 0.1 0\n 0.5 0.1 0\n 0.5 0.15 0\n 0 0.2 0\n 0.5 0.2 0\n 0.5 0.25 0\n 0 0.3 0\n 0.5 0.3 0\n 0.5 0.35 0\n 0 0.4 0\n 0.5 0.4 0\n 0.5 0.45 0\n 0 0.5 0\n 0.5 0.5 0\n 0.5 0.55 0\n 0 0.6 0\n 0.5 0.6 0\n 0.5 0.65 0\n 0 0.7 0\n 0.5 0.7 0\n 0.5 0.75 0\n 0 0.8 0\n 0.5 0.8 0\n 0.5 0.85 0\n 0 0.9 0\n 0.5 0.9 0\n 0.5 0.95 0\n 1 0 0\n 1 0.1 0\n 1 0.2 0\n 1 0.3 0\n 1 0.4 0\n 1 0.5 0\n 1 0.6 0\n 1 0.7 0\n 1 0.8 0\n 1 0.9 0\n 1 0.05 0\n 1 0.15 0\n 1 0.25 0\n 1 0.35 0\n 1 0.45 0\n 1 0.55 0\n 1 0.65 0\n 1 0.75 0\n 1 0.85 0\n 1 0.95 0\n 1 1 0\n 0 1 0\n 0.5 1 0\n 0 0.05 0\n 0 0.15 0\n 0 0.25 0\n 0 0.35 0\n 0 0.45 0\n 0 0.55 0\n 0 0.65 0\n 0 0.75 0\n 0 0.85 0\n 0 0.95 0\n2.9268\n3.1707\n3\n1\n 0 0 0\n 0 0 0\n 0 0 0.05\n 0 0 0.1\n 0 0 0.1\n 0 0 0.15\n 0 0 0.2\n 0 0 0.2\n 0 0 0.25\n 0 0 0.3\n 0 0 0.3\n 0 0 0.35\n 0 0 0.4\n 0 0 0.4\n 0 0 0.45\n 0 0 0.5\n 0 0 0.5\n 0 0 0.55\n 0 0 0.6\n 0 0 0.6\n 0 0 0.65\n 0 0 0.7\n 0 0 0.7\n 0 0 0.75\n 0 0 0.8\n 0 0 0.8\n 0 0 0.85\n 0 0 0.9\n 0 0 0.9\n 0 0 0.95\n 0 0 0\n 0 0 0.1\n 0 0 0.2\n 0 0 0.3\n 0 0 0.4\n 0 0 0.5\n 0 0 0.6\n 0 0 0.7\n 0 0 0.8\n 0 0 0.9\n 0 0 0.05\n 0 0 0.15\n 0 0 0.25\n 0 0 0.35\n 0 0 0.45\n 0 0 0.55\n 0 0 0.65\n 0 0 0.75\n 0 0 0.85\n 0 0 0.95\n 0 0 1\n 0 0 1\n 0 0 1\n 0 0 0.05\n 0 0 0.15\n 0 0 0.25\n 0 0 0.35\n 0 0 0.45\n 0 0 0.55\n 0 0 0.65\n 0 0 0.75\n 0 0 0.85\n 0 0 0.95\n1\n";
254
255   const char save_vit_out_0_2[]="vitesse_in_chaude\n0\n2\n3\n22\n4\n10\n-1081737852\n 0 1 3 2\n 2 3 5 4\n 4 5 7 6\n 6 7 9 8\n 8 9 11 10\n 10 11 13 12\n 12 13 15 14\n 14 15 17 16\n 16 17 19 18\n 18 19 21 20\n 0 0 0\n 1 0 0\n 0 0.1 0\n 1 0.1 0\n 0 0.2 0\n 1 0.2 0\n 0 0.3 0\n 1 0.3 0\n 0 0.4 0\n 1 0.4 0\n 0 0.5 0\n 1 0.5 0\n 0 0.6 0\n 1 0.6 0\n 0 0.7 0\n 1 0.7 0\n 0 0.8 0\n 1 0.8 0\n 0 0.9 0\n 1 0.9 0\n 0 1 0\n 1 1 0\n2.9268\n3.1707\n3\n1\n 0 0 0.05\n 0 0 0.15\n 0 0 0.25\n 0 0 0.35\n 0 0 0.45\n 0 0 0.55\n 0 0 0.65\n 0 0 0.75\n 0 0 0.85\n 0 0 0.95\n0\n";
256   const char save_vit_out_1_2[]="vitesse_in_chaude\n1\n2\n3\n22\n4\n10\n-1081737852\n 0 1 3 2\n 2 3 5 4\n 4 5 7 6\n 6 7 9 8\n 8 9 11 10\n 10 11 13 12\n 12 13 15 14\n 14 15 17 16\n 16 17 19 18\n 18 19 21 20\n 0 0 0\n 1 0 0\n 0 0.1 0\n 1 0.1 0\n 0 0.2 0\n 1 0.2 0\n 0 0.3 0\n 1 0.3 0\n 0 0.4 0\n 1 0.4 0\n 0 0.5 0\n 1 0.5 0\n 0 0.6 0\n 1 0.6 0\n 0 0.7 0\n 1 0.7 0\n 0 0.8 0\n 1 0.8 0\n 0 0.9 0\n 1 0.9 0\n 0 1 0\n 1 1 0\n2.9268\n3.1707\n3\n1\n 0 0 0.029375\n 0 0 0.029375\n 0 0 0.1\n 0 0 0.1\n 0 0 0.2\n 0 0 0.2\n 0 0 0.3\n 0 0 0.3\n 0 0 0.4\n 0 0 0.4\n 0 0 0.5\n 0 0 0.5\n 0 0 0.6\n 0 0 0.6\n 0 0 0.7\n 0 0 0.7\n 0 0 0.8\n 0 0 0.8\n 0 0 0.9\n 0 0 0.9\n 0 0 0.970625\n 0 0 0.970625\n0\n";
257
258   const char *save_vit_outs[2]={save_vit_out_1_2,save_vit_out_0_2};
259
260   const char save_vit_out_1_0[]="vitesse_in_chaude\n1\n2\n3\n22\n4\n10\n-1081737852\n 0 1 3 2\n 2 3 5 4\n 4 5 7 6\n 6 7 9 8\n 8 9 11 10\n 10 11 13 12\n 12 13 15 14\n 14 15 17 16\n 16 17 19 18\n 18 19 21 20\n 0 0 0\n 1 0 0\n 0 0.1 0\n 1 0.1 0\n 0 0.2 0\n 1 0.2 0\n 0 0.3 0\n 1 0.3 0\n 0 0.4 0\n 1 0.4 0\n 0 0.5 0\n 1 0.5 0\n 0 0.6 0\n 1 0.6 0\n 0 0.7 0\n 1 0.7 0\n 0 0.8 0\n 1 0.8 0\n 0 0.9 0\n 1 0.9 0\n 0 1 0\n 1 1 0\n2.9268\n3.1707\n3\n1\n 0 0 0.029375\n 0 0 0.029375\n 0 0 0.1\n 0 0 0.1\n 0 0 0.2\n 0 0 0.2\n 0 0 0.3\n 0 0 0.3\n 0 0 0.4\n 0 0 0.4\n 0 0 0.5\n 0 0 0.5\n 0 0 0.6\n 0 0 0.6\n 0 0 0.7\n 0 0 0.7\n 0 0 0.8\n 0 0 0.8\n 0 0 0.9\n 0 0 0.9\n 0 0 0.970625\n 0 0 0.970625\n0\n";
261   
262   const char save_vit_in[]="VITESSE_P1_OUT\n1\n2\n3\n63\n3\n80\n0\n 0 1 2\n 3 4 5\n 6 7 8\n 9 10 11\n 12 13 14\n 15 16 17\n 18 19 20\n 21 22 23\n 24 25 26\n 27 28 29\n 30 2 1\n 31 5 4\n 32 8 7\n 33 11 10\n 34 14 13\n 35 17 16\n 36 20 19\n 37 23 22\n 38 26 25\n 39 29 28\n 30 40 2\n 31 41 5\n 32 42 8\n 33 43 11\n 34 44 14\n 35 45 17\n 36 46 20\n 37 47 23\n 38 48 26\n 39 49 29\n 31 2 40\n 32 5 41\n 33 8 42\n 34 11 43\n 35 14 44\n 36 17 45\n 37 20 46\n 38 23 47\n 39 26 48\n 50 29 49\n 3 2 4\n 6 5 7\n 9 8 10\n 12 11 13\n 15 14 16\n 18 17 19\n 21 20 22\n 24 23 25\n 27 26 28\n 51 29 52\n 31 4 2\n 32 7 5\n 33 10 8\n 34 13 11\n 35 16 14\n 36 19 17\n 37 22 20\n 38 25 23\n 39 28 26\n 50 52 29\n 0 2 53\n 3 5 54\n 6 8 55\n 9 11 56\n 12 14 57\n 15 17 58\n 18 20 59\n 21 23 60\n 24 26 61\n 27 29 62\n 3 53 2\n 6 54 5\n 9 55 8\n 12 56 11\n 15 57 14\n 18 58 17\n 21 59 20\n 24 60 23\n 27 61 26\n 51 62 29\n 0 0 0\n 0.5 0 0\n 0.5 0.05 0\n 0 0.1 0\n 0.5 0.1 0\n 0.5 0.15 0\n 0 0.2 0\n 0.5 0.2 0\n 0.5 0.25 0\n 0 0.3 0\n 0.5 0.3 0\n 0.5 0.35 0\n 0 0.4 0\n 0.5 0.4 0\n 0.5 0.45 0\n 0 0.5 0\n 0.5 0.5 0\n 0.5 0.55 0\n 0 0.6 0\n 0.5 0.6 0\n 0.5 0.65 0\n 0 0.7 0\n 0.5 0.7 0\n 0.5 0.75 0\n 0 0.8 0\n 0.5 0.8 0\n 0.5 0.85 0\n 0 0.9 0\n 0.5 0.9 0\n 0.5 0.95 0\n 1 0 0\n 1 0.1 0\n 1 0.2 0\n 1 0.3 0\n 1 0.4 0\n 1 0.5 0\n 1 0.6 0\n 1 0.7 0\n 1 0.8 0\n 1 0.9 0\n 1 0.05 0\n 1 0.15 0\n 1 0.25 0\n 1 0.35 0\n 1 0.45 0\n 1 0.55 0\n 1 0.65 0\n 1 0.75 0\n 1 0.85 0\n 1 0.95 0\n 1 1 0\n 0 1 0\n 0.5 1 0\n 0 0.05 0\n 0 0.15 0\n 0 0.25 0\n 0 0.35 0\n 0 0.45 0\n 0 0.55 0\n 0 0.65 0\n 0 0.75 0\n 0 0.85 0\n 0 0.95 0\n2.9268\n3.1707\n3\n1\n 0 0 0\n 0 0 0\n 0 0 0.05\n 0 0 0.1\n 0 0 0.1\n 0 0 0.15\n 0 0 0.2\n 0 0 0.2\n 0 0 0.25\n 0 0 0.3\n 0 0 0.3\n 0 0 0.35\n 0 0 0.4\n 0 0 0.4\n 0 0 0.45\n 0 0 0.5\n 0 0 0.5\n 0 0 0.55\n 0 0 0.6\n 0 0 0.6\n 0 0 0.65\n 0 0 0.7\n 0 0 0.7\n 0 0 0.75\n 0 0 0.8\n 0 0 0.8\n 0 0 0.85\n 0 0 0.9\n 0 0 0.9\n 0 0 0.95\n 0 0 0\n 0 0 0.1\n 0 0 0.2\n 0 0 0.3\n 0 0 0.4\n 0 0 0.5\n 0 0 0.6\n 0 0 0.7\n 0 0 0.8\n 0 0 0.9\n 0 0 0.05\n 0 0 0.15\n 0 0 0.25\n 0 0 0.35\n 0 0 0.45\n 0 0 0.55\n 0 0 0.65\n 0 0 0.75\n 0 0 0.85\n 0 0 0.95\n 0 0 1\n 0 0 1\n 0 0 1\n 0 0 0.05\n 0 0 0.15\n 0 0 0.25\n 0 0 0.35\n 0 0 0.45\n 0 0 0.55\n 0 0 0.65\n 0 0 0.75\n 0 0 0.85\n 0 0 0.95\n1\n";
263
264   double valuesExpected1[2]={0.,0.};
265   double valuesExpected2[2]={0.95,0.970625};
266   
267   double valuesExpected30[]={0., 0., 0.05, 0., 0., 0.15, 0., 0., 0.25, 0., 0., 0.35, 0., 0., 0.45, 0., 0., 0.55, 0., 0., 0.65, 0., 0., 0.75, 0., 0., 0.85, 0., 0., 0.95};
268   double valuesExpected31[]={0.,  0.,  0.029375,  0.,  0.,  0.029375,  0.,  0.,  0.1,  0.,  0.,  0.1,  0.,  0.,  0.2,  0.,  0.,  0.2,  0.,  0.,  0.3,  0.,  0.,  0.3,  0.,  0.,  0.4,  0.,  0.,  0.4,  0.,  0.,  0.5,  0.,  0.,  0.5,  0.,  0.,  0.6,  0.,  0.,  0.6,  0.,  0.,  0.7,  0.,  0.,  0.7,  0.,  0.,  0.8,  0.,  0.,  0.8,  0.,  0.,  0.9,  0.,  0.,  0.9,  0.,  0.,  0.970625,  0.,  0.,  0.970625 };
269
270   double *valuesExpected3[2]={valuesExpected30,valuesExpected31};
271
272   int rank, size;
273   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
274   MPI_Comm_size(MPI_COMM_WORLD,&size);
275   if (size <2)
276     return ;
277   CommInterface comm;
278   set<int> Genepi_ids;
279   set<int> entree_chaude_ids;
280   Genepi_ids.insert(0);
281   for (int i=1;i<size;i++)
282     entree_chaude_ids.insert(i);
283   for (int type=0;type<2;type++)
284     {
285       MPIProcessorGroup entree_chaude_group(comm,entree_chaude_ids);
286       MPIProcessorGroup Genepi_group(comm,Genepi_ids);
287
288       TrioField vitesse;
289       InterpKernelDEC dec_vit_in_chaude(entree_chaude_group, Genepi_group);
290
291       if ( entree_chaude_group.containsMyRank())
292         {
293           istringstream save_vit(save_vit_in);
294           vitesse.restore(save_vit);
295         }
296       else
297         {
298           istringstream save_vit(save_vit_out_1_0);
299           vitesse.restore(save_vit);
300           vitesse._has_field_ownership=false;
301       
302           if (vitesse._field)
303             {
304               delete [] vitesse._field;
305               // cette ligne est super importante sinon c'est tout faux !!!!!!!
306               vitesse._field=0;
307             }
308           // pour tester P1->P0
309           vitesse._type=type;  
310         }
311   
312       if (vitesse._type==1)
313         dec_vit_in_chaude.setMethod("P1");
314   
315   
316
317       dec_vit_in_chaude.attachLocalField((ICoCo::Field*) &vitesse);
318       
319       dec_vit_in_chaude.synchronize();
320   
321   
322       // Envois - receptions
323       if (entree_chaude_group.containsMyRank())
324         {
325           dec_vit_in_chaude.sendData();
326         }
327       else
328         {
329           dec_vit_in_chaude.recvData(); 
330         }
331       if (entree_chaude_group.containsMyRank() )
332         {
333           if (1)
334             {
335               ostringstream save_vit(save_vit_in_2);
336               vitesse.save(save_vit);
337             }
338         }
339       else
340         {
341       
342           double pmin=1e38, pmax=-1e38;
343       
344           for(int i=0;i<vitesse.nb_values()*vitesse._nb_field_components;i++)
345             {
346               double p=*(vitesse._field+i);
347               if (p<pmin) pmin=p;
348               if (p>pmax) pmax=p;
349             }
350           CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesExpected1[type],pmin,1e-12);
351           CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesExpected2[type],pmax,1e-12);
352       
353           ostringstream save_vit(save_vit_outs[type]);
354           vitesse.save(save_vit);
355
356           for(int i=0;i<vitesse.nb_values();i++)
357             {
358               for(int c=0;c<vitesse._nb_field_components;c++)
359                 {
360                   double p=vitesse._field[i*vitesse._nb_field_components+c];
361                   CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesExpected3[type][i*vitesse._nb_field_components+c],p,1e-12);
362                 }
363             }
364       
365         }
366     }
367 }
368
369 /*!
370  * Non regression test testing copy constructor of InterpKernelDEC. 
371  */
372 void ParaMEDMEMTest::testGauthier3()
373 {
374   int num_cas=0;
375   int rank, size;
376   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
377   MPI_Comm_size(MPI_COMM_WORLD,&size);
378   
379   int is_master=0;
380
381   CommInterface comm;
382   set<int> emetteur_ids;
383   set<int> recepteur_ids;
384   emetteur_ids.insert(0);
385   if(size!=4)
386     return;
387   recepteur_ids.insert(1);
388   if (size >2) 
389     recepteur_ids.insert(2);
390   if (size >2) 
391     emetteur_ids.insert(3);
392   if ((rank==0)||(rank==1)) 
393     is_master=1;
394   
395   MPIProcessorGroup recepteur_group(comm,recepteur_ids);
396   MPIProcessorGroup emetteur_group(comm,emetteur_ids);
397
398   string cas;
399   if (recepteur_group.containsMyRank())
400     {
401       cas="recepteur";
402       //freopen("recpeteur.out","w",stdout);
403       //freopen("recepteur.err","w",stderr);
404     }
405   else
406     {
407       cas="emetteur";
408       // freopen("emetteur.out","w",stdout);
409       //freopen("emetteur.err","w",stderr);
410     }
411   double expected[8][4]={
412     {1.,1.,1.,1.},
413     {40., 40., 1., 1.},
414     {1.,1.,1e200,1e200},
415     {40.,1.,1e200,1e200},
416     {1.,1.,1.,1.},
417     {40.,1.,1.,1.},
418     {1.,1.,1e200,1e200},
419     {20.5,1.,1e200,1e200}
420   };
421   int expectedLgth[8]={4,4,2,2,4,4,2,2};
422   
423   for (int send=0;send<2;send++)
424     for (int rec=0;rec<2;rec++)
425       {
426         std::vector<InterpKernelDEC> decu(1);
427         decu[0]=InterpKernelDEC(emetteur_group,recepteur_group);
428         InterpKernelDEC& dec_emetteur=decu[0];
429         ParaMEDMEM::ParaFIELD *champ_emetteur(0),*champ_recepteur(0);
430         ParaMEDMEM::ParaMESH *paramesh(0);
431         MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh;
432         dec_emetteur.setOrientation(2);
433         if (send==0)
434           {
435             mesh=init_quadGauthier1(is_master);
436           }
437         else
438           {
439             mesh=init_triangleGauthier1(is_master);
440           }
441         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"emetteur mesh");
442         ParaMEDMEM::ComponentTopology comptopo;
443         champ_emetteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
444         champ_emetteur->getField()->setNature(ConservativeVolumic);
445         champ_emetteur->setOwnSupport(true);
446         if (rec==0)
447           {
448             mesh=init_triangleGauthier1(is_master);
449           }
450         else
451           {
452             mesh=init_quadGauthier1(is_master);
453           }
454         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"recepteur mesh");
455         champ_recepteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
456         champ_recepteur->getField()->setNature(ConservativeVolumic);
457         champ_recepteur->setOwnSupport(true);
458         if (cas=="emetteur") 
459           {
460             champ_emetteur->getField()->getArray()->fillWithValue(1.);
461           }
462   
463   
464         MPI_Barrier(MPI_COMM_WORLD);
465
466         //clock_t clock0= clock ();
467         int compti=0;
468
469         bool init=true; // first time step ??
470         bool stop=false;
471         //boucle sur les pas de quads
472         while (!stop) {
473   
474           compti++;
475           //clock_t clocki= clock ();
476           //cout << compti << " CLOCK " << (clocki-clock0)*1.e-6 << endl; 
477           for (int non_unif=0;non_unif<2;non_unif++)
478             {
479               if (cas=="emetteur") 
480                 {
481                   if (non_unif)
482                     if(rank!=3)
483                       champ_emetteur->getField()->getArray()->setIJ(0,0,40);
484                 }
485               //bool ok=false; // Is the time interval successfully solved ?
486     
487               // Loop on the time interval tries
488               if(1) {
489       
490
491                 if (cas=="emetteur")
492                   dec_emetteur.attachLocalField(champ_emetteur);
493                 else
494                   dec_emetteur.attachLocalField(champ_recepteur);
495
496
497                 if(init) dec_emetteur.synchronize();
498                 init=false;
499
500                 if (cas=="emetteur") {
501                   //    affiche(champ_emetteur);
502                   dec_emetteur.sendData();
503                 }
504                 else if (cas=="recepteur")
505                   {
506                     dec_emetteur.recvData();
507                     if (is_master)
508                       afficheGauthier1(*champ_recepteur,expected[num_cas],expectedLgth[num_cas]);
509                   }
510                 else
511                   throw 0;
512                 MPI_Barrier(MPI_COMM_WORLD);
513               }
514               stop=true;
515               num_cas++;
516             }
517         }
518         delete champ_emetteur;
519         delete champ_recepteur;
520       }
521 }
522
523 /*!
524  * This test is the parallel version of MEDCouplingBasicsTest.test3D1DOnP1P0_1 test.
525  */
526 void ParaMEDMEMTest::testGauthier4()
527 {
528   //
529   const double sourceCoords[19*3]={0.5,0.5,0.1,0.5,0.5,1.2,0.5,0.5,1.6,0.5,0.5,1.8,0.5,0.5,2.43,0.5,0.5,2.55,0.5,0.5,4.1,0.5,0.5,4.4,0.5,0.5,4.9,0.5,0.5,5.1,0.5,0.5,7.6,0.5,0.5,7.7,0.5,0.5,8.2,0.5,0.5,8.4,0.5,0.5,8.6,0.5,0.5,8.8,0.5,0.5,9.2,0.5,0.5,9.6,0.5,0.5,11.5};
530   const int sourceConn[18*2]={0,1,1,2,2,3,3,4,4,5,5,6,6,7,7,8,8,9,9,10,10,11,11,12,12,13,13,14,14,15,15,16,16,17,17,18};
531   const double sourceVals[19]={0.49,2.8899999999999997,7.29,13.69,22.09,32.49,44.89,59.29,75.69,94.09, 114.49,136.89,161.29,187.69,216.09,246.49,278.89,313.29,349.69};
532   const double targetCoords0[20*3]={0.,0.,0.,1.,0.,0.,0.,1.,0.,1.,1.,0.,0.,0.,1.,1.,0.,1.,0.,1.,1.,1.,1.,1.,0.,0.,2.,1.,0.,2.,0.,1.,2.,1.,1.,2.,0.,0.,3.,1.,0.,3.,0.,1.,3.,1.,1.,3.,0.,0.,4.,1.,0.,4.,0.,1.,4.,1.,1.,4.};
533   const int targetConn0[8*4]={1,0,2,3,5,4,6,7,5,4,6,7,9,8,10,11,9,8,10,11,13,12,14,15,13,12,14,15,17,16,18,19};
534   const double targetCoords1[28*3]={0.,0.,4.,1.,0.,4.,0.,1.,4.,1.,1.,4.,0.,0.,5.,1.,0.,5.,0.,1.,5.,1.,1.,5.,0.,0.,6.,1.,0.,6.,0.,1.,6.,1.,1.,6.,0.,0.,7.,1.,0.,7.,0.,1.,7.,1.,1.,7.,0.,0.,8.,1.,0.,8.,0.,1.,8.,1.,1.,8.,0.,0.,9.,1.,0.,9.,0.,1.,9.,1.,1.,9.,0.,0.,10.,1.,0.,10.,0.,1.,10.,1.,1.,10.};
535   const int targetConn1[8*6]={1,0,2,3,5,4,6,7,5,4,6,7,9,8,10,11,9,8,10,11,13,12,14,15,13,12,14,15,17,16,18,19,17,16,18,19,21,20,22,23,21,20,22,23,25,24,26,27};
536   //
537   int size;
538   int rank;
539   MPI_Comm_size(MPI_COMM_WORLD,&size);
540   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
541   //
542   if(size!=3)
543     return ;
544   int nproc_source = 1;
545   set<int> self_procs;
546   set<int> procs_source;
547   set<int> procs_target;
548   
549   for (int i=0; i<nproc_source; i++)
550     procs_source.insert(i);
551   for (int i=nproc_source; i<size; i++)
552     procs_target.insert(i);
553   self_procs.insert(rank);
554   //
555   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
556   ParaMEDMEM::ParaMESH *paramesh=0;
557   ParaMEDMEM::ParaFIELD* parafield=0;
558   //
559   ParaMEDMEM::CommInterface interface;
560   //
561   ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
562   ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
563   ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
564   //
565   MPI_Barrier(MPI_COMM_WORLD);
566   if(source_group->containsMyRank())
567     {
568       std::ostringstream stream; stream << "sourcemesh2D proc " << rank;
569       mesh=MEDCouplingUMesh::New(stream.str().c_str(),1);
570       mesh->allocateCells();
571       for(int i=0;i<18;i++)
572         mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,sourceConn+2*i);
573       mesh->finishInsertingCells();
574       DataArrayDouble *myCoords=DataArrayDouble::New();
575       myCoords->alloc(19,3);
576       std::copy(sourceCoords,sourceCoords+19*3,myCoords->getPointer());
577       mesh->setCoords(myCoords);
578       myCoords->decrRef();
579       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
580       ParaMEDMEM::ComponentTopology comptopo;
581       parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh,comptopo);
582       double *value=parafield->getField()->getArray()->getPointer();
583       std::copy(sourceVals,sourceVals+19,value);
584     }
585   else
586     {
587       if(rank==1)
588         {
589           std::ostringstream stream; stream << "targetmesh2D proc " << rank-nproc_source;
590           mesh=MEDCouplingUMesh::New(stream.str().c_str(),3);
591           mesh->allocateCells();
592           for(int i=0;i<4;i++)
593             mesh->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,targetConn0+8*i);
594           mesh->finishInsertingCells();
595           DataArrayDouble *myCoords=DataArrayDouble::New();
596           myCoords->alloc(20,3);
597           std::copy(targetCoords0,targetCoords0+20*3,myCoords->getPointer());
598           mesh->setCoords(myCoords);
599           myCoords->decrRef();
600           paramesh=new ParaMESH (mesh,*target_group,"target mesh");
601           ParaMEDMEM::ComponentTopology comptopo;
602           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
603         }
604       else if(rank==2)
605         {
606           std::ostringstream stream; stream << "targetmesh2D proc " << rank-nproc_source;
607           mesh=MEDCouplingUMesh::New(stream.str().c_str(),3);
608           mesh->allocateCells();
609           for(int i=0;i<6;i++)
610             mesh->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,targetConn1+8*i);
611           mesh->finishInsertingCells();
612           DataArrayDouble *myCoords=DataArrayDouble::New();
613           myCoords->alloc(28,3);
614           std::copy(targetCoords1,targetCoords1+28*3,myCoords->getPointer());
615           mesh->setCoords(myCoords);
616           myCoords->decrRef();
617           paramesh=new ParaMESH (mesh,*target_group,"target mesh");
618           ParaMEDMEM::ComponentTopology comptopo;
619           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
620         }
621     }
622   //test 1 - primaire -> secondaire
623   ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
624   dec.setIntersectionType(INTERP_KERNEL::PointLocator);
625   parafield->getField()->setNature(ConservativeVolumic);//very important
626   if (source_group->containsMyRank())
627     { 
628       dec.setMethod("P1");
629       dec.attachLocalField(parafield);
630       dec.synchronize();
631       dec.setForcedRenormalization(false);
632       dec.sendData();
633     }
634   else
635     {
636       dec.setMethod("P0");
637       dec.attachLocalField(parafield);
638       dec.synchronize();
639       dec.setForcedRenormalization(false);
640       dec.recvData();
641       const double *res(parafield->getField()->getArray()->getConstPointer());
642       if(rank==1)
643         {
644           const double expected0[4]={0.49,7.956666666666667,27.29,0.};
645           for(int i=0;i<4;i++)
646             CPPUNIT_ASSERT_DOUBLES_EQUAL(expected0[i],res[i],1e-13);
647         }
648       else
649         {
650           const double expected1[6]={59.95666666666667,94.09,0.,125.69,202.89,296.09};
651           for(int i=0;i<6;i++)
652             CPPUNIT_ASSERT_DOUBLES_EQUAL(expected1[i],res[i],1e-13);
653         }
654     }
655   MPI_Barrier(MPI_COMM_WORLD);
656   if (source_group->containsMyRank())
657     {
658       dec.recvData();
659       const double expected2[19]={0.49,7.956666666666667,7.956666666666667,7.956666666666667,27.29,27.29,59.95666666666667,59.95666666666667,59.95666666666667,94.09,125.69,125.69,202.89,202.89,202.89,202.89,296.09,296.09,0.};
660       const double *res(parafield->getField()->getArray()->getConstPointer());
661       for(int i=0;i<19;i++)
662         CPPUNIT_ASSERT_DOUBLES_EQUAL(expected2[i],res[i],1e-13);
663     }
664   else
665     {
666       dec.sendData();
667     }
668   delete parafield;
669   mesh->decrRef();
670   delete paramesh;
671   delete self_group;
672   delete target_group;
673   delete source_group;
674   //
675   MPI_Barrier(MPI_COMM_WORLD);
676 }