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