2 // Copyright (C) 2020 CEA/DEN, EDF R&D
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 // Author : Anthony Geay (EDF R&D)
22 #include "ParaUMesh.hxx"
23 #include "ProcessorGroup.hxx"
24 #include "MPIProcessorGroup.hxx"
25 #include "Topology.hxx"
26 #include "BlockTopology.hxx"
27 #include "CommInterface.hxx"
28 #include "MEDCouplingMemArray.hxx"
38 using namespace MEDCoupling;
40 ParaUMesh *ParaUMesh::New(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
42 return new ParaUMesh(mesh,globalCellIds,globalNodeIds);
45 ParaUMesh::ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
48 _cell_global.takeRef(globalCellIds);
49 _node_global.takeRef(globalNodeIds);
51 _cell_global.checkNotNull();
52 _node_global.checkNotNull();
53 _mesh->checkConsistencyLight();
54 if(_mesh->getNumberOfNodes() != _node_global->getNumberOfTuples())
55 throw INTERP_KERNEL::Exception("ParaUMesh constructor : mismatch between # nodes and len of global # nodes.");
56 if(_mesh->getNumberOfCells() != _cell_global->getNumberOfTuples())
57 throw INTERP_KERNEL::Exception("ParaUMesh constructor : mismatch between # cells and len of global # cells.");
60 std::size_t ParaUMesh::getHeapMemorySizeWithoutChildren() const
65 std::vector<const BigMemoryObject *> ParaUMesh::getDirectChildrenWithNull() const
67 return {_mesh,_cell_global,_node_global};
71 * This method computes the cells part of distributed mesh lying on \a globalNodeIds nodes.
72 * The input \a globalNodeIds are not supposed to reside on the current process.
74 MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
77 throw INTERP_KERNEL::Exception("ParaUMesh::getCellIdsLyingOnNodes : not implemented yet for fullyIn == True !");
78 MPI_Comm comm(MPI_COMM_WORLD);
81 ci.commSize(comm,&size);
82 std::unique_ptr<mcIdType[]> nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[size]),nbOfElems3(new mcIdType[size]);
83 mcIdType nbOfNodeIdsLoc(globalNodeIds->getNumberOfTuples());
84 ci.allGather(&nbOfNodeIdsLoc,1,MPI_ID_TYPE,nbOfElems.get(),1,MPI_ID_TYPE,comm);
85 std::vector< MCAuto<DataArrayIdType> > tabs(size);
86 // loop to avoid to all procs to have all the nodes per proc
87 for(int subDiv = 0 ; subDiv < size ; ++subDiv)
89 std::unique_ptr<mcIdType[]> nbOfElemsSp(CommInterface::SplitArrayOfLength(nbOfElems,size,subDiv,size));
90 mcIdType nbOfNodeIdsSum(std::accumulate(nbOfElemsSp.get(),nbOfElemsSp.get()+size,0));
91 std::unique_ptr<mcIdType[]> allGlobalNodeIds(new mcIdType[nbOfNodeIdsSum]);
92 std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElemsSp,size) );
93 std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) );
94 mcIdType startGlobalNodeIds,endGlobalNodeIds;
95 DataArray::GetSlice(0,globalNodeIds->getNumberOfTuples(),1,subDiv,size,startGlobalNodeIds,endGlobalNodeIds);
96 ci.allGatherV(globalNodeIds->begin()+startGlobalNodeIds,endGlobalNodeIds-startGlobalNodeIds,MPI_ID_TYPE,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
98 for(int curRk = 0 ; curRk < size ; ++curRk)
100 MCAuto<DataArrayIdType> globalNodeIdsOfCurProc(DataArrayIdType::New());
101 globalNodeIdsOfCurProc->useArray(allGlobalNodeIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElemsSp[curRk],1);
102 offset += nbOfElemsSp[curRk];
103 MCAuto<DataArrayIdType> globalNodeIdsCaptured(_node_global->buildIntersection(globalNodeIdsOfCurProc));
104 MCAuto<DataArrayIdType> localNodeIdsToLocate(_node_global->findIdForEach(globalNodeIdsCaptured->begin(),globalNodeIdsCaptured->end()));
105 MCAuto<DataArrayIdType> localCellCaptured(_mesh->getCellIdsLyingOnNodes(localNodeIdsToLocate->begin(),localNodeIdsToLocate->end(),false));
106 MCAuto<DataArrayIdType> localCellCapturedGlob(_cell_global->selectByTupleIdSafe(localCellCaptured->begin(),localCellCaptured->end()));
107 if(tabs[curRk].isNull())
108 tabs[curRk] = localCellCapturedGlob;
110 tabs[curRk]->insertAtTheEnd(localCellCapturedGlob->begin(),localCellCapturedGlob->end());
113 for(int curRk = 0 ; curRk < size ; ++curRk)
115 tabs[curRk] = tabs[curRk]->buildUniqueNotSorted();
116 nbOfElems3[curRk] = tabs[curRk]->getNumberOfTuples();
118 std::vector<const DataArrayIdType *> tabss(tabs.begin(),tabs.end());
119 MCAuto<DataArrayIdType> cells(DataArrayIdType::Aggregate(tabss));
120 ci.allToAll(nbOfElems3.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm);
121 mcIdType nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0));
122 MCAuto<DataArrayIdType> cellIdsFromProcs(DataArrayIdType::New());
123 cellIdsFromProcs->alloc(nbOfCellIdsSum,1);
125 std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElems3,size) ),nbOfElemsOutInt( CommInterface::ToIntArray<mcIdType>(nbOfElems2,size) );
126 std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(nbOfElemsOutInt,size) );
127 ci.allToAllV(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,
128 cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_ID_TYPE,comm);
130 cellIdsFromProcs->sort();
131 return cellIdsFromProcs;
134 DataArrayIdType *ParaUMesh::redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const
136 return this->redistributeCellFieldT<mcIdType>(globalCellIds,fieldValueToRed);
139 DataArrayDouble *ParaUMesh::redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const
141 return this->redistributeCellFieldT<double>(globalCellIds,fieldValueToRed);
144 DataArrayIdType *ParaUMesh::redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const
146 return this->redistributeNodeFieldT<mcIdType>(globalCellIds,fieldValueToRed);
149 DataArrayDouble *ParaUMesh::redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const
151 return this->redistributeNodeFieldT<double>(globalCellIds,fieldValueToRed);
155 * Return part of \a this mesh split over COMM_WORLD. Part is defined by global cell ids array \a globaCellIds.
157 ParaUMesh *ParaUMesh::redistributeCells(const DataArrayIdType *globalCellIds) const
159 MPI_Comm comm(MPI_COMM_WORLD);
161 std::unique_ptr<mcIdType[]> allGlobalCellIds,allGlobalCellIdsIndex;
162 int size(ci.allGatherArrays(comm,globalCellIds,allGlobalCellIds,allGlobalCellIdsIndex));
163 // Prepare ParaUMesh parts to be sent : compute for each proc the contribution of current rank.
164 std::vector< MCAuto<DataArrayIdType> > globalCellIdsToBeSent(size),globalNodeIdsToBeSent(size);
165 std::vector< MCAuto<MEDCouplingUMesh> > meshPartsToBeSent(size);
166 for(int curRk = 0 ; curRk < size ; ++curRk)
168 mcIdType offset(allGlobalCellIdsIndex[curRk]);
169 MCAuto<DataArrayIdType> globalCellIdsOfCurProc(DataArrayIdType::New());
170 globalCellIdsOfCurProc->useArray(allGlobalCellIds.get()+offset,false,DeallocType::CPP_DEALLOC,allGlobalCellIdsIndex[curRk+1]-offset,1);
171 // the key call is here : compute for rank curRk the cells to be sent
172 MCAuto<DataArrayIdType> globalCellIdsCaptured(_cell_global->buildIntersection(globalCellIdsOfCurProc));// OK for the global cellIds
173 MCAuto<DataArrayIdType> localCellIdsCaptured(_cell_global->findIdForEach(globalCellIdsCaptured->begin(),globalCellIdsCaptured->end()));
174 MCAuto<MEDCouplingUMesh> meshPart(_mesh->buildPartOfMySelf(localCellIdsCaptured->begin(),localCellIdsCaptured->end(),true));
175 MCAuto<DataArrayIdType> o2n(meshPart->zipCoordsTraducer());// OK for the mesh
176 MCAuto<DataArrayIdType> n2o(o2n->invertArrayO2N2N2O(meshPart->getNumberOfNodes()));
177 MCAuto<DataArrayIdType> globalNodeIdsPart(_node_global->selectByTupleIdSafe(n2o->begin(),n2o->end())); // OK for the global nodeIds
178 meshPartsToBeSent[curRk] = meshPart;
179 globalCellIdsToBeSent[curRk] = globalCellIdsCaptured;
180 globalNodeIdsToBeSent[curRk] = globalNodeIdsPart;
183 std::vector< MCAuto<DataArrayIdType> > globalCellIdsReceived,globalNodeIdsReceived;
184 ci.allToAllArrays(comm,globalCellIdsToBeSent,globalCellIdsReceived);
185 ci.allToAllArrays(comm,globalNodeIdsToBeSent,globalNodeIdsReceived);
186 //now exchange the 3 arrays for the umesh : connectivity, connectivityindex and coordinates
187 std::vector<const MEDCouplingUMesh *> meshPartsToBeSent2(FromVecAutoToVecOfConst<MEDCouplingUMesh>(meshPartsToBeSent));
189 std::vector< MCAuto<DataArrayIdType> > connectivityIndexReceived,connectivityReceived;
191 std::vector<const DataArrayIdType *> connectivityIndexToBeSent(UMeshConnectivityIndexIterator(0,&meshPartsToBeSent2),UMeshConnectivityIndexIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
192 ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayIdType>(connectivityIndexToBeSent),connectivityIndexReceived);
196 std::vector<const DataArrayIdType *> connectivityToBeSent(UMeshConnectivityIterator(0,&meshPartsToBeSent2),UMeshConnectivityIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
197 ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayIdType>(connectivityToBeSent),connectivityReceived);
200 MCAuto<DataArrayDouble> coords;
202 std::vector<const DataArrayDouble *> coordsToBeSent(UMeshCoordsIterator(0,&meshPartsToBeSent2),UMeshCoordsIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
203 ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayDouble>(coordsToBeSent),coords);
205 /////// Sort it all !
206 // firstly deal with nodes.
207 MCAuto<DataArrayIdType> aggregatedNodeIds( DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(globalNodeIdsReceived)) );
208 MCAuto<DataArrayIdType> aggregatedNodeIdsSorted(aggregatedNodeIds->copySorted());
209 MCAuto<DataArrayIdType> nodeIdsIntoAggregatedIds(DataArrayIdType::FindPermutationFromFirstToSecondDuplicate(aggregatedNodeIdsSorted,aggregatedNodeIds));
210 MCAuto<DataArrayIdType> idxOfSameNodeIds(aggregatedNodeIdsSorted->indexOfSameConsecutiveValueGroups());
211 MCAuto<DataArrayIdType> n2o_nodes(nodeIdsIntoAggregatedIds->selectByTupleIdSafe(idxOfSameNodeIds->begin(),idxOfSameNodeIds->end()-1));//new == new ordering so that global node ids are sorted . old == coarse ordering implied by the aggregation
212 MCAuto<DataArrayIdType> finalGlobalNodeIds(aggregatedNodeIdsSorted->selectByTupleIdSafe(idxOfSameNodeIds->begin(),idxOfSameNodeIds->end()-1));
213 MCAuto<DataArrayDouble> finalCoords(coords->selectByTupleIdSafe(n2o_nodes->begin(),n2o_nodes->end()));
214 finalCoords->copyStringInfoFrom(*_mesh->getCoords());
215 // secondly renumbering of node ids in connectivityReceived
216 for(int curRk = 0 ; curRk < size ; ++curRk)
218 auto current(globalNodeIdsReceived[curRk]);
219 MCAuto<DataArrayIdType> aa(finalGlobalNodeIds->findIdForEach(current->begin(),current->end()));
220 // work on connectivityReceived[curRk] with transformWithIndArr but do not forget type of cells that should be excluded !
221 auto connectivityToModify(connectivityReceived[curRk]);
222 auto connectivityIndex(connectivityIndexReceived[curRk]);
223 MCAuto<DataArrayIdType> types(connectivityToModify->selectByTupleIdSafe(connectivityIndex->begin(),connectivityIndex->end()-1));
224 connectivityToModify->setPartOfValuesSimple3(0,connectivityIndex->begin(),connectivityIndex->end()-1,0,1,1);
225 connectivityToModify->transformWithIndArr(aa->begin(),aa->end());
226 connectivityToModify->setPartOfValues3(types,connectivityIndex->begin(),connectivityIndex->end()-1,0,1,1,true);
228 // thirdly renumber cells
229 MCAuto<DataArrayIdType> aggregatedCellIds( DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(globalCellIdsReceived)) );
230 MCAuto<DataArrayIdType> aggregatedCellIdsSorted(aggregatedCellIds->copySorted());
231 MCAuto<DataArrayIdType> idsIntoAggregatedIds(DataArrayIdType::FindPermutationFromFirstToSecondDuplicate(aggregatedCellIdsSorted,aggregatedCellIds));
232 MCAuto<DataArrayIdType> cellIdsOfSameNodeIds(aggregatedCellIdsSorted->indexOfSameConsecutiveValueGroups());
233 MCAuto<DataArrayIdType> n2o_cells(idsIntoAggregatedIds->selectByTupleIdSafe(cellIdsOfSameNodeIds->begin(),cellIdsOfSameNodeIds->end()-1));//new == new ordering so that global cell ids are sorted . old == coarse ordering implied by the aggregation
234 // TODO : check coordsReceived==globalCellIds
235 MCAuto<DataArrayIdType> connSorted,indicesSorted;
237 MCAuto<DataArrayIdType> conn(DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(connectivityReceived)));
238 MCAuto<DataArrayIdType> connIndex(DataArrayIdType::AggregateIndexes(FromVecAutoToVecOfConst<DataArrayIdType>(connectivityIndexReceived)));
240 DataArrayIdType *indicesSortedTmp(nullptr),*valuesSortedTmp(nullptr);
241 DataArrayIdType::ExtractFromIndexedArrays(n2o_cells->begin(),n2o_cells->end(),conn,connIndex,valuesSortedTmp,indicesSortedTmp);
242 indicesSorted = indicesSortedTmp; connSorted=valuesSortedTmp;
246 MCAuto<MEDCouplingUMesh> mesh(MEDCouplingUMesh::New(_mesh->getName(),_mesh->getMeshDimension()));
247 mesh->setConnectivity(connSorted,indicesSorted,true);
248 mesh->setCoords(finalCoords);
249 mesh->setDescription(_mesh->getDescription());
250 MCAuto<ParaUMesh> ret(ParaUMesh::New(mesh,aggregatedCellIdsSorted,finalGlobalNodeIds));