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 "TetrahedralMesh.hpp"
00030
00031 #include <iostream>
00032 #include <cassert>
00033 #include <sstream>
00034 #include <map>
00035
00036 #include "BoundaryElement.hpp"
00037 #include "Element.hpp"
00038 #include "Exception.hpp"
00039 #include "Node.hpp"
00040 #include "OutputFileHandler.hpp"
00041 #include "PetscTools.hpp"
00042 #include "RandomNumberGenerator.hpp"
00043
00044
00045 #define REAL double
00046 #define VOID void
00047 #include "triangle.h"
00048 #include "tetgen.h"
00049 #undef REAL
00050 #undef VOID
00051
00053
00055
00056
00057 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00058 TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::TetrahedralMesh()
00059 {
00060 Clear();
00061 }
00062
00063 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00064 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ConstructFromMeshReader(
00065 AbstractMeshReader<ELEMENT_DIM, SPACE_DIM>& rMeshReader)
00066 {
00067 this->mMeshFileBaseName = rMeshReader.GetMeshFileBaseName();
00068
00069
00070 unsigned num_nodes = rMeshReader.GetNumNodes();
00071
00072
00073
00074 this->mNodes.reserve(num_nodes);
00075
00076 rMeshReader.Reset();
00077
00078
00079
00080
00081
00082 std::vector<double> coords;
00083 for (unsigned i=0; i < num_nodes; i++)
00084 {
00085 coords = rMeshReader.GetNextNode();
00086 Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(i, coords, false);
00087
00088 for (unsigned i = 0; i < rMeshReader.GetNodeAttributes().size(); i++)
00089 {
00090 double attribute = rMeshReader.GetNodeAttributes()[i];
00091 p_node->AddNodeAttribute(attribute);
00092 }
00093
00094 this->mNodes.push_back(p_node);
00095 }
00096
00097
00098
00099 rMeshReader.Reset();
00100
00101
00102 this->mElements.reserve(rMeshReader.GetNumElements());
00103
00104 for (unsigned element_index=0; element_index < (unsigned) rMeshReader.GetNumElements(); element_index++)
00105 {
00106 ElementData element_data = rMeshReader.GetNextElementData();
00107 std::vector<Node<SPACE_DIM>*> nodes;
00108
00109
00110
00111
00112
00113
00114 for (unsigned j=0; j<ELEMENT_DIM+1; j++)
00115 {
00116 assert(element_data.NodeIndices[j] < this->mNodes.size());
00117 nodes.push_back(this->mNodes[element_data.NodeIndices[j]]);
00118 }
00119
00120 Element<ELEMENT_DIM,SPACE_DIM>* p_element = new Element<ELEMENT_DIM,SPACE_DIM>(element_index, nodes);
00121 this->mElements.push_back(p_element);
00122
00123 if (rMeshReader.GetNumElementAttributes() > 0)
00124 {
00125 assert(rMeshReader.GetNumElementAttributes() == 1);
00126 unsigned attribute_value = element_data.AttributeValue;
00127 p_element->SetRegion(attribute_value);
00128 }
00129 }
00130
00131
00132 for (unsigned face_index=0; face_index<(unsigned)rMeshReader.GetNumFaces(); face_index++)
00133 {
00134 ElementData face_data = rMeshReader.GetNextFaceData();
00135 std::vector<unsigned> node_indices = face_data.NodeIndices;
00136
00137
00138
00139
00140
00141
00142 std::vector<Node<SPACE_DIM>*> nodes;
00143 for (unsigned node_index=0; node_index<ELEMENT_DIM; node_index++)
00144 {
00145 assert(node_indices[node_index] < this->mNodes.size());
00146
00147 nodes.push_back(this->mNodes[node_indices[node_index]]);
00148 }
00149
00150
00151
00152
00153 assert(nodes.size()==ELEMENT_DIM);
00154 for (unsigned j=0; j<nodes.size(); j++)
00155 {
00156 if (!nodes[j]->IsBoundaryNode())
00157 {
00158 nodes[j]->SetAsBoundaryNode();
00159 this->mBoundaryNodes.push_back(nodes[j]);
00160 }
00161
00162
00163 nodes[j]->AddBoundaryElement(face_index);
00164 }
00165
00166
00167 BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* p_boundary_element = new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(face_index, nodes);
00168 this->mBoundaryElements.push_back(p_boundary_element);
00169
00170 if (rMeshReader.GetNumFaceAttributes() > 0)
00171 {
00172 assert(rMeshReader.GetNumFaceAttributes() == 1);
00173 unsigned attribute_value = face_data.AttributeValue;
00174 p_boundary_element->SetRegion(attribute_value);
00175 }
00176 }
00177
00178 RefreshJacobianCachedData();
00179 rMeshReader.Reset();
00180 }
00181
00182 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00183 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ReadNodesPerProcessorFile(const std::string& rNodesPerProcessorFile)
00184 {
00185 std::vector<unsigned> nodes_per_processor_vec;
00186
00187 std::ifstream file_stream(rNodesPerProcessorFile.c_str());
00188 if (file_stream.is_open())
00189 {
00190 while (file_stream)
00191 {
00192 unsigned nodes_per_processor;
00193 file_stream >> nodes_per_processor;
00194
00195 if (file_stream)
00196 {
00197 nodes_per_processor_vec.push_back(nodes_per_processor);
00198 }
00199 }
00200 }
00201 else
00202 {
00203 EXCEPTION("Unable to read nodes per processor file " + rNodesPerProcessorFile);
00204 }
00205
00206 unsigned sum = 0;
00207 for (unsigned i=0; i<nodes_per_processor_vec.size(); i++)
00208 {
00209 sum += nodes_per_processor_vec[i];
00210 }
00211
00212 if (sum != this->GetNumNodes())
00213 {
00214 std::stringstream string_stream;
00215 string_stream << "Sum of nodes per processor, " << sum
00216 << ", not equal to number of nodes in mesh, " << this->GetNumNodes();
00217 EXCEPTION(string_stream.str());
00218 }
00219
00220 unsigned num_owned=nodes_per_processor_vec[PetscTools::GetMyRank()];
00221
00222 if (nodes_per_processor_vec.size() != PetscTools::GetNumProcs())
00223 {
00224 EXCEPTION("Number of processes doesn't match the size of the nodes-per-processor file");
00225 }
00226 delete this->mpDistributedVectorFactory;
00227 this->mpDistributedVectorFactory=new DistributedVectorFactory(this->GetNumNodes(), num_owned);
00228 }
00229 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00230 bool TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::CheckIsConforming()
00231 {
00232
00233
00234
00235
00236
00237
00238 std::set< std::set<unsigned> > odd_parity_faces;
00239
00240 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00241 iter != this->GetElementIteratorEnd();
00242 ++iter)
00243 {
00244 for (unsigned face_index=0; face_index<=ELEMENT_DIM; face_index++)
00245 {
00246 std::set<unsigned> face_info;
00247 for (unsigned node_index=0; node_index<=ELEMENT_DIM; node_index++)
00248 {
00249
00250 if (node_index != face_index)
00251 {
00252 face_info.insert(iter->GetNodeGlobalIndex(node_index));
00253 }
00254 }
00255
00256 std::set< std::set<unsigned> >::iterator find_face=odd_parity_faces.find(face_info);
00257 if( find_face != odd_parity_faces.end())
00258 {
00259
00260
00261 odd_parity_faces.erase(find_face);
00262 }
00263 else
00264 {
00265
00266 odd_parity_faces.insert(face_info);
00267 }
00268
00269 }
00270 }
00271
00272
00273
00274
00275 return( odd_parity_faces.size() == this->GetNumBoundaryElements());
00276 }
00277
00278 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00279 double TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetVolume()
00280 {
00281 double mesh_volume = 0.0;
00282
00283 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00284 iter != this->GetElementIteratorEnd();
00285 ++iter)
00286 {
00287 mesh_volume += iter->GetVolume(mElementJacobianDeterminants[iter->GetIndex()]);
00288 }
00289
00290 return mesh_volume;
00291 }
00292
00293 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00294 double TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetSurfaceArea()
00295 {
00296
00297 assert(ELEMENT_DIM >= 1);
00298 const unsigned bound_element_dim = ELEMENT_DIM-1;
00299 assert(bound_element_dim < 3);
00300 if ( bound_element_dim == 0)
00301 {
00302 return 0.0;
00303 }
00304
00305 double mesh_surface = 0.0;
00306 typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator it = this->GetBoundaryElementIteratorBegin();
00307
00308 while (it != this->GetBoundaryElementIteratorEnd())
00309 {
00310 mesh_surface += mBoundaryElementJacobianDeterminants[(*it)->GetIndex()];
00311 it++;
00312 }
00313
00314 if ( bound_element_dim == 2)
00315 {
00316 mesh_surface /= 2.0;
00317 }
00318
00319 return mesh_surface;
00320 }
00321
00322 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00323 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::PermuteNodes()
00324 {
00325
00326 RandomNumberGenerator* p_rng = RandomNumberGenerator::Instance();
00327 std::vector<unsigned> perm;
00328 p_rng->Shuffle(this->mNodes.size(), perm);
00329
00330 PermuteNodes(perm);
00331 }
00332
00333 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00334 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::PermuteNodes(const std::vector<unsigned>& perm)
00335 {
00336
00337 assert( this->GetNumAllNodes() == this->GetNumNodes());
00338
00339 assert(perm.size() == this->mNodes.size());
00340
00341
00342 std::vector< Node<SPACE_DIM>* > copy_m_nodes;
00343 copy_m_nodes.assign(this->mNodes.begin(), this->mNodes.end());
00344
00345 for (unsigned original_index=0; original_index<this->mNodes.size(); original_index++)
00346 {
00347 assert(perm[original_index] < this->mNodes.size());
00348
00349 this->mNodes[ perm[original_index] ] = copy_m_nodes[original_index];
00350 }
00351
00352
00353 for (unsigned index=0; index<this->mNodes.size(); index++)
00354 {
00355 this->mNodes[index]->SetIndex(index);
00356 }
00357
00358 this->mNodesPermutation = perm;
00359 }
00360
00361
00362 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00363 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetContainingElementIndex(const ChastePoint<SPACE_DIM>& rTestPoint,
00364 bool strict,
00365 std::set<unsigned> testElements,
00366 bool onlyTryWithTestElements)
00367 {
00368 for (std::set<unsigned>::iterator iter=testElements.begin(); iter!=testElements.end(); iter++)
00369 {
00370 assert(*iter<this->GetNumElements());
00371 if (this->mElements[*iter]->IncludesPoint(rTestPoint, strict))
00372 {
00373 assert(!this->mElements[*iter]->IsDeleted());
00374 return *iter;
00375 }
00376 }
00377
00378 if(!onlyTryWithTestElements)
00379 {
00380 for (unsigned i=0; i<this->mElements.size(); i++)
00381 {
00382 if (this->mElements[i]->IncludesPoint(rTestPoint, strict))
00383 {
00384 assert(!this->mElements[i]->IsDeleted());
00385 return i;
00386 }
00387 }
00388 }
00389
00390
00391 std::stringstream ss;
00392 ss << "Point [";
00393 for(unsigned j=0; (int)j<(int)SPACE_DIM-1; j++)
00394 {
00395 ss << rTestPoint[j] << ",";
00396 }
00397 ss << rTestPoint[SPACE_DIM-1] << "] is not in ";
00398 if(!onlyTryWithTestElements)
00399 {
00400 ss << "mesh - all elements tested";
00401 }
00402 else
00403 {
00404 ss << "set of elements given";
00405 }
00406 EXCEPTION(ss.str());
00407 }
00408
00409 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00410 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetContainingElementIndexWithInitialGuess(const ChastePoint<SPACE_DIM>& rTestPoint, unsigned startingElementGuess, bool strict)
00411 {
00412 assert(startingElementGuess<this->GetNumElements());
00413
00414
00415
00416
00417 unsigned i = startingElementGuess;
00418 bool reached_end = false;
00419
00420 while(!reached_end)
00421 {
00422 if (this->mElements[i]->IncludesPoint(rTestPoint, strict))
00423 {
00424 assert(!this->mElements[i]->IsDeleted());
00425 return i;
00426 }
00427
00428
00429 i++;
00430 if(i==this->GetNumElements())
00431 {
00432 i=0;
00433 }
00434
00435
00436 if(i==startingElementGuess)
00437 {
00438 reached_end = true;
00439 }
00440 }
00441
00442
00443 std::stringstream ss;
00444 ss << "Point [";
00445 for(unsigned j=0; (int)j<(int)SPACE_DIM-1; j++)
00446 {
00447 ss << rTestPoint[j] << ",";
00448 }
00449 ss << rTestPoint[SPACE_DIM-1] << "] is not in mesh - all elements tested";
00450 EXCEPTION(ss.str());
00451 }
00452
00453
00454 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00455 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNearestElementIndex(const ChastePoint<SPACE_DIM>& rTestPoint)
00456 {
00457 double max_min_weight = -INFINITY;
00458 unsigned closest_index = 0;
00459 for (unsigned i=0; i<this->mElements.size(); i++)
00460 {
00461 c_vector<double, ELEMENT_DIM+1> weight=this->mElements[i]->CalculateInterpolationWeights(rTestPoint);
00462 double neg_weight_sum=0.0;
00463 for (unsigned j=0; j<=ELEMENT_DIM; j++)
00464 {
00465 if (weight[j] < 0.0)
00466 {
00467 neg_weight_sum += weight[j];
00468 }
00469 }
00470 if (neg_weight_sum > max_min_weight)
00471 {
00472 max_min_weight = neg_weight_sum;
00473 closest_index = i;
00474 }
00475
00476 }
00477 assert(!this->mElements[closest_index]->IsDeleted());
00478 return closest_index;
00479 }
00480
00481 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00482 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNearestElementIndexFromTestElements(const ChastePoint<SPACE_DIM>& rTestPoint,
00483 std::set<unsigned> testElements)
00484 {
00485 assert(testElements.size()>0);
00486
00487 double max_min_weight = -INFINITY;
00488 unsigned closest_index = 0;
00489 for(std::set<unsigned>::iterator iter = testElements.begin();
00490 iter != testElements.end();
00491 iter++)
00492 {
00493 c_vector<double, ELEMENT_DIM+1> weight=this->mElements[*iter]->CalculateInterpolationWeights(rTestPoint);
00494 double neg_weight_sum=0.0;
00495 for (unsigned j=0; j<=ELEMENT_DIM; j++)
00496 {
00497 if (weight[j] < 0.0)
00498 {
00499 neg_weight_sum += weight[j];
00500 }
00501 }
00502 if (neg_weight_sum > max_min_weight)
00503 {
00504 max_min_weight = neg_weight_sum;
00505 closest_index = *iter;
00506 }
00507
00508 }
00509 assert(!this->mElements[closest_index]->IsDeleted());
00510 return closest_index;
00511 }
00512
00513
00514 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00515 std::vector<unsigned> TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetContainingElementIndices(const ChastePoint<SPACE_DIM> &rTestPoint)
00516 {
00517 std::vector<unsigned> element_indices;
00518 for (unsigned i=0; i<this->mElements.size(); i++)
00519 {
00520 if (this->mElements[i]->IncludesPoint(rTestPoint))
00521 {
00522 assert(!this->mElements[i]->IsDeleted());
00523 element_indices.push_back(i);
00524 }
00525 }
00526 return element_indices;
00527 }
00528
00529
00530 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00531 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::Clear()
00532 {
00533
00534 for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00535 {
00536 delete this->mBoundaryElements[i];
00537 }
00538 for (unsigned i=0; i<this->mElements.size(); i++)
00539 {
00540 delete this->mElements[i];
00541 }
00542 for (unsigned i=0; i<this->mNodes.size(); i++)
00543 {
00544 delete this->mNodes[i];
00545 }
00546
00547 this->mNodes.clear();
00548 this->mElements.clear();
00549 this->mBoundaryElements.clear();
00550 this->mBoundaryNodes.clear();
00551 }
00552
00553 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00554 std::set<unsigned> TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::CalculateBoundaryOfFlaggedRegion()
00555 {
00556
00557 typedef std::set<unsigned> FaceNodes;
00558
00559
00560
00561 std::map<FaceNodes,bool> face_on_boundary;
00562
00563
00564 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00565 iter != this->GetElementIteratorEnd();
00566 ++iter)
00567 {
00568 if (iter->IsFlagged())
00569 {
00570
00571 std::set<unsigned> all_nodes;
00572 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
00573 {
00574 all_nodes.insert(iter->GetNodeGlobalIndex(i));
00575 }
00576
00577
00578 for (unsigned i=0; i<ELEMENT_DIM+1; i++)
00579 {
00580 FaceNodes face_nodes = all_nodes;
00581 face_nodes.erase(iter->GetNodeGlobalIndex(i));
00582
00583
00584 std::map<FaceNodes,bool>::iterator it = face_on_boundary.find(face_nodes);
00585
00586 if (it == face_on_boundary.end())
00587 {
00588
00589 face_on_boundary[face_nodes]=true;
00590 }
00591 else
00592 {
00593
00594 it->second = false;
00595 }
00596 }
00597 }
00598 }
00599
00600
00601 std::set<unsigned> boundary_of_flagged_region;
00602
00603
00604 std::map<FaceNodes,bool>::iterator it=face_on_boundary.begin();
00605 while (it!=face_on_boundary.end())
00606 {
00607
00608 if (it->second==true)
00609 {
00610
00611 boundary_of_flagged_region.insert(it->first.begin(),it->first.end());
00612 }
00613 it++;
00614 }
00615
00616 return boundary_of_flagged_region;
00617 }
00618
00619 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00620 double TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetAngleBetweenNodes(unsigned indexA, unsigned indexB)
00621 {
00622 assert(SPACE_DIM == 2);
00623 assert(SPACE_DIM == ELEMENT_DIM);
00624
00625 double x_diff = this->mNodes[indexB]->rGetLocation()[0] - this->mNodes[indexA]->rGetLocation()[0];
00626 double y_diff = this->mNodes[indexB]->rGetLocation()[1] - this->mNodes[indexA]->rGetLocation()[1];
00627
00628 if (x_diff==0)
00629 {
00630 if (y_diff>0)
00631 {
00632 return M_PI/2.0;
00633 }
00634 else if (y_diff<0)
00635 {
00636 return -M_PI/2.0;
00637 }
00638 else
00639 {
00640 EXCEPTION("Tried to compute polar angle of (0,0)");
00641 }
00642 }
00643
00644 double angle = atan2(y_diff,x_diff);
00645 return angle;
00646 }
00647
00648 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00649 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::UnflagAllElements()
00650 {
00651 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00652 iter != this->GetElementIteratorEnd();
00653 ++iter)
00654 {
00655 iter->Unflag();
00656 }
00657 }
00658
00659 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00660 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::FlagElementsNotContainingNodes(std::set<unsigned> nodesList)
00661 {
00662 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00663 iter != this->GetElementIteratorEnd();
00664 ++iter)
00665 {
00666 bool found_node = false;
00667
00668 for (unsigned i=0; i<iter->GetNumNodes(); i++)
00669 {
00670 unsigned node_index = iter->GetNodeGlobalIndex(i);
00671
00672 std::set<unsigned>::iterator set_iter = nodesList.find(node_index);
00673 if (set_iter != nodesList.end())
00674 {
00675 found_node = true;
00676 }
00677 }
00678
00679 if (!found_node)
00680 {
00681 iter->Flag();
00682 }
00683 }
00684 }
00685
00686
00688
00690
00691 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00692 Node<SPACE_DIM>* TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator::GetNodeA()
00693 {
00694 assert((*this) != mrMesh.EdgesEnd());
00695 Element<ELEMENT_DIM,SPACE_DIM>* p_element = mrMesh.GetElement(mElemIndex);
00696 return p_element->GetNode(mNodeALocalIndex);
00697 }
00698
00699 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00700 Node<SPACE_DIM>* TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator::GetNodeB()
00701 {
00702 assert((*this) != mrMesh.EdgesEnd());
00703 Element<ELEMENT_DIM,SPACE_DIM>* p_element = mrMesh.GetElement(mElemIndex);
00704 return p_element->GetNode(mNodeBLocalIndex);
00705 }
00706
00707
00708 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00709 bool TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator::operator!=(const TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator& rOther)
00710 {
00711 return (mElemIndex != rOther.mElemIndex ||
00712 mNodeALocalIndex != rOther.mNodeALocalIndex ||
00713 mNodeBLocalIndex != rOther.mNodeBLocalIndex);
00714 }
00715
00716 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00717 typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator& TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator::operator++()
00718 {
00719 bool already_seen_this_edge;
00720
00721 unsigned num_elements = mrMesh.GetNumAllElements();
00722 std::pair<unsigned, unsigned> current_node_pair;
00723 do
00724 {
00725
00726
00727 mNodeBLocalIndex = (mNodeBLocalIndex + 1) % (ELEMENT_DIM+1);
00728 if (mNodeBLocalIndex == mNodeALocalIndex)
00729 {
00730 mNodeALocalIndex = (mNodeALocalIndex + 1) % (ELEMENT_DIM+1);
00731 mNodeBLocalIndex = (mNodeALocalIndex + 1) % (ELEMENT_DIM+1);
00732 }
00733
00734 if (mNodeALocalIndex == 0 && mNodeBLocalIndex == 1)
00735 {
00736
00737
00738 do
00739 {
00740 mElemIndex++;
00741 }
00742 while (mElemIndex!=num_elements && mrMesh.GetElement(mElemIndex)->IsDeleted());
00743 }
00744
00745 if (mElemIndex != num_elements)
00746 {
00747 Element<ELEMENT_DIM, SPACE_DIM>* p_element = mrMesh.GetElement(mElemIndex);
00748 unsigned node_a_global_index = p_element->GetNodeGlobalIndex(mNodeALocalIndex);
00749 unsigned node_b_global_index = p_element->GetNodeGlobalIndex(mNodeBLocalIndex);
00750 if (node_b_global_index < node_a_global_index)
00751 {
00752
00753 unsigned temp = node_a_global_index;
00754 node_a_global_index = node_b_global_index;
00755 node_b_global_index = temp;
00756 }
00757
00758
00759 current_node_pair = std::pair<unsigned, unsigned>(node_a_global_index, node_b_global_index);
00760 already_seen_this_edge = (mEdgesVisited.count(current_node_pair) != 0);
00761 }
00762 else
00763 {
00764 already_seen_this_edge = false;
00765 }
00766 }
00767 while (already_seen_this_edge);
00768 mEdgesVisited.insert(current_node_pair);
00769
00770 return (*this);
00771 }
00772
00773
00774 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00775 TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator::EdgeIterator(TetrahedralMesh& rMesh, unsigned elemIndex)
00776 : mrMesh(rMesh),
00777 mElemIndex(elemIndex),
00778 mNodeALocalIndex(0),
00779 mNodeBLocalIndex(1)
00780 {
00781 if (elemIndex==mrMesh.GetNumAllElements())
00782 {
00783 return;
00784 }
00785
00786 mEdgesVisited.clear();
00787
00788
00789
00790 unsigned node_a_global_index = mrMesh.GetElement(mElemIndex)->GetNodeGlobalIndex(mNodeALocalIndex);
00791 unsigned node_b_global_index = mrMesh.GetElement(mElemIndex)->GetNodeGlobalIndex(mNodeBLocalIndex);
00792 if (node_b_global_index < node_a_global_index)
00793 {
00794
00795 unsigned temp = node_a_global_index;
00796 node_a_global_index = node_b_global_index;
00797 node_b_global_index = temp;
00798 }
00799
00800
00801 std::pair<unsigned, unsigned> current_node_pair=std::pair<unsigned, unsigned>(node_a_global_index, node_b_global_index);
00802
00803 mEdgesVisited.insert(current_node_pair);
00804 }
00805
00806 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00807 typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgesBegin()
00808 {
00809 unsigned first_element_index=0;
00810 while (first_element_index!=this->GetNumAllElements() && this->GetElement(first_element_index)->IsDeleted())
00811 {
00812 first_element_index++;
00813 }
00814 return EdgeIterator(*this, first_element_index);
00815 }
00816
00817 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00818 typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgeIterator TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::EdgesEnd()
00819 {
00820 return EdgeIterator(*this, this->GetNumAllElements());
00821 }
00822
00823 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00824 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::RefreshMesh()
00825 {
00826 RefreshJacobianCachedData();
00827 }
00828
00829 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00830 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SolveNodeMapping(unsigned index) const
00831 {
00832 assert(index < this->mNodes.size() );
00833 return index;
00834 }
00835
00836 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00837 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SolveElementMapping(unsigned index) const
00838 {
00839 assert(index < this->mElements.size() );
00840 return index;
00841 }
00842
00843 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00844 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::SolveBoundaryElementMapping(unsigned index) const
00845 {
00846 assert(index < this->mBoundaryElements.size() );
00847 return index;
00848 }
00849
00850 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00851 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::RefreshJacobianCachedData()
00852 {
00853 unsigned num_elements = this->GetNumAllElements();
00854 unsigned num_boundary_elements = this->GetNumAllBoundaryElements();
00855
00856
00857 this->mElementJacobians.resize(num_elements);
00858 this->mElementInverseJacobians.resize(num_elements);
00859
00860 if (ELEMENT_DIM < SPACE_DIM)
00861 {
00862 this->mElementWeightedDirections.resize(num_elements);
00863 }
00864
00865 this->mBoundaryElementWeightedDirections.resize(num_boundary_elements);
00866
00867 this->mElementJacobianDeterminants.resize(num_elements);
00868 this->mBoundaryElementJacobianDeterminants.resize(num_boundary_elements);
00869
00870
00871 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00872 iter != this->GetElementIteratorEnd();
00873 ++iter)
00874 {
00875 unsigned index = iter->GetIndex();
00876 iter->CalculateInverseJacobian(this->mElementJacobians[index], this->mElementJacobianDeterminants[index], this->mElementInverseJacobians[index]);
00877 }
00878
00879 if (ELEMENT_DIM < SPACE_DIM)
00880 {
00881 for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
00882 iter != this->GetElementIteratorEnd();
00883 ++iter)
00884 {
00885 unsigned index = iter->GetIndex();
00886 iter->CalculateWeightedDirection(this->mElementWeightedDirections[index], this->mElementJacobianDeterminants[index]);
00887 }
00888 }
00889
00890 for ( typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator itb = this->GetBoundaryElementIteratorBegin();
00891 itb != this->GetBoundaryElementIteratorEnd();
00892 itb++)
00893 {
00894 unsigned index = (*itb)->GetIndex();
00895 (*itb)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[index], this->mBoundaryElementJacobianDeterminants[index]);
00896 }
00897 }
00898
00899 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00900 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetJacobianForElement(unsigned elementIndex, c_matrix<double, SPACE_DIM, SPACE_DIM>& rJacobian, double& rJacobianDeterminant) const
00901 {
00902 assert(ELEMENT_DIM <= SPACE_DIM);
00903 assert(elementIndex < this->mElementJacobians.size());
00904 rJacobian = this->mElementJacobians[elementIndex];
00905 rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
00906 }
00907
00908 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00909 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetInverseJacobianForElement(unsigned elementIndex, c_matrix<double, SPACE_DIM, ELEMENT_DIM>& rJacobian, double& rJacobianDeterminant, c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian) const
00910 {
00911 assert(ELEMENT_DIM <= SPACE_DIM);
00912 assert(elementIndex < this->mElementInverseJacobians.size());
00913 rInverseJacobian = this->mElementInverseJacobians[elementIndex];
00914 rJacobian = this->mElementJacobians[elementIndex];
00915 rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
00916 }
00917
00918 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00919 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetWeightedDirectionForElement(unsigned elementIndex, c_vector<double, SPACE_DIM>& rWeightedDirection, double& rJacobianDeterminant) const
00920 {
00921 assert(ELEMENT_DIM < SPACE_DIM);
00922 assert(elementIndex < this->mElementWeightedDirections.size());
00923 rWeightedDirection = this->mElementWeightedDirections[elementIndex];
00924 rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
00925 }
00926
00927 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00928 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetWeightedDirectionForBoundaryElement(unsigned elementIndex, c_vector<double, SPACE_DIM>& rWeightedDirection, double& rJacobianDeterminant) const
00929 {
00930 assert(elementIndex < this->mBoundaryElementWeightedDirections.size());
00931 rWeightedDirection = this->mBoundaryElementWeightedDirections[elementIndex];
00932 rJacobianDeterminant = this->mBoundaryElementJacobianDeterminants[elementIndex];
00933 }
00934
00935
00936 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00937 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::InitialiseTriangulateIo(triangulateio& mesherIo)
00938 {
00939 mesherIo.numberofpoints = 0;
00940 mesherIo.pointlist = NULL;
00941 mesherIo.numberofpointattributes = 0;
00942 mesherIo.pointattributelist = (double *) NULL;
00943 mesherIo.pointmarkerlist = (int *) NULL;
00944 mesherIo.numberofsegments = 0;
00945 mesherIo.numberofholes = 0;
00946 mesherIo.numberofregions = 0;
00947 mesherIo.trianglelist = (int *) NULL;
00948 mesherIo.triangleattributelist = (double *) NULL;
00949 mesherIo.numberoftriangleattributes = 0;
00950 mesherIo.edgelist = (int *) NULL;
00951 mesherIo.edgemarkerlist = (int *) NULL;
00952 }
00953 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00954 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::FreeTriangulateIo(triangulateio& mesherIo)
00955 {
00956 if (mesherIo.numberofpoints != 0)
00957 {
00958 mesherIo.numberofpoints=0;
00959 free(mesherIo.pointlist);
00960 }
00961
00962
00963 free(mesherIo.pointattributelist);
00964 free(mesherIo.pointmarkerlist);
00965 free(mesherIo.trianglelist);
00966 free(mesherIo.triangleattributelist);
00967 free(mesherIo.edgelist);
00968 free(mesherIo.edgemarkerlist);
00969 }
00970 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00971 template <class MESHER_IO>
00972 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ExportToMesher(NodeMap& map, MESHER_IO& mesherInput, int *elementList)
00973 {
00974 if (SPACE_DIM == 2)
00975 {
00976 mesherInput.pointlist = (double *) malloc(this->GetNumNodes() * SPACE_DIM * sizeof(double));
00977 }
00978 else
00979 {
00980 mesherInput.pointlist = new double[this->GetNumNodes() * SPACE_DIM];
00981 }
00982
00983 mesherInput.numberofpoints = this->GetNumNodes();
00984 unsigned new_index = 0;
00985 for (unsigned i=0; i<this->GetNumAllNodes(); i++)
00986 {
00987 if (this->mNodes[i]->IsDeleted())
00988 {
00989 map.SetDeleted(i);
00990 }
00991 else
00992 {
00993 map.SetNewIndex(i, new_index);
00994 for (unsigned j=0; j<SPACE_DIM; j++)
00995 {
00996 mesherInput.pointlist[SPACE_DIM*new_index + j] = this->mNodes[i]->rGetLocation()[j];
00997 }
00998 new_index++;
00999 }
01000 }
01001 if(elementList != NULL)
01002 {
01003 unsigned element_index=0u;
01004
01005 mesherInput.numberofcorners=ELEMENT_DIM+1;
01006 for (typename TetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::ElementIterator elem_iter = this->GetElementIteratorBegin();
01007 elem_iter != this->GetElementIteratorEnd();
01008 ++elem_iter)
01009 {
01010
01011 for (unsigned j=0; j<=ELEMENT_DIM; j++)
01012 {
01013 elementList[element_index*(ELEMENT_DIM+1) + j] = (*elem_iter).GetNodeGlobalIndex(j);
01014
01015 }
01016 element_index++;
01017 }
01018 }
01019 }
01020
01021 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
01022 template <class MESHER_IO>
01023 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ImportFromMesher(MESHER_IO& mesherOutput, unsigned numberOfElements, int *elementList, unsigned numberOfFaces, int *faceList, int *edgeMarkerList)
01024 {
01025 unsigned nodes_per_element = mesherOutput.numberofcorners;
01026
01027 assert( nodes_per_element == ELEMENT_DIM+1 || nodes_per_element == (ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2 );
01028
01029 Clear();
01030
01031 for (unsigned node_index=0; node_index<(unsigned)mesherOutput.numberofpoints; node_index++)
01032 {
01033 this->mNodes.push_back(new Node<SPACE_DIM>(node_index, &mesherOutput.pointlist[node_index * SPACE_DIM], false));
01034 }
01035
01036
01037 this->mElements.reserve(numberOfElements);
01038
01039 unsigned real_element_index=0;
01040 for (unsigned element_index=0; element_index<numberOfElements; element_index++)
01041 {
01042 std::vector<Node<SPACE_DIM>*> nodes;
01043 for (unsigned j=0; j<ELEMENT_DIM+1; j++)
01044 {
01045 unsigned global_node_index = elementList[element_index*(nodes_per_element) + j];
01046 assert(global_node_index < this->mNodes.size());
01047 nodes.push_back(this->mNodes[global_node_index]);
01048
01049 }
01050
01051
01052 Element<ELEMENT_DIM, SPACE_DIM>* p_element;
01053 try
01054 {
01055 p_element = new Element<ELEMENT_DIM, SPACE_DIM>(real_element_index, nodes);
01056
01057 this->mElements.push_back(p_element);
01058
01059
01060 for (unsigned j=ELEMENT_DIM+1; j<nodes_per_element; j++)
01061 {
01062 unsigned global_node_index = elementList[element_index*nodes_per_element + j];
01063 assert(global_node_index < this->mNodes.size());
01064 this->mElements[real_element_index]->AddNode( this->mNodes[global_node_index] );
01065 this->mNodes[global_node_index]->AddElement(real_element_index);
01066 this->mNodes[global_node_index]->MarkAsInternal();
01067 }
01068 real_element_index++;
01069 }
01070 catch (Exception &e)
01071 {
01072
01073 assert(SPACE_DIM == 3);
01074
01075
01076
01077
01078 }
01079 }
01080
01081
01082 unsigned next_boundary_element_index = 0;
01083 for (unsigned boundary_element_index=0; boundary_element_index<numberOfFaces; boundary_element_index++)
01084 {
01085
01086
01087 if (edgeMarkerList == NULL || edgeMarkerList[boundary_element_index] == 1)
01088 {
01089 std::vector<Node<SPACE_DIM>*> nodes;
01090 for (unsigned j=0; j<ELEMENT_DIM; j++)
01091 {
01092 unsigned global_node_index = faceList[boundary_element_index*ELEMENT_DIM + j];
01093 assert(global_node_index < this->mNodes.size());
01094 nodes.push_back(this->mNodes[global_node_index]);
01095 if (!nodes[j]->IsBoundaryNode())
01096 {
01097 nodes[j]->SetAsBoundaryNode();
01098 this->mBoundaryNodes.push_back(nodes[j]);
01099 }
01100 }
01101
01102
01103 BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>* p_b_element;
01104 try
01105 {
01106 p_b_element = new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(next_boundary_element_index, nodes);
01107 this->mBoundaryElements.push_back(p_b_element);
01108 next_boundary_element_index++;
01109 }
01110 catch (Exception &e)
01111 {
01112
01113 assert(SPACE_DIM == 3);
01114 }
01115 }
01116 }
01117
01118 this->RefreshJacobianCachedData();
01119 }
01120
01121
01122
01124
01126
01127 template class TetrahedralMesh<1,1>;
01128 template class TetrahedralMesh<1,2>;
01129 template class TetrahedralMesh<1,3>;
01130 template class TetrahedralMesh<2,2>;
01131 template class TetrahedralMesh<2,3>;
01132 template class TetrahedralMesh<3,3>;
01133
01138 template void TetrahedralMesh<2,2>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
01139 template void TetrahedralMesh<2,2>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
01140
01141 template void TetrahedralMesh<3,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
01142 template void TetrahedralMesh<3,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
01143
01144
01145 template void TetrahedralMesh<1,2>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
01146 template void TetrahedralMesh<1,2>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
01147
01148 template void TetrahedralMesh<1,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
01149 template void TetrahedralMesh<1,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
01150 template void TetrahedralMesh<2,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
01151 template void TetrahedralMesh<2,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
01152
01153
01154
01155 template void TetrahedralMesh<3u, 3u>::ImportFromMesher<triangulateio>(triangulateio&, unsigned int, int*, unsigned int, int*, int*);
01156 template void TetrahedralMesh<1u, 1u>::ImportFromMesher<triangulateio>(triangulateio&, unsigned int, int*, unsigned int, int*, int*);
01157 template void TetrahedralMesh<1u, 1u>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned int, int*, unsigned int, int*, int*);
01158 template void TetrahedralMesh<2u, 2u>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned int, int*, unsigned int, int*, int*);
01159 template void TetrahedralMesh<1u, 1u>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
01160
01161
01162
01163
01164
01165
01166
01167
01168
01169
01170
01176
01177 #define CHASTE_SERIALIZATION_CPP
01178 #include "SerializationExportWrapper.hpp"
01179 EXPORT_TEMPLATE_CLASS_ALL_DIMS(TetrahedralMesh)