Salome HOME
Various fixes for test runs - introducing MEDCOUPLING_RESOURCE_DIR env variable
[tools/medcoupling.git] / src / ParaMEDMEM / ParaUMesh.cxx
1 //
2 // Copyright (C) 2020  CEA/DEN, EDF R&D
3 //
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.
8 //
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.
13 //
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
17 //
18 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 //
20 // Author : Anthony Geay (EDF R&D)
21
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"
29
30 #include "mpi.h"
31
32 #include <fstream>
33 #include <sstream>
34 #include <numeric>
35 #include <memory>
36 #include <vector>
37
38 using namespace MEDCoupling;
39
40 ParaUMesh *ParaUMesh::New(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
41 {
42   return new ParaUMesh(mesh,globalCellIds,globalNodeIds);
43 }
44
45 ParaUMesh::ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
46 {
47   _mesh.takeRef(mesh);
48   _cell_global.takeRef(globalCellIds);
49   _node_global.takeRef(globalNodeIds);
50   _mesh.checkNotNull();
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.");
58 }
59
60 std::size_t ParaUMesh::getHeapMemorySizeWithoutChildren() const
61 {
62   return 0;
63 }
64
65 std::vector<const BigMemoryObject *> ParaUMesh::getDirectChildrenWithNull() const
66 {
67   return {_mesh,_cell_global,_node_global};
68 }
69
70 /*!
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.
73 */
74 MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
75 {
76   if(fullyIn)
77     return this->getCellIdsLyingOnNodesTrue(globalNodeIds);
78   else
79     return this->getCellIdsLyingOnNodesFalse(globalNodeIds);
80 }
81
82 MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodesTrue(const DataArrayIdType *globalNodeIds) const
83 {
84   MPI_Comm comm(MPI_COMM_WORLD);
85   CommInterface ci;
86   int size;
87   ci.commSize(comm,&size);
88   std::unique_ptr<mcIdType[]> nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[size]),nbOfElems3(new mcIdType[size]);
89   mcIdType nbOfNodeIdsLoc(globalNodeIds->getNumberOfTuples());
90   ci.allGather(&nbOfNodeIdsLoc,1,MPI_ID_TYPE,nbOfElems.get(),1,MPI_ID_TYPE,comm);
91   //store for each proc the local nodeids intercepted by current proc
92   std::vector< MCAuto<DataArrayIdType> > tabs(size);
93   // loop to avoid to all procs to have all the nodes per proc
94   for(int subDiv = 0 ; subDiv < size ; ++subDiv)
95   {
96     std::unique_ptr<mcIdType[]> nbOfElemsSp(CommInterface::SplitArrayOfLength(nbOfElems,size,subDiv,size));
97     mcIdType nbOfNodeIdsSum(std::accumulate(nbOfElemsSp.get(),nbOfElemsSp.get()+size,0));
98     std::unique_ptr<mcIdType[]> allGlobalNodeIds(new mcIdType[nbOfNodeIdsSum]);
99     std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElemsSp,size) );
100     std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) );
101     mcIdType startGlobalNodeIds,endGlobalNodeIds;
102     DataArray::GetSlice(0,globalNodeIds->getNumberOfTuples(),1,subDiv,size,startGlobalNodeIds,endGlobalNodeIds);
103     ci.allGatherV(globalNodeIds->begin()+startGlobalNodeIds,endGlobalNodeIds-startGlobalNodeIds,MPI_ID_TYPE,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
104     mcIdType offset(0);
105     for(int curRk = 0 ; curRk < size ; ++curRk)
106     {
107       MCAuto<DataArrayIdType> globalNodeIdsOfCurProc(DataArrayIdType::New());
108       globalNodeIdsOfCurProc->useArray(allGlobalNodeIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElemsSp[curRk],1);
109       offset += nbOfElemsSp[curRk];
110       MCAuto<DataArrayIdType> globalNodeIdsCaptured(_node_global->buildIntersection(globalNodeIdsOfCurProc));
111       MCAuto<DataArrayIdType> localNodeIdsToLocate(_node_global->findIdForEach(globalNodeIdsCaptured->begin(),globalNodeIdsCaptured->end()));
112       if(tabs[curRk].isNull())
113         tabs[curRk] = localNodeIdsToLocate;
114       else
115         tabs[curRk]->insertAtTheEnd(localNodeIdsToLocate->begin(),localNodeIdsToLocate->end());
116     }
117   }
118
119   for(int curRk = 0 ; curRk < size ; ++curRk)
120   {
121     MCAuto<DataArrayIdType> localNodeIds(tabs[curRk]);
122     localNodeIds->sort();
123     MCAuto<DataArrayIdType> localNodeIdsUnique(localNodeIds->buildUnique());
124     MCAuto<DataArrayIdType> localCellCaptured(_mesh->getCellIdsLyingOnNodes(localNodeIdsUnique->begin(),localNodeIdsUnique->end(),true));
125     MCAuto<DataArrayIdType> localCellCapturedGlob(_cell_global->selectByTupleIdSafe(localCellCaptured->begin(),localCellCaptured->end()));
126     tabs[curRk] = localCellCapturedGlob;
127   }
128   
129   for(int curRk = 0 ; curRk < size ; ++curRk)
130   {
131     tabs[curRk] = tabs[curRk]->buildUniqueNotSorted();
132     nbOfElems3[curRk] = tabs[curRk]->getNumberOfTuples();
133   }
134   std::vector<const DataArrayIdType *> tabss(tabs.begin(),tabs.end());
135   MCAuto<DataArrayIdType> cells(DataArrayIdType::Aggregate(tabss));
136   ci.allToAll(nbOfElems3.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm);
137   mcIdType nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0));
138   MCAuto<DataArrayIdType> cellIdsFromProcs(DataArrayIdType::New());
139   cellIdsFromProcs->alloc(nbOfCellIdsSum,1);
140   {
141     std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElems3,size) ),nbOfElemsOutInt( CommInterface::ToIntArray<mcIdType>(nbOfElems2,size) );
142     std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(nbOfElemsOutInt,size) );
143     ci.allToAllV(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,
144                  cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_ID_TYPE,comm);
145   }
146   cellIdsFromProcs->sort();
147   return cellIdsFromProcs;
148 }
149
150 MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodesFalse(const DataArrayIdType *globalNodeIds) const
151 {
152   MPI_Comm comm(MPI_COMM_WORLD);
153   CommInterface ci;
154   int size;
155   ci.commSize(comm,&size);
156   std::unique_ptr<mcIdType[]> nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[size]),nbOfElems3(new mcIdType[size]);
157   mcIdType nbOfNodeIdsLoc(globalNodeIds->getNumberOfTuples());
158   ci.allGather(&nbOfNodeIdsLoc,1,MPI_ID_TYPE,nbOfElems.get(),1,MPI_ID_TYPE,comm);
159   std::vector< MCAuto<DataArrayIdType> > tabs(size);
160   // loop to avoid to all procs to have all the nodes per proc
161   for(int subDiv = 0 ; subDiv < size ; ++subDiv)
162   {
163     std::unique_ptr<mcIdType[]> nbOfElemsSp(CommInterface::SplitArrayOfLength(nbOfElems,size,subDiv,size));
164     mcIdType nbOfNodeIdsSum(std::accumulate(nbOfElemsSp.get(),nbOfElemsSp.get()+size,0));
165     std::unique_ptr<mcIdType[]> allGlobalNodeIds(new mcIdType[nbOfNodeIdsSum]);
166     std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElemsSp,size) );
167     std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) );
168     mcIdType startGlobalNodeIds,endGlobalNodeIds;
169     DataArray::GetSlice(0,globalNodeIds->getNumberOfTuples(),1,subDiv,size,startGlobalNodeIds,endGlobalNodeIds);
170     ci.allGatherV(globalNodeIds->begin()+startGlobalNodeIds,endGlobalNodeIds-startGlobalNodeIds,MPI_ID_TYPE,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
171     mcIdType offset(0);
172     for(int curRk = 0 ; curRk < size ; ++curRk)
173     {
174       MCAuto<DataArrayIdType> globalNodeIdsOfCurProc(DataArrayIdType::New());
175       globalNodeIdsOfCurProc->useArray(allGlobalNodeIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElemsSp[curRk],1);
176       offset += nbOfElemsSp[curRk];
177       MCAuto<DataArrayIdType> globalNodeIdsCaptured(_node_global->buildIntersection(globalNodeIdsOfCurProc));
178       MCAuto<DataArrayIdType> localNodeIdsToLocate(_node_global->findIdForEach(globalNodeIdsCaptured->begin(),globalNodeIdsCaptured->end()));
179       MCAuto<DataArrayIdType> localCellCaptured(_mesh->getCellIdsLyingOnNodes(localNodeIdsToLocate->begin(),localNodeIdsToLocate->end(),false));
180       MCAuto<DataArrayIdType> localCellCapturedGlob(_cell_global->selectByTupleIdSafe(localCellCaptured->begin(),localCellCaptured->end()));
181       if(tabs[curRk].isNull())
182         tabs[curRk] = localCellCapturedGlob;
183       else
184         tabs[curRk]->insertAtTheEnd(localCellCapturedGlob->begin(),localCellCapturedGlob->end());
185     }
186   }
187   for(int curRk = 0 ; curRk < size ; ++curRk)
188   {
189     tabs[curRk] = tabs[curRk]->buildUniqueNotSorted();
190     nbOfElems3[curRk] = tabs[curRk]->getNumberOfTuples();
191   }
192   std::vector<const DataArrayIdType *> tabss(tabs.begin(),tabs.end());
193   MCAuto<DataArrayIdType> cells(DataArrayIdType::Aggregate(tabss));
194   ci.allToAll(nbOfElems3.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm);
195   mcIdType nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0));
196   MCAuto<DataArrayIdType> cellIdsFromProcs(DataArrayIdType::New());
197   cellIdsFromProcs->alloc(nbOfCellIdsSum,1);
198   {
199     std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElems3,size) ),nbOfElemsOutInt( CommInterface::ToIntArray<mcIdType>(nbOfElems2,size) );
200     std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(nbOfElemsOutInt,size) );
201     ci.allToAllV(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,
202                  cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_ID_TYPE,comm);
203   }
204   cellIdsFromProcs->sort();
205   return cellIdsFromProcs;
206 }
207
208 DataArrayIdType *ParaUMesh::redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const
209 {
210   return this->redistributeCellFieldT<mcIdType>(globalCellIds,fieldValueToRed);
211 }
212
213 DataArrayDouble *ParaUMesh::redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const
214 {
215   return this->redistributeCellFieldT<double>(globalCellIds,fieldValueToRed);
216 }
217
218 DataArrayIdType *ParaUMesh::redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const
219 {
220   return this->redistributeNodeFieldT<mcIdType>(globalCellIds,fieldValueToRed);
221 }
222
223 DataArrayDouble *ParaUMesh::redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const
224 {
225   return this->redistributeNodeFieldT<double>(globalCellIds,fieldValueToRed);
226 }
227
228 /*!
229  * Return part of \a this mesh split over COMM_WORLD. Part is defined by global cell ids array \a globaCellIds.
230  */
231 ParaUMesh *ParaUMesh::redistributeCells(const DataArrayIdType *globalCellIds) const
232 {
233   MPI_Comm comm(MPI_COMM_WORLD);
234   CommInterface ci;
235   std::unique_ptr<mcIdType[]> allGlobalCellIds,allGlobalCellIdsIndex;
236   int size(ci.allGatherArrays(comm,globalCellIds,allGlobalCellIds,allGlobalCellIdsIndex));
237   // Prepare ParaUMesh parts to be sent : compute for each proc the contribution of current rank.
238   std::vector< MCAuto<DataArrayIdType> > globalCellIdsToBeSent(size),globalNodeIdsToBeSent(size);
239   std::vector< MCAuto<MEDCouplingUMesh> > meshPartsToBeSent(size);
240   for(int curRk = 0 ; curRk < size ; ++curRk)
241   {
242     mcIdType offset(allGlobalCellIdsIndex[curRk]);
243     MCAuto<DataArrayIdType> globalCellIdsOfCurProc(DataArrayIdType::New());
244     globalCellIdsOfCurProc->useArray(allGlobalCellIds.get()+offset,false,DeallocType::CPP_DEALLOC,allGlobalCellIdsIndex[curRk+1]-offset,1);
245     // the key call is here : compute for rank curRk the cells to be sent
246     MCAuto<DataArrayIdType> globalCellIdsCaptured(_cell_global->buildIntersection(globalCellIdsOfCurProc));// OK for the global cellIds
247     MCAuto<DataArrayIdType> localCellIdsCaptured(_cell_global->findIdForEach(globalCellIdsCaptured->begin(),globalCellIdsCaptured->end()));
248     MCAuto<MEDCouplingUMesh> meshPart(_mesh->buildPartOfMySelf(localCellIdsCaptured->begin(),localCellIdsCaptured->end(),true));
249     MCAuto<DataArrayIdType> o2n(meshPart->zipCoordsTraducer());// OK for the mesh
250     MCAuto<DataArrayIdType> n2o(o2n->invertArrayO2N2N2O(meshPart->getNumberOfNodes()));
251     MCAuto<DataArrayIdType> globalNodeIdsPart(_node_global->selectByTupleIdSafe(n2o->begin(),n2o->end())); // OK for the global nodeIds
252     meshPartsToBeSent[curRk] = meshPart;
253     globalCellIdsToBeSent[curRk] = globalCellIdsCaptured;
254     globalNodeIdsToBeSent[curRk] = globalNodeIdsPart;
255   }
256   // Receive
257   std::vector< MCAuto<DataArrayIdType> > globalCellIdsReceived,globalNodeIdsReceived;
258   ci.allToAllArrays(comm,globalCellIdsToBeSent,globalCellIdsReceived);
259   ci.allToAllArrays(comm,globalNodeIdsToBeSent,globalNodeIdsReceived);
260   //now exchange the 3 arrays for the umesh : connectivity, connectivityindex and coordinates
261   std::vector<const MEDCouplingUMesh *> meshPartsToBeSent2(FromVecAutoToVecOfConst<MEDCouplingUMesh>(meshPartsToBeSent));
262   //connectivityindex
263   std::vector< MCAuto<DataArrayIdType> > connectivityIndexReceived,connectivityReceived;
264   {
265     std::vector<const DataArrayIdType *> connectivityIndexToBeSent(UMeshConnectivityIndexIterator(0,&meshPartsToBeSent2),UMeshConnectivityIndexIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
266     ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayIdType>(connectivityIndexToBeSent),connectivityIndexReceived);
267   }
268   //connectivity
269   {
270     std::vector<const DataArrayIdType *> connectivityToBeSent(UMeshConnectivityIterator(0,&meshPartsToBeSent2),UMeshConnectivityIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
271     ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayIdType>(connectivityToBeSent),connectivityReceived);
272   }
273   //coordinates
274   MCAuto<DataArrayDouble> coords;
275   {
276     std::vector<const DataArrayDouble *> coordsToBeSent(UMeshCoordsIterator(0,&meshPartsToBeSent2),UMeshCoordsIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
277     ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayDouble>(coordsToBeSent),coords);
278   }
279   /////// Sort it all !
280   // firstly deal with nodes.
281   MCAuto<DataArrayIdType> aggregatedNodeIds( DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(globalNodeIdsReceived)) );
282   MCAuto<DataArrayIdType> aggregatedNodeIdsSorted(aggregatedNodeIds->copySorted());
283   MCAuto<DataArrayIdType> nodeIdsIntoAggregatedIds(DataArrayIdType::FindPermutationFromFirstToSecondDuplicate(aggregatedNodeIdsSorted,aggregatedNodeIds));
284   MCAuto<DataArrayIdType> idxOfSameNodeIds(aggregatedNodeIdsSorted->indexOfSameConsecutiveValueGroups());
285   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
286   MCAuto<DataArrayIdType> finalGlobalNodeIds(aggregatedNodeIdsSorted->selectByTupleIdSafe(idxOfSameNodeIds->begin(),idxOfSameNodeIds->end()-1));
287   MCAuto<DataArrayDouble> finalCoords(coords->selectByTupleIdSafe(n2o_nodes->begin(),n2o_nodes->end()));
288   finalCoords->copyStringInfoFrom(*_mesh->getCoords());
289   // secondly renumbering of node ids in connectivityReceived
290   for(int curRk = 0 ; curRk < size ; ++curRk)
291   {
292     auto current(globalNodeIdsReceived[curRk]);
293     MCAuto<DataArrayIdType> aa(finalGlobalNodeIds->findIdForEach(current->begin(),current->end()));
294     // work on connectivityReceived[curRk] with transformWithIndArr but do not forget type of cells that should be excluded !
295     auto connectivityToModify(connectivityReceived[curRk]);
296     auto connectivityIndex(connectivityIndexReceived[curRk]);
297     MCAuto<DataArrayIdType> types(connectivityToModify->selectByTupleIdSafe(connectivityIndex->begin(),connectivityIndex->end()-1));
298     connectivityToModify->setPartOfValuesSimple3(0,connectivityIndex->begin(),connectivityIndex->end()-1,0,1,1);
299     connectivityToModify->transformWithIndArr(aa->begin(),aa->end());
300     connectivityToModify->setPartOfValues3(types,connectivityIndex->begin(),connectivityIndex->end()-1,0,1,1,true);
301   }
302   // thirdly renumber cells
303   MCAuto<DataArrayIdType> aggregatedCellIds( DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(globalCellIdsReceived)) );
304   MCAuto<DataArrayIdType> aggregatedCellIdsSorted(aggregatedCellIds->copySorted());
305   MCAuto<DataArrayIdType> idsIntoAggregatedIds(DataArrayIdType::FindPermutationFromFirstToSecondDuplicate(aggregatedCellIdsSorted,aggregatedCellIds));
306   MCAuto<DataArrayIdType> cellIdsOfSameNodeIds(aggregatedCellIdsSorted->indexOfSameConsecutiveValueGroups());
307   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
308   // TODO : check coordsReceived==globalCellIds
309   MCAuto<DataArrayIdType> connSorted,indicesSorted;
310   {
311     MCAuto<DataArrayIdType> conn(DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(connectivityReceived)));
312     MCAuto<DataArrayIdType> connIndex(DataArrayIdType::AggregateIndexes(FromVecAutoToVecOfConst<DataArrayIdType>(connectivityIndexReceived))); 
313     {
314       DataArrayIdType *indicesSortedTmp(nullptr),*valuesSortedTmp(nullptr);
315       DataArrayIdType::ExtractFromIndexedArrays(n2o_cells->begin(),n2o_cells->end(),conn,connIndex,valuesSortedTmp,indicesSortedTmp);
316       indicesSorted = indicesSortedTmp; connSorted=valuesSortedTmp;
317     }
318   }
319   // finalize all
320   MCAuto<MEDCouplingUMesh> mesh(MEDCouplingUMesh::New(_mesh->getName(),_mesh->getMeshDimension()));
321   mesh->setConnectivity(connSorted,indicesSorted,true);
322   mesh->setCoords(finalCoords);
323   mesh->setDescription(_mesh->getDescription());
324   MCAuto<ParaUMesh> ret(ParaUMesh::New(mesh,aggregatedCellIdsSorted,finalGlobalNodeIds));
325   return ret.retn();
326 }