00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include "DistributedTetrahedralMesh.hpp"
00030
00031 #include <cassert>
00032 #include <sstream>
00033 #include <string>
00034
00035 #include "Exception.hpp"
00036 #include "Element.hpp"
00037 #include "BoundaryElement.hpp"
00038
00039 #include "PetscTools.hpp"
00040 #include "DistributedVectorFactory.hpp"
00041 #include "OutputFileHandler.hpp"
00042
00043 #include "Timer.hpp"
00044
00046
00048
00049 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00050 DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::DistributedTetrahedralMesh(PartitionType metisPartitioning)
00051 :
00052 mTotalNumElements(0u),
00053 mTotalNumBoundaryElements(0u),
00054 mTotalNumNodes(0u),
00055 mMetisPartitioning(metisPartitioning)
00056 {
00057 if (ELEMENT_DIM == 1)
00058 {
00059
00060 mMetisPartitioning=DUMB;
00061 }
00062 }
00063
00064 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00065 DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::~DistributedTetrahedralMesh()
00066 {
00067 for (unsigned i=0; i<this->mHaloNodes.size(); i++)
00068 {
00069 delete this->mHaloNodes[i];
00070 }
00071 }
00072
00073 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00074 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SetDistributedVectorFactory(DistributedVectorFactory* pFactory)
00075 {
00076 AbstractMesh<ELEMENT_DIM,SPACE_DIM>::SetDistributedVectorFactory(pFactory);
00077 mMetisPartitioning = DUMB;
00078 }
00079
00080 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00081 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ComputeMeshPartitioning(
00082 AbstractMeshReader<ELEMENT_DIM, SPACE_DIM>& rMeshReader,
00083 std::set<unsigned>& rNodesOwned,
00084 std::set<unsigned>& rHaloNodesOwned,
00085 std::set<unsigned>& rElementsOwned,
00086 std::vector<unsigned>& rProcessorsOffset)
00087 {
00089
00090 if (mMetisPartitioning==PARMETIS_LIBRARY && !PetscTools::IsSequential())
00091 {
00092
00093
00094
00095
00096 ParMetisLibraryNodePartitioning(rMeshReader, rElementsOwned, rNodesOwned, rHaloNodesOwned, rProcessorsOffset);
00097 }
00098 else
00099 {
00100
00101
00102
00103 if (mMetisPartitioning==METIS_LIBRARY && !PetscTools::IsSequential())
00104 {
00105 MetisLibraryNodePartitioning(rMeshReader, rNodesOwned, rProcessorsOffset);
00106 }
00107 else
00108 {
00109 DumbNodePartitioning(rMeshReader, rNodesOwned);
00110 }
00111
00112 for (unsigned element_number = 0; element_number < mTotalNumElements; element_number++)
00113 {
00114 ElementData element_data = rMeshReader.GetNextElementData();
00115
00116 bool element_owned = false;
00117 std::set<unsigned> temp_halo_nodes;
00118
00119 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
00120 {
00121 if (rNodesOwned.find(element_data.NodeIndices[i]) != rNodesOwned.end())
00122 {
00123 element_owned = true;
00124 rElementsOwned.insert(element_number);
00125 }
00126 else
00127 {
00128 temp_halo_nodes.insert(element_data.NodeIndices[i]);
00129 }
00130 }
00131
00132 if (element_owned)
00133 {
00134 rHaloNodesOwned.insert(temp_halo_nodes.begin(), temp_halo_nodes.end());
00135 }
00136 }
00137 }
00138 rMeshReader.Reset();
00139 }
00140
00141 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00142 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ConstructFromMeshReader(
00143 AbstractMeshReader<ELEMENT_DIM, SPACE_DIM>& rMeshReader)
00144 {
00145 std::set<unsigned> nodes_owned;
00146 std::set<unsigned> halo_nodes_owned;
00147 std::set<unsigned> elements_owned;
00148 std::vector<unsigned> proc_offsets;
00149
00150 this->mMeshFileBaseName = rMeshReader.GetMeshFileBaseName();
00151 mTotalNumElements = rMeshReader.GetNumElements();
00152 mTotalNumBoundaryElements = rMeshReader.GetNumFaces();
00153 mTotalNumNodes = rMeshReader.GetNumNodes();
00154
00155
00156 PetscTools::Barrier();
00157 Timer::Reset();
00158 ComputeMeshPartitioning(rMeshReader, nodes_owned, halo_nodes_owned, elements_owned, proc_offsets);
00159 PetscTools::Barrier();
00160
00161
00162
00163 this->mElements.reserve(elements_owned.size());
00164 this->mNodes.reserve(nodes_owned.size());
00165
00166 if ( rMeshReader.IsFileFormatBinary() )
00167 {
00168 std::vector<double> coords;
00169
00170 for (std::set<unsigned>::const_iterator it=nodes_owned.begin(); it!=nodes_owned.end(); it++)
00171 {
00172
00173 unsigned global_node_index=*it;
00174 coords = rMeshReader.GetNode(global_node_index);
00175 RegisterNode(global_node_index);
00176 this->mNodes.push_back(new Node<SPACE_DIM>(global_node_index, coords, false));
00177 }
00178 for (std::set<unsigned>::const_iterator it=halo_nodes_owned.begin(); it!=halo_nodes_owned.end(); it++)
00179 {
00180
00181 unsigned global_node_index=*it;
00182 coords = rMeshReader.GetNode(global_node_index);
00183 RegisterHaloNode(global_node_index);
00184 mHaloNodes.push_back(new Node<SPACE_DIM>(global_node_index, coords, false));
00185 }
00186 }
00187 else
00188 {
00189
00190 for (unsigned node_index=0; node_index < mTotalNumNodes; node_index++)
00191 {
00192 std::vector<double> coords;
00194 coords = rMeshReader.GetNextNode();
00195
00196
00197 if (nodes_owned.find(node_index) != nodes_owned.end())
00198 {
00199 RegisterNode(node_index);
00200 this->mNodes.push_back(new Node<SPACE_DIM>(node_index, coords, false));
00201 }
00202
00203
00204 if (halo_nodes_owned.find(node_index) != halo_nodes_owned.end())
00205 {
00206 RegisterHaloNode(node_index);
00207 mHaloNodes.push_back(new Node<SPACE_DIM>(node_index, coords, false));
00208 }
00209 }
00210 }
00211
00212 if ( rMeshReader.IsFileFormatBinary() )
00213 {
00214
00215 for (std::set<unsigned>::const_iterator elem_it=elements_owned.begin(); elem_it!=elements_owned.end(); elem_it++)
00216 {
00217 unsigned global_element_index=*elem_it;
00218 ElementData element_data;
00219 element_data = rMeshReader.GetElementData(global_element_index);
00220
00221 std::vector<Node<SPACE_DIM>*> nodes;
00222 for (unsigned j=0; j<ELEMENT_DIM+1; j++)
00223 {
00224
00225 nodes.push_back(this->GetAnyNode(element_data.NodeIndices[j]));
00226 }
00227
00228 RegisterElement(global_element_index);
00229 Element<ELEMENT_DIM,SPACE_DIM>* p_element = new Element<ELEMENT_DIM,SPACE_DIM>(global_element_index, nodes);
00230 this->mElements.push_back(p_element);
00231
00232 if (rMeshReader.GetNumElementAttributes() > 0)
00233 {
00234 assert(rMeshReader.GetNumElementAttributes() == 1);
00235 unsigned attribute_value = element_data.AttributeValue;
00236 p_element->SetRegion(attribute_value);
00237 }
00238 }
00239 }
00240 else
00241 {
00242
00243 for (unsigned element_index=0; element_index < mTotalNumElements; element_index++)
00244 {
00245 ElementData element_data;
00246
00247 element_data = rMeshReader.GetNextElementData();
00248
00249
00250 if (elements_owned.find(element_index) != elements_owned.end())
00251 {
00252 std::vector<Node<SPACE_DIM>*> nodes;
00253 for (unsigned j=0; j<ELEMENT_DIM+1; j++)
00254 {
00255
00256 nodes.push_back(this->GetAnyNode(element_data.NodeIndices[j]));
00257 }
00258
00259 RegisterElement(element_index);
00260
00261 Element<ELEMENT_DIM,SPACE_DIM>* p_element = new Element<ELEMENT_DIM,SPACE_DIM>(element_index, nodes);
00262 this->mElements.push_back(p_element);
00263
00264 if (rMeshReader.GetNumElementAttributes() > 0)
00265 {
00266 assert(rMeshReader.GetNumElementAttributes() == 1);
00267 unsigned attribute_value = element_data.AttributeValue;
00268 p_element->SetRegion(attribute_value);
00269 }
00270 }
00271 }
00272 }
00273
00274
00275 try
00276 {
00277 for (unsigned face_index=0; face_index<mTotalNumBoundaryElements; face_index++)
00278 {
00279 ElementData face_data = rMeshReader.GetNextFaceData();
00280 std::vector<unsigned> node_indices = face_data.NodeIndices;
00281
00282 bool own = false;
00283
00284 for (unsigned node_index=0; node_index<node_indices.size(); node_index++)
00285 {
00286
00287 if (mNodesMapping.find(node_indices[node_index]) != mNodesMapping.end())
00288 {
00289 own = true;
00290 break;
00291 }
00292 }
00293
00294 if (!own)
00295 {
00296 continue;
00297 }
00298
00299
00300 std::set<unsigned> containing_element_indices;
00301 std::vector<Node<SPACE_DIM>*> nodes;
00302
00303 for (unsigned node_index=0; node_index<node_indices.size(); node_index++)
00304 {
00305
00306
00307 try
00308 {
00309 nodes.push_back(this->GetAnyNode(node_indices[node_index]));
00310 }
00311 catch (Exception &e)
00312 {
00313 std::stringstream message;
00314 message << "Face does not appear in element file (Face " << face_index << " in "<<this->mMeshFileBaseName<< ")";
00315 EXCEPTION(message.str().c_str());
00316 }
00317 }
00318
00319
00320
00321 for (unsigned j=0; j<nodes.size(); j++)
00322 {
00323 if (!nodes[j]->IsBoundaryNode())
00324 {
00325 nodes[j]->SetAsBoundaryNode();
00326 this->mBoundaryNodes.push_back(nodes[j]);
00327 }
00328
00329 nodes[j]->AddBoundaryElement(face_index);
00330 }
00331
00332 RegisterBoundaryElement(face_index);
00333 BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* p_boundary_element = new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(face_index, nodes);
00334 this->mBoundaryElements.push_back(p_boundary_element);
00335
00336 if (rMeshReader.GetNumFaceAttributes() > 0)
00337 {
00338 assert(rMeshReader.GetNumFaceAttributes() == 1);
00339 unsigned attribute_value = face_data.AttributeValue;
00340 p_boundary_element->SetRegion(attribute_value);
00341 }
00342 }
00343 }
00344 catch (Exception &e)
00345 {
00346 PetscTools::ReplicateException(true);
00347 throw e;
00348 }
00349 PetscTools::ReplicateException(false);
00350
00351 if (mMetisPartitioning != DUMB && !PetscTools::IsSequential())
00352 {
00353 assert(this->mNodesPermutation.size() != 0);
00354
00355 ReorderNodes();
00356
00357 unsigned num_owned;
00358 unsigned rank = PetscTools::GetMyRank();
00359 if ( !PetscTools::AmTopMost() )
00360 {
00361 num_owned = proc_offsets[rank+1]-proc_offsets[rank];
00362 }
00363 else
00364 {
00365 num_owned = mTotalNumNodes - proc_offsets[rank];
00366 }
00367
00368 assert(!this->mpDistributedVectorFactory);
00369 this->mpDistributedVectorFactory = new DistributedVectorFactory(this->GetNumNodes(), num_owned);
00370 }
00371 else
00372 {
00373
00374 assert(this->mpDistributedVectorFactory);
00375 }
00376 rMeshReader.Reset();
00377 }
00378
00379 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00380 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumLocalNodes() const
00381 {
00382 return this->mNodes.size();
00383 }
00384
00385 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00386 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumLocalElements() const
00387 {
00388 return this->mElements.size();
00389 }
00390
00391 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00392 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumLocalBoundaryElements() const
00393 {
00394 return this->mBoundaryElements.size();
00395 }
00396
00397 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00398 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumNodes() const
00399 {
00400 return mTotalNumNodes;
00401 }
00402
00403 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00404 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumAllNodes() const
00405 {
00406 return mTotalNumNodes;
00407 }
00408
00409 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00410 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumElements() const
00411 {
00412 return mTotalNumElements;
00413 }
00414
00415 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00416 typename DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::PartitionType DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetPartitionType() const
00417 {
00418 return mMetisPartitioning;
00419 }
00420
00421 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00422 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumBoundaryElements() const
00423 {
00424 return mTotalNumBoundaryElements;
00425 }
00426
00427 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00428 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetHaloNodeIndices(std::vector<unsigned>& rHaloIndices) const
00429 {
00430
00431 rHaloIndices.clear();
00432 for (unsigned i=0; i<mHaloNodes.size(); i++)
00433 {
00434 rHaloIndices.push_back(mHaloNodes[i]->GetIndex());
00435 }
00436 }
00437
00438 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00439 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SetElementOwnerships(unsigned lo, unsigned hi)
00440 {
00441
00442
00443
00444
00445 assert(hi>=lo);
00446 for (unsigned element_index=0; element_index<this->mElements.size(); element_index++)
00447 {
00448 Element<ELEMENT_DIM, SPACE_DIM>* p_element=this->mElements[element_index];
00449 p_element->SetOwnership(true);
00450 }
00451 }
00452
00453 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00454 bool DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::CalculateDesignatedOwnershipOfElement( unsigned elementIndex )
00455 {
00456 try
00457 {
00458 return(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::CalculateDesignatedOwnershipOfElement(elementIndex));
00459 }
00460 catch(Exception& e)
00461 {
00462 return false;
00463 }
00464 }
00465
00466 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00467 bool DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::CalculateDesignatedOwnershipOfBoundaryElement( unsigned faceIndex )
00468 {
00469 try
00470 {
00471 return(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::CalculateDesignatedOwnershipOfBoundaryElement(faceIndex));
00472 }
00473 catch(Exception& e)
00474 {
00475 return false;
00476 }
00477 }
00478
00479 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00480 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::RegisterNode(unsigned index)
00481 {
00482 mNodesMapping[index] = this->mNodes.size();
00483 }
00484
00485 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00486 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::RegisterHaloNode(unsigned index)
00487 {
00488 mHaloNodesMapping[index] = mHaloNodes.size();
00489 }
00490
00491 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00492 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::RegisterElement(unsigned index)
00493 {
00494 mElementsMapping[index] = this->mElements.size();
00495 }
00496
00497 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00498 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::RegisterBoundaryElement(unsigned index)
00499 {
00500 mBoundaryElementsMapping[index] = this->mBoundaryElements.size();
00501 }
00502
00503 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00504 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SolveNodeMapping(unsigned index) const
00505 {
00506 std::map<unsigned, unsigned>::const_iterator node_position = mNodesMapping.find(index);
00507
00508 if (node_position == mNodesMapping.end())
00509 {
00510 std::stringstream message;
00511 message << "Requested node " << index << " does not belong to processor " << PetscTools::GetMyRank();
00512 EXCEPTION(message.str().c_str());
00513 }
00514 return node_position->second;
00515 }
00516
00517
00518
00519
00520
00521
00522
00523
00524 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00525 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SolveElementMapping(unsigned index) const
00526 {
00527 std::map<unsigned, unsigned>::const_iterator element_position = mElementsMapping.find(index);
00528
00529 if (element_position == mElementsMapping.end())
00530 {
00531 std::stringstream message;
00532 message << "Requested element " << index << " does not belong to processor " << PetscTools::GetMyRank();
00533 EXCEPTION(message.str().c_str());
00534 }
00535
00536 return element_position->second;
00537 }
00538
00539 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00540 unsigned DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SolveBoundaryElementMapping(unsigned index) const
00541 {
00542 std::map<unsigned, unsigned>::const_iterator boundary_element_position = mBoundaryElementsMapping.find(index);
00543
00544 if (boundary_element_position == mBoundaryElementsMapping.end())
00545 {
00546 std::stringstream message;
00547 message << "Requested boundary element " << index << " does not belong to processor " << PetscTools::GetMyRank();
00548 EXCEPTION(message.str().c_str());
00549 }
00550
00551 return boundary_element_position->second;
00552 }
00553
00554 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00555 Node<SPACE_DIM> * DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetAnyNode(unsigned index) const
00556 {
00557 std::map<unsigned, unsigned>::const_iterator node_position;
00558
00559 if ((node_position=mHaloNodesMapping.find(index)) != mHaloNodesMapping.end())
00560 {
00561 return mHaloNodes[node_position->second];
00562 }
00563
00564 if ((node_position=mNodesMapping.find(index)) != mNodesMapping.end())
00565 {
00566
00567 return this->mNodes[node_position->second];
00568 }
00569
00570 std::stringstream message;
00571 message << "Requested node/halo " << index << " does not belong to processor " << PetscTools::GetMyRank();
00572 EXCEPTION(message.str().c_str());
00573 }
00574
00575 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00576 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::DumbNodePartitioning(AbstractMeshReader<ELEMENT_DIM, SPACE_DIM>& rMeshReader,
00577 std::set<unsigned>& rNodesOwned)
00578 {
00579 if (this->mpDistributedVectorFactory)
00580 {
00581
00582 if (this->mpDistributedVectorFactory->GetProblemSize() != mTotalNumNodes)
00583 {
00584
00585 this->mpDistributedVectorFactory = NULL;
00586 this->mTotalNumNodes = 0u;
00587 this->mTotalNumElements = 0u;
00588 this->mTotalNumBoundaryElements = 0u;
00589 EXCEPTION("The distributed vector factory size in the mesh doesn't match the total number of nodes.");
00590 }
00591 }
00592 else
00593 {
00594 this->mpDistributedVectorFactory = new DistributedVectorFactory(mTotalNumNodes);
00595 }
00596 for (unsigned node_index = this->mpDistributedVectorFactory->GetLow();
00597 node_index < this->mpDistributedVectorFactory->GetHigh();
00598 node_index++)
00599 {
00600 rNodesOwned.insert(node_index);
00601 }
00602 }
00603
00604
00605 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00606 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::MetisLibraryNodePartitioning(AbstractMeshReader<ELEMENT_DIM, SPACE_DIM>& rMeshReader,
00607 std::set<unsigned>& rNodesOwned,
00608 std::vector<unsigned>& rProcessorsOffset)
00609 {
00610 assert(!PetscTools::IsSequential());
00611
00612 assert(ELEMENT_DIM==2 || ELEMENT_DIM==3);
00613
00614 int nn = rMeshReader.GetNumNodes();
00615 idxtype* npart = new idxtype[nn];
00616 assert(npart != NULL);
00617
00618
00619 if (PetscTools::AmMaster())
00620 {
00621 int ne = rMeshReader.GetNumElements();
00622 idxtype* elmnts = new idxtype[ne * (ELEMENT_DIM+1)];
00623 assert(elmnts != NULL);
00624
00625 unsigned counter=0;
00626 for (unsigned element_number = 0; element_number < mTotalNumElements; element_number++)
00627 {
00628 ElementData element_data = rMeshReader.GetNextElementData();
00629
00630 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
00631 {
00632 elmnts[counter++] = element_data.NodeIndices[i];
00633 }
00634 }
00635 rMeshReader.Reset();
00636
00637 int etype;
00638
00639 switch (ELEMENT_DIM)
00640 {
00641 case 2:
00642 etype = 1;
00643 break;
00644 case 3:
00645 etype = 2;
00646 break;
00647 default:
00648 NEVER_REACHED;
00649 }
00650
00651 int numflag = 0;
00652 int nparts = PetscTools::GetNumProcs();
00653 int edgecut;
00654 idxtype* epart = new idxtype[ne];
00655 assert(epart != NULL);
00656
00657 Timer::Reset();
00658 METIS_PartMeshNodal(&ne, &nn, elmnts, &etype, &numflag, &nparts, &edgecut, epart, npart);
00659
00660
00661 delete[] elmnts;
00662 delete[] epart;
00663 }
00664
00665
00666
00667 assert(sizeof(idxtype) == sizeof(int));
00668 MPI_Bcast(npart , nn , MPI_INT, 0 , PETSC_COMM_WORLD);
00669
00670 assert(rProcessorsOffset.size() == 0);
00671 rProcessorsOffset.resize(PetscTools::GetNumProcs(), 0);
00672
00673 for (unsigned node_index=0; node_index<this->GetNumNodes(); node_index++)
00674 {
00675 unsigned part_read = npart[node_index];
00676
00677
00678 if (part_read == PetscTools::GetMyRank())
00679 {
00680 rNodesOwned.insert(node_index);
00681 }
00682
00683
00684
00685
00686 for (unsigned proc=part_read+1; proc<PetscTools::GetNumProcs(); proc++)
00687 {
00688 rProcessorsOffset[proc]++;
00689 }
00690 }
00691
00692
00693
00694
00695 std::vector<unsigned> local_index(PetscTools::GetNumProcs(), 0);
00696
00697 this->mNodesPermutation.resize(this->GetNumNodes());
00698
00699 for (unsigned node_index=0; node_index<this->GetNumNodes(); node_index++)
00700 {
00701 unsigned part_read = npart[node_index];
00702
00703 this->mNodesPermutation[node_index] = rProcessorsOffset[part_read] + local_index[part_read];
00704
00705 local_index[part_read]++;
00706 }
00707
00708 delete[] npart;
00709 }
00710
00711 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00712 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ReorderNodes()
00713 {
00714 assert(!PetscTools::IsSequential());
00715
00716
00717 mNodesMapping.clear();
00718 mHaloNodesMapping.clear();
00719
00720
00721 for (unsigned index=0; index<this->mNodes.size(); index++)
00722 {
00723 unsigned old_index = this->mNodes[index]->GetIndex();
00724 unsigned new_index = this->mNodesPermutation[old_index];
00725
00726 this->mNodes[index]->SetIndex(new_index);
00727 mNodesMapping[new_index] = index;
00728 }
00729
00730 for (unsigned index=0; index<mHaloNodes.size(); index++)
00731 {
00732 unsigned old_index = mHaloNodes[index]->GetIndex();
00733 unsigned new_index = this->mNodesPermutation[old_index];
00734
00735 mHaloNodes[index]->SetIndex(new_index);
00736 mHaloNodesMapping[new_index] = index;
00737 }
00738 }
00739
00740 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00741 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ConstructLinearMesh(unsigned width)
00742 {
00743 assert(ELEMENT_DIM == 1);
00744
00745
00746 if (width==0)
00747 {
00748 EXCEPTION("There aren't enough nodes to make parallelisation worthwhile");
00749 }
00750
00751 mMetisPartitioning=DUMB;
00752 mTotalNumNodes=width+1;
00753 mTotalNumBoundaryElements=2u;
00754 mTotalNumElements=width;
00755
00756
00757 assert(!this->mpDistributedVectorFactory);
00758 this->mpDistributedVectorFactory = new DistributedVectorFactory(mTotalNumNodes);
00759 if (this->mpDistributedVectorFactory->GetLocalOwnership() == 0)
00760 {
00761
00762 return;
00763 }
00764
00765
00766
00767
00768
00769 bool am_top_most = (this->mpDistributedVectorFactory->GetHigh() == mTotalNumNodes);
00770
00771 unsigned lo_node=this->mpDistributedVectorFactory->GetLow();
00772 unsigned hi_node=this->mpDistributedVectorFactory->GetHigh();
00773 if (!PetscTools::AmMaster())
00774 {
00775
00776 lo_node--;
00777 }
00778 if (!am_top_most)
00779 {
00780
00781 hi_node++;
00782 }
00783 Node<SPACE_DIM>* p_old_node=NULL;
00784 for (unsigned node_index=lo_node; node_index<hi_node; node_index++)
00785 {
00786
00787 Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(node_index, node_index==0 || node_index==width, node_index);
00788 if (node_index<this->mpDistributedVectorFactory->GetLow() ||
00789 node_index==this->mpDistributedVectorFactory->GetHigh() )
00790 {
00791
00792 RegisterHaloNode(node_index);
00793 mHaloNodes.push_back(p_node);
00794 }
00795 else
00796 {
00797 RegisterNode(node_index);
00798 this->mNodes.push_back(p_node);
00799
00800
00801
00802 if (node_index==0)
00803 {
00804 this->mBoundaryNodes.push_back(p_node);
00805 RegisterBoundaryElement(0);
00806 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(0, p_node) );
00807 }
00808 if (node_index==width)
00809 {
00810 this->mBoundaryNodes.push_back(p_node);
00811 RegisterBoundaryElement(1);
00812 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(1, p_node) );
00813 }
00814 }
00815 if (node_index>lo_node)
00816 {
00817 std::vector<Node<SPACE_DIM>*> nodes;
00818 nodes.push_back(p_old_node);
00819 nodes.push_back(p_node);
00820 RegisterElement(node_index-1);
00821 this->mElements.push_back(new Element<ELEMENT_DIM,SPACE_DIM>(node_index-1, nodes) );
00822 }
00823
00824 p_old_node=p_node;
00825 }
00826 }
00827
00828 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00829 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ConstructRectangularMesh(unsigned width, unsigned height, bool stagger)
00830 {
00831 assert(SPACE_DIM == 2);
00832 assert(ELEMENT_DIM == 2);
00833
00834 if (height==0)
00835 {
00836 EXCEPTION("There aren't enough nodes to make parallelisation worthwhile");
00837 }
00838
00839 mMetisPartitioning=DUMB;
00840
00841 mTotalNumNodes=(width+1)*(height+1);
00842 mTotalNumBoundaryElements=(width+height)*2;
00843 mTotalNumElements=width*height*2;
00844
00845
00846 DistributedVectorFactory y_partition(height+1);
00847 unsigned lo_y = y_partition.GetLow();
00848 unsigned hi_y = y_partition.GetHigh();
00849
00850 assert(!this->mpDistributedVectorFactory);
00851 this->mpDistributedVectorFactory = new DistributedVectorFactory(mTotalNumNodes, (width+1)*y_partition.GetLocalOwnership());
00852 if (this->mpDistributedVectorFactory->GetLocalOwnership() == 0)
00853 {
00854
00855 return;
00856 }
00857
00858
00859
00860
00861 bool am_top_most = (this->mpDistributedVectorFactory->GetHigh() == mTotalNumNodes);
00862
00863
00864 if (!PetscTools::AmMaster())
00865 {
00866
00867 lo_y--;
00868 }
00869 if (!am_top_most)
00870 {
00871
00872 hi_y++;
00873 }
00874
00875
00876 for (unsigned j=lo_y; j<hi_y; j++)
00877 {
00878 for (unsigned i=0; i<width+1; i++)
00879 {
00880 bool is_boundary=false;
00881 if (i==0 || j==0 || i==width || j==height)
00882 {
00883 is_boundary=true;
00884 }
00885 unsigned global_node_index=((width+1)*(j) + i);
00886 Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(global_node_index, is_boundary, i, j);
00887 if (j<y_partition.GetLow() || j==y_partition.GetHigh() )
00888 {
00889
00890 RegisterHaloNode(global_node_index);
00891 mHaloNodes.push_back(p_node);
00892 }
00893 else
00894 {
00895 RegisterNode(global_node_index);
00896 this->mNodes.push_back(p_node);
00897 }
00898 if (is_boundary)
00899 {
00900 this->mBoundaryNodes.push_back(p_node);
00901 }
00902 }
00903 }
00904
00905
00906 unsigned belem_index;
00907
00908 if (am_top_most)
00909 {
00910 for (unsigned i=0; i<width; i++)
00911 {
00912 std::vector<Node<SPACE_DIM>*> nodes;
00913 nodes.push_back(GetAnyNode( height*(width+1)+i ));
00914 nodes.push_back(GetAnyNode( height*(width+1)+i+1 ));
00915 belem_index=i;
00916 RegisterBoundaryElement(belem_index);
00917 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index,nodes));
00918 }
00919 }
00920
00921
00922 for (unsigned j=lo_y+1; j<hi_y; j++)
00923 {
00924 std::vector<Node<SPACE_DIM>*> nodes;
00925 nodes.push_back(GetAnyNode( (width+1)*(j+1)-1 ));
00926 nodes.push_back(GetAnyNode( (width+1)*j-1 ));
00927 belem_index=width+j-1;
00928 RegisterBoundaryElement(belem_index);
00929 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index,nodes));
00930 }
00931
00932
00933 if (PetscTools::AmMaster())
00934 {
00935 for (unsigned i=0; i<width; i++)
00936 {
00937 std::vector<Node<SPACE_DIM>*> nodes;
00938 nodes.push_back(GetAnyNode( i+1 ));
00939 nodes.push_back(GetAnyNode( i ));
00940 belem_index=width+height+i;
00941 RegisterBoundaryElement(belem_index);
00942 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index,nodes));
00943 }
00944 }
00945
00946
00947 for (unsigned j=lo_y; j<hi_y-1; j++)
00948 {
00949 std::vector<Node<SPACE_DIM>*> nodes;
00950 nodes.push_back(GetAnyNode( (width+1)*(j+1) ));
00951 nodes.push_back(GetAnyNode( (width+1)*(j) ));
00952 belem_index=2*width+height+j;
00953 RegisterBoundaryElement(belem_index);
00954 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index,nodes));
00955 }
00956
00957
00958
00959 unsigned elem_index;
00960 for (unsigned j=lo_y; j<hi_y-1; j++)
00961 {
00962 for (unsigned i=0; i<width; i++)
00963 {
00964 unsigned parity=(i+(height-j))%2;
00965 unsigned nw=(j+1)*(width+1)+i;
00966 unsigned sw=(j)*(width+1)+i;
00967 std::vector<Node<SPACE_DIM>*> upper_nodes;
00968 upper_nodes.push_back(GetAnyNode( nw ));
00969 upper_nodes.push_back(GetAnyNode( nw+1 ));
00970 if (stagger==false || parity == 1)
00971 {
00972 upper_nodes.push_back(GetAnyNode( sw+1 ));
00973 }
00974 else
00975 {
00976 upper_nodes.push_back(GetAnyNode( sw ));
00977 }
00978 elem_index=2*(j*width+i);
00979 RegisterElement(elem_index);
00980 this->mElements.push_back(new Element<ELEMENT_DIM,SPACE_DIM>(elem_index,upper_nodes));
00981 std::vector<Node<SPACE_DIM>*> lower_nodes;
00982 lower_nodes.push_back(GetAnyNode( sw+1 ));
00983 lower_nodes.push_back(GetAnyNode( sw ));
00984 if (stagger==false ||parity == 1)
00985 {
00986 lower_nodes.push_back(GetAnyNode( nw ));
00987 }
00988 else
00989 {
00990 lower_nodes.push_back(GetAnyNode( nw+1 ));
00991 }
00992 elem_index++;
00993 RegisterElement(elem_index);
00994 this->mElements.push_back(new Element<ELEMENT_DIM,SPACE_DIM>(elem_index,lower_nodes));
00995 }
00996 }
00997
00998 }
00999
01000
01001 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
01002 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ConstructCuboid(unsigned width,
01003 unsigned height,
01004 unsigned depth)
01005 {
01006 assert(SPACE_DIM == 3);
01007 assert(ELEMENT_DIM == 3);
01008
01009 if (depth==0)
01010 {
01011 EXCEPTION("There aren't enough nodes to make parallelisation worthwhile");
01012 }
01013
01014
01015 mMetisPartitioning=DUMB;
01016
01017 mTotalNumNodes=(width+1)*(height+1)*(depth+1);
01018 mTotalNumBoundaryElements=((width*height)+(width*depth)+(height*depth))*4;
01019 mTotalNumElements=width*height*depth*6;
01020
01021
01022 DistributedVectorFactory z_partition(depth+1);
01023 unsigned lo_z = z_partition.GetLow();
01024 unsigned hi_z = z_partition.GetHigh();
01025
01026
01027 assert(!this->mpDistributedVectorFactory);
01028 this->mpDistributedVectorFactory = new DistributedVectorFactory(mTotalNumNodes, (width+1)*(height+1)*z_partition.GetLocalOwnership());
01029 if (this->mpDistributedVectorFactory->GetLocalOwnership() == 0)
01030 {
01031
01032 return;
01033 }
01034
01035
01036
01037
01038 bool am_top_most = (this->mpDistributedVectorFactory->GetHigh() == mTotalNumNodes);
01039
01040
01041
01042 if (!PetscTools::AmMaster())
01043 {
01044
01045 lo_z--;
01046 }
01047 if (!am_top_most)
01048 {
01049
01050 hi_z++;
01051 }
01052
01053
01054 unsigned global_node_index;
01055 for (unsigned k=lo_z; k<hi_z; k++)
01056 {
01057 for (unsigned j=0; j<height+1; j++)
01058 {
01059 for (unsigned i=0; i<width+1; i++)
01060 {
01061 bool is_boundary = false;
01062 if (i==0 || j==0 || k==0 || i==width || j==height || k==depth)
01063 {
01064 is_boundary = true;
01065 }
01066 global_node_index = (k*(height+1)+j)*(width+1)+i;
01067
01068 Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(global_node_index, is_boundary, i, j, k);
01069
01070 if (k<z_partition.GetLow() || k==z_partition.GetHigh() )
01071 {
01072
01073 RegisterHaloNode(global_node_index);
01074 mHaloNodes.push_back(p_node);
01075 }
01076 else
01077 {
01078 RegisterNode(global_node_index);
01079 this->mNodes.push_back(p_node);
01080 }
01081
01082 if (is_boundary)
01083 {
01084 this->mBoundaryNodes.push_back(p_node);
01085 }
01086 }
01087 }
01088 }
01089
01090
01091
01092 unsigned element_nodes[6][4] = {{0, 1, 5, 7}, {0, 1, 3, 7},
01093 {0, 2, 3, 7}, {0, 2, 6, 7},
01094 {0, 4, 6, 7}, {0, 4, 5, 7}};
01095 std::vector<Node<SPACE_DIM>*> tetrahedra_nodes;
01096
01097 for (unsigned k=lo_z; k<hi_z-1; k++)
01098 {
01099 unsigned belem_index = 0;
01100 if (k != 0)
01101 {
01102
01103 belem_index = 2*(height*width+k*2*(height+width));
01104 }
01105
01106 for (unsigned j=0; j<height; j++)
01107 {
01108 for (unsigned i=0; i<width; i++)
01109 {
01110
01111 unsigned global_node_indices[8];
01112 unsigned local_node_index = 0;
01113
01114 for (unsigned z = 0; z < 2; z++)
01115 {
01116 for (unsigned y = 0; y < 2; y++)
01117 {
01118 for (unsigned x = 0; x < 2; x++)
01119 {
01120 global_node_indices[local_node_index] = i+x+(width+1)*(j+y+(height+1)*(k+z));
01121
01122 local_node_index++;
01123 }
01124 }
01125 }
01126
01127 for (unsigned m = 0; m < 6; m++)
01128 {
01129
01130
01131 tetrahedra_nodes.clear();
01132
01133 for (unsigned n = 0; n < 4; n++)
01134 {
01135 tetrahedra_nodes.push_back(GetAnyNode( global_node_indices[element_nodes[m][n]] ));
01136 }
01137 unsigned elem_index = 6 * ((k*height+j)*width+i)+m;
01138 RegisterElement(elem_index);
01139 this->mElements.push_back(new Element<ELEMENT_DIM,SPACE_DIM>(elem_index, tetrahedra_nodes));
01140 }
01141
01142
01143 std::vector<Node<SPACE_DIM>*> triangle_nodes;
01144
01145 if (i == 0)
01146 {
01147 triangle_nodes.clear();
01148 triangle_nodes.push_back(GetAnyNode( global_node_indices[0] ));
01149 triangle_nodes.push_back(GetAnyNode( global_node_indices[2] ));
01150 triangle_nodes.push_back(GetAnyNode( global_node_indices[6] ));
01151 RegisterBoundaryElement(belem_index);
01152 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01153 triangle_nodes.clear();
01154 triangle_nodes.push_back(GetAnyNode( global_node_indices[0] ));
01155 triangle_nodes.push_back(GetAnyNode( global_node_indices[6] ));
01156 triangle_nodes.push_back(GetAnyNode( global_node_indices[4] ));
01157 RegisterBoundaryElement(belem_index);
01158 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01159 }
01160 if (i == width-1)
01161 {
01162 triangle_nodes.clear();
01163 triangle_nodes.push_back(GetAnyNode( global_node_indices[1] ));
01164 triangle_nodes.push_back(GetAnyNode( global_node_indices[5] ));
01165 triangle_nodes.push_back(GetAnyNode( global_node_indices[7] ));
01166 RegisterBoundaryElement(belem_index);
01167 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01168 triangle_nodes.clear();
01169 triangle_nodes.push_back(GetAnyNode( global_node_indices[1] ));
01170 triangle_nodes.push_back(GetAnyNode( global_node_indices[7] ));
01171 triangle_nodes.push_back(GetAnyNode( global_node_indices[3] ));
01172 RegisterBoundaryElement(belem_index);
01173 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01174 }
01175 if (j == 0)
01176 {
01177 triangle_nodes.clear();
01178 triangle_nodes.push_back(GetAnyNode( global_node_indices[0] ));
01179 triangle_nodes.push_back(GetAnyNode( global_node_indices[5] ));
01180 triangle_nodes.push_back(GetAnyNode( global_node_indices[1] ));
01181 RegisterBoundaryElement(belem_index);
01182 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01183 triangle_nodes.clear();
01184 triangle_nodes.push_back(GetAnyNode( global_node_indices[0] ));
01185 triangle_nodes.push_back(GetAnyNode( global_node_indices[4] ));
01186 triangle_nodes.push_back(GetAnyNode( global_node_indices[5] ));
01187 RegisterBoundaryElement(belem_index);
01188 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01189 }
01190 if (j == height-1)
01191 {
01192 triangle_nodes.clear();
01193 triangle_nodes.push_back(GetAnyNode( global_node_indices[2] ));
01194 triangle_nodes.push_back(GetAnyNode( global_node_indices[3] ));
01195 triangle_nodes.push_back(GetAnyNode( global_node_indices[7] ));
01196 RegisterBoundaryElement(belem_index);
01197 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01198 triangle_nodes.clear();
01199 triangle_nodes.push_back(GetAnyNode( global_node_indices[2] ));
01200 triangle_nodes.push_back(GetAnyNode( global_node_indices[7] ));
01201 triangle_nodes.push_back(GetAnyNode( global_node_indices[6] ));
01202 RegisterBoundaryElement(belem_index);
01203 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01204 }
01205 if (k == 0)
01206 {
01207 triangle_nodes.clear();
01208 triangle_nodes.push_back(GetAnyNode( global_node_indices[0] ));
01209 triangle_nodes.push_back(GetAnyNode( global_node_indices[3] ));
01210 triangle_nodes.push_back(GetAnyNode( global_node_indices[2] ));
01211 RegisterBoundaryElement(belem_index);
01212 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01213 triangle_nodes.clear();
01214 triangle_nodes.push_back(GetAnyNode( global_node_indices[0] ));
01215 triangle_nodes.push_back(GetAnyNode( global_node_indices[1] ));
01216 triangle_nodes.push_back(GetAnyNode( global_node_indices[3] ));
01217 RegisterBoundaryElement(belem_index);
01218 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01219 }
01220 if (k == depth-1)
01221 {
01222 triangle_nodes.clear();
01223 triangle_nodes.push_back(GetAnyNode( global_node_indices[4] ));
01224 triangle_nodes.push_back(GetAnyNode( global_node_indices[7] ));
01225 triangle_nodes.push_back(GetAnyNode( global_node_indices[5] ));
01226 RegisterBoundaryElement(belem_index);
01227 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01228 triangle_nodes.clear();
01229 triangle_nodes.push_back(GetAnyNode( global_node_indices[4] ));
01230 triangle_nodes.push_back(GetAnyNode( global_node_indices[6] ));
01231 triangle_nodes.push_back(GetAnyNode( global_node_indices[7] ));
01232 RegisterBoundaryElement(belem_index);
01233 this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(belem_index++,triangle_nodes));
01234 }
01235 }
01236 }
01237 }
01238 }
01239
01240 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
01241 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::Scale(const double xFactor, const double yFactor, const double zFactor)
01242 {
01243
01244 AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::Scale(xFactor, yFactor, zFactor);
01245
01246 for (unsigned i=0; i<mHaloNodes.size(); i++)
01247 {
01248 c_vector<double, SPACE_DIM>& r_location = mHaloNodes[i]->rGetModifiableLocation();
01249 if (SPACE_DIM>=3)
01250 {
01251 r_location[2] *= zFactor;
01252 }
01253 if (SPACE_DIM>=2)
01254 {
01255 r_location[1] *= yFactor;
01256 }
01257 r_location[0] *= xFactor;
01258 }
01259
01260 }
01261
01262
01263 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
01264 void DistributedTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ParMetisLibraryNodePartitioning(AbstractMeshReader<ELEMENT_DIM, SPACE_DIM>& rMeshReader,
01265 std::set<unsigned>& rElementsOwned,
01266 std::set<unsigned>& rNodesOwned,
01267 std::set<unsigned>& rHaloNodesOwned,
01268 std::vector<unsigned>& rProcessorsOffset)
01269 {
01270 assert(!PetscTools::IsSequential());
01271 assert(ELEMENT_DIM==2 || ELEMENT_DIM==3);
01272
01273 unsigned num_elements = rMeshReader.GetNumElements();
01274 unsigned num_procs = PetscTools::GetNumProcs();
01275 unsigned local_proc_index = PetscTools::GetMyRank();
01276
01277
01278
01279
01280 idxtype element_distribution[num_procs+1];
01281 idxtype element_count[num_procs];
01282
01283 element_distribution[0]=0;
01284
01285 for (unsigned proc_index=1; proc_index<num_procs; proc_index++)
01286 {
01287 element_distribution[proc_index] = element_distribution[proc_index-1] + num_elements/num_procs;
01288 element_count[proc_index-1] = element_distribution[proc_index] - element_distribution[proc_index-1];
01289 }
01290
01291 element_distribution[num_procs] = num_elements;
01292 element_count[num_procs-1] = element_distribution[num_procs] - element_distribution[num_procs-1];
01293
01294
01295
01296
01297 unsigned first_local_element = element_distribution[local_proc_index];
01298 unsigned last_plus_one_element = element_distribution[local_proc_index+1];
01299 unsigned num_local_elements = last_plus_one_element - first_local_element;
01300
01301 idxtype* eind = new idxtype[num_local_elements*(ELEMENT_DIM+1)];
01302 idxtype* eptr = new idxtype[num_local_elements+1];
01303
01304
01305 for (unsigned element_index = 0; element_index < first_local_element; element_index++)
01306 {
01307 ElementData element_data = rMeshReader.GetNextElementData();
01308 }
01309
01310 unsigned counter=0;
01311 for (unsigned element_index = 0; element_index < num_local_elements; element_index++)
01312 {
01313 ElementData element_data = rMeshReader.GetNextElementData();
01314
01315 eptr[element_index] = counter;
01316 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
01317 {
01318 eind[counter++] = element_data.NodeIndices[i];
01319 }
01320
01321 }
01322 eptr[num_local_elements] = counter;
01323
01324 rMeshReader.Reset();
01325
01326 int numflag = 0;
01327 int ncommonnodes = 3;
01328 MPI_Comm communicator = PETSC_COMM_WORLD;
01329
01330 idxtype* xadj;
01331 idxtype* adjncy;
01332
01333 Timer::Reset();
01334 ParMETIS_V3_Mesh2Dual(element_distribution, eptr, eind,
01335 &numflag, &ncommonnodes, &xadj, &adjncy, &communicator);
01336
01337
01338 delete[] eind;
01339 delete[] eptr;
01340
01341 int weight_flag = 0;
01342 int n_constrains = 0;
01343 int n_subdomains = PetscTools::GetNumProcs();
01344 int options[3];
01345 options[0] = 0;
01346 int edgecut;
01347
01348 idxtype* local_partition = new idxtype[num_local_elements];
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358
01359
01360
01361
01362
01363 Timer::Reset();
01364 ParMETIS_V3_PartKway(element_distribution, xadj, adjncy, NULL, NULL, &weight_flag, &numflag,
01365 &n_constrains, &n_subdomains, NULL, NULL,
01366 options, &edgecut, local_partition, &communicator);
01367
01368
01369
01370 idxtype* global_element_partition = new idxtype[num_elements];
01371
01372 MPI_Allgatherv(local_partition, num_local_elements, MPI_INT,
01373 global_element_partition, element_count, element_distribution, MPI_INT, PETSC_COMM_WORLD);
01374
01375 delete[] local_partition;
01376
01377 for(unsigned elem_index=0; elem_index<num_elements; elem_index++)
01378 {
01379 if ((unsigned) global_element_partition[elem_index] == local_proc_index)
01380 {
01381 rElementsOwned.insert(elem_index);
01382 }
01383 }
01384
01385 rMeshReader.Reset();
01386 free(xadj);
01387 free(adjncy);
01388
01389 unsigned num_nodes = rMeshReader.GetNumNodes();
01390
01391
01392 std::vector<unsigned> global_node_partition;
01393 global_node_partition.resize(num_nodes, UNASSIGNED_NODE);
01394
01395 assert(rProcessorsOffset.size() == 0);
01396 rProcessorsOffset.resize(PetscTools::GetNumProcs(), 0);
01397
01398
01399
01400
01401
01402
01403
01404
01405
01406 for (unsigned element_number = 0; element_number < mTotalNumElements; element_number++)
01407 {
01408 unsigned element_owner = global_element_partition[element_number];
01409 ElementData element_data = rMeshReader.GetNextElementData();
01410
01411 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
01412 {
01413
01414
01415
01416
01417 if( global_node_partition[element_data.NodeIndices[i]] == UNASSIGNED_NODE )
01418 {
01419 if (element_owner == local_proc_index)
01420 {
01421 rNodesOwned.insert(element_data.NodeIndices[i]);
01422 }
01423
01424 global_node_partition[element_data.NodeIndices[i]] = element_owner;
01425
01426
01427
01428
01429 for (unsigned proc=element_owner+1; proc<PetscTools::GetNumProcs(); proc++)
01430 {
01431 rProcessorsOffset[proc]++;
01432 }
01433 }
01434 else
01435 {
01436 if (element_owner == local_proc_index)
01437 {
01438
01439 if (global_node_partition[element_data.NodeIndices[i]] != local_proc_index)
01440 {
01441 rHaloNodesOwned.insert(element_data.NodeIndices[i]);
01442 }
01443 }
01444 }
01445 }
01446 }
01447
01448 delete[] global_element_partition;
01449
01450
01451
01452
01453
01454
01455 rMeshReader.Reset();
01456
01457 for (unsigned element_number = 0; element_number < mTotalNumElements; element_number++)
01458 {
01459 ElementData element_data = rMeshReader.GetNextElementData();
01460
01461 bool element_owned = false;
01462 std::set<unsigned> temp_halo_nodes;
01463
01464 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
01465 {
01466 if (rNodesOwned.find(element_data.NodeIndices[i]) != rNodesOwned.end())
01467 {
01468 element_owned = true;
01469 rElementsOwned.insert(element_number);
01470 }
01471 else
01472 {
01473 temp_halo_nodes.insert(element_data.NodeIndices[i]);
01474 }
01475 }
01476
01477 if (element_owned)
01478 {
01479 rHaloNodesOwned.insert(temp_halo_nodes.begin(), temp_halo_nodes.end());
01480 }
01481 }
01482
01483 rMeshReader.Reset();
01484
01485
01486
01487
01488 std::vector<unsigned> local_index(PetscTools::GetNumProcs(), 0);
01489
01490 this->mNodesPermutation.resize(this->GetNumNodes());
01491
01492 for (unsigned node_index=0; node_index<this->GetNumNodes(); node_index++)
01493 {
01494 unsigned partition = global_node_partition[node_index];
01495 assert(partition!=UNASSIGNED_NODE);
01496
01497 this->mNodesPermutation[node_index] = rProcessorsOffset[partition] + local_index[partition];
01498
01499 local_index[partition]++;
01500 }
01501
01502 }
01503
01505
01507
01508 template class DistributedTetrahedralMesh<1,1>;
01509 template class DistributedTetrahedralMesh<1,2>;
01510 template class DistributedTetrahedralMesh<1,3>;
01511 template class DistributedTetrahedralMesh<2,2>;
01512 template class DistributedTetrahedralMesh<2,3>;
01513 template class DistributedTetrahedralMesh<3,3>;
01514
01515
01516
01517 #include "SerializationExportWrapperForCpp.hpp"
01518 EXPORT_TEMPLATE_CLASS_ALL_DIMS(DistributedTetrahedralMesh);