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