Chaste Release::3.1
MutableMesh.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include <map>
00037 #include <cstring>
00038 
00039 #include "MutableMesh.hpp"
00040 #include "OutputFileHandler.hpp"
00041 #include "TrianglesMeshReader.hpp"
00042 
00043 //Jonathan Shewchuk's triangle and Hang Si's tetgen
00044 #define REAL double
00045 #define VOID void
00046 #include "triangle.h"
00047 #include "tetgen.h"
00048 #undef REAL
00049 #undef VOID
00050 
00051 
00052 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00053 MutableMesh<ELEMENT_DIM, SPACE_DIM>::MutableMesh()
00054     : mAddedNodes(false)
00055 {
00056     this->mMeshChangesDuringSimulation = true;
00057 }
00058 
00059 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00060 MutableMesh<ELEMENT_DIM, SPACE_DIM>::MutableMesh(std::vector<Node<SPACE_DIM> *> nodes)
00061 {
00062     this->mMeshChangesDuringSimulation = true;
00063     Clear();
00064     for (unsigned index=0; index<nodes.size(); index++)
00065     {
00066         Node<SPACE_DIM>* p_temp_node = nodes[index];
00067         this->mNodes.push_back(p_temp_node);
00068     }
00069     mAddedNodes = true;
00070     NodeMap node_map(nodes.size());
00071     ReMesh(node_map);
00072 }
00073 
00074 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00075 MutableMesh<ELEMENT_DIM, SPACE_DIM>::~MutableMesh()
00076 {
00077     Clear();
00078 }
00079 
00080 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00081 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::AddNode(Node<SPACE_DIM>* pNewNode)
00082 {
00083     if (mDeletedNodeIndices.empty())
00084     {
00085         pNewNode->SetIndex(this->mNodes.size());
00086         this->mNodes.push_back(pNewNode);
00087     }
00088     else
00089     {
00090         unsigned index = mDeletedNodeIndices.back();
00091         pNewNode->SetIndex(index);
00092         mDeletedNodeIndices.pop_back();
00093         delete this->mNodes[index];
00094         this->mNodes[index] = pNewNode;
00095     }
00096     mAddedNodes = true;
00097     return pNewNode->GetIndex();
00098 }
00099 
00100 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00101 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::Clear()
00102 {
00103     mDeletedElementIndices.clear();
00104     mDeletedBoundaryElementIndices.clear();
00105     mDeletedNodeIndices.clear();
00106     mAddedNodes = false;
00107 
00108     TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::Clear();
00109 }
00110 
00111 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00112 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::GetNumBoundaryElements() const
00113 {
00114     return this->mBoundaryElements.size() - mDeletedBoundaryElementIndices.size();
00115 }
00116 
00117 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00118 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::GetNumElements() const
00119 {
00120     return this->mElements.size() - mDeletedElementIndices.size();
00121 }
00122 
00123 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00124 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::GetNumNodes() const
00125 {
00126     return this->mNodes.size() - mDeletedNodeIndices.size();
00127 }
00128 
00135 template<>
00136 void MutableMesh<1, 1>::RescaleMeshFromBoundaryNode(ChastePoint<1> updatedPoint, unsigned boundaryNodeIndex)
00137 {
00138     assert(this->GetNode(boundaryNodeIndex)->IsBoundaryNode());
00139     double scaleFactor = updatedPoint[0] / this->GetNode(boundaryNodeIndex)->GetPoint()[0];
00140     double temp;
00141     for (unsigned i=0; i < boundaryNodeIndex+1; i++)
00142     {
00143         temp = scaleFactor * this->mNodes[i]->GetPoint()[0];
00144         ChastePoint<1> newPoint(temp);
00145         this->mNodes[i]->SetPoint(newPoint);
00146     }
00147     this->RefreshMesh();
00148 }
00149 
00150 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00151 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::SetNode(unsigned index,
00152         ChastePoint<SPACE_DIM> point,
00153         bool concreteMove)
00154 {
00155     this->mNodes[index]->SetPoint(point);
00156 
00157     if (concreteMove)
00158     {
00159         for (typename Node<SPACE_DIM>::ContainingElementIterator it = this->mNodes[index]->ContainingElementsBegin();
00160              it != this->mNodes[index]->ContainingElementsEnd();
00161              ++it)
00162         {
00163             if (ELEMENT_DIM == SPACE_DIM)
00164             {
00165                 try
00166                 {
00167                     this->GetElement(*it)->CalculateInverseJacobian(this->mElementJacobians[ (*it) ],
00168                                                               this->mElementJacobianDeterminants[ (*it) ],
00169                                                               this->mElementInverseJacobians[ (*it) ]);
00170                 }
00171                 catch (Exception e)
00172                 {
00173                         EXCEPTION("Moving node caused an element to have a non-positive Jacobian determinant");
00174                 }
00175             }
00176             else
00177             {
00178                 c_vector<double,SPACE_DIM> previous_direction = this->mElementWeightedDirections[ (*it) ];
00179 
00180                 this->GetElement(*it)->CalculateWeightedDirection(this->mElementWeightedDirections[ (*it) ],
00181                                                             this->mElementJacobianDeterminants[ (*it) ]);
00182 
00183                 if ( inner_prod(previous_direction, this->mElementWeightedDirections[ (*it) ]) < 0)
00184                 {
00185                     EXCEPTION("Moving node caused an subspace element to change direction");
00186                 }
00187 
00188             }
00189         }
00190         for (typename Node<SPACE_DIM>::ContainingBoundaryElementIterator it = this->mNodes[index]->ContainingBoundaryElementsBegin();
00191              it != this->mNodes[index]->ContainingBoundaryElementsEnd();
00192              ++it)
00193         {
00194             try
00195             {
00196                 this->GetBoundaryElement(*it)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[ (*it) ],
00197                                                                     this->mBoundaryElementJacobianDeterminants[ (*it) ]);
00198             }
00199             catch (Exception e)
00200             {
00201                 EXCEPTION("Moving node caused a boundary element to have a non-positive Jacobian determinant");
00202             }
00203         }
00204     }
00205 }
00206 
00207 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00208 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteNode(unsigned index)
00209 {
00210     if (this->mNodes[index]->IsDeleted())
00211     {
00212         EXCEPTION("Trying to delete a deleted node");
00213     }
00214     unsigned target_index = (unsigned)(-1);
00215     bool found_target=false;
00216     for (typename Node<SPACE_DIM>::ContainingElementIterator it = this->mNodes[index]->ContainingElementsBegin();
00217          !found_target && it != this->mNodes[index]->ContainingElementsEnd();
00218          ++it)
00219     {
00220         Element <ELEMENT_DIM,SPACE_DIM>* p_element = this->GetElement(*it);
00221         for (unsigned i=0; i<=ELEMENT_DIM && !found_target; i++)
00222         {
00223             target_index = p_element->GetNodeGlobalIndex(i);
00224             try
00225             {
00226                 MoveMergeNode(index, target_index, false);
00227                 found_target = true;
00228             }
00229             catch (Exception e)
00230             {
00231                 // Just try the next node
00232             }
00233         }
00234     }
00235     if (!found_target)
00236     {
00237         EXCEPTION("Failure to delete node");
00238     }
00239 
00240     MoveMergeNode(index, target_index);
00241 }
00242 
00243 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00244 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteNodePriorToReMesh(unsigned index)
00245 {
00246     this->mNodes[index]->MarkAsDeleted();
00247     mDeletedNodeIndices.push_back(index);
00248 }
00249 
00250 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00251 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::MoveMergeNode(unsigned index,
00252         unsigned targetIndex,
00253         bool concreteMove)
00254 {
00255     if (this->mNodes[index]->IsDeleted())
00256     {
00257         EXCEPTION("Trying to move a deleted node");
00258     }
00259 
00260     if (index == targetIndex)
00261     {
00262         EXCEPTION("Trying to merge a node with itself");
00263     }
00264     if (this->mNodes[index]->IsBoundaryNode())
00265     {
00266         if (!this->mNodes[targetIndex]->IsBoundaryNode())
00267         {
00268             EXCEPTION("A boundary node can only be moved on to another boundary node");
00269         }
00270     }
00271     std::set<unsigned> unshared_element_indices;
00272     std::set_difference(this->mNodes[index]->rGetContainingElementIndices().begin(),
00273                         this->mNodes[index]->rGetContainingElementIndices().end(),
00274                         this->mNodes[targetIndex]->rGetContainingElementIndices().begin(),
00275                         this->mNodes[targetIndex]->rGetContainingElementIndices().end(),
00276                         std::inserter(unshared_element_indices, unshared_element_indices.begin()));
00277 
00278 
00279     if (unshared_element_indices.size() == this->mNodes[index]->rGetContainingElementIndices().size())
00280     {
00281         EXCEPTION("These nodes cannot be merged since they are not neighbours");
00282     }
00283 
00284     std::set<unsigned> unshared_boundary_element_indices;
00285     std::set_difference(this->mNodes[index]->rGetContainingBoundaryElementIndices().begin(),
00286                         this->mNodes[index]->rGetContainingBoundaryElementIndices().end(),
00287                         this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().begin(),
00288                         this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().end(),
00289                         std::inserter(unshared_boundary_element_indices, unshared_boundary_element_indices.begin()));
00290 
00291 
00292     if (this->mNodes[index]->IsBoundaryNode())
00293     {
00294         if (unshared_boundary_element_indices.size()
00295             == this->mNodes[index]->rGetContainingBoundaryElementIndices().size())
00296         {
00297             //May be redundant (only thrown in 1D tests)
00298             EXCEPTION("These nodes cannot be merged since they are not neighbours on the boundary");
00299         }
00300     }
00301 
00302     this->mNodes[index]->rGetModifiableLocation() = this->mNodes[targetIndex]->rGetLocation();
00303 
00304     for (std::set<unsigned>::const_iterator element_iter=unshared_element_indices.begin();
00305              element_iter != unshared_element_indices.end();
00306              element_iter++)
00307     {
00308         try
00309         {
00310             if (SPACE_DIM == ELEMENT_DIM)
00311             {
00312                 this->GetElement(*element_iter)->CalculateInverseJacobian(this->mElementJacobians[(*element_iter)],
00313                                                                           this->mElementJacobianDeterminants[(*element_iter)],
00314                                                                           this->mElementInverseJacobians[ (*element_iter) ]);
00315             }
00316             else
00317             {
00318                 this->GetElement(*element_iter)->CalculateWeightedDirection(this->mElementWeightedDirections[(*element_iter)],
00319                                                                             this->mElementJacobianDeterminants[(*element_iter)]);
00320             }
00321 
00322             if (concreteMove)
00323             {
00324                 this->GetElement(*element_iter)->ReplaceNode(this->mNodes[index], this->mNodes[targetIndex]);
00325             }
00326 
00327         }
00328         catch (Exception e)
00329         {
00330             EXCEPTION("Moving node caused an element to have a non-positive Jacobian determinant");
00331         }
00332     }
00333 
00334     for (std::set<unsigned>::const_iterator boundary_element_iter=
00335                  unshared_boundary_element_indices.begin();
00336              boundary_element_iter != unshared_boundary_element_indices.end();
00337              boundary_element_iter++)
00338     {
00339 
00340         this->GetBoundaryElement(*boundary_element_iter)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[(*boundary_element_iter)],
00341                                                                                      this->mBoundaryElementJacobianDeterminants[(*boundary_element_iter)]);
00342 
00343         if (concreteMove)
00344         {
00345             this->GetBoundaryElement(*boundary_element_iter)->ReplaceNode(this->mNodes[index], this->mNodes[targetIndex]);
00346         }
00347     }
00348 
00349     std::set<unsigned> shared_element_indices;
00350     std::set_intersection(this->mNodes[index]->rGetContainingElementIndices().begin(),
00351                           this->mNodes[index]->rGetContainingElementIndices().end(),
00352                           this->mNodes[targetIndex]->rGetContainingElementIndices().begin(),
00353                           this->mNodes[targetIndex]->rGetContainingElementIndices().end(),
00354                           std::inserter(shared_element_indices, shared_element_indices.begin()));
00355 
00356     for (std::set<unsigned>::const_iterator element_iter=shared_element_indices.begin();
00357              element_iter != shared_element_indices.end();
00358              element_iter++)
00359     {
00360         if (concreteMove)
00361         {
00362             this->GetElement(*element_iter)->MarkAsDeleted();
00363             this->mElementJacobianDeterminants[ (*element_iter) ] = 0.0; //This used to be done in MarkAsDeleted
00364             mDeletedElementIndices.push_back(*element_iter);
00365         }
00366         else
00367         {
00368             this->mElementJacobianDeterminants[ (*element_iter) ] = 0.0;
00369         }
00370     }
00371 
00372 
00373     std::set<unsigned> shared_boundary_element_indices;
00374     std::set_intersection(this->mNodes[index]->rGetContainingBoundaryElementIndices().begin(),
00375                           this->mNodes[index]->rGetContainingBoundaryElementIndices().end(),
00376                           this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().begin(),
00377                           this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().end(),
00378                           std::inserter(shared_boundary_element_indices, shared_boundary_element_indices.begin()));
00379 
00380     for (std::set<unsigned>::const_iterator boundary_element_iter=shared_boundary_element_indices.begin();
00381              boundary_element_iter != shared_boundary_element_indices.end();
00382              boundary_element_iter++)
00383     {
00384         if (concreteMove)
00385         {
00386             this->GetBoundaryElement(*boundary_element_iter)->MarkAsDeleted();
00387             this->mBoundaryElementJacobianDeterminants[ (*boundary_element_iter) ] = 0.0; //This used to be done in MarkAsDeleted
00388             mDeletedBoundaryElementIndices.push_back(*boundary_element_iter);
00389         }
00390         else
00391         {
00392             this->mBoundaryElementJacobianDeterminants[ (*boundary_element_iter) ] = 0.0;
00393             this->mBoundaryElementWeightedDirections[ (*boundary_element_iter) ] = zero_vector<double>(SPACE_DIM);
00394         }
00395     }
00396 
00397     if (concreteMove)
00398     {
00399         this->mNodes[index]->MarkAsDeleted();
00400         mDeletedNodeIndices.push_back(index);
00401     }
00402 }
00403 
00404 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00405 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::RefineElement(
00406     Element<ELEMENT_DIM,SPACE_DIM>* pElement,
00407     ChastePoint<SPACE_DIM> point)
00408 {
00409     //Check that the point is in the element
00410     if (pElement->IncludesPoint(point, true) == false)
00411     {
00412         EXCEPTION("RefineElement could not be started (point is not in element)");
00413     }
00414 
00415     // Add a new node from the point that is passed to RefineElement
00416     unsigned new_node_index = AddNode(new Node<SPACE_DIM>(0, point.rGetLocation()));
00417     // Note: the first argument is the index of the node, which is going to be
00418     //       overriden by AddNode, so it can safely be ignored
00419 
00420     // This loop constructs the extra elements which are going to fill the space
00421     for (unsigned i = 0; i < ELEMENT_DIM; i++)
00422     {
00423 
00424         // First, make a copy of the current element making sure we update its index
00425         unsigned new_elt_index;
00426         if (mDeletedElementIndices.empty())
00427         {
00428             new_elt_index = this->mElements.size();
00429         }
00430         else
00431         {
00432             new_elt_index = mDeletedElementIndices.back();
00433             mDeletedElementIndices.pop_back();
00434         }
00435 
00436         Element<ELEMENT_DIM,SPACE_DIM>* p_new_element=
00437             new Element<ELEMENT_DIM,SPACE_DIM>(*pElement, new_elt_index);
00438 
00439         // Second, update the node in the element with the new one
00440         p_new_element->UpdateNode(ELEMENT_DIM-1-i, this->mNodes[new_node_index]);
00441 
00442         // Third, add the new element to the set
00443         if ((unsigned) new_elt_index == this->mElements.size())
00444         {
00445             this->mElements.push_back(p_new_element);
00446         }
00447         else
00448         {
00449             delete this->mElements[new_elt_index];
00450             this->mElements[new_elt_index] = p_new_element;
00451         }
00452     }
00453 
00454     // Lastly, update the last node in the element to be refined
00455     pElement->UpdateNode(ELEMENT_DIM, this->mNodes[new_node_index]);
00456 
00457     return new_node_index;
00458 }
00459 
00460 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00461 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteBoundaryNodeAt(unsigned index)
00462 {
00463     if (!this->mNodes[index]->IsBoundaryNode() )
00464     {
00465         EXCEPTION(" You may only delete a boundary node ");
00466     }
00467 
00468     this->mNodes[index]->MarkAsDeleted();
00469     mDeletedNodeIndices.push_back(index);
00470 
00471     // Update the boundary node vector
00472     typename std::vector<Node<SPACE_DIM>*>::iterator b_node_iter
00473     = std::find(this->mBoundaryNodes.begin(), this->mBoundaryNodes.end(), this->mNodes[index]);
00474     this->mBoundaryNodes.erase(b_node_iter);
00475 
00476     // Remove boundary elements containing this node
00477     std::set<unsigned> boundary_element_indices = this->mNodes[index]->rGetContainingBoundaryElementIndices();
00478     std::set<unsigned>::const_iterator boundary_element_indices_iterator = boundary_element_indices.begin();
00479     while (boundary_element_indices_iterator != boundary_element_indices.end())
00480     {
00481         BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>* p_boundary_element = this->GetBoundaryElement(*boundary_element_indices_iterator);
00482         p_boundary_element->MarkAsDeleted();
00483         mDeletedBoundaryElementIndices.push_back(*boundary_element_indices_iterator);
00484         boundary_element_indices_iterator++;
00485     }
00486 
00487     // Remove elements containing this node
00488     std::set<unsigned> element_indices = this->mNodes[index]->rGetContainingElementIndices();
00489     std::set<unsigned>::const_iterator element_indices_iterator = element_indices.begin();
00490     while (element_indices_iterator != element_indices.end())
00491     {
00492         Element<ELEMENT_DIM, SPACE_DIM>* p_element = this->GetElement(*element_indices_iterator);
00493         for (unsigned i=0; i<p_element->GetNumNodes(); i++)
00494         {
00495             Node<SPACE_DIM>* p_node = p_element->GetNode(i);
00496             if (!p_node->IsDeleted())
00497             {
00498                 p_node->SetAsBoundaryNode();
00499                 // Update the boundary node vector
00500                 this->mBoundaryNodes.push_back(p_node);
00501             }
00502         }
00503         p_element->MarkAsDeleted();
00504         mDeletedElementIndices.push_back(p_element->GetIndex());
00505         element_indices_iterator++;
00506     }
00507 }
00508 
00509 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00510 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::ReIndex(NodeMap& map)
00511 {
00512     assert(!mAddedNodes);
00513     map.Resize(this->GetNumAllNodes());
00514 
00515     std::vector<Element<ELEMENT_DIM, SPACE_DIM> *> live_elements;
00516 
00517     for (unsigned i=0; i<this->mElements.size(); i++)
00518     {
00519         assert(i==this->mElements[i]->GetIndex()); // We need this to be true to be able to reindex the Jacobian cache
00520         if (this->mElements[i]->IsDeleted())
00521         {
00522             delete this->mElements[i];
00523         }
00524         else
00525         {
00526             live_elements.push_back(this->mElements[i]);
00527 
00528             unsigned this_element_index = this->mElements[i]->GetIndex();
00529             if (SPACE_DIM == ELEMENT_DIM)
00530             {
00531                 this->mElementJacobians[live_elements.size()-1] = this->mElementJacobians[this_element_index];
00532                 this->mElementInverseJacobians[live_elements.size()-1] = this->mElementInverseJacobians[this_element_index];
00533             }
00534             else
00535             {
00536                 this->mElementWeightedDirections[live_elements.size()-1] = this->mElementWeightedDirections[this_element_index];
00537             }
00538             this->mElementJacobianDeterminants[live_elements.size()-1] = this->mElementJacobianDeterminants[this_element_index];
00539         }
00540     }
00541 
00542     assert(mDeletedElementIndices.size() == this->mElements.size()-live_elements.size());
00543     mDeletedElementIndices.clear();
00544     this->mElements = live_elements;
00545     unsigned num_elements = this->mElements.size();
00546 
00547     if (SPACE_DIM == ELEMENT_DIM)
00548     {
00549         this->mElementJacobians.resize(num_elements);
00550         this->mElementInverseJacobians.resize(num_elements);
00551     }
00552     else
00553     {
00554         this->mElementWeightedDirections.resize(num_elements);
00555     }
00556     this->mElementJacobianDeterminants.resize(num_elements);
00557 
00558     std::vector<Node<SPACE_DIM> *> live_nodes;
00559     for (unsigned i=0; i<this->mNodes.size(); i++)
00560     {
00561         if (this->mNodes[i]->IsDeleted())
00562         {
00563             delete this->mNodes[i];
00564             map.SetDeleted(i);
00565         }
00566         else
00567         {
00568             live_nodes.push_back(this->mNodes[i]);
00569             // the nodes will have their index set to be the index into the live_nodes
00570             // vector further down
00571             map.SetNewIndex(i, (unsigned)(live_nodes.size()-1));
00572         }
00573     }
00574 
00575     assert(mDeletedNodeIndices.size() == this->mNodes.size()-live_nodes.size());
00576     this->mNodes = live_nodes;
00577     mDeletedNodeIndices.clear();
00578 
00579     std::vector<BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *> live_boundary_elements;
00580     for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00581     {
00582         if (this->mBoundaryElements[i]->IsDeleted())
00583         {
00584             delete this->mBoundaryElements[i];
00585         }
00586         else
00587         {
00588             live_boundary_elements.push_back(this->mBoundaryElements[i]);
00589 
00590             this->mBoundaryElementWeightedDirections[live_boundary_elements.size()-1] = this->mBoundaryElementWeightedDirections[this->mBoundaryElements[i]->GetIndex()];
00591             this->mBoundaryElementJacobianDeterminants[live_boundary_elements.size()-1] = this->mBoundaryElementJacobianDeterminants[this->mBoundaryElements[i]->GetIndex()];
00592         }
00593     }
00594 
00595     assert(mDeletedBoundaryElementIndices.size() == this->mBoundaryElements.size()-live_boundary_elements.size());
00596     this->mBoundaryElements = live_boundary_elements;
00597     mDeletedBoundaryElementIndices.clear();
00598 
00599     unsigned num_boundary_elements = this->mBoundaryElements.size();
00600 
00601     this->mBoundaryElementWeightedDirections.resize(num_boundary_elements);
00602     this->mBoundaryElementJacobianDeterminants.resize(num_boundary_elements);
00603 
00604     for (unsigned i=0; i<this->mNodes.size(); i++)
00605     {
00606         this->mNodes[i]->SetIndex(i);
00607     }
00608     for (unsigned i=0; i<this->mElements.size(); i++)
00609     {
00610 
00611         this->mElements[i]->ResetIndex(i);
00612     }
00613     for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00614     {
00615         this->mBoundaryElements[i]->ResetIndex(i);
00616     }
00617 }
00618 
00619 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00620 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::ReMesh(NodeMap& map)
00621 {
00622     // Make sure that we are in the correct dimension - this code will be eliminated at compile time
00623     #define COVERAGE_IGNORE
00624     assert( ELEMENT_DIM == SPACE_DIM );
00625     #undef COVERAGE_IGNORE
00626 
00627     // Avoid some triangle/tetgen errors: need at least four
00628     // nodes for tetgen, and at least three for triangle
00629     if (GetNumNodes() <= SPACE_DIM)
00630     {
00631         EXCEPTION("The number of nodes must exceed the spatial dimension.");
00632     }
00633 
00634     // Make sure the map is big enough
00635     map.Resize(this->GetNumAllNodes());
00636     if (mAddedNodes || !mDeletedNodeIndices.empty())
00637     {
00638         // Size of mesh is about to change
00639         if (this->mpDistributedVectorFactory)
00640         {
00641             delete this->mpDistributedVectorFactory;
00642             this->mpDistributedVectorFactory = new DistributedVectorFactory(this->GetNumNodes());
00643         }
00644     }
00645     if (SPACE_DIM==1)
00646     {
00647         // Store the node locations
00648         std::vector<c_vector<double, SPACE_DIM> > old_node_locations;
00649         unsigned new_index = 0;
00650         for (unsigned i=0; i<this->GetNumAllNodes(); i++)
00651         {
00652             if (this->mNodes[i]->IsDeleted())
00653             {
00654                 map.SetDeleted(i);
00655             }
00656             else
00657             {
00658                 map.SetNewIndex(i, new_index);
00659                 old_node_locations.push_back(this->mNodes[i]->rGetLocation());
00660                 new_index++;
00661             }
00662         }
00663 
00664         // Remove current data
00665         Clear();
00666 
00667         // Construct the nodes and boundary nodes
00668         for (unsigned node_index=0; node_index<old_node_locations.size(); node_index++)
00669         {
00670             Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(node_index, old_node_locations[node_index], false);
00671             this->mNodes.push_back(p_node);
00672 
00673             // As we're in 1D, the boundary nodes are simply at either end of the mesh
00674             if ( node_index==0 || node_index==old_node_locations.size()-1 )
00675             {
00676                 this->mBoundaryNodes.push_back(p_node);
00677             }
00678         }
00679 
00680         // Create a map between node indices and node locations
00681         std::map<double, unsigned> location_index_map;
00682         for (unsigned i=0; i<this->mNodes.size(); i++)
00683         {
00684             location_index_map[this->mNodes[i]->rGetLocation()[0]] = this->mNodes[i]->GetIndex();
00685         }
00686 
00687         // Use this map to generate a vector of node indices that are ordered spatially
00688         std::vector<unsigned> node_indices_ordered_spatially;
00689         for (std::map<double, unsigned>::iterator iter = location_index_map.begin();
00690              iter != location_index_map.end();
00691              ++iter)
00692         {
00693             node_indices_ordered_spatially.push_back(iter->second);
00694         }
00695 
00696         // Construct the elements
00697         this->mElements.reserve(old_node_locations.size()-1);
00698         for (unsigned element_index=0; element_index<old_node_locations.size()-1; element_index++)
00699         {
00700             std::vector<Node<SPACE_DIM>*> nodes;
00701             for (unsigned j=0; j<2; j++)
00702             {
00703                 unsigned global_node_index = node_indices_ordered_spatially[element_index + j];
00704                 assert(global_node_index < this->mNodes.size());
00705                 nodes.push_back(this->mNodes[global_node_index]);
00706             }
00707             this->mElements.push_back(new Element<ELEMENT_DIM, SPACE_DIM>(element_index, nodes));
00708         }
00709 
00710         // Construct the two boundary elements - as we're in 1D, these are simply at either end of the mesh
00711         std::vector<Node<SPACE_DIM>*> nodes;
00712         nodes.push_back(this->mNodes[0]);
00713         this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(0, nodes));
00714 
00715         nodes.clear();
00716         nodes.push_back(this->mNodes[old_node_locations.size()-1]);
00717         this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(1, nodes));
00718 
00719         this->RefreshJacobianCachedData();
00720     }
00721     else if (SPACE_DIM==2)  // In 2D, remesh using triangle via library calls
00722     {
00723         struct triangulateio mesher_input, mesher_output;
00724         this->InitialiseTriangulateIo(mesher_input);
00725         this->InitialiseTriangulateIo(mesher_output);
00726 
00727         this->ExportToMesher(map, mesher_input);
00728 
00729         // Library call
00730         triangulate((char*)"Qze", &mesher_input, &mesher_output, NULL);
00731 
00732         this->ImportFromMesher(mesher_output, mesher_output.numberoftriangles, mesher_output.trianglelist, mesher_output.numberofedges, mesher_output.edgelist, mesher_output.edgemarkerlist);
00733 
00734         //Tidy up triangle
00735         this->FreeTriangulateIo(mesher_input);
00736         this->FreeTriangulateIo(mesher_output);
00737     }
00738     else // in 3D, remesh using tetgen
00739     {
00740 
00741         struct tetgen::tetgenio mesher_input, mesher_output;
00742 
00743         this->ExportToMesher(map, mesher_input);
00744 
00745         // Library call
00746         tetgen::tetrahedralize((char*)"Qz", &mesher_input, &mesher_output);
00747 
00748         this->ImportFromMesher(mesher_output, mesher_output.numberoftetrahedra, mesher_output.tetrahedronlist, mesher_output.numberoftrifaces, mesher_output.trifacelist, NULL);
00749     }
00750 }
00751 
00752 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00753 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::ReMesh()
00754 {
00755     NodeMap map(GetNumNodes());
00756     ReMesh(map);
00757 }
00758 
00759 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00760 bool MutableMesh<ELEMENT_DIM, SPACE_DIM>::CheckIsVoronoi(Element<ELEMENT_DIM, SPACE_DIM>* pElement, double maxPenetration)
00761 {
00762     assert(ELEMENT_DIM == SPACE_DIM);
00763     unsigned num_nodes = pElement->GetNumNodes();
00764     std::set<unsigned> neighbouring_elements_indices;
00765     std::set< Element<ELEMENT_DIM,SPACE_DIM> *> neighbouring_elements;
00766     std::set<unsigned> neighbouring_nodes_indices;
00767 
00768     // Form a set of neighbouring elements via the nodes
00769     for (unsigned i=0; i<num_nodes; i++)
00770     {
00771         Node<SPACE_DIM>* p_node = pElement->GetNode(i);
00772         neighbouring_elements_indices = p_node->rGetContainingElementIndices();
00773         for (std::set<unsigned>::const_iterator it = neighbouring_elements_indices.begin();
00774              it != neighbouring_elements_indices.end();
00775              ++it)
00776         {
00777             neighbouring_elements.insert(this->GetElement(*it));
00778         }
00779     }
00780     neighbouring_elements.erase(pElement);
00781 
00782     // For each neighbouring element find the supporting nodes
00783     typedef typename std::set<Element<ELEMENT_DIM,SPACE_DIM> *>::const_iterator ElementIterator;
00784 
00785     for (ElementIterator it = neighbouring_elements.begin();
00786          it != neighbouring_elements.end();
00787          ++it)
00788     {
00789         for (unsigned i=0; i<num_nodes; i++)
00790         {
00791             neighbouring_nodes_indices.insert((*it)->GetNodeGlobalIndex(i));
00792         }
00793     }
00794 
00795     // Remove the nodes that support this element
00796     for (unsigned i = 0; i < num_nodes; i++)
00797     {
00798         neighbouring_nodes_indices.erase(pElement->GetNodeGlobalIndex(i));
00799     }
00800 
00801     // Get the circumsphere information
00802     c_vector<double, SPACE_DIM+1> this_circum_centre;
00803 
00804     this_circum_centre = pElement->CalculateCircumsphere(this->mElementJacobians[pElement->GetIndex()], this->mElementInverseJacobians[pElement->GetIndex()]);
00805 
00806     // Copy the actualy circumcentre into a smaller vector
00807     c_vector<double, ELEMENT_DIM> circum_centre;
00808     for (unsigned i=0; i<ELEMENT_DIM; i++)
00809     {
00810         circum_centre[i] = this_circum_centre[i];
00811     }
00812 
00813     for (std::set<unsigned>::const_iterator it = neighbouring_nodes_indices.begin();
00814          it != neighbouring_nodes_indices.end();
00815          ++it)
00816     {
00817         c_vector<double, ELEMENT_DIM> node_location = this->GetNode(*it)->rGetLocation();
00818 
00819         // Calculate vector from circumcenter to node
00820         node_location -= circum_centre;
00821 
00822         // This is to calculate the squared distance betweeen them
00823         double squared_distance = inner_prod(node_location, node_location);
00824 
00825         // If the squared idstance is less than the elements circum-radius(squared),
00826         // then the Voronoi property is violated.
00827         if (squared_distance < this_circum_centre[ELEMENT_DIM])
00828         {
00829             // We know the node is inside the circumsphere, but we don't know how far
00830             double radius = sqrt(this_circum_centre[ELEMENT_DIM]);
00831             double distance = radius - sqrt(squared_distance);
00832 
00833             // If the node penetration is greater than supplied maximum penetration factor
00834             if (distance/radius > maxPenetration)
00835             {
00836                 return false;
00837             }
00838         }
00839     }
00840     return true;
00841 }
00842 
00843 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00844 bool MutableMesh<ELEMENT_DIM, SPACE_DIM>::CheckIsVoronoi(double maxPenetration)
00845 {
00846     // Looping through all the elements in the mesh
00848     for (unsigned i=0; i<this->mElements.size(); i++)
00849     {
00850         // Check if the element is not deleted
00851         if (!this->mElements[i]->IsDeleted())
00852         {
00853             // Checking the Voronoi of the Element
00854             if (CheckIsVoronoi(this->mElements[i], maxPenetration) == false)
00855             {
00856                 return false;
00857             }
00858         }
00859     }
00860     return true;
00861 }
00862 
00863 
00865 // Explicit instantiation
00867 
00868 template class MutableMesh<1,1>;
00869 template class MutableMesh<1,2>;
00870 template class MutableMesh<1,3>;
00871 template class MutableMesh<2,2>;
00872 template class MutableMesh<2,3>;
00873 template class MutableMesh<3,3>;
00874 
00875 
00876 // Serialization for Boost >= 1.36
00877 #include "SerializationExportWrapperForCpp.hpp"
00878 EXPORT_TEMPLATE_CLASS_ALL_DIMS(MutableMesh)