]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
[Partial load] Adding test cases
authorAnida Khizar <anida.khizar@cea.fr>
Mon, 30 Jan 2023 15:34:18 +0000 (16:34 +0100)
committerAnida Khizar <anida.khizar@cea.fr>
Mon, 30 Jan 2023 15:49:27 +0000 (16:49 +0100)
resources/CMakeLists.txt
resources/SimpleTest2D.med [new file with mode: 0644]
src/ParaMEDLoader/ParaMEDFileMesh.cxx
src/ParaMEDMEMTest/CMakeLists.txt
src/ParaMEDMEMTest/ParaMEDMEMTest.hxx
src/ParaMEDMEMTest/ParaMEDMEMTest_MEDLoader.cxx

index 85482604ea174614619404fe574a885a47559fa6..12135b962242fe5994c761322056246408e93876 100644 (file)
@@ -39,6 +39,7 @@ SET(MED_other_FILES
   Test2Dpoly.med
   Test3D.med
   Test3Dpoly.med
+  SimpleTest2D.med
   UnitTetraDegenT.med
   DegenEdgeXY.med
   DegenFaceXYZ.med
diff --git a/resources/SimpleTest2D.med b/resources/SimpleTest2D.med
new file mode 100644 (file)
index 0000000..d4cfcd7
Binary files /dev/null and b/resources/SimpleTest2D.med differ
index 364141bd0777a6bd784977a05ff47885cee44df9..013079d059546b307f2224139cdc3e11cd672f5c 100644 (file)
@@ -150,6 +150,7 @@ MEDFileUMesh *ParaMEDFileUMesh::NewPrivate(med_idt fid, const MPI_Comm& com, con
   med_geometry_type geoMedType(typmai3[geoType]);
   med_bool changement,transformation;
   med_int totalNumberOfElements(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_CELL,geoMedType,MED_CONNECTIVITY,MED_NODAL,&changement,&transformation));
+
   mcIdType nbEltsInDistribLoc = distrib.size();
   mcIdType nbEltsInDistribTot = -1;
 #ifdef HAVE_MPI
index b6813e2a519a32f6d131e46d48b0baab5499403f..83c1b7648bc9a622508e98ab6ef9cf131e6ad010 100644 (file)
@@ -23,6 +23,7 @@ ADD_DEFINITIONS(${MPI_DEFINITIONS} ${CPPUNIT_DEFINITIONS})
 
 INCLUDE_DIRECTORIES(
   ${MPI_INCLUDE_DIRS}
+  ${MEDFILE_INCLUDE_DIRS}
   ${CPPUNIT_INCLUDE_DIRS}
   ${CMAKE_CURRENT_SOURCE_DIR}/../ParaMEDLoader
   ${CMAKE_CURRENT_SOURCE_DIR}/../ParaMEDMEM
@@ -46,6 +47,7 @@ SET(ParaMEDMEMTest_SOURCES
   ParaMEDMEMTest_FabienAPI.cxx
   ParaMEDMEMTest_NonCoincidentDEC.cxx
   ParaMEDMEMTest_OverlapDEC.cxx
+  ParaMEDMEMTest_MEDLoader.cxx
 )
 
 ADD_LIBRARY(ParaMEDMEMTest ${ParaMEDMEMTest_SOURCES})
index 240f948300cb62d617eab32e6604d3ac82381cc6..b26be735664873eac9f30fa68096e0bfcebe9493 100644 (file)
@@ -85,6 +85,8 @@ class ParaMEDMEMTest : public CppUnit::TestFixture
   CPPUNIT_TEST(testFabienAPI1);       // 3 procs
   CPPUNIT_TEST(testFabienAPI2);       // 3 procs
 
+  CPPUNIT_TEST(testParallelLoad1);   // 2 procs
+  CPPUNIT_TEST(testParallelLoad2);   // 2 procs
   CPPUNIT_TEST_SUITE_END();
 
 public:
@@ -148,6 +150,9 @@ public:
   void testFabienAPI1();
   void testFabienAPI2();
 
+  void testParallelLoad1();
+  void testParallelLoad2();
+
   std::string getTmpDirectory();
   std::string makeTmpFile( const std::string&, const std::string& = "" );
 
@@ -166,6 +171,8 @@ private:
   void testInterpKernelDEC2_2D_(const char *srcMeth, const char *targetMeth);
   void testInterpKernelDEC_3D_(const char *srcMeth, const char *targetMeth);
   void testGauthier3_GEN(bool, int);
+
+
 };
 
 // to automatically remove temporary files from disk
index 485186b09ec2ad797abb3b2a13e53a5d85f0c59c..379846bd0c55993002e450d1fd342089204766e4 100644 (file)
 
 #include "ParaMEDMEMTest.hxx"
 #include "MEDLoader.hxx"
-#include "MEDCouplingUMesh.hxx"
+
+#include "ParaMEDFileMesh.hxx"
+#include "MEDFileMesh.hxx"
+#include "MEDFileField1TS.hxx"
+#include "TestInterpKernelUtils.hxx"
 #include "MEDCouplingFieldDouble.hxx"
 
 #include <cppunit/TestAssert.h>
 #include <iostream>
 #include <iterator>
 
-using namespace std;
-using namespace INTERP_KERNEL;
 using namespace MEDCoupling;
 
+MEDCouplingUMesh* genLocMesh(int rk)
+{
+  int nxTot=4,nyTot=2;
+  int nx=2,ny=2;
+  MCAuto<MEDCouplingCMesh> msh = MEDCouplingCMesh::New("mesh");
+  MCAuto<DataArrayDouble> dax = DataArrayDouble::New(); dax->alloc(nx+1,1);
+  MCAuto<DataArrayDouble> day = DataArrayDouble::New(); day->alloc(ny+1,1);
+  dax->iota(); day->iota();
+  if (rk == 0)
+    {
+      std::transform(dax->begin(), dax->end(),
+                     dax->rwBegin(),
+                     [nxTot](const int& c){return c/(float)nxTot;});
+      std::transform(day->begin(), day->end(),
+                     day->rwBegin(),
+                     [nyTot](const int& c){return c/(float)nyTot;});
+    }
+  else
+    {
+      std::transform(dax->begin(), dax->end(),
+                     dax->rwBegin(),
+                     [nxTot](const int& c){return c/(float)nxTot+0.5; });
+      std::transform(day->begin(), day->end(),
+                     day->rwBegin(),
+                     [nyTot](const int& c){return c/(float)nyTot;});
+    }
+  msh->setCoords(dax, day);
+  MCAuto<MEDCouplingUMesh> ret = msh->buildUnstructured();
+  return ret.retn();
+}
+
+MEDCouplingFieldDouble *genLocField(int rank)
+{
+  MCAuto<MEDCouplingUMesh> mesh = genLocMesh(rank);
+  MCAuto<MEDCouplingFieldDouble> f1=MEDCouplingFieldDouble::New(ON_CELLS,ONE_TIME);
+  f1->setName("field");
+  f1->setMesh(mesh);
+
+  MCAuto<DataArrayDouble> array(DataArrayDouble::New());
+  double* values = new double[4];
+  if(rank == 0)
+    values[0]= 0.,values[1]=10.,values[2]=40.,values[3]=50.;
+  else
+    values[0]=20.,values[1]=30.,values[2]=60.,values[3]=70.;
+  array->useArray(values,true, DeallocType::CPP_DEALLOC,4,1);
+  array->setInfoOnComponent(0,"s");
+  f1->setArray(array);
+  return f1.retn();
+}
+
+
+void ParaMEDMEMTest::testParallelLoad1()
+{
+  int size;
+  int rank;
+  MPI_Comm_size(MPI_COMM_WORLD,&size);
+  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
+  //
+  if(size!=2)
+    return ;
+
+  std::vector<mcIdType> distrib;
+  if (rank == 0)
+      distrib = {0,1,4,5};    //c++ type of indexing: index starts from zero!
+  else
+      distrib = {2,3,6,7};
+
+  std::string filename=INTERP_TEST::getResourceFile("SimpleTest2D.med");
+  MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
+  MCAuto<MEDCouplingUMesh> meshRef = genLocMesh(rank);
+  CPPUNIT_ASSERT(mu->getMeshAtLevel(0)->isEqual(meshRef,1e-12));
+  MPI_Barrier(MPI_COMM_WORLD);
+
+}
+
+
+void ParaMEDMEMTest::testParallelLoad2()
+{
+  int size;
+  int rank;
+  MPI_Comm_size(MPI_COMM_WORLD,&size);
+  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
+  //
+  if(size!=2)
+    return ;
+
+  std::vector<mcIdType> distrib;
+  if (rank == 0)
+      distrib = {0,1,4,5};    //c++ type of indexing: index starts from zero!
+  else
+      distrib = {2,3,6,7};
+
+  std::string filename=INTERP_TEST::getResourceFile("SimpleTest2D.med");
+  MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
+  MCAuto<MEDFileField1TS> f1TS = ParaMEDFileField1TS::ParaNew(mu, filename, "mesh", "field");
+  MCAuto<MEDCouplingFieldDouble> partField = f1TS->field(mu);
+  MCAuto<MEDCouplingFieldDouble> fieldRef = genLocField(rank);
+  CPPUNIT_ASSERT_EQUAL(4,(int)partField->getNumberOfTuples());
+  CPPUNIT_ASSERT(partField->getArray()->isEqual(*fieldRef->getArray(),1e-12));
+  MPI_Barrier(MPI_COMM_WORLD);
+}
+