MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, mcIdType strt, mcIdType stp, mcIdType step, MEDFileMeshReadSelector *mrs)
{
- int geoElt2i((int)geoElt2);
- if(geoElt2i<0 || geoElt2i>=INTERP_KERNEL::NORM_MAXTYPE)
- throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : Not recognized MEDCoupling/MEDLoader geometric type !");
- med_geometry_type geoElt(typmai3[geoElt2]);
+ med_geometry_type geoElt;
med_entity_type whichEntity;
- if(!isExisting(fid,mName,dt,it,geoElt,whichEntity))
- throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : The specified geo type is not present in the specified mesh !");
+ getMedTypes(fid, mName, dt, it, geoElt2, geoElt, whichEntity);
+
MCAuto<MEDFileUMeshPerType> ret(new MEDFileUMeshPerType);
ret->loadPart(fid,mName,dt,it,mdim,geoElt,geoElt2,whichEntity,strt,stp,step,mrs);
return ret.retn();
}
-MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs)
+
+MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, med_geometry_type geoElt, med_entity_type whichEntity, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs)
{
- int geoElt2i((int)geoElt2);
- if(geoElt2i<0 || geoElt2i>=INTERP_KERNEL::NORM_MAXTYPE)
- throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : Not recognized MEDCoupling/MEDLoader geometric type !");
- med_geometry_type geoElt(typmai3[geoElt2]);
- med_entity_type whichEntity;
- if(!isExisting(fid,mName,dt,it,geoElt,whichEntity))
- throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : The specified geo type is not present in the specified mesh !");
- MCAuto<MEDFileUMeshPerType> ret(new MEDFileUMeshPerType);
- ret->loadPart(fid,mName,dt,it,mdim,geoElt,geoElt2,whichEntity,distrib,mrs);
+ MCAuto<MEDFileUMeshPerType> ret(new MEDFileUMeshPerType);
+ ret->loadPart(fid,mName,dt,it,mdim,geoElt,geoElt2,whichEntity,distrib,mrs);
return ret.retn();
}
return ret;
}
+void MEDFileUMeshPerType::getMedTypes(med_idt fid, const char *mName, int dt, int it, INTERP_KERNEL::NormalizedCellType MCgeoElt, med_geometry_type& geoElt, med_entity_type& whichEntity)
+{
+ int geoElt2i((int)MCgeoElt);
+ if(geoElt2i<0 || geoElt2i>=INTERP_KERNEL::NORM_MAXTYPE)
+ throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::getMedTypes : Not recognized MEDCoupling/MEDLoader geometric type !");
+ geoElt = typmai3[MCgeoElt];
+ if(!isExisting(fid,mName,dt,it,geoElt,whichEntity))
+ throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::getMedTypes : The specified geo type is not present in the specified mesh !");
+}
+
bool MEDFileUMeshPerType::isExisting(med_idt fid, const char *mName, int dt, int it, med_geometry_type geoElt, med_entity_type& whichEntity)
{
static const med_entity_type entities[3]={MED_CELL,MED_DESCENDING_FACE,MED_DESCENDING_EDGE};
public:
static MEDFileUMeshPerType *New(med_idt fid, const char *mName, int dt, int it, int mdim, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType geoElt2, MEDFileMeshReadSelector *mrs);
static MEDFileUMeshPerType *NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, mcIdType strt, mcIdType stp, mcIdType step, MEDFileMeshReadSelector *mrs);
- static MEDFileUMeshPerType *NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs);
+ static MEDFileUMeshPerType *NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, med_geometry_type geoElt, med_entity_type whichEntity, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs);
std::string getClassName() const override { return std::string("MEDFileUMeshPerType"); }
static bool isExisting(med_idt fid, const char *mName, int dt, int it, med_geometry_type geoElt, med_entity_type& whichEntity);
+ static void getMedTypes(med_idt fid, const char *mName, int dt, int it, INTERP_KERNEL::NormalizedCellType MCgeoElt, med_geometry_type& geoElt, med_entity_type& whichEntity);
std::size_t getHeapMemorySizeWithoutChildren() const;
std::vector<const BigMemoryObject *> getDirectChildrenWithNull() const;
int getDim() const;
MCAuto<DataArrayIdType> fni(DataArrayIdType::BuildListOfSwitchedOn(fetchedNodeIds));
MCAuto< MapKeyVal<mcIdType, mcIdType> > o2n(fni->invertArrayN2O2O2NOptimized());
for(std::vector< std::vector< MCAuto<MEDFileUMeshPerType> > >::const_iterator it0=_per_type_mesh.begin();it0!=_per_type_mesh.end();it0++)
- for(std::vector< MCAuto<MEDFileUMeshPerType> >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
- (*it1)->getMesh()->renumberNodesInConn(o2n->data());
+ for(std::vector< MCAuto<MEDFileUMeshPerType> >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
+ (*it1)->getMesh()->renumberNodesInConn(o2n->data());
+
// loading coordinates of fetched nodes
std::vector<mcIdType> distribNodes;
for(std::map<mcIdType,mcIdType>::const_iterator mapIter = o2n->data().begin(); mapIter != o2n->data().end(); ++mapIter)
- distribNodes.push_back(mapIter->first);
- this->loadPartCoords(fid,infosOnComp,mName,dt,it,distribNodes);
+ distribNodes.push_back(mapIter->first);
+ if(distribNodes.size() != 0)
+ this->loadPartCoords(fid,infosOnComp,mName,dt,it,distribNodes);
}
else
throw INTERP_KERNEL::Exception("MEDFileUMeshL2::loadPartFromUserDistrib: multiple load sessions not handled!");
std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>::const_iterator iter;
for (iter = distrib.begin(); iter != distrib.end(); iter++)
{
- MCAuto<MEDFileUMeshPerType> tmp(MEDFileUMeshPerType::NewPart(fid,mName.c_str(),dt,it,mdim,iter->first/*type*/,iter->second/*distrib over the current type*/,mrs));
- _per_type_mesh[0].push_back(tmp);
+ INTERP_KERNEL::NormalizedCellType MCGeoType = iter->first;
+ med_geometry_type geoElt;
+ med_entity_type whichEntity;
+ MEDFileUMeshPerType::getMedTypes(fid, mName.c_str(), dt, it, MCGeoType, geoElt, whichEntity);
+ if(iter->second.size() !=0)
+ {
+ MCAuto<MEDFileUMeshPerType> tmp(MEDFileUMeshPerType::NewPart(fid,mName.c_str(),dt,it,mdim,MCGeoType,geoElt,whichEntity,iter->second/*distrib over the current type*/,mrs));
+ _per_type_mesh[0].push_back(tmp);
+ }
}
sortTypes();
}
using namespace MEDCoupling;
-void checkDistribution(const MPI_Comm& com, mcIdType totalNumberOfElements, const std::vector<mcIdType>& distrib)
-{
- mcIdType nbEltsInDistribLoc = distrib.size();
- mcIdType nbEltsInDistribTot = -1;
-#ifdef HAVE_MPI
- MPI_Allreduce(&nbEltsInDistribLoc, &nbEltsInDistribTot, 1, MPI_LONG, MPI_SUM, com);
-#else
- throw INTERP_KERNEL::Exception("not(HAVE_MPI) incompatible with MPI_World_Size>1");
-#endif
- if(nbEltsInDistribTot != totalNumberOfElements)
- {
- if(nbEltsInDistribTot > totalNumberOfElements)
- throw INTERP_KERNEL::Exception("ParaMEDFileMesh : Some of your partitions overlap each other ! Each element in your distribution vector must appear only once ! ");
- else
- throw INTERP_KERNEL::Exception("ParaMEDFileMesh : The distribution does not cover the whole mesh ! Each element of the mesh must appear once in your distribution vector ");
- }
-}
-
-
MEDFileMesh *ParaMEDFileMesh::New(int iPart, int nbOfParts, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
{
MEDFileUtilities::CheckFileForRead(fileName);
MEDFileUMesh *ParaMEDFileUMesh::NewPrivate(med_idt fid, const MPI_Comm& com, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
{
MCAuto<MEDFileUMesh> ret;
- for(std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>::const_iterator iter=distrib.begin(); iter!= distrib.end(); iter++)
- {
- med_geometry_type geoMedType(typmai3[iter->first /*current geometric type*/]);
- med_bool changement,transformation;
- med_int totalNumberOfElements(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_CELL,geoMedType,MED_CONNECTIVITY,MED_NODAL,&changement,&transformation));
- checkDistribution(com,totalNumberOfElements,iter->second /*distrib over this geometric type*/);
- }
ret=MEDFileUMesh::LoadPartOfFromUserDistrib(fid,mName,distrib,dt,it,mrs);
return ret.retn();
}
CPPUNIT_TEST(testParallelLoad3); // 2 procs
CPPUNIT_TEST(testParallelLoad4); // 2 procs
CPPUNIT_TEST(testParallelLoad5); // 2 procs
+ CPPUNIT_TEST(testParallelLoad6); // 3 procs
CPPUNIT_TEST_SUITE_END();
public:
void testParallelLoad3();
void testParallelLoad4();
void testParallelLoad5();
+ void testParallelLoad6();
std::string getTmpDirectory();
std::string makeTmpFile( const std::string&, const std::string& = "" );
return ret.retn();
}
+/*
+ * Generate a 2D mesh that is supposed to match the part that will be loaded by proc0 in testParallelLoad6
+ */
+MEDCouplingUMesh* genPartialLocMeshMultipleTypes1()
+{
+ MCAuto<MEDCouplingUMesh> ret= MEDCouplingUMesh::New("mesh",2);
+ double coords[8] = {0.,2., 1.,2., 0.,3., 1.,3.};
+ DataArrayDouble *myCoords=DataArrayDouble::New();
+ myCoords->alloc(4,2);
+ std::copy(coords,coords+8,myCoords->getPointer());
+ ret->setCoords(myCoords);
+ myCoords->decrRef();
+ mcIdType conn[4]={0,1,3,2};
+ ret->allocateCells(1);
+ ret->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
+ ret->finishInsertingCells();
+ return ret.retn();
+}
+
+/*
+ * Generate a 2D mesh that is supposed to match the part that will be loaded by proc1 in testParallelLoad6
+ */
+MEDCouplingUMesh* genPartialLocMeshMultipleTypes2()
+{
+ MCAuto<MEDCouplingUMesh> ret= MEDCouplingUMesh::New("mesh",2);
+ double coords[6] = {0.,1., 1.,1., 1.,2.};
+ DataArrayDouble *myCoords=DataArrayDouble::New();
+ myCoords->alloc(3,2);
+ std::copy(coords,coords+6,myCoords->getPointer());
+ ret->setCoords(myCoords);
+ myCoords->decrRef();
+ mcIdType conn[3]={0,1,2};
+ ret->allocateCells(1);
+ ret->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
+ ret->finishInsertingCells();
+ return ret.retn();
+}
+
+
+/*
+ * Generate a 2D mesh that is supposed to match the part that will be loaded by proc2 in testParallelLoad6
+ */
+MEDCouplingUMesh* genPartialLocMeshMultipleTypes3()
+{
+ MCAuto<MEDCouplingUMesh> ret= MEDCouplingUMesh::New("mesh",2);
+ double coords[12] = {1.,0., 2.,0., 1.,1., 2.,1., 1.,2., 2.,2. };
+ DataArrayDouble *myCoords=DataArrayDouble::New();
+ myCoords->alloc(6,2);
+ std::copy(coords,coords+12,myCoords->getPointer());
+ ret->setCoords(myCoords);
+ myCoords->decrRef();
+ mcIdType conn[7]={0,1,3, 2,3,5,4};
+ ret->allocateCells(2);
+ ret->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
+ ret->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+3);
+ ret->finishInsertingCells();
+ return ret.retn();
+}
+
+
+
/*
* Generate a 2D field that is supposed to match the local field loaded by each proc in testParallelLoad4
*/
distrib= { {INTERP_KERNEL::NORM_TRI3,{0,1}} , {INTERP_KERNEL::NORM_QUAD4,{1,3}} };
std::string filename=INTERP_TEST::getResourceFile("Test2DMultiGeoType.med");
+ //partial loading
MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
MCAuto<MEDCouplingUMesh> mesh = mu->getMeshAtLevel(0);
+
MEDCouplingUMesh *meshRef;
if(rank==0)
meshRef=genLocMeshMultipleTypes1();
int allEqual = -1;
MPI_Allreduce(&equal, &allEqual, 1, MPI_INT,MPI_SUM,MPI_COMM_WORLD);
CPPUNIT_ASSERT(allEqual==3);
-
meshRef->decrRef();
+
MPI_Barrier(MPI_COMM_WORLD);
}
}
+/*!
+ * Test case to load a 2D mesh with multiple geometric types in parallel on 3 procs.
+ * Some procs may have empty partition to load for some types
+ */
+void ParaMEDMEMTest::testParallelLoad6()
+{
+ int size;
+ int rank;
+ MPI_Comm_size(MPI_COMM_WORLD,&size);
+ MPI_Comm_rank(MPI_COMM_WORLD,&rank);
+ //
+ if(size!=3)
+ return ;
+
+ std::map<INTERP_KERNEL::NormalizedCellType, std::vector<mcIdType>> distrib;
+ // independant numerotation for each geometric type!
+ if (rank == 0)
+ distrib = { {INTERP_KERNEL::NORM_TRI3,{}} , {INTERP_KERNEL::NORM_QUAD4,{2}} };
+ else if(rank == 1)
+ distrib = { {INTERP_KERNEL::NORM_TRI3,{2}} , {INTERP_KERNEL::NORM_QUAD4,{}} };
+ else
+ distrib= { {INTERP_KERNEL::NORM_TRI3,{0}} , {INTERP_KERNEL::NORM_QUAD4,{1}} };
+
+ std::string filename=INTERP_TEST::getResourceFile("Test2DMultiGeoType.med");
+ MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
+ MCAuto<MEDCouplingUMesh> mesh = mu->getMeshAtLevel(0);
+
+ MEDCouplingUMesh *meshRef;
+ if(rank==0)
+ meshRef=genPartialLocMeshMultipleTypes1();
+ else if(rank==1)
+ meshRef=genPartialLocMeshMultipleTypes2();
+ else
+ meshRef=genPartialLocMeshMultipleTypes3();
+ //checking that all 3 procs have correctly loaded their part
+ int equal = (int)mesh->isEqual(meshRef,1e-12);
+ int allEqual = -1;
+ MPI_Allreduce(&equal, &allEqual, 1, MPI_INT,MPI_SUM,MPI_COMM_WORLD);
+ CPPUNIT_ASSERT(allEqual==3);
+ meshRef->decrRef();
+
+ MPI_Barrier(MPI_COMM_WORLD);
+}
+