1 // Copyright (C) 2007-2020 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay (CEA/DEN)
21 #include "MEDLoaderTest.hxx"
22 #include "MEDLoader.hxx"
23 #include "MEDLoaderBase.hxx"
24 #include "MEDCouplingUMesh.hxx"
25 #include "MEDCouplingFieldDouble.hxx"
26 #include "MEDCouplingFieldInt.hxx"
27 #include "MEDCouplingFieldFloat.hxx"
28 #include "MEDCouplingMemArray.hxx"
29 #include "TestInterpKernelUtils.hxx" // getResourceFile()
34 using namespace MEDCoupling;
36 void MEDLoaderTest::testMesh1DRW()
38 MEDCouplingUMesh *mesh=build1DMesh_1();
39 mesh->checkConsistencyLight();
40 WriteUMesh("file1.med",mesh,true);
41 MEDCouplingUMesh *mesh_rw=ReadUMeshFromFile("file1.med",mesh->getName().c_str(),0);
42 CPPUNIT_ASSERT(mesh->isEqual(mesh_rw,1e-12));
47 void MEDLoaderTest::testMesh2DCurveRW()
49 MEDCouplingUMesh *mesh=build2DCurveMesh_1();
50 mesh->checkConsistencyLight();
51 WriteUMesh("file2.med",mesh,true);
52 MEDCouplingUMesh *mesh_rw=ReadUMeshFromFile("file2.med",mesh->getName().c_str(),0);
53 CPPUNIT_ASSERT(mesh->isEqual(mesh_rw,1e-12));
58 void MEDLoaderTest::testMesh2DRW()
60 MEDCouplingUMesh *mesh=build2DMesh_1();
61 mesh->checkConsistencyLight();
62 WriteUMesh("file3.med",mesh,true);
63 MEDCouplingUMesh *mesh_rw=ReadUMeshFromFile("file3.med",mesh->getName().c_str(),0);
64 CPPUNIT_ASSERT(mesh->isEqual(mesh_rw,1e-12));
69 void MEDLoaderTest::testMesh3DSurfRW()
71 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
72 mesh->checkConsistencyLight();
73 WriteUMesh("file4.med",mesh,true);
74 MEDCouplingUMesh *mesh_rw=ReadUMeshFromFile("file4.med",mesh->getName().c_str(),0);
75 CPPUNIT_ASSERT(mesh->isEqual(mesh_rw,1e-12));
80 void MEDLoaderTest::testMesh3DRW()
82 MEDCouplingUMesh *mesh=build3DMesh_1();
83 mesh->checkConsistencyLight();
84 WriteUMesh("file5.med",mesh,true);
85 MEDCouplingUMesh *mesh_rw=ReadUMeshFromFile("file5.med",mesh->getName().c_str(),0);
86 CPPUNIT_ASSERT(mesh->isEqual(mesh_rw,1e-12));
92 * Most basic test : one and only one MEDCoupling field in a new file.
94 void MEDLoaderTest::testFieldRW1()
96 MEDCouplingFieldDouble *f1=buildVecFieldOnCells_1();
97 WriteField("file6.med",f1,true);
98 MEDCouplingFieldDouble *f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell("file6.med",f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),0,1));
99 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
103 f1=buildVecFieldOnNodes_1();
104 WriteField("file7.med",f1,true);
105 f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode("file7.med",f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,3));
106 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
107 // testing kind message on error of field type.
108 CPPUNIT_ASSERT_THROW(ReadFieldCell("file7.med",f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,3),INTERP_KERNEL::Exception);
115 * Multi field writing in a same file.
117 void MEDLoaderTest::testFieldRW2()
119 const char fileName[]="file8.med";
120 static const double VAL1=12345.67890314;
121 static const double VAL2=-1111111111111.;
122 MEDCouplingFieldDouble *f1=buildVecFieldOnCells_1();
123 MEDCouplingFieldInt *f1_int=buildIntVecFieldOnCells_1();
124 MEDCouplingFieldFloat *f1_fl=buildFloatVecFieldOnCells_1();
125 WriteField(fileName,f1,true);
126 f1->setTime(10.,8,9);
127 f1_int->setTime(10.,8,9);
128 f1_fl->setTime(10.,8,9);
129 double *tmp=f1->getArray()->getPointer();
131 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
132 f1->setTime(10.14,18,19);
134 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
135 // Write int and float fields:
136 WriteFieldUsingAlreadyWrittenMesh(fileName,f1_int);
137 WriteFieldUsingAlreadyWrittenMesh(fileName,f1_fl);
138 //retrieving time steps...
139 MEDCouplingFieldDouble *f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),8,9));
140 f1->setTime(10.,8,9);
142 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
144 f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),0,1));
145 MEDCouplingFieldDouble *f3=buildVecFieldOnCells_1();
146 CPPUNIT_ASSERT(f3->isEqual(f2,1e-12,1e-12));
149 f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),18,19));
150 f1->setTime(10.14,18,19);
152 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
153 //test of throw on invalid (dt,it)
154 CPPUNIT_ASSERT_THROW(ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),28,19),INTERP_KERNEL::Exception);
157 // Reading Int and Float fields:
158 MEDCouplingFieldInt *f2_int=dynamic_cast<MEDCouplingFieldInt *>(ReadFieldCell(fileName,f1_int->getMesh()->getName().c_str(),0,f1_int->getName().c_str(),8,9));
159 CPPUNIT_ASSERT(f1_int->isEqual(f2_int,1e-12,0)); // exact equality for int values
161 MEDCouplingFieldFloat *f2_fl=dynamic_cast<MEDCouplingFieldFloat *>(ReadFieldCell(fileName,f1_fl->getMesh()->getName().c_str(),0,f1_fl->getName().c_str(),8,9));
162 CPPUNIT_ASSERT(f1_fl->isEqual(f2_fl,1e-12,1e-12));
165 f1=buildVecFieldOnNodes_1();
166 const char fileName2[]="file9.med";
167 WriteField(fileName2,f1,true);
168 f1->setTime(110.,108,109);
169 tmp=f1->getArray()->getPointer();
171 WriteFieldUsingAlreadyWrittenMesh(fileName2,f1);
172 f1->setTime(210.,208,209);
174 WriteFieldUsingAlreadyWrittenMesh(fileName2,f1);
175 f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName2,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),108,109));
176 f1->setTime(110.,108,109);
178 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
180 f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName2,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,3));
181 f3=buildVecFieldOnNodes_1();
182 CPPUNIT_ASSERT(f3->isEqual(f2,1e-12,1e-12));
185 f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName2,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),208,209));
186 f1->setTime(210.,208,209);
188 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
194 * Multi field in a same file, but this field has several timesteps
196 void MEDLoaderTest::testFieldRW3()
198 const char fileName[]="file11.med";
199 static const double VAL1=12345.67890314;
200 static const double VAL2=-1111111111111.;
201 const char name1[]="AField";
202 const char name3[]="AMesh1";
203 MEDCouplingFieldDouble *f1=buildVecFieldOnCells_1();
204 (const_cast<MEDCouplingMesh *>(f1->getMesh()))->setName(name3);
206 f1->setTime(10.,8,9);
207 double *tmp=f1->getArray()->getPointer();
209 WriteField(fileName,f1,true);
210 f1->setTime(10.14,18,19);
212 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
213 f1->setTime(10.55,28,29);
215 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
216 f1->setTime(10.66,38,39);
218 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
219 f1->setTime(10.77,48,49);
221 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
224 f1=buildVecFieldOnNodes_1();
226 (const_cast<MEDCouplingMesh *>(f1->getMesh()))->setName(name3);
227 f1->setTime(110.,8,9);
228 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
229 f1->setTime(110.,108,109);
230 tmp=f1->getArray()->getPointer();
232 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
233 f1->setTime(210.,208,209);
235 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
237 std::vector< std::pair<int,int> > it1=GetCellFieldIterations(fileName,name3,name1);
238 CPPUNIT_ASSERT_EQUAL(5,(int)it1.size());
239 CPPUNIT_ASSERT_EQUAL(8,it1[0].first); CPPUNIT_ASSERT_EQUAL(9,it1[0].second);
240 CPPUNIT_ASSERT_EQUAL(18,it1[1].first); CPPUNIT_ASSERT_EQUAL(19,it1[1].second);
241 CPPUNIT_ASSERT_EQUAL(28,it1[2].first); CPPUNIT_ASSERT_EQUAL(29,it1[2].second);
242 CPPUNIT_ASSERT_EQUAL(38,it1[3].first); CPPUNIT_ASSERT_EQUAL(39,it1[3].second);
243 CPPUNIT_ASSERT_EQUAL(48,it1[4].first); CPPUNIT_ASSERT_EQUAL(49,it1[4].second);
244 std::vector< std::pair<int,int> > it3=GetNodeFieldIterations(fileName,name3,name1);
245 CPPUNIT_ASSERT_EQUAL(3,(int)it3.size());
246 CPPUNIT_ASSERT_EQUAL(8,it3[0].first); CPPUNIT_ASSERT_EQUAL(9,it3[0].second);
247 CPPUNIT_ASSERT_EQUAL(108,it3[1].first); CPPUNIT_ASSERT_EQUAL(109,it3[1].second);
248 CPPUNIT_ASSERT_EQUAL(208,it3[2].first); CPPUNIT_ASSERT_EQUAL(209,it3[2].second);
252 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,name3,0,name1,8,9));
253 CPPUNIT_ASSERT_DOUBLES_EQUAL(VAL1,f1->getArray()->getConstPointer()[0],1e-13);
255 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,name3,0,name1,18,19));
256 CPPUNIT_ASSERT_DOUBLES_EQUAL(VAL2,f1->getArray()->getConstPointer()[0],1e-13);
258 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,name3,0,name1,28,29));
259 CPPUNIT_ASSERT_DOUBLES_EQUAL(3*VAL1,f1->getArray()->getConstPointer()[0],1e-13);
261 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,name3,0,name1,38,39));
262 CPPUNIT_ASSERT_DOUBLES_EQUAL(3*VAL2,f1->getArray()->getConstPointer()[0],1e-13);
264 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,name3,0,name1,48,49));
265 CPPUNIT_ASSERT_DOUBLES_EQUAL(4*VAL2,f1->getArray()->getConstPointer()[0],1e-13);
268 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName,name3,0,name1,8,9));
269 CPPUNIT_ASSERT_DOUBLES_EQUAL(71.,f1->getArray()->getConstPointer()[3],1e-13);
271 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName,name3,0,name1,108,109));
272 CPPUNIT_ASSERT_DOUBLES_EQUAL(VAL1,f1->getArray()->getConstPointer()[3],1e-13);
274 f1=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName,name3,0,name1,208,209));
275 CPPUNIT_ASSERT_DOUBLES_EQUAL(VAL2,f1->getArray()->getConstPointer()[3],1e-13);
279 void MEDLoaderTest::testMultiMeshRW1()
281 const char fileName[]="file10.med";
282 MEDCouplingUMesh *mesh1=build3DMesh_1();
283 const mcIdType part1[5]={1,2,4,13,15};
284 MEDCouplingUMesh *mesh2=(MEDCouplingUMesh *)mesh1->buildPartOfMySelf(part1,part1+5,true);
285 mesh2->setName("mesh2");
286 const mcIdType part2[4]={3,4,13,14};
287 MEDCouplingUMesh *mesh3=(MEDCouplingUMesh *)mesh1->buildPartOfMySelf(part2,part2+4,true);
288 mesh3->setName("mesh3");
289 MEDCouplingUMesh *mesh4=MEDCouplingUMesh::New();
290 mesh4->setName("mesh4");
291 mesh4->setMeshDimension(3);
292 mesh4->allocateCells(1);
293 mcIdType conn[4]={0,11,1,3};
294 mesh4->insertNextCell(INTERP_KERNEL::NORM_TETRA4,4,conn);
295 mesh4->finishInsertingCells();
296 mesh4->setCoords(mesh1->getCoords());
297 std::vector<const MEDCouplingUMesh *> meshes;
298 meshes.push_back(mesh1);
299 meshes.push_back(mesh2);
300 meshes.push_back(mesh3);
301 meshes.push_back(mesh4);
302 const char mnane[]="3DToto";
303 WriteUMeshesPartition(fileName,mnane,meshes,true);
305 MEDCouplingUMesh *mesh5=ReadUMeshFromFile(fileName,mnane);
306 mesh1->setName(mnane);
307 const mcIdType part3[18]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17};
308 MEDCouplingUMesh *mesh6=(MEDCouplingUMesh *)mesh5->buildPartOfMySelf(part3,part3+18,true);
309 mesh6->setName(mnane);
311 CPPUNIT_ASSERT(mesh6->isEqual(mesh1,1e-12));
313 std::vector<std::string> grps=GetMeshGroupsNames(fileName,mnane);
314 CPPUNIT_ASSERT_EQUAL(4,(int)grps.size());
315 CPPUNIT_ASSERT(std::find(grps.begin(),grps.end(),std::string("mesh2"))!=grps.end());
316 CPPUNIT_ASSERT(std::find(grps.begin(),grps.end(),std::string("mesh3"))!=grps.end());
317 CPPUNIT_ASSERT(std::find(grps.begin(),grps.end(),std::string("mesh4"))!=grps.end());
318 CPPUNIT_ASSERT(std::find(grps.begin(),grps.end(),std::string("3DMesh_1"))!=grps.end());
320 std::vector<std::string> vec;
321 vec.push_back(std::string("mesh2"));
322 MEDCouplingUMesh *mesh2_2=ReadUMeshFromGroups(fileName,mnane,0,vec);
323 CPPUNIT_ASSERT(mesh2_2->isEqual(mesh2,1e-12));
325 vec.clear(); vec.push_back(std::string("mesh3"));
326 MEDCouplingUMesh *mesh3_2=ReadUMeshFromGroups(fileName,mnane,0,vec);
327 CPPUNIT_ASSERT(mesh3_2->isEqual(mesh3,1e-12));
329 vec.clear(); vec.push_back(std::string("mesh4"));
330 MEDCouplingUMesh *mesh4_2=ReadUMeshFromGroups(fileName,mnane,0,vec);
331 CPPUNIT_ASSERT(mesh4_2->isEqual(mesh4,1e-12));
333 vec.clear(); vec.push_back(std::string("3DMesh_1"));
334 MEDCouplingUMesh *mesh1_2=ReadUMeshFromGroups(fileName,mnane,0,vec);
335 mesh1->setName("3DMesh_1");
336 CPPUNIT_ASSERT(mesh1_2->isEqual(mesh1,1e-12));
339 vec.clear(); vec.push_back(std::string("Family_-3")); vec.push_back(std::string("Family_-5"));
340 mesh2_2=ReadUMeshFromFamilies(fileName,mnane,0,vec);
341 mesh2_2->setName("mesh2");
342 CPPUNIT_ASSERT(mesh2_2->isEqual(mesh2,1e-12));
345 std::vector<std::string> ret(GetMeshFamiliesNamesOnGroup(fileName,"3DToto","3DMesh_1"));
346 std::set<std::string> s(ret.begin(),ret.end());
347 std::set<std::string> ref_s;
348 ref_s.insert("Family_-2");
349 ref_s.insert("Family_-3");
350 ref_s.insert("Family_-4");
351 ref_s.insert("Family_-5");
352 CPPUNIT_ASSERT_EQUAL(4,(int)ret.size());
353 CPPUNIT_ASSERT(s==ref_s);
355 std::vector<std::string> ret1=GetMeshGroupsNamesOnFamily(fileName,"3DToto","Family_-3");
356 CPPUNIT_ASSERT_EQUAL(2,(int)ret1.size());
357 CPPUNIT_ASSERT(ret1[0]=="3DMesh_1");
358 CPPUNIT_ASSERT(ret1[1]=="mesh2");
366 void MEDLoaderTest::testFieldProfilRW1()
368 const char fileName[]="file12.med";
369 MEDCouplingUMesh *mesh1=build3DMesh_1();
371 mcIdType newNbOfNodes;
372 DataArrayIdType *da=mesh1->mergeNodes(1e-12,b,newNbOfNodes);
374 WriteUMesh(fileName,mesh1,true);
375 const mcIdType part1[5]={1,2,4,13,15};
376 MEDCouplingUMesh *mesh2=(MEDCouplingUMesh *)mesh1->buildPartOfMySelf(part1,part1+5,true);
377 mesh2->setName(mesh1->getName().c_str());//<- important for the test
379 mcIdType nbOfCells=mesh2->getNumberOfCells();
380 CPPUNIT_ASSERT_EQUAL(ToIdType(5),nbOfCells);
381 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
382 f1->setName("VectorFieldOnCells");
384 DataArrayDouble *array=DataArrayDouble::New();
385 array->alloc(nbOfCells,2);
388 double *tmp=array->getPointer();
389 const double arr1[10]={71.,171.,10.,110.,20.,120.,30.,130.,40.,140.};
390 std::copy(arr1,arr1+10,tmp);
391 f1->setTime(3.14,2,7);
392 f1->checkConsistencyLight();
394 WriteField(fileName,f1,false);//<- false important for the test
396 MEDCouplingFieldDouble *f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,7));
397 std::vector<MEDCoupling::TypeOfField> types=GetTypesOfField(fileName,f1->getMesh()->getName().c_str(),f1->getName().c_str());
398 CPPUNIT_ASSERT_EQUAL(1,(int)types.size());
399 CPPUNIT_ASSERT(types[0]==ON_CELLS);
400 f2->checkConsistencyLight();
401 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
410 * Test MED file profiles.
412 void MEDLoaderTest::testFieldNodeProfilRW1()
414 const char fileName[]="file19.med";
415 const char fileName2[]="file20.med";
416 MEDCouplingUMesh *m=build2DMesh_1();
417 mcIdType nbOfNodes=m->getNumberOfNodes();
418 WriteUMesh(fileName,m,true);
419 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_NODES,ONE_TIME);
420 f1->setName("VFieldOnNodes");
422 DataArrayDouble *array=DataArrayDouble::New();
423 const double arr1[24]={1.,101.,2.,102.,3.,103.,4.,104.,5.,105.,6.,106.,7.,107.,8.,108.,9.,109.,10.,110.,11.,111.,12.,112.};
424 array->alloc(nbOfNodes,2);
425 std::copy(arr1,arr1+24,array->getPointer());
427 array->setInfoOnComponent(0,"tyty [mm]");
428 array->setInfoOnComponent(1,"uiop [MW]");
430 f1->setTime(3.14,2,7);
431 f1->checkConsistencyLight();
432 const mcIdType arr2[2]={1,4};//node ids are 2,4,5,3,6,7
433 MEDCouplingFieldDouble *f2=f1->buildSubPart(arr2,arr2+2);
434 (const_cast<MEDCouplingMesh *>(f2->getMesh()))->setName(f1->getMesh()->getName().c_str());
435 WriteField(fileName,f2,false);//<- false important for the test
437 MEDCouplingFieldDouble *f3=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName,f2->getMesh()->getName().c_str(),0,f2->getName().c_str(),2,7));
438 f3->checkConsistencyLight();
439 CPPUNIT_ASSERT(f3->isEqual(f2,1e-12,1e-12));
442 const mcIdType arr3[6]={1,3,0,5,2,4};
443 f2->renumberNodes(arr3);
444 WriteUMesh(fileName2,m,true);
445 WriteField(fileName2,f2,false);//<- false important for the test
446 f3=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName2,f2->getMesh()->getName().c_str(),0,f2->getName().c_str(),2,7));
447 f3->checkConsistencyLight();
448 CPPUNIT_ASSERT(f3->isEqual(f2,1e-12,1e-12));
456 void MEDLoaderTest::testFieldNodeProfilRW2()
458 const char fileName[]="file23.med";
459 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
460 WriteUMesh(fileName,mesh,true);
462 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_NODES,ONE_TIME);
463 f1->setName("FieldMix");
465 const double arr2[24]={
466 1071.,1171.,1010.,1110.,1020.,1120.,1030.,1130.,1040.,1140.,1050.,1150.,
467 1060.,1160.,1070.,1170.,1080.,1180.,1090.,1190.,1091.,1191.,1092.,1192.
469 DataArrayDouble *array=DataArrayDouble::New();
472 array->setInfoOnComponent(0,"plkj [mm]");
473 array->setInfoOnComponent(1,"pqqqss [mm]");
475 double *tmp=array->getPointer();
476 std::copy(arr2,arr2+24,tmp);
477 f1->setTime(3.17,2,7);
479 const mcIdType renumArr[12]={3,7,2,1,5,11,10,0,9,6,8,4};
480 f1->renumberNodes(renumArr);
481 f1->checkConsistencyLight();
482 WriteField(fileName,f1,false);//<- false important for the test
483 MEDCouplingFieldDouble *f2=dynamic_cast<MEDCouplingFieldDouble *>(ReadFieldNode(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,7));
484 CPPUNIT_ASSERT(f2->isEqual(f1,1e-12,1e-12));
491 void MEDLoaderTest::testFieldGaussRW1()
493 const char fileName[]="file13.med";
494 MEDCouplingFieldDouble *f1=buildVecFieldOnGauss_1();
495 WriteField(fileName,f1,true);
496 MCAuto<MEDCouplingField> f2Tmp(ReadField(ON_GAUSS_PT,fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),1,5));
497 MCAuto<MEDCouplingFieldDouble> f2(MEDCoupling::DynamicCast<MEDCouplingField,MEDCouplingFieldDouble>(f2Tmp));
498 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
502 void MEDLoaderTest::testFieldGaussNERW1()
504 const char fileName[]="file14.med";
505 MEDCouplingFieldDouble *f1=buildVecFieldOnGaussNE_1();
506 WriteField(fileName,f1,true);
507 std::vector<MEDCoupling::TypeOfField> tof(GetTypesOfField(fileName,"2DMesh_2","MyFieldOnGaussNE"));
508 CPPUNIT_ASSERT_EQUAL(1,(int)tof.size());
509 CPPUNIT_ASSERT(ON_GAUSS_NE==tof[0]);
510 MCAuto<MEDCouplingField> f2Tmp(ReadField(ON_GAUSS_NE,fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),1,5));
511 MCAuto<MEDCouplingFieldDouble> f2(MEDCoupling::DynamicCast<MEDCouplingField,MEDCouplingFieldDouble>(f2Tmp));
512 CPPUNIT_ASSERT(f1->isEqual(f2,1e-12,1e-12));
516 void MEDLoaderTest::testLittleStrings1()
518 std::string s("azeeeerrrtty");
519 MEDLoaderBase::zipEqualConsChar(s,3);
520 CPPUNIT_ASSERT(s=="azertty");
523 void MEDLoaderTest::testSplitIntoNameAndUnit1()
525 std::string s(" []");
527 MEDLoaderBase::splitIntoNameAndUnit(s,c,u);
528 CPPUNIT_ASSERT(c.empty());
529 CPPUNIT_ASSERT(u.empty());
531 MEDLoaderBase::strip(s);
532 CPPUNIT_ASSERT(s=="lmmm kki jjj");
534 MEDLoaderBase::strip(s);
535 CPPUNIT_ASSERT(s.empty());
537 MEDLoaderBase::strip(s);
538 CPPUNIT_ASSERT(s.empty());
540 MEDLoaderBase::strip(s);
541 CPPUNIT_ASSERT(s.empty());
543 MEDLoaderBase::strip(s);
544 CPPUNIT_ASSERT(s=="pp");
547 void MEDLoaderTest::testMesh3DSurfShuffleRW()
549 const char fileName[]="file15.med";
550 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
551 const mcIdType renumber1[6]={2,5,1,0,3,4};
552 mesh->renumberCells(renumber1,false);
553 mesh->checkConsistencyLight();
554 WriteUMesh(fileName,mesh,true);
555 MEDCouplingUMesh *mesh_rw=ReadUMeshFromFile(fileName,mesh->getName().c_str(),0);
556 CPPUNIT_ASSERT(mesh->isEqual(mesh_rw,1e-12));
561 void MEDLoaderTest::testFieldShuffleRW1()
563 const char fileName[]="file16.med";
564 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
565 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
566 f1->setName("FieldOnCellsShuffle");
568 DataArrayDouble *array=DataArrayDouble::New();
572 double *tmp=array->getPointer();
573 const double arr1[12]={71.,171.,10.,110.,20.,120.,30.,130.,40.,140.,50.,150.};
574 std::copy(arr1,arr1+12,tmp);
575 f1->setTime(3.14,2,7);
576 f1->checkConsistencyLight();
578 const mcIdType renumber1[6]={2,1,5,0,3,4};
579 f1->renumberCells(renumber1,false);
580 WriteField(fileName,f1,true);
581 MEDCouplingFieldDouble *f2=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldCell(fileName,mesh->getName().c_str(),0,f1->getName().c_str(),2,7));
582 CPPUNIT_ASSERT(f2->isEqual(f1,1e-12,1e-12));
590 * Shuffle de cells but no profile. Like pointe.med
592 void MEDLoaderTest::testMultiFieldShuffleRW1()
594 const char fileName[]="file17.med";
595 MEDCouplingUMesh *m=build3DMesh_2();
596 CPPUNIT_ASSERT_EQUAL(20,(int)m->getNumberOfCells());
597 CPPUNIT_ASSERT_EQUAL(ToIdType(45),m->getNumberOfNodes());
598 const mcIdType polys[3]={1,4,6};
599 std::vector<mcIdType> poly2(polys,polys+3);
600 m->convertToPolyTypes(&poly2[0],&poly2[0]+poly2.size());
601 const mcIdType renum[20]={1,3,2,8,9,12,13,16,19,0,4,7,5,15,14,17,10,18,6,11};
602 m->renumberCells(renum,false);
603 m->orientCorrectlyPolyhedrons();
605 WriteUMesh(fileName,m,true);
606 MEDCouplingFieldDouble *f1Tmp=m->getMeasureField(false);
607 MEDCouplingFieldDouble *f1=f1Tmp->buildNewTimeReprFromThis(ONE_TIME,false);
610 MEDCouplingFieldDouble *f_1=f1->cloneWithMesh(true);
611 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
612 f1->applyFunc("2*x");
613 f1->setTime(0.01,3,4);
614 MEDCouplingFieldDouble *f_2=f1->cloneWithMesh(true);
615 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
616 f1->applyFunc("2*x/3");
617 f1->setTime(0.02,5,6);
618 MEDCouplingFieldDouble *f_3=f1->cloneWithMesh(true);
619 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
622 std::vector<std::pair<int,int> > its;
623 its.push_back(std::pair<int,int>(1,2));
624 its.push_back(std::pair<int,int>(3,4));
625 its.push_back(std::pair<int,int>(5,6));
626 std::vector<MEDCouplingFieldDouble *> fs=ReadFieldsOnSameMesh(ON_CELLS,fileName,f_1->getMesh()->getName().c_str(),0,f_1->getName().c_str(),its);
627 CPPUNIT_ASSERT_EQUAL(3,(int)fs.size());
628 const MEDCouplingMesh *mm=fs[0]->getMesh();
629 CPPUNIT_ASSERT(fs[0]->isEqual(f_1,1e-12,1e-12));
630 CPPUNIT_ASSERT(fs[1]->isEqual(f_2,1e-12,1e-12));
631 CPPUNIT_ASSERT(fs[2]->isEqual(f_3,1e-12,1e-12));
632 CPPUNIT_ASSERT(mm==fs[1]->getMesh());// <- important for the test
633 CPPUNIT_ASSERT(mm==fs[2]->getMesh());// <- important for the test
634 for(std::vector<MEDCouplingFieldDouble *>::iterator iter=fs.begin();iter!=fs.end();iter++)
644 void MEDLoaderTest::testWriteUMeshesRW1()
646 const char fileName[]="file18.med";
647 MEDCouplingUMesh *m3d=build3DMesh_2();
648 const double pt[3]={0.,0.,-0.3};
649 const double vec[3]={0.,0.,1.};
650 std::vector<mcIdType> nodes;
651 m3d->findNodesOnPlane(pt,vec,1e-12,nodes);
652 MEDCouplingUMesh *m2d=(MEDCouplingUMesh *)m3d->buildFacePartOfMySelfNode(&nodes[0],&nodes[0]+nodes.size(),true);
653 const mcIdType renumber[5]={1,2,0,4,3};
654 m2d->renumberCells(renumber,false);
655 m2d->setName("ExampleOfMultiDimW");
656 std::vector<const MEDCouplingUMesh *> meshes;
657 meshes.push_back(m2d);
658 meshes.push_back(m3d);
659 WriteUMeshes(fileName,meshes,true);
660 MEDCouplingUMesh *m3d_bis=ReadUMeshFromFile(fileName,m2d->getName().c_str(),0);
661 CPPUNIT_ASSERT(!m3d_bis->isEqual(m3d,1e-12));
662 m3d_bis->setName(m3d->getName().c_str());
663 CPPUNIT_ASSERT(m3d_bis->isEqual(m3d,1e-12));
664 MEDCouplingUMesh *m2d_bis=ReadUMeshFromFile(fileName,m2d->getName().c_str(),-1);//-1 for faces
665 CPPUNIT_ASSERT(m2d_bis->isEqual(m2d,1e-12));
666 // Creation of a field on faces.
667 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
668 f1->setName("FieldOnFacesShuffle");
670 DataArrayDouble *array=DataArrayDouble::New();
671 array->alloc(m2d->getNumberOfCells(),2);
672 array->setInfoOnComponent(0,"plkj [mm]");
673 array->setInfoOnComponent(1,"pqqqss [mm]");
676 double *tmp=array->getPointer();
677 const double arr1[10]={71.,171.,10.,110.,20.,120.,30.,130.,40.,140.};
678 std::copy(arr1,arr1+10,tmp);
679 f1->setTime(3.14,2,7);
680 f1->checkConsistencyLight();
681 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
682 MEDCouplingFieldDouble *f2=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),-1,f1->getName().c_str(),2,7));
683 CPPUNIT_ASSERT(f2->isEqual(f1,1e-12,1e-12));
693 void MEDLoaderTest::testMixCellAndNodesFieldRW1()
695 const char fileName[]="file21.med";
696 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
697 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
698 f1->setName("FieldMix");
700 DataArrayDouble *array=DataArrayDouble::New();
703 array->setInfoOnComponent(0,"plkj [mm]");
704 array->setInfoOnComponent(1,"pqqqss [mm]");
706 double *tmp=array->getPointer();
707 const double arr1[12]={71.,171.,10.,110.,20.,120.,30.,130.,40.,140.,50.,150.};
708 std::copy(arr1,arr1+12,tmp);
709 f1->setTime(3.14,2,7);
710 f1->checkConsistencyLight();
712 MEDCouplingFieldDouble *f2=MEDCouplingFieldDouble::New(ON_NODES,ONE_TIME);
713 f2->setName("FieldMix");
715 array=DataArrayDouble::New();
718 array->setInfoOnComponent(0,"plkj [mm]");
719 array->setInfoOnComponent(1,"pqqqss [mm]");
721 tmp=array->getPointer();
722 const double arr2[24]={
723 1071.,1171.,1010.,1110.,1020.,1120.,1030.,1130.,1040.,1140.,1050.,1150.,
724 1060.,1160.,1070.,1170.,1080.,1180.,1090.,1190.,1091.,1191.,1092.,1192.
726 std::copy(arr2,arr2+24,tmp);
727 f2->setTime(3.14,2,7);
728 f2->checkConsistencyLight();
730 WriteField(fileName,f1,true);
731 std::vector<MEDCoupling::TypeOfField> ts=GetTypesOfField(fileName,f1->getMesh()->getName().c_str(),f1->getName().c_str());
732 CPPUNIT_ASSERT_EQUAL(1,(int)ts.size());
733 CPPUNIT_ASSERT_EQUAL(ON_CELLS,ts[0]);
734 std::vector<std::string> fs=GetAllFieldNamesOnMesh(fileName,f1->getMesh()->getName().c_str());
735 CPPUNIT_ASSERT_EQUAL(1,(int)fs.size());
736 CPPUNIT_ASSERT(fs[0]=="FieldMix");
737 WriteFieldUsingAlreadyWrittenMesh(fileName,f2);
738 fs=GetAllFieldNamesOnMesh(fileName,f1->getMesh()->getName().c_str());
739 CPPUNIT_ASSERT_EQUAL(1,(int)fs.size());
740 CPPUNIT_ASSERT(fs[0]=="FieldMix");
742 ts=GetTypesOfField(fileName,f1->getMesh()->getName().c_str(),f1->getName().c_str());
743 CPPUNIT_ASSERT_EQUAL(2,(int)ts.size());
744 CPPUNIT_ASSERT_EQUAL(ON_NODES,ts[0]);
745 CPPUNIT_ASSERT_EQUAL(ON_CELLS,ts[1]);
747 MEDCouplingFieldDouble *f3=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldNode(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,7));
748 CPPUNIT_ASSERT(f3->isEqual(f2,1e-12,1e-12));
750 f3=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldCell(fileName,f1->getMesh()->getName().c_str(),0,f1->getName().c_str(),2,7));
751 CPPUNIT_ASSERT(f3->isEqual(f1,1e-12,1e-12));
759 void MEDLoaderTest::testGetAllFieldNamesRW1()
761 const char fileName[]="file22.med";
762 MEDCouplingUMesh *mesh=build2DMesh_2();
763 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_NODES,ONE_TIME);
764 f1->setName("Field1");
765 f1->setTime(3.44,5,6);
767 f1->fillFromAnalytic(2,"x+y");
768 WriteField(fileName,f1,true);
769 f1->setTime(1002.3,7,8);
770 f1->fillFromAnalytic(2,"x+77.*y");
771 WriteFieldUsingAlreadyWrittenMesh(fileName,f1);
772 f1->setName("Field2");
773 WriteField(fileName,f1,false);
774 f1->setName("Field3");
775 mesh->setName("2DMesh_2Bis");
776 WriteField(fileName,f1,false);
778 f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
779 f1->setName("Field8");
780 f1->setTime(8.99,7,9);
782 f1->fillFromAnalytic(3,"3*x+y");
783 WriteField(fileName,f1,false);
785 std::vector<std::string> fs=GetAllFieldNames(fileName);
786 CPPUNIT_ASSERT_EQUAL(4,(int)fs.size());
787 CPPUNIT_ASSERT(fs[0]=="Field1");
788 CPPUNIT_ASSERT(fs[1]=="Field2");
789 CPPUNIT_ASSERT(fs[2]=="Field3");
790 CPPUNIT_ASSERT(fs[3]=="Field8");
795 void MEDLoaderTest::testMEDLoaderRead1()
798 using namespace INTERP_KERNEL;
800 string fileName= INTERP_TEST::getResourceFile("pointe.med", 3);
801 vector<string> meshNames=GetMeshNames(fileName.c_str());
802 CPPUNIT_ASSERT_EQUAL(1,(int)meshNames.size());
803 MEDCouplingUMesh *mesh=ReadUMeshFromFile(fileName.c_str(),meshNames[0].c_str(),0);
804 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
805 CPPUNIT_ASSERT_EQUAL(3,mesh->getMeshDimension());
806 CPPUNIT_ASSERT_EQUAL(16,(int)mesh->getNumberOfCells());
807 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
808 CPPUNIT_ASSERT_EQUAL(3,(int)mesh->getAllGeoTypes().size());
809 for(int i=0;i<12;i++)
810 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,mesh->getTypeOfCell(i));
811 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,mesh->getTypeOfCell(12));
812 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,mesh->getTypeOfCell(13));
813 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,mesh->getTypeOfCell(14));
814 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,mesh->getTypeOfCell(15));
815 CPPUNIT_ASSERT_EQUAL(ToIdType(90),mesh->getNodalConnectivity()->getNbOfElems());
816 CPPUNIT_ASSERT_EQUAL(701,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+90,0));
817 CPPUNIT_ASSERT_EQUAL(705,std::accumulate(mesh->getNodalConnectivityIndex()->getPointer(),mesh->getNodalConnectivityIndex()->getPointer()+17,0));
818 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+57,0),1e-12);
821 vector<string> families=GetMeshFamiliesNames(fileName.c_str(),meshNames[0].c_str());
822 CPPUNIT_ASSERT_EQUAL(8,(int)families.size());
823 CPPUNIT_ASSERT(families[2]=="FAMILLE_ELEMENT_3");
825 vector<string> families2;
826 families2.push_back(families[2]);
827 mesh=ReadUMeshFromFamilies(fileName.c_str(),meshNames[0].c_str(),0,families2);
828 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
829 CPPUNIT_ASSERT_EQUAL(3,mesh->getMeshDimension());
830 CPPUNIT_ASSERT_EQUAL(2,(int)mesh->getNumberOfCells());
831 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
832 CPPUNIT_ASSERT_EQUAL(2,(int)mesh->getAllGeoTypes().size());
833 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,mesh->getTypeOfCell(0));
834 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,mesh->getTypeOfCell(1));
835 CPPUNIT_ASSERT_EQUAL(ToIdType(11),mesh->getNodalConnectivity()->getNbOfElems());
836 CPPUNIT_ASSERT_EQUAL(132,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+11,0));
837 CPPUNIT_ASSERT_EQUAL(16,std::accumulate(mesh->getNodalConnectivityIndex()->getPointer(),mesh->getNodalConnectivityIndex()->getPointer()+3,0));
838 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+57,0),1e-12);
841 vector<string> groups=GetMeshGroupsNames(fileName.c_str(),meshNames[0].c_str());
842 CPPUNIT_ASSERT_EQUAL(5,(int)groups.size());
843 CPPUNIT_ASSERT(groups[0]=="groupe1");
844 CPPUNIT_ASSERT(groups[1]=="groupe2");
845 CPPUNIT_ASSERT(groups[2]=="groupe3");
846 CPPUNIT_ASSERT(groups[3]=="groupe4");
847 CPPUNIT_ASSERT(groups[4]=="groupe5");
848 vector<string> groups2;
849 groups2.push_back(groups[0]);
850 mesh=ReadUMeshFromGroups(fileName.c_str(),meshNames[0].c_str(),0,groups2);
851 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
852 CPPUNIT_ASSERT_EQUAL(3,mesh->getMeshDimension());
853 CPPUNIT_ASSERT_EQUAL(7,(int)mesh->getNumberOfCells());
854 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
855 CPPUNIT_ASSERT_EQUAL(2,(int)mesh->getAllGeoTypes().size());
857 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,mesh->getTypeOfCell(i));
858 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,mesh->getTypeOfCell(6));
859 CPPUNIT_ASSERT_EQUAL(ToIdType(36),mesh->getNodalConnectivity()->getNbOfElems());
860 CPPUNIT_ASSERT_EQUAL(254,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+36,0));
861 CPPUNIT_ASSERT_EQUAL(141,std::accumulate(mesh->getNodalConnectivityIndex()->getPointer(),mesh->getNodalConnectivityIndex()->getPointer()+8,0));
862 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+57,0),1e-12);
865 std::vector<std::string> fieldsName=GetCellFieldNamesOnMesh(fileName.c_str(),meshNames[0].c_str());
866 CPPUNIT_ASSERT_EQUAL(2,(int)fieldsName.size());
867 CPPUNIT_ASSERT(fieldsName[0]=="fieldcelldoublescalar");
868 CPPUNIT_ASSERT(fieldsName[1]=="fieldcelldoublevector");
869 std::vector<std::pair<int,int> > its0=GetCellFieldIterations(fileName.c_str(),meshNames[0].c_str(),fieldsName[0].c_str());
870 CPPUNIT_ASSERT_EQUAL(1,(int)its0.size());
871 CPPUNIT_ASSERT_EQUAL(-1,its0[0].first);
872 CPPUNIT_ASSERT_EQUAL(-1,its0[0].second);
873 std::vector<std::pair<int,int> > its1=GetCellFieldIterations(fileName.c_str(),meshNames[0].c_str(),fieldsName[1].c_str());
874 CPPUNIT_ASSERT_EQUAL(1,(int)its1.size());
875 CPPUNIT_ASSERT_EQUAL(-1,its1[0].first);
876 CPPUNIT_ASSERT_EQUAL(-1,its1[0].second);
878 MEDCouplingFieldDouble *field0=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldCell(fileName.c_str(),meshNames[0].c_str(),0,fieldsName[0].c_str(),its0[0].first,its0[0].second));
879 field0->checkConsistencyLight();
880 CPPUNIT_ASSERT(field0->getName()==fieldsName[0]);
881 CPPUNIT_ASSERT_EQUAL(1,(int)field0->getNumberOfComponents());
882 CPPUNIT_ASSERT_EQUAL(16,(int)field0->getNumberOfTuples());
883 const double expectedValues[16]={1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,2.,3.,3.,2.};
884 double diffValue[16];
885 std::transform(field0->getArray()->getPointer(),field0->getArray()->getPointer()+16,expectedValues,diffValue,std::minus<double>());
886 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue,diffValue+16),1e-12);
887 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue,diffValue+16),1e-12);
888 const MEDCouplingUMesh *constMesh=dynamic_cast<const MEDCouplingUMesh *>(field0->getMesh());
889 CPPUNIT_ASSERT(constMesh);
890 CPPUNIT_ASSERT_EQUAL(3,constMesh->getSpaceDimension());
891 CPPUNIT_ASSERT_EQUAL(3,constMesh->getMeshDimension());
892 CPPUNIT_ASSERT_EQUAL(16,(int)constMesh->getNumberOfCells());
893 CPPUNIT_ASSERT_EQUAL(ToIdType(19),constMesh->getNumberOfNodes());
894 CPPUNIT_ASSERT_EQUAL(3,(int)constMesh->getAllGeoTypes().size());
895 for(int i=0;i<12;i++)
896 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,constMesh->getTypeOfCell(i));
897 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(12));
898 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(13));
899 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(14));
900 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(15));
901 CPPUNIT_ASSERT_EQUAL(ToIdType(90),constMesh->getNodalConnectivity()->getNbOfElems());
902 CPPUNIT_ASSERT_EQUAL(701,std::accumulate(constMesh->getNodalConnectivity()->getConstPointer(),constMesh->getNodalConnectivity()->getConstPointer()+90,0));
903 CPPUNIT_ASSERT_EQUAL(705,std::accumulate(constMesh->getNodalConnectivityIndex()->getConstPointer(),constMesh->getNodalConnectivityIndex()->getConstPointer()+17,0));
904 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(constMesh->getCoords()->getConstPointer(),constMesh->getCoords()->getConstPointer()+57,0),1e-12);
907 MEDCouplingFieldDouble *field1=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldCell(fileName.c_str(),meshNames[0].c_str(),0,fieldsName[1].c_str(),its1[0].first,its1[0].second));
908 field1->checkConsistencyLight();
909 CPPUNIT_ASSERT(field1->getName()==fieldsName[1]);
910 CPPUNIT_ASSERT_EQUAL(3,(int)field1->getNumberOfComponents());
911 CPPUNIT_ASSERT_EQUAL(16,(int)field1->getNumberOfTuples());
912 const double expectedValues2[48]={1.,0.,1.,1.,0.,1.,1.,0.,1.,2.,1.,0.,2.,1.,0.,2.,1.,0.,3.,0.,1.,3.,0.,1.,3.,0.,1.,4.,1.,0.,4.,1.,0.,4.,1.,0.,5.,0.,0.,6.,1.,1.,6.,0.,0.,5.,1.,1.};
913 double diffValue2[48];
914 std::transform(field1->getArray()->getPointer(),field1->getArray()->getPointer()+48,expectedValues2,diffValue2,std::minus<double>());
915 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue2,diffValue2+48),1e-12);
916 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue2,diffValue2+48),1e-12);
917 constMesh=dynamic_cast<const MEDCouplingUMesh *>(field1->getMesh());
918 CPPUNIT_ASSERT(constMesh);
919 CPPUNIT_ASSERT_EQUAL(3,constMesh->getSpaceDimension());
920 CPPUNIT_ASSERT_EQUAL(3,constMesh->getMeshDimension());
921 CPPUNIT_ASSERT_EQUAL(16,(int)constMesh->getNumberOfCells());
922 CPPUNIT_ASSERT_EQUAL(ToIdType(19),constMesh->getNumberOfNodes());
923 CPPUNIT_ASSERT_EQUAL(3,(int)constMesh->getAllGeoTypes().size());
924 for(int i=0;i<12;i++)
925 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,constMesh->getTypeOfCell(i));
926 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(12));
927 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(13));
928 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(14));
929 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(15));
930 CPPUNIT_ASSERT_EQUAL(ToIdType(90),constMesh->getNodalConnectivity()->getNbOfElems());
931 CPPUNIT_ASSERT_EQUAL(701,std::accumulate(constMesh->getNodalConnectivity()->getConstPointer(),constMesh->getNodalConnectivity()->getConstPointer()+90,0));
932 CPPUNIT_ASSERT_EQUAL(705,std::accumulate(constMesh->getNodalConnectivityIndex()->getConstPointer(),constMesh->getNodalConnectivityIndex()->getConstPointer()+17,0));
933 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(constMesh->getCoords()->getConstPointer(),constMesh->getCoords()->getConstPointer()+57,0),1e-12);
936 std::vector<std::string> fieldsNameNode=GetNodeFieldNamesOnMesh(fileName.c_str(),meshNames[0].c_str());
937 CPPUNIT_ASSERT_EQUAL(2,(int)fieldsNameNode.size());
938 CPPUNIT_ASSERT(fieldsNameNode[0]=="fieldnodedouble");
939 CPPUNIT_ASSERT(fieldsNameNode[1]=="fieldnodeint");
940 std::vector<std::pair<int,int> > its0Node=GetNodeFieldIterations(fileName.c_str(),meshNames[0].c_str(),fieldsNameNode[0].c_str());
941 CPPUNIT_ASSERT_EQUAL(3,(int)its0Node.size());
942 CPPUNIT_ASSERT_EQUAL(-1,its0Node[0].first);
943 CPPUNIT_ASSERT_EQUAL(-1,its0Node[0].second);
944 CPPUNIT_ASSERT_EQUAL(1,its0Node[1].first);
945 CPPUNIT_ASSERT_EQUAL(-1,its0Node[1].second);
946 CPPUNIT_ASSERT_EQUAL(2,its0Node[2].first);
947 CPPUNIT_ASSERT_EQUAL(-1,its0Node[2].second);
948 MEDCouplingFieldDouble *field0Nodes=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldNode(fileName.c_str(),meshNames[0].c_str(),0,fieldsNameNode[0].c_str(),its0Node[0].first,its0Node[0].second));
949 field0Nodes->checkConsistencyLight();
950 CPPUNIT_ASSERT(field0Nodes->getName()==fieldsNameNode[0]);
951 CPPUNIT_ASSERT_EQUAL(1,(int)field0Nodes->getNumberOfComponents());
952 CPPUNIT_ASSERT_EQUAL(19,(int)field0Nodes->getNumberOfTuples());
953 const double expectedValues3[19]={1.,1.,1.,2.,2.,2.,3.,3.,3.,4.,4.,4.,5.,5.,5.,6.,6.,6.,7.};
954 double diffValue3[19];
955 std::transform(field0Nodes->getArray()->getPointer(),field0Nodes->getArray()->getPointer()+19,expectedValues3,diffValue3,std::minus<double>());
956 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue3,diffValue3+19),1e-12);
957 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue3,diffValue3+19),1e-12);
958 constMesh=dynamic_cast<const MEDCouplingUMesh *>(field0Nodes->getMesh());
959 CPPUNIT_ASSERT(constMesh);
960 field0Nodes->decrRef();
962 field0Nodes=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldNode(fileName.c_str(),meshNames[0].c_str(),0,fieldsNameNode[0].c_str(),its0Node[2].first,its0Node[2].second));
963 field0Nodes->checkConsistencyLight();
964 CPPUNIT_ASSERT(field0Nodes->getName()==fieldsNameNode[0]);
965 CPPUNIT_ASSERT_EQUAL(1,(int)field0Nodes->getNumberOfComponents());
966 CPPUNIT_ASSERT_EQUAL(19,(int)field0Nodes->getNumberOfTuples());
967 const double expectedValues4[19]={1.,2.,2.,2.,3.,3.,3.,4.,4.,4.,5.,5.,5.,6.,6.,6.,7.,7.,7.};
968 std::transform(field0Nodes->getArray()->getPointer(),field0Nodes->getArray()->getPointer()+19,expectedValues4,diffValue3,std::minus<double>());
969 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue3,diffValue3+19),1e-12);
970 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue3,diffValue3+19),1e-12);
971 constMesh=dynamic_cast<const MEDCouplingUMesh *>(field0Nodes->getMesh());
972 CPPUNIT_ASSERT(constMesh);
973 CPPUNIT_ASSERT_EQUAL(3,constMesh->getSpaceDimension());
974 CPPUNIT_ASSERT_EQUAL(3,constMesh->getMeshDimension());
975 CPPUNIT_ASSERT_EQUAL(16,(int)constMesh->getNumberOfCells());
976 CPPUNIT_ASSERT_EQUAL(ToIdType(19),constMesh->getNumberOfNodes());
977 CPPUNIT_ASSERT_EQUAL(3,(int)constMesh->getAllGeoTypes().size());
978 for(int i=0;i<12;i++)
979 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,constMesh->getTypeOfCell(i));
980 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(12));
981 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(13));
982 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(14));
983 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(15));
984 CPPUNIT_ASSERT_EQUAL(ToIdType(90),constMesh->getNodalConnectivity()->getNbOfElems());
985 CPPUNIT_ASSERT_EQUAL(701,std::accumulate(constMesh->getNodalConnectivity()->getConstPointer(),constMesh->getNodalConnectivity()->getConstPointer()+90,0));
986 CPPUNIT_ASSERT_EQUAL(705,std::accumulate(constMesh->getNodalConnectivityIndex()->getConstPointer(),constMesh->getNodalConnectivityIndex()->getConstPointer()+17,0));
987 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(constMesh->getCoords()->getConstPointer(),constMesh->getCoords()->getConstPointer()+57,0),1e-12);
988 field0Nodes->decrRef();
990 field0Nodes=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldNode(fileName.c_str(),meshNames[0].c_str(),0,fieldsNameNode[0].c_str(),its0Node[0].first,its0Node[0].second));
991 field0Nodes->checkConsistencyLight();
992 CPPUNIT_ASSERT(field0Nodes->getName()==fieldsNameNode[0]);
993 CPPUNIT_ASSERT_EQUAL(1,(int)field0Nodes->getNumberOfComponents());
994 CPPUNIT_ASSERT_EQUAL(19,(int)field0Nodes->getNumberOfTuples());
995 const double expectedValues5[19]={1.,1.,1.,2.,2.,2.,3.,3.,3.,4.,4.,4.,5.,5.,5.,6.,6.,6.,7.};
996 std::transform(field0Nodes->getArray()->getPointer(),field0Nodes->getArray()->getPointer()+19,expectedValues5,diffValue3,std::minus<double>());
997 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue3,diffValue3+19),1e-12);
998 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue3,diffValue3+19),1e-12);
999 constMesh=dynamic_cast<const MEDCouplingUMesh *>(field0Nodes->getMesh());
1000 CPPUNIT_ASSERT(constMesh);
1001 CPPUNIT_ASSERT_EQUAL(3,constMesh->getSpaceDimension());
1002 CPPUNIT_ASSERT_EQUAL(3,constMesh->getMeshDimension());
1003 CPPUNIT_ASSERT_EQUAL(16,(int)constMesh->getNumberOfCells());
1004 CPPUNIT_ASSERT_EQUAL(ToIdType(19),constMesh->getNumberOfNodes());
1005 CPPUNIT_ASSERT_EQUAL(3,(int)constMesh->getAllGeoTypes().size());
1006 for(int i=0;i<12;i++)
1007 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,constMesh->getTypeOfCell(i));
1008 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(12));
1009 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(13));
1010 CPPUNIT_ASSERT_EQUAL(NORM_HEXA8,constMesh->getTypeOfCell(14));
1011 CPPUNIT_ASSERT_EQUAL(NORM_PYRA5,constMesh->getTypeOfCell(15));
1012 CPPUNIT_ASSERT_EQUAL(ToIdType(90),constMesh->getNodalConnectivity()->getNbOfElems());
1013 CPPUNIT_ASSERT_EQUAL(701,std::accumulate(constMesh->getNodalConnectivity()->getConstPointer(),constMesh->getNodalConnectivity()->getConstPointer()+90,0));
1014 CPPUNIT_ASSERT_EQUAL(705,std::accumulate(constMesh->getNodalConnectivityIndex()->getConstPointer(),constMesh->getNodalConnectivityIndex()->getConstPointer()+17,0));
1015 CPPUNIT_ASSERT_DOUBLES_EQUAL(46.,std::accumulate(constMesh->getCoords()->getConstPointer(),constMesh->getCoords()->getConstPointer()+57,0),1e-12);
1016 field0Nodes->decrRef();
1019 void MEDLoaderTest::testMEDLoaderPolygonRead()
1021 using namespace std;
1022 using namespace INTERP_KERNEL;
1024 string fileName=INTERP_TEST::getResourceFile("polygones.med", 3);
1025 vector<string> meshNames=GetMeshNames(fileName.c_str());
1026 CPPUNIT_ASSERT_EQUAL(1,(int)meshNames.size());
1027 CPPUNIT_ASSERT(meshNames[0]=="Bord");
1028 MEDCouplingUMesh *mesh=ReadUMeshFromFile(fileName.c_str(),meshNames[0].c_str(),0);
1029 mesh->checkConsistencyLight();
1030 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
1031 CPPUNIT_ASSERT_EQUAL(2,mesh->getMeshDimension());
1032 CPPUNIT_ASSERT_EQUAL(538,(int)mesh->getNumberOfCells());
1033 CPPUNIT_ASSERT_EQUAL(ToIdType(579),mesh->getNumberOfNodes());
1034 CPPUNIT_ASSERT_EQUAL(2,(int)mesh->getAllGeoTypes().size());
1035 for(int i=0;i<514;i++)
1036 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(i));
1037 for(int i=514;i<538;i++)
1038 CPPUNIT_ASSERT_EQUAL(NORM_POLYGON,mesh->getTypeOfCell(i));
1039 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,std::accumulate(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+1737,0),1e-12);
1040 const double expectedVals1[12]={1.4851585216522212,-0.5,0.,1.4851585216522212,-0.4,0.,1.4851585216522212,-0.3,0., 1.5741585216522211, -0.5, 0. };
1041 double diffValue1[12];
1042 std::transform(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+12,expectedVals1,diffValue1,std::minus<double>());
1043 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue1,diffValue1+12),1e-12);
1044 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue1,diffValue1+12),1e-12);
1045 CPPUNIT_ASSERT_EQUAL(ToIdType(2768),mesh->getNodalConnectivity()->getNbOfElems());
1046 CPPUNIT_ASSERT_EQUAL(651050,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+2768,0));
1047 CPPUNIT_ASSERT_EQUAL(725943,std::accumulate(mesh->getNodalConnectivityIndex()->getPointer(),mesh->getNodalConnectivityIndex()->getPointer()+539,0));
1050 std::vector<std::string> fieldsName=GetCellFieldNamesOnMesh(fileName.c_str(),meshNames[0].c_str());
1051 CPPUNIT_ASSERT_EQUAL(3,(int)fieldsName.size());
1052 CPPUNIT_ASSERT(fieldsName[0]=="bord_:_distorsion");
1053 CPPUNIT_ASSERT(fieldsName[1]=="bord_:_familles");
1054 CPPUNIT_ASSERT(fieldsName[2]=="bord_:_non-ortho");
1055 std::vector<std::pair<int,int> > its0=GetCellFieldIterations(fileName.c_str(),meshNames[0].c_str(),fieldsName[0].c_str());
1056 CPPUNIT_ASSERT_EQUAL(1,(int)its0.size());
1057 MEDCouplingFieldDouble *field=dynamic_cast<MEDCouplingFieldDouble *>((MEDCouplingField *)ReadFieldCell(fileName.c_str(),meshNames[0].c_str(),0,fieldsName[0].c_str(),its0[0].first,its0[0].second));
1058 field->checkConsistencyLight();
1059 CPPUNIT_ASSERT(field->getName()==fieldsName[0]);
1060 CPPUNIT_ASSERT_EQUAL(1,(int)field->getNumberOfComponents());
1061 CPPUNIT_ASSERT_EQUAL(538,(int)field->getNumberOfTuples());
1062 const MEDCouplingUMesh *constMesh=dynamic_cast<const MEDCouplingUMesh *>(field->getMesh());
1063 CPPUNIT_ASSERT(constMesh);
1064 CPPUNIT_ASSERT_EQUAL(3,constMesh->getSpaceDimension());
1065 CPPUNIT_ASSERT_EQUAL(2,constMesh->getMeshDimension());
1066 CPPUNIT_ASSERT_EQUAL(538,(int)constMesh->getNumberOfCells());
1067 CPPUNIT_ASSERT_EQUAL(ToIdType(579),constMesh->getNumberOfNodes());
1068 CPPUNIT_ASSERT_EQUAL(2,(int)constMesh->getAllGeoTypes().size());
1069 for(int i=0;i<514;i++)
1070 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,constMesh->getTypeOfCell(i));
1071 for(int i=514;i<538;i++)
1072 CPPUNIT_ASSERT_EQUAL(NORM_POLYGON,constMesh->getTypeOfCell(i));
1073 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,std::accumulate(constMesh->getCoords()->getConstPointer(),constMesh->getCoords()->getConstPointer()+1737,0),1e-12);
1074 std::transform(constMesh->getCoords()->getConstPointer(),constMesh->getCoords()->getConstPointer()+12,expectedVals1,diffValue1,std::minus<double>());
1075 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::max_element(diffValue1,diffValue1+12),1e-12);
1076 CPPUNIT_ASSERT_DOUBLES_EQUAL(0.,*std::min_element(diffValue1,diffValue1+12),1e-12);
1077 CPPUNIT_ASSERT_EQUAL(ToIdType(2768),constMesh->getNodalConnectivity()->getNbOfElems());
1078 CPPUNIT_ASSERT_EQUAL(651050,std::accumulate(constMesh->getNodalConnectivity()->getConstPointer(),constMesh->getNodalConnectivity()->getConstPointer()+2768,0));
1079 CPPUNIT_ASSERT_EQUAL(725943,std::accumulate(constMesh->getNodalConnectivityIndex()->getConstPointer(),constMesh->getNodalConnectivityIndex()->getConstPointer()+539,0));
1080 const double *values=field->getArray()->getPointer();
1081 CPPUNIT_ASSERT_DOUBLES_EQUAL(2.87214203182918,std::accumulate(values,values+538,0.),1e-12);
1085 void MEDLoaderTest::testMEDLoaderPolyhedronRead()
1087 using namespace std;
1088 using namespace INTERP_KERNEL;
1090 string fileName=INTERP_TEST::getResourceFile("poly3D.med", 3);
1091 vector<string> meshNames=GetMeshNames(fileName.c_str());
1092 CPPUNIT_ASSERT_EQUAL(1,(int)meshNames.size());
1093 CPPUNIT_ASSERT(meshNames[0]=="poly3D");
1094 MEDCouplingUMesh *mesh=ReadUMeshFromFile(fileName.c_str(),meshNames[0].c_str(),0);
1095 mesh->checkConsistencyLight();
1096 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
1097 CPPUNIT_ASSERT_EQUAL(3,mesh->getMeshDimension());
1098 CPPUNIT_ASSERT_EQUAL(3,(int)mesh->getNumberOfCells());
1099 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
1100 CPPUNIT_ASSERT_EQUAL(2,(int)mesh->getAllGeoTypes().size());
1101 CPPUNIT_ASSERT_EQUAL(NORM_TETRA4,mesh->getTypeOfCell(0));
1102 CPPUNIT_ASSERT_EQUAL(NORM_POLYHED,mesh->getTypeOfCell(1));
1103 CPPUNIT_ASSERT_EQUAL(NORM_POLYHED,mesh->getTypeOfCell(2));
1104 CPPUNIT_ASSERT_EQUAL(ToIdType(98),mesh->getNodalConnectivity()->getNbOfElems());
1105 CPPUNIT_ASSERT_EQUAL(725,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+98,0));
1106 CPPUNIT_ASSERT_DOUBLES_EQUAL(110.,std::accumulate(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+57,0),1e-12);
1107 CPPUNIT_ASSERT_EQUAL(155,std::accumulate(mesh->getNodalConnectivityIndex()->getPointer(),mesh->getNodalConnectivityIndex()->getPointer()+4,0));
1110 mesh=ReadUMeshFromFile(fileName.c_str(),meshNames[0].c_str(),-1);
1111 mesh->checkConsistencyLight();
1112 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
1113 CPPUNIT_ASSERT_EQUAL(2,mesh->getMeshDimension());
1114 CPPUNIT_ASSERT_EQUAL(17,(int)mesh->getNumberOfCells());
1115 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
1116 CPPUNIT_ASSERT_EQUAL(3,(int)mesh->getAllGeoTypes().size());
1117 CPPUNIT_ASSERT_EQUAL(NORM_POLYGON,mesh->getTypeOfCell(0));
1118 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(1));
1119 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(2));
1120 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(3));
1121 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(4));
1122 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(5));
1123 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(6));
1124 CPPUNIT_ASSERT_EQUAL(NORM_TRI3,mesh->getTypeOfCell(7));
1125 CPPUNIT_ASSERT_EQUAL(NORM_POLYGON,mesh->getTypeOfCell(8));
1126 CPPUNIT_ASSERT_EQUAL(NORM_POLYGON,mesh->getTypeOfCell(9));
1127 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(10));
1128 CPPUNIT_ASSERT_EQUAL(NORM_TRI3,mesh->getTypeOfCell(11));
1129 CPPUNIT_ASSERT_EQUAL(NORM_TRI3,mesh->getTypeOfCell(12));
1130 CPPUNIT_ASSERT_EQUAL(NORM_TRI3,mesh->getTypeOfCell(13));
1131 CPPUNIT_ASSERT_EQUAL(NORM_TRI3,mesh->getTypeOfCell(14));
1132 CPPUNIT_ASSERT_EQUAL(NORM_QUAD4,mesh->getTypeOfCell(15));
1133 CPPUNIT_ASSERT_EQUAL(NORM_TRI3,mesh->getTypeOfCell(16));
1134 CPPUNIT_ASSERT_DOUBLES_EQUAL(110.,std::accumulate(mesh->getCoords()->getPointer(),mesh->getCoords()->getPointer()+57,0),1e-12);
1135 CPPUNIT_ASSERT_EQUAL(ToIdType(83),mesh->getNodalConnectivity()->getNbOfElems());
1136 CPPUNIT_ASSERT_EQUAL(619,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+83,0));
1139 vector<string> families=GetMeshFamiliesNames(fileName.c_str(),meshNames[0].c_str());
1140 CPPUNIT_ASSERT_EQUAL(4,(int)families.size());
1141 CPPUNIT_ASSERT(families[0]=="FAMILLE_FACE_POLYGONS3");
1142 CPPUNIT_ASSERT(families[1]=="FAMILLE_FACE_QUAD41");
1143 CPPUNIT_ASSERT(families[2]=="FAMILLE_FACE_TRIA32");
1144 CPPUNIT_ASSERT(families[3]=="FAMILLE_ZERO");
1145 vector<string> families2;
1146 families2.push_back(families[0]);
1147 mesh=ReadUMeshFromFamilies(fileName.c_str(),meshNames[0].c_str(),-1,families2);
1148 mesh->checkConsistencyLight();
1149 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
1150 CPPUNIT_ASSERT_EQUAL(2,mesh->getMeshDimension());
1151 CPPUNIT_ASSERT_EQUAL(3,(int)mesh->getNumberOfCells());
1152 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
1153 CPPUNIT_ASSERT_EQUAL(1,(int)mesh->getAllGeoTypes().size());
1154 for(int i=0;i<3;i++)
1155 CPPUNIT_ASSERT_EQUAL(NORM_POLYGON,mesh->getTypeOfCell(i));
1156 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNodalConnectivity()->getNbOfElems());
1157 CPPUNIT_ASSERT_EQUAL(117,std::accumulate(mesh->getNodalConnectivity()->getPointer(),mesh->getNodalConnectivity()->getPointer()+19,0));
1160 mesh=ReadUMeshFromFamilies(fileName.c_str(),meshNames[0].c_str(),0,families2);
1161 CPPUNIT_ASSERT_EQUAL(3,mesh->getSpaceDimension());
1162 CPPUNIT_ASSERT_EQUAL(0,(int)mesh->getNumberOfCells());
1163 CPPUNIT_ASSERT_EQUAL(ToIdType(19),mesh->getNumberOfNodes());
1164 CPPUNIT_ASSERT_EQUAL(3,mesh->getMeshDimension());
1165 CPPUNIT_ASSERT_EQUAL(0,(int)mesh->getAllGeoTypes().size());
1169 MEDCouplingUMesh *MEDLoaderTest::build1DMesh_1()
1171 double coords[6]={ 0.0, 0.3, 0.75, 1.0, 1.4, 1.3 };
1172 mcIdType conn[9]={ 0,1, 1,2, 2,3 , 3,4,5};
1173 MEDCouplingUMesh *mesh=MEDCouplingUMesh::New();
1174 mesh->setName("1DMesh_1");
1175 mesh->setMeshDimension(1);
1176 mesh->allocateCells(4);
1177 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
1178 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+2);
1179 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+4);
1180 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG3,3,conn+6);
1181 mesh->finishInsertingCells();
1182 DataArrayDouble *myCoords=DataArrayDouble::New();
1183 myCoords->alloc(6,1);
1184 myCoords->setInfoOnComponent(0,"tototototototot [m*m*m*m*m*m*m*m]");
1185 std::copy(coords,coords+6,myCoords->getPointer());
1186 mesh->setCoords(myCoords);
1187 myCoords->decrRef();
1191 MEDCouplingUMesh *MEDLoaderTest::build2DCurveMesh_1()
1193 double coords[12]={ 0.0,0.0, 0.3,0.3, 0.75,0.75, 1.0,1.0, 1.4,1.4, 1.3,1.3 };
1194 mcIdType conn[9]={ 0,1, 1,2, 2,3 , 3,4,5};
1195 MEDCouplingUMesh *mesh=MEDCouplingUMesh::New();
1196 mesh->setName("2DCurveMesh_1");
1197 mesh->setMeshDimension(1);
1198 mesh->allocateCells(4);
1199 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
1200 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+2);
1201 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+4);
1202 mesh->insertNextCell(INTERP_KERNEL::NORM_SEG3,3,conn+6);
1203 mesh->finishInsertingCells();
1204 DataArrayDouble *myCoords=DataArrayDouble::New();
1205 myCoords->alloc(6,2);
1206 std::copy(coords,coords+12,myCoords->getPointer());
1207 mesh->setCoords(myCoords);
1208 myCoords->decrRef();
1212 MEDCouplingUMesh *MEDLoaderTest::build2DMesh_1()
1214 double targetCoords[24]={-0.3,-0.3, 0.2,-0.3, 0.7,-0.3, -0.3,0.2, 0.2,0.2, 0.7,0.2, -0.3,0.7, 0.2,0.7, 0.7,0.7, -0.05,0.95, 0.2,1.2, 0.45,0.95 };
1215 mcIdType targetConn[24]={1,4,2, 4,5,2, 6,10,8,9,11,7, 0,3,4,1, 6,7,4,3, 7,8,5,4};
1216 MEDCouplingUMesh *targetMesh=MEDCouplingUMesh::New();
1217 targetMesh->setMeshDimension(2);
1218 targetMesh->allocateCells(6);
1219 targetMesh->setName("2DMesh_1");
1220 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1221 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
1222 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI6,6,targetConn+6);
1223 targetMesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+12);
1224 targetMesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+16);
1225 targetMesh->insertNextCell(INTERP_KERNEL::NORM_POLYGON,4,targetConn+20);
1226 targetMesh->finishInsertingCells();
1227 DataArrayDouble *myCoords=DataArrayDouble::New();
1228 myCoords->alloc(12,2);
1229 myCoords->setInfoOnComponent(0,"tototototototot [m]");
1230 myCoords->setInfoOnComponent(1,"energie [kW]");
1231 std::copy(targetCoords,targetCoords+24,myCoords->getPointer());
1232 targetMesh->setCoords(myCoords);
1233 myCoords->decrRef();
1237 MEDCouplingUMesh *MEDLoaderTest::build2DMesh_2()
1239 double targetCoords[24]={
1240 -0.3,-0.3, 0.2,-0.3, 0.7,-0.3, -0.3,0.2, 0.2,0.2, 0.7,0.2, -0.3,0.7, 0.2,0.7, 0.7,0.7,
1241 -0.05,0.95, 0.2,1.2, 0.45,0.95
1243 mcIdType targetConn[24]={1,4,2, 4,5,2, 6,10,8,9,11,7, 0,3,4,1, 6,7,4,3, 7,8,5,4};
1244 MEDCouplingUMesh *targetMesh=MEDCouplingUMesh::New();
1245 targetMesh->setMeshDimension(2);
1246 targetMesh->allocateCells(5);
1247 targetMesh->setName("2DMesh_2");
1248 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1249 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
1250 targetMesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+12);
1251 targetMesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+16);
1252 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI6,6,targetConn+6);
1253 targetMesh->finishInsertingCells();
1254 DataArrayDouble *myCoords=DataArrayDouble::New();
1255 myCoords->alloc(12,2);
1256 myCoords->setInfoOnComponent(0,"toto [m]");
1257 myCoords->setInfoOnComponent(1,"energie [kW]");
1258 std::copy(targetCoords,targetCoords+24,myCoords->getPointer());
1259 targetMesh->setCoords(myCoords);
1260 myCoords->decrRef();
1264 MEDCouplingUMesh *MEDLoaderTest::build3DSurfMesh_1()
1266 double targetCoords[36]={
1267 -0.3,-0.3,-0.3, 0.2,-0.3,-0.3, 0.7,-0.3,-0.3, -0.3,0.2,-0.3, 0.2,0.2,-0.3, 0.7,0.2,-0.3, -0.3,0.7,-0.3, 0.2,0.7,-0.3, 0.7,0.7,-0.3
1268 ,-0.05,0.95,-0.3, 0.2,1.2,-0.3, 0.45,0.95,-0.3
1270 mcIdType targetConn[24]={1,4,2, 4,5,2, 6,10,8,9,11,7, 0,3,4,1, 6,7,4,3, 7,8,5,4};
1271 MEDCouplingUMesh *targetMesh=MEDCouplingUMesh::New();
1272 targetMesh->setMeshDimension(2);
1273 targetMesh->allocateCells(6);
1274 targetMesh->setName("3DSurfMesh_1");
1275 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1276 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
1277 targetMesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+12);
1278 targetMesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+16);
1279 targetMesh->insertNextCell(INTERP_KERNEL::NORM_TRI6,6,targetConn+6);
1280 targetMesh->insertNextCell(INTERP_KERNEL::NORM_POLYGON,4,targetConn+20);
1281 targetMesh->finishInsertingCells();
1282 DataArrayDouble *myCoords=DataArrayDouble::New();
1283 myCoords->alloc(12,3);
1284 myCoords->setInfoOnComponent(0,"toto [m]");
1285 myCoords->setInfoOnComponent(2,"ff [km]");//component 1 is not set for test
1286 std::copy(targetCoords,targetCoords+36,myCoords->getPointer());
1287 targetMesh->setCoords(myCoords);
1288 myCoords->decrRef();
1292 MEDCouplingUMesh *MEDLoaderTest::build3DMesh_1()
1294 double coords[180]={
1295 0.,0.,0., 1.,1.,0., 1.,1.25,0., 0.,1.,0., 1.,1.5,0., 2.,0.,0., 2.,1.,0., 1.,2.,0., 0.,2.,0., 3.,1.,0.,
1296 3.,2.,0., 0.,1.,0., 1.,3.,0., 2.,2.,0., 2.,3.,0.,
1297 0.,0.,1., 1.,1.,1., 1.,1.25,1., 0.,1.,1., 1.,1.5,1., 2.,0.,1., 2.,1.,1., 1.,2.,1., 0.,2.,1., 3.,1.,1.,
1298 3.,2.,1., 0.,1.,1., 1.,3.,1., 2.,2.,1., 2.,3.,1.,
1299 0.,0.,2., 1.,1.,2., 1.,1.25,2., 0.,1.,2., 1.,1.5,2., 2.,0.,2., 2.,1.,2., 1.,2.,2., 0.,2.,2., 3.,1.,2.,
1300 3.,2.,2., 0.,1.,2., 1.,3.,2., 2.,2.,2., 2.,3.,2.,
1301 0.,0.,3., 1.,1.,3., 1.,1.25,3., 0.,1.,3., 1.,1.5,3., 2.,0.,3., 2.,1.,3., 1.,2.,3., 0.,2.,3., 3.,1.,3.,
1302 3.,2.,3., 0.,1.,3., 1.,3.,3., 2.,2.,3., 2.,3.,3.};
1304 mcIdType conn[354]={
1306 0,11,1,3,15,26,16,18, 1,2,4,7,13,6,-1,1,16,21,6,-1,6,21,28,13,-1,13,7,22,28,-1,7,4,19,22,-1,4,2,17,19,-1,2,1,16,17,-1,16,21,28,22,19,17,
1307 1,6,5,3,16,21,20,18, 13,10,9,6,28,25,24,21,
1308 11,8,7,4,2,1,-1,11,26,16,1,-1,1,16,17,2,-1,2,17,19,4,-1,4,19,22,7,-1,7,8,23,22,-1,8,11,26,23,-1,26,16,17,19,22,23,
1309 7,12,14,13,22,27,29,28,
1311 15,26,16,18,30,41,31,33, 16,17,19,22,28,21,-1,16,31,36,21,-1,21,36,43,28,-1,28,22,37,43,-1,22,19,34,37,-1,19,17,32,34,-1,17,16,31,32,-1,31,36,43,37,34,32,
1312 16,21,20,18,31,36,35,33, 28,25,24,21,43,40,39,36,
1313 26,23,22,19,17,16,-1,26,41,31,16,-1,16,31,32,17,-1,17,32,34,19,-1,19,34,37,22,-1,22,23,38,37,-1,23,26,41,38,-1,41,31,32,34,37,38,
1314 22,27,29,28,37,42,44,43,
1316 30,41,31,33,45,56,46,48, 31,32,34,37,43,36,-1,31,46,51,36,-1,36,51,58,43,-1,43,37,52,58,-1,37,34,49,52,-1,34,32,47,49,-1,32,31,46,47,-1,46,51,58,52,49,47,
1317 31,36,35,33,46,51,50,48, 43,40,39,36,58,55,54,51,
1318 41,38,37,34,32,31,-1,41,56,46,31,-1,31,46,47,32,-1,32,47,49,34,-1,34,49,52,37,-1,37,38,53,52,-1,38,41,56,53,-1,56,46,47,49,52,53,
1319 37,42,44,43,52,57,59,58
1322 MEDCouplingUMesh *ret=MEDCouplingUMesh::New();
1323 ret->setName("3DMesh_1");
1324 ret->setMeshDimension(3);
1325 ret->allocateCells(18);
1327 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn);
1328 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+51);
1329 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+59);
1330 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+110);
1332 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+118);
1333 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+169);
1334 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+177);
1335 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+228);
1337 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+236);
1338 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+287);
1339 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+295);
1340 ret->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,conn+346);
1342 ret->insertNextCell(INTERP_KERNEL::NORM_POLYHED,43,conn+8);
1343 ret->insertNextCell(INTERP_KERNEL::NORM_POLYHED,43,conn+67);
1344 ret->insertNextCell(INTERP_KERNEL::NORM_POLYHED,43,conn+126);
1345 ret->insertNextCell(INTERP_KERNEL::NORM_POLYHED,43,conn+185);
1346 ret->insertNextCell(INTERP_KERNEL::NORM_POLYHED,43,conn+244);
1347 ret->insertNextCell(INTERP_KERNEL::NORM_POLYHED,43,conn+303);
1349 ret->finishInsertingCells();
1350 DataArrayDouble *myCoords=DataArrayDouble::New();
1351 myCoords->alloc(60,3);
1352 myCoords->setInfoOnComponent(0,"titi [m]");
1353 myCoords->setInfoOnComponent(1,"density power [MW/m^3]");
1354 myCoords->setInfoOnComponent(2,"t [kW]");
1355 std::copy(coords,coords+180,myCoords->getPointer());
1356 ret->setCoords(myCoords);
1357 myCoords->decrRef();
1361 MEDCouplingUMesh *MEDLoaderTest::build3DMesh_2()
1363 MEDCouplingUMesh *m3dsurfBase=build3DSurfMesh_1();
1364 mcIdType numbers[5]={0,1,2,3,5};
1365 MEDCouplingUMesh *m3dsurf=(MEDCouplingUMesh *)m3dsurfBase->buildPartOfMySelf(numbers,numbers+5,false);
1366 m3dsurfBase->decrRef();
1367 MEDCouplingUMesh *m1dBase=build1DMesh_1();
1368 mcIdType numbers2[4]={0,1,2,3};
1369 MEDCouplingUMesh *m1d=(MEDCouplingUMesh *)m1dBase->buildPartOfMySelf(numbers2,numbers2+4,false);
1371 m1d->changeSpaceDimension(3);
1372 const double vec[3]={0.,1.,0.};
1373 const double pt[3]={0.,0.,0.};
1374 m1d->rotate(pt,vec,-M_PI/2.);
1375 MEDCouplingUMesh *ret=m3dsurf->buildExtrudedMesh(m1d,0);
1381 MEDCouplingFieldDouble *MEDLoaderTest::buildVecFieldOnCells_1()
1383 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
1384 mcIdType nbOfCells=mesh->getNumberOfCells();
1385 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
1386 f1->setName("VectorFieldOnCells");
1388 DataArrayDouble *array=DataArrayDouble::New();
1389 array->alloc(nbOfCells,3);
1390 array->setInfoOnComponent(0,"power [MW/m^3]");
1391 array->setInfoOnComponent(1,"density [g/cm^3]");
1392 array->setInfoOnComponent(2,"temperature [K]");
1393 f1->setArray(array);
1395 double *tmp=array->getPointer();
1396 const double arr1[18]={0.,10.,20.,1.,11.,21.,2.,12.,22.,3.,13.,23.,4.,14.,24.,5.,15.,25.};
1397 std::copy(arr1,arr1+18,tmp);
1398 f1->setTime(2.,0,1);
1399 f1->checkConsistencyLight();
1404 MEDCouplingFieldInt *MEDLoaderTest::buildIntVecFieldOnCells_1()
1406 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
1407 mcIdType nbOfCells=mesh->getNumberOfCells();
1408 MEDCouplingFieldInt *f1=MEDCouplingFieldInt::New(ON_CELLS,ONE_TIME);
1409 f1->setName("IntVectorFieldOnCells");
1411 DataArrayInt *array=DataArrayInt::New();
1412 array->alloc(nbOfCells,3);
1413 array->setInfoOnComponent(0,"val1 [MW/m^3]");
1414 array->setInfoOnComponent(1,"va2 [g/cm^3]");
1415 array->setInfoOnComponent(2,"val3 [K]");
1416 f1->setArray(array);
1418 mcIdType *tmp=array->getPointer();
1419 const mcIdType arr1[18]={0,10,20,1,11,21,2,12,22,3,13,23,4,14,24,5,15,25};
1420 std::copy(arr1,arr1+18,tmp);
1421 f1->setTime(2.,0,1);
1422 f1->checkConsistencyLight();
1427 MEDCouplingFieldFloat *MEDLoaderTest::buildFloatVecFieldOnCells_1()
1429 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
1430 mcIdType nbOfCells=mesh->getNumberOfCells();
1431 MEDCouplingFieldFloat *f1=MEDCouplingFieldFloat::New(ON_CELLS,ONE_TIME);
1432 f1->setName("FloatVectorFieldOnCells");
1434 DataArrayFloat *array=DataArrayFloat::New();
1435 array->alloc(nbOfCells,3);
1436 array->setInfoOnComponent(0,"power [MW/m^3]");
1437 array->setInfoOnComponent(1,"density [g/cm^3]");
1438 array->setInfoOnComponent(2,"temperature [K]");
1439 f1->setArray(array);
1441 float *tmp=array->getPointer();
1442 const float arr1[18]={0.,10.,20.,1.,11.,21.,2.,12.,22.,3.,13.,23.,4.,14.,24.,5.,15.,25.};
1443 std::copy(arr1,arr1+18,tmp);
1444 f1->setTime(2.,0,1);
1445 f1->checkConsistencyLight();
1451 MEDCouplingFieldDouble *MEDLoaderTest::buildVecFieldOnNodes_1()
1453 MEDCouplingUMesh *mesh=build3DSurfMesh_1();
1454 mcIdType nbOfNodes=mesh->getNumberOfNodes();
1455 MEDCouplingFieldDouble *f1=MEDCouplingFieldDouble::New(ON_NODES,ONE_TIME);
1456 f1->setName("VectorFieldOnNodes");
1458 DataArrayDouble *array=DataArrayDouble::New();
1459 array->alloc(nbOfNodes,3);
1460 f1->setArray(array);
1461 array->setInfoOnComponent(0,"power [MW/m^3]");
1462 array->setInfoOnComponent(1,"density [g/cm^3]");
1463 array->setInfoOnComponent(2,"temperature [K]");
1465 double *tmp=array->getPointer();
1466 const double arr1[36]={
1467 70.,80.,90.,71.,81.,91.,72.,82.,92.,73.,83.,93.,74.,84.,94.,75.,85.,95.,
1468 1000.,10010.,10020.,1001.,10011.,10021.,1002.,10012.,10022.,1003.,10013.,10023.,1004.,10014.,10024.,1005.,10015.,10025.,
1470 std::copy(arr1,arr1+36,tmp);
1471 f1->setTime(2.12,2,3);
1472 f1->checkConsistencyLight();
1477 MEDCouplingFieldDouble *MEDLoaderTest::buildVecFieldOnGauss_1()
1479 const double _a=0.446948490915965;
1480 const double _b=0.091576213509771;
1481 const double _p1=0.11169079483905;
1482 const double _p2=0.0549758718227661;
1483 const double refCoo1[6]={ 0.,0., 1.,0., 0.,1. };
1484 const double gsCoo1[12]={ 2*_b-1, 1-4*_b, 2*_b-1, 2.07*_b-1, 1-4*_b,
1485 2*_b-1, 1-4*_a, 2*_a-1, 2*_a-1, 1-4*_a, 2*_a-1, 2*_a-1 };
1486 const double wg1[6]={ 4*_p2, 4*_p2, 4*_p2, 4*_p1, 4*_p1, 4*_p1 };
1487 std::vector<double> _refCoo1(refCoo1,refCoo1+6);
1488 std::vector<double> _gsCoo1(gsCoo1,gsCoo1+12);
1489 std::vector<double> _wg1(wg1,wg1+6);
1490 MEDCouplingUMesh *m=build2DMesh_2();
1491 MEDCouplingFieldDouble *f=MEDCouplingFieldDouble::New(ON_GAUSS_PT,ONE_TIME);
1492 f->setTime(3.14,1,5);
1494 f->setGaussLocalizationOnType(INTERP_KERNEL::NORM_TRI3,_refCoo1,_gsCoo1,_wg1);
1495 const double refCoo2[12]={-1.0,1.0, -1.0,-1.0, 1.0,-1.0, -1.0,0.0, 0.0,-1.0, 0.0,0.0 };
1496 std::vector<double> _refCoo2(refCoo2,refCoo2+12);
1497 std::vector<double> _gsCoo2(_gsCoo1);
1498 std::vector<double> _wg2(_wg1);
1499 _gsCoo2.resize(6); _wg2.resize(3);
1500 const double refCoo3[8]={ 0.,0., 1.,0., 1.,1., 0.,1. };
1501 std::vector<double> _refCoo3(refCoo3,refCoo3+8);
1502 _gsCoo1.resize(4); _wg1.resize(2);
1503 f->setGaussLocalizationOnType(INTERP_KERNEL::NORM_QUAD4,_refCoo3,_gsCoo1,_wg1);
1504 f->setGaussLocalizationOnType(INTERP_KERNEL::NORM_TRI6,_refCoo2,_gsCoo2,_wg2);
1505 DataArrayDouble *array=DataArrayDouble::New();
1507 double *ptr=array->getPointer();
1508 for(int i=0;i<19*2;i++)
1509 ptr[i]=(double)(i+7);
1511 f->setName("MyFirstFieldOnGaussPoint");
1512 array->setInfoOnComponent(0,"power [MW/m^3]");
1513 array->setInfoOnComponent(1,"density");
1515 f->checkConsistencyLight();
1520 MEDCouplingFieldDouble *MEDLoaderTest::buildVecFieldOnGaussNE_1()
1522 MEDCouplingUMesh *m=build2DMesh_2();
1523 MEDCouplingFieldDouble *f=MEDCouplingFieldDouble::New(ON_GAUSS_NE,ONE_TIME);
1524 f->setTime(3.14,1,5);
1526 DataArrayDouble *array=DataArrayDouble::New();
1528 double *ptr=array->getPointer();
1529 for(int i=0;i<20*2;i++)
1530 ptr[i]=(double)(i+8);
1532 array->setInfoOnComponent(0,"power [W]");
1533 array->setInfoOnComponent(1,"temperature");
1534 f->setName("MyFieldOnGaussNE");
1536 f->checkConsistencyLight();