MutableMesh.cpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, 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 
00042 //Jonathan Shewchuk's triangle and Hang Si's tetgen
00043 #define REAL double
00044 #define VOID void
00045 #include "triangle.h"
00046 #include "tetgen.h"
00047 #undef REAL
00048 #undef VOID
00049 
00050 
00051 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00052 MutableMesh<ELEMENT_DIM, SPACE_DIM>::MutableMesh()
00053     : mAddedNodes(false)
00054 {
00055     this->mMeshChangesDuringSimulation = true;
00056 }
00057 
00058 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00059 MutableMesh<ELEMENT_DIM, SPACE_DIM>::MutableMesh(std::vector<Node<SPACE_DIM> *> nodes)
00060 {
00061     this->mMeshChangesDuringSimulation = true;
00062     Clear();
00063     for (unsigned index=0; index<nodes.size(); index++)
00064     {
00065         Node<SPACE_DIM>* p_temp_node = nodes[index];
00066         this->mNodes.push_back(p_temp_node);
00067     }
00068     mAddedNodes = true;
00069     NodeMap node_map(nodes.size());
00070     ReMesh(node_map);
00071 }
00072 
00073 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00074 MutableMesh<ELEMENT_DIM, SPACE_DIM>::~MutableMesh()
00075 {
00076     Clear();
00077 }
00078 
00079 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00080 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::AddNode(Node<SPACE_DIM>* pNewNode)
00081 {
00082     if (mDeletedNodeIndices.empty())
00083     {
00084         pNewNode->SetIndex(this->mNodes.size());
00085         this->mNodes.push_back(pNewNode);
00086     }
00087     else
00088     {
00089         unsigned index = mDeletedNodeIndices.back();
00090         pNewNode->SetIndex(index);
00091         mDeletedNodeIndices.pop_back();
00092         delete this->mNodes[index];
00093         this->mNodes[index] = pNewNode;
00094     }
00095     mAddedNodes = true;
00096     return pNewNode->GetIndex();
00097 }
00098 
00099 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00100 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::AddElement(Element<ELEMENT_DIM,SPACE_DIM>* pNewElement)
00101 {
00102     unsigned new_elt_index;
00103 
00104     if (mDeletedElementIndices.empty())
00105     {
00106         new_elt_index = this->mElements.size();
00107         this->mElements.push_back(pNewElement);
00108         pNewElement->ResetIndex(new_elt_index);
00109     }
00110     else
00111     {
00112         unsigned index = mDeletedElementIndices.back();
00113         pNewElement->ResetIndex(index);
00114         mDeletedElementIndices.pop_back();
00115         delete this->mElements[index];
00116         this->mElements[index] = pNewElement;
00117     }
00118 
00119     return pNewElement->GetIndex();
00120 }
00121 
00122 
00123 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00124 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::Clear()
00125 {
00126     mDeletedElementIndices.clear();
00127     mDeletedBoundaryElementIndices.clear();
00128     mDeletedNodeIndices.clear();
00129     mAddedNodes = false;
00130 
00131     TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::Clear();
00132 }
00133 
00134 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00135 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::GetNumBoundaryElements() const
00136 {
00137     return this->mBoundaryElements.size() - mDeletedBoundaryElementIndices.size();
00138 }
00139 
00140 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00141 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::GetNumElements() const
00142 {
00143     return this->mElements.size() - mDeletedElementIndices.size();
00144 }
00145 
00146 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00147 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::GetNumNodes() const
00148 {
00149     return this->mNodes.size() - mDeletedNodeIndices.size();
00150 }
00151 
00158 template<>
00159 void MutableMesh<1, 1>::RescaleMeshFromBoundaryNode(ChastePoint<1> updatedPoint, unsigned boundaryNodeIndex)
00160 {
00161     assert(this->GetNode(boundaryNodeIndex)->IsBoundaryNode());
00162     double scaleFactor = updatedPoint[0] / this->GetNode(boundaryNodeIndex)->GetPoint()[0];
00163     double temp;
00164     for (unsigned i=0; i < boundaryNodeIndex+1; i++)
00165     {
00166         temp = scaleFactor * this->mNodes[i]->GetPoint()[0];
00167         ChastePoint<1> newPoint(temp);
00168         this->mNodes[i]->SetPoint(newPoint);
00169     }
00170     this->RefreshMesh();
00171 }
00172 
00173 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00174 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::SetNode(unsigned index,
00175         ChastePoint<SPACE_DIM> point,
00176         bool concreteMove)
00177 {
00178     this->mNodes[index]->SetPoint(point);
00179 
00180     if (concreteMove)
00181     {
00182         for (typename Node<SPACE_DIM>::ContainingBoundaryElementIterator it = this->mNodes[index]->ContainingBoundaryElementsBegin();
00183              it != this->mNodes[index]->ContainingBoundaryElementsEnd();
00184              ++it)
00185         {
00186             try
00187             {
00188                 this->GetBoundaryElement(*it)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[ (*it) ],
00189                                                                     this->mBoundaryElementJacobianDeterminants[ (*it) ]);
00190             }
00191             catch (Exception&)
00192             {
00193                 EXCEPTION("Moving node caused a boundary element to have a non-positive Jacobian determinant");
00194             }
00195         }
00196         for (typename Node<SPACE_DIM>::ContainingElementIterator it = this->mNodes[index]->ContainingElementsBegin();
00197              it != this->mNodes[index]->ContainingElementsEnd();
00198              ++it)
00199         {
00200             if (ELEMENT_DIM == SPACE_DIM)
00201             {
00202                 try
00203                 {
00204                     this->GetElement(*it)->CalculateInverseJacobian(this->mElementJacobians[ (*it) ],
00205                                                               this->mElementJacobianDeterminants[ (*it) ],
00206                                                               this->mElementInverseJacobians[ (*it) ]);
00207                 }
00208                 catch (Exception&)
00209                 {
00210                         EXCEPTION("Moving node caused an element to have a non-positive Jacobian determinant");
00211                 }
00212             }
00213             else
00214             {
00215                 c_vector<double,SPACE_DIM> previous_direction = this->mElementWeightedDirections[ (*it) ];
00216 
00217                 this->GetElement(*it)->CalculateWeightedDirection(this->mElementWeightedDirections[ (*it) ],
00218                                                             this->mElementJacobianDeterminants[ (*it) ]);
00219 
00220                 if ( inner_prod(previous_direction, this->mElementWeightedDirections[ (*it) ]) < 0)
00221                 {
00222                     EXCEPTION("Moving node caused an subspace element to change direction");
00223                 }
00224 
00225             }
00226         }
00227     }
00228 }
00229 
00230 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00231 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteNode(unsigned index)
00232 {
00233     if (this->mNodes[index]->IsDeleted())
00234     {
00235         EXCEPTION("Trying to delete a deleted node");
00236     }
00237     unsigned target_index = (unsigned)(-1);
00238     bool found_target=false;
00239     for (typename Node<SPACE_DIM>::ContainingElementIterator it = this->mNodes[index]->ContainingElementsBegin();
00240          !found_target && it != this->mNodes[index]->ContainingElementsEnd();
00241          ++it)
00242     {
00243         Element <ELEMENT_DIM,SPACE_DIM>* p_element = this->GetElement(*it);
00244         for (unsigned i=0; i<=ELEMENT_DIM && !found_target; i++)
00245         {
00246             target_index = p_element->GetNodeGlobalIndex(i);
00247             try
00248             {
00249                 MoveMergeNode(index, target_index, false);
00250                 found_target = true;
00251             }
00252             catch (Exception&)
00253             {
00254                 // Just try the next node
00255             }
00256         }
00257     }
00258     if (!found_target)
00259     {
00260         EXCEPTION("Failure to delete node");
00261     }
00262 
00263     MoveMergeNode(index, target_index);
00264 }
00265 
00266 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00267 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteElement(unsigned index)
00268 {
00269     assert(!this->mElements[index]->IsDeleted());
00270     this->mElements[index]->MarkAsDeleted();
00271     mDeletedElementIndices.push_back(index);
00272 
00273     //Delete any nodes that are no longer attached to mesh.
00274     for (unsigned node_index = 0; node_index < this->mElements[index]->GetNumNodes(); ++node_index)
00275     {
00276         if (this->mElements[index]->GetNode(node_index)->GetNumContainingElements() == 0u)
00277         {
00278             if (this->mElements[index]->GetNode(node_index)->GetNumBoundaryElements() == 0u)
00279             {
00280                 this->mElements[index]->GetNode(node_index)->MarkAsDeleted();
00281                 mDeletedNodeIndices.push_back(this->mElements[index]->GetNode(node_index)->GetIndex());
00282             }
00283             else if (this->mElements[index]->GetNode(node_index)->GetNumBoundaryElements() == 1u && ELEMENT_DIM == 1)
00284             {
00285                 std::set<unsigned> indices = this->mElements[index]->GetNode(node_index)->rGetContainingBoundaryElementIndices();
00286                 assert(indices.size() == 1u);
00287                 this->mBoundaryElements[*indices.begin()]->MarkAsDeleted();
00288                 mDeletedBoundaryElementIndices.push_back(*indices.begin());
00289 
00290                 this->mElements[index]->GetNode(node_index)->MarkAsDeleted();
00291                 mDeletedNodeIndices.push_back(this->mElements[index]->GetNode(node_index)->GetIndex());
00292             }
00293         }
00294     }
00295 }
00296 
00297 
00298 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00299 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteNodePriorToReMesh(unsigned index)
00300 {
00301     this->mNodes[index]->MarkAsDeleted();
00302     mDeletedNodeIndices.push_back(index);
00303 }
00304 
00305 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00306 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::MoveMergeNode(unsigned index,
00307         unsigned targetIndex,
00308         bool concreteMove)
00309 {
00310     if (this->mNodes[index]->IsDeleted())
00311     {
00312         EXCEPTION("Trying to move a deleted node");
00313     }
00314 
00315     if (index == targetIndex)
00316     {
00317         EXCEPTION("Trying to merge a node with itself");
00318     }
00319     if (this->mNodes[index]->IsBoundaryNode())
00320     {
00321         if (!this->mNodes[targetIndex]->IsBoundaryNode())
00322         {
00323             EXCEPTION("A boundary node can only be moved on to another boundary node");
00324         }
00325     }
00326     std::set<unsigned> unshared_element_indices;
00327     std::set_difference(this->mNodes[index]->rGetContainingElementIndices().begin(),
00328                         this->mNodes[index]->rGetContainingElementIndices().end(),
00329                         this->mNodes[targetIndex]->rGetContainingElementIndices().begin(),
00330                         this->mNodes[targetIndex]->rGetContainingElementIndices().end(),
00331                         std::inserter(unshared_element_indices, unshared_element_indices.begin()));
00332 
00333 
00334     if (unshared_element_indices.size() == this->mNodes[index]->rGetContainingElementIndices().size())
00335     {
00336         EXCEPTION("These nodes cannot be merged since they are not neighbours");
00337     }
00338 
00339     std::set<unsigned> unshared_boundary_element_indices;
00340     std::set_difference(this->mNodes[index]->rGetContainingBoundaryElementIndices().begin(),
00341                         this->mNodes[index]->rGetContainingBoundaryElementIndices().end(),
00342                         this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().begin(),
00343                         this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().end(),
00344                         std::inserter(unshared_boundary_element_indices, unshared_boundary_element_indices.begin()));
00345 
00346 
00347     if (this->mNodes[index]->IsBoundaryNode())
00348     {
00349         if (unshared_boundary_element_indices.size()
00350             == this->mNodes[index]->rGetContainingBoundaryElementIndices().size())
00351         {
00352             //May be redundant (only thrown in 1D tests)
00353             EXCEPTION("These nodes cannot be merged since they are not neighbours on the boundary");
00354         }
00355     }
00356 
00357     this->mNodes[index]->rGetModifiableLocation() = this->mNodes[targetIndex]->rGetLocation();
00358 
00359     for (std::set<unsigned>::const_iterator element_iter=unshared_element_indices.begin();
00360              element_iter != unshared_element_indices.end();
00361              element_iter++)
00362     {
00363         try
00364         {
00365             if (SPACE_DIM == ELEMENT_DIM)
00366             {
00367                 this->GetElement(*element_iter)->CalculateInverseJacobian(this->mElementJacobians[(*element_iter)],
00368                                                                           this->mElementJacobianDeterminants[(*element_iter)],
00369                                                                           this->mElementInverseJacobians[ (*element_iter) ]);
00370             }
00371             else
00372             {
00373                 this->GetElement(*element_iter)->CalculateWeightedDirection(this->mElementWeightedDirections[(*element_iter)],
00374                                                                             this->mElementJacobianDeterminants[(*element_iter)]);
00375             }
00376 
00377             if (concreteMove)
00378             {
00379                 this->GetElement(*element_iter)->ReplaceNode(this->mNodes[index], this->mNodes[targetIndex]);
00380             }
00381 
00382         }
00383         catch (Exception&)
00384         {
00385             EXCEPTION("Moving node caused an element to have a non-positive Jacobian determinant");
00386         }
00387     }
00388 
00389     for (std::set<unsigned>::const_iterator boundary_element_iter=
00390                  unshared_boundary_element_indices.begin();
00391              boundary_element_iter != unshared_boundary_element_indices.end();
00392              boundary_element_iter++)
00393     {
00394 
00395         this->GetBoundaryElement(*boundary_element_iter)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[(*boundary_element_iter)],
00396                                                                                      this->mBoundaryElementJacobianDeterminants[(*boundary_element_iter)]);
00397 
00398         if (concreteMove)
00399         {
00400             this->GetBoundaryElement(*boundary_element_iter)->ReplaceNode(this->mNodes[index], this->mNodes[targetIndex]);
00401         }
00402     }
00403 
00404     std::set<unsigned> shared_element_indices;
00405     std::set_intersection(this->mNodes[index]->rGetContainingElementIndices().begin(),
00406                           this->mNodes[index]->rGetContainingElementIndices().end(),
00407                           this->mNodes[targetIndex]->rGetContainingElementIndices().begin(),
00408                           this->mNodes[targetIndex]->rGetContainingElementIndices().end(),
00409                           std::inserter(shared_element_indices, shared_element_indices.begin()));
00410 
00411     for (std::set<unsigned>::const_iterator element_iter=shared_element_indices.begin();
00412              element_iter != shared_element_indices.end();
00413              element_iter++)
00414     {
00415         if (concreteMove)
00416         {
00417             this->GetElement(*element_iter)->MarkAsDeleted();
00418             this->mElementJacobianDeterminants[ (*element_iter) ] = 0.0; //This used to be done in MarkAsDeleted
00419             mDeletedElementIndices.push_back(*element_iter);
00420         }
00421         else
00422         {
00423             this->mElementJacobianDeterminants[ (*element_iter) ] = 0.0;
00424         }
00425     }
00426 
00427 
00428     std::set<unsigned> shared_boundary_element_indices;
00429     std::set_intersection(this->mNodes[index]->rGetContainingBoundaryElementIndices().begin(),
00430                           this->mNodes[index]->rGetContainingBoundaryElementIndices().end(),
00431                           this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().begin(),
00432                           this->mNodes[targetIndex]->rGetContainingBoundaryElementIndices().end(),
00433                           std::inserter(shared_boundary_element_indices, shared_boundary_element_indices.begin()));
00434 
00435     for (std::set<unsigned>::const_iterator boundary_element_iter=shared_boundary_element_indices.begin();
00436              boundary_element_iter != shared_boundary_element_indices.end();
00437              boundary_element_iter++)
00438     {
00439         if (concreteMove)
00440         {
00441             this->GetBoundaryElement(*boundary_element_iter)->MarkAsDeleted();
00442             this->mBoundaryElementJacobianDeterminants[ (*boundary_element_iter) ] = 0.0; //This used to be done in MarkAsDeleted
00443             mDeletedBoundaryElementIndices.push_back(*boundary_element_iter);
00444         }
00445         else
00446         {
00447             this->mBoundaryElementJacobianDeterminants[ (*boundary_element_iter) ] = 0.0;
00448             this->mBoundaryElementWeightedDirections[ (*boundary_element_iter) ] = zero_vector<double>(SPACE_DIM);
00449         }
00450     }
00451 
00452     if (concreteMove)
00453     {
00454         this->mNodes[index]->MarkAsDeleted();
00455         mDeletedNodeIndices.push_back(index);
00456     }
00457 }
00458 
00459 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00460 unsigned MutableMesh<ELEMENT_DIM, SPACE_DIM>::RefineElement(
00461     Element<ELEMENT_DIM,SPACE_DIM>* pElement,
00462     ChastePoint<SPACE_DIM> point)
00463 {
00464     //Check that the point is in the element
00465     if (pElement->IncludesPoint(point, true) == false)
00466     {
00467         EXCEPTION("RefineElement could not be started (point is not in element)");
00468     }
00469 
00470     // Add a new node from the point that is passed to RefineElement
00471     unsigned new_node_index = AddNode(new Node<SPACE_DIM>(0, point.rGetLocation()));
00472     // Note: the first argument is the index of the node, which is going to be
00473     //       overridden by AddNode, so it can safely be ignored
00474 
00475     // This loop constructs the extra elements which are going to fill the space
00476     for (unsigned i = 0; i < ELEMENT_DIM; i++)
00477     {
00478 
00479         // First, make a copy of the current element making sure we update its index
00480         unsigned new_elt_index;
00481         if (mDeletedElementIndices.empty())
00482         {
00483             new_elt_index = this->mElements.size();
00484         }
00485         else
00486         {
00487             new_elt_index = mDeletedElementIndices.back();
00488             mDeletedElementIndices.pop_back();
00489         }
00490 
00491         Element<ELEMENT_DIM,SPACE_DIM>* p_new_element=
00492             new Element<ELEMENT_DIM,SPACE_DIM>(*pElement, new_elt_index);
00493 
00494         // Second, update the node in the element with the new one
00495         p_new_element->UpdateNode(ELEMENT_DIM-1-i, this->mNodes[new_node_index]);
00496 
00497         // Third, add the new element to the set
00498         if ((unsigned) new_elt_index == this->mElements.size())
00499         {
00500             this->mElements.push_back(p_new_element);
00501         }
00502         else
00503         {
00504             delete this->mElements[new_elt_index];
00505             this->mElements[new_elt_index] = p_new_element;
00506         }
00507     }
00508 
00509     // Lastly, update the last node in the element to be refined
00510     pElement->UpdateNode(ELEMENT_DIM, this->mNodes[new_node_index]);
00511 
00512     return new_node_index;
00513 }
00514 
00515 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00516 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::DeleteBoundaryNodeAt(unsigned index)
00517 {
00518     if (!this->mNodes[index]->IsBoundaryNode() )
00519     {
00520         EXCEPTION(" You may only delete a boundary node ");
00521     }
00522 
00523     this->mNodes[index]->MarkAsDeleted();
00524     mDeletedNodeIndices.push_back(index);
00525 
00526     // Update the boundary node vector
00527     typename std::vector<Node<SPACE_DIM>*>::iterator b_node_iter
00528     = std::find(this->mBoundaryNodes.begin(), this->mBoundaryNodes.end(), this->mNodes[index]);
00529     this->mBoundaryNodes.erase(b_node_iter);
00530 
00531     // Remove boundary elements containing this node
00532     std::set<unsigned> boundary_element_indices = this->mNodes[index]->rGetContainingBoundaryElementIndices();
00533     std::set<unsigned>::const_iterator boundary_element_indices_iterator = boundary_element_indices.begin();
00534     while (boundary_element_indices_iterator != boundary_element_indices.end())
00535     {
00536         BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>* p_boundary_element = this->GetBoundaryElement(*boundary_element_indices_iterator);
00537         p_boundary_element->MarkAsDeleted();
00538         mDeletedBoundaryElementIndices.push_back(*boundary_element_indices_iterator);
00539         boundary_element_indices_iterator++;
00540     }
00541 
00542     // Remove elements containing this node
00543     std::set<unsigned> element_indices = this->mNodes[index]->rGetContainingElementIndices();
00544     std::set<unsigned>::const_iterator element_indices_iterator = element_indices.begin();
00545     while (element_indices_iterator != element_indices.end())
00546     {
00547         Element<ELEMENT_DIM, SPACE_DIM>* p_element = this->GetElement(*element_indices_iterator);
00548         for (unsigned i=0; i<p_element->GetNumNodes(); i++)
00549         {
00550             Node<SPACE_DIM>* p_node = p_element->GetNode(i);
00551             if (!p_node->IsDeleted())
00552             {
00553                 p_node->SetAsBoundaryNode();
00554                 // Update the boundary node vector
00555                 this->mBoundaryNodes.push_back(p_node);
00556             }
00557         }
00558         p_element->MarkAsDeleted();
00559         mDeletedElementIndices.push_back(p_element->GetIndex());
00560         element_indices_iterator++;
00561     }
00562 }
00563 
00564 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00565 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::ReIndex(NodeMap& map)
00566 {
00567     assert(!mAddedNodes);
00568     map.Resize(this->GetNumAllNodes());
00569 
00570     std::vector<Element<ELEMENT_DIM, SPACE_DIM> *> live_elements;
00571 
00572     for (unsigned i=0; i<this->mElements.size(); i++)
00573     {
00574         assert(i==this->mElements[i]->GetIndex()); // We need this to be true to be able to reindex the Jacobian cache
00575         if (this->mElements[i]->IsDeleted())
00576         {
00577             delete this->mElements[i];
00578         }
00579         else
00580         {
00581             live_elements.push_back(this->mElements[i]);
00582 
00583             unsigned this_element_index = this->mElements[i]->GetIndex();
00584             if (SPACE_DIM == ELEMENT_DIM)
00585             {
00586                 this->mElementJacobians[live_elements.size()-1] = this->mElementJacobians[this_element_index];
00587                 this->mElementInverseJacobians[live_elements.size()-1] = this->mElementInverseJacobians[this_element_index];
00588             }
00589             else
00590             {
00591                 this->mElementWeightedDirections[live_elements.size()-1] = this->mElementWeightedDirections[this_element_index];
00592             }
00593             this->mElementJacobianDeterminants[live_elements.size()-1] = this->mElementJacobianDeterminants[this_element_index];
00594         }
00595     }
00596 
00597     assert(mDeletedElementIndices.size() == this->mElements.size()-live_elements.size());
00598     mDeletedElementIndices.clear();
00599     this->mElements = live_elements;
00600     unsigned num_elements = this->mElements.size();
00601 
00602     if (SPACE_DIM == ELEMENT_DIM)
00603     {
00604         this->mElementJacobians.resize(num_elements);
00605         this->mElementInverseJacobians.resize(num_elements);
00606     }
00607     else
00608     {
00609         this->mElementWeightedDirections.resize(num_elements);
00610     }
00611     this->mElementJacobianDeterminants.resize(num_elements);
00612 
00613     std::vector<Node<SPACE_DIM> *> live_nodes;
00614     for (unsigned i=0; i<this->mNodes.size(); i++)
00615     {
00616         if (this->mNodes[i]->IsDeleted())
00617         {
00618             delete this->mNodes[i];
00619             map.SetDeleted(i);
00620         }
00621         else
00622         {
00623             live_nodes.push_back(this->mNodes[i]);
00624             // the nodes will have their index set to be the index into the live_nodes
00625             // vector further down
00626             map.SetNewIndex(i, (unsigned)(live_nodes.size()-1));
00627         }
00628     }
00629 
00630     assert(mDeletedNodeIndices.size() == this->mNodes.size()-live_nodes.size());
00631     this->mNodes = live_nodes;
00632     mDeletedNodeIndices.clear();
00633 
00634     std::vector<BoundaryElement<ELEMENT_DIM-1, SPACE_DIM> *> live_boundary_elements;
00635     for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00636     {
00637         if (this->mBoundaryElements[i]->IsDeleted())
00638         {
00639             delete this->mBoundaryElements[i];
00640         }
00641         else
00642         {
00643             live_boundary_elements.push_back(this->mBoundaryElements[i]);
00644 
00645             this->mBoundaryElementWeightedDirections[live_boundary_elements.size()-1] = this->mBoundaryElementWeightedDirections[this->mBoundaryElements[i]->GetIndex()];
00646             this->mBoundaryElementJacobianDeterminants[live_boundary_elements.size()-1] = this->mBoundaryElementJacobianDeterminants[this->mBoundaryElements[i]->GetIndex()];
00647         }
00648     }
00649 
00650     assert(mDeletedBoundaryElementIndices.size() == this->mBoundaryElements.size()-live_boundary_elements.size());
00651     this->mBoundaryElements = live_boundary_elements;
00652     mDeletedBoundaryElementIndices.clear();
00653 
00654     unsigned num_boundary_elements = this->mBoundaryElements.size();
00655 
00656     this->mBoundaryElementWeightedDirections.resize(num_boundary_elements);
00657     this->mBoundaryElementJacobianDeterminants.resize(num_boundary_elements);
00658 
00659     for (unsigned i=0; i<this->mNodes.size(); i++)
00660     {
00661         this->mNodes[i]->SetIndex(i);
00662     }
00663 
00664     for (unsigned i=0; i<this->mElements.size(); i++)
00665     {
00666 
00667         this->mElements[i]->ResetIndex(i);
00668     }
00669 
00670     for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00671     {
00672         this->mBoundaryElements[i]->ResetIndex(i);
00673     }
00674 }
00675 
00676 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00677 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::ReMesh(NodeMap& map)
00678 {
00679     // Make sure that we are in the correct dimension - this code will be eliminated at compile time
00680     #define COVERAGE_IGNORE
00681     assert( ELEMENT_DIM == SPACE_DIM );
00682     #undef COVERAGE_IGNORE
00683 
00684     // Avoid some triangle/tetgen errors: need at least four
00685     // nodes for tetgen, and at least three for triangle
00686     if (GetNumNodes() <= SPACE_DIM)
00687     {
00688         EXCEPTION("The number of nodes must exceed the spatial dimension.");
00689     }
00690 
00691     // Make sure the map is big enough
00692     map.Resize(this->GetNumAllNodes());
00693     if (mAddedNodes || !mDeletedNodeIndices.empty())
00694     {
00695         // Size of mesh is about to change
00696         if (this->mpDistributedVectorFactory)
00697         {
00698             delete this->mpDistributedVectorFactory;
00699             this->mpDistributedVectorFactory = new DistributedVectorFactory(this->GetNumNodes());
00700         }
00701     }
00702     if (SPACE_DIM==1)
00703     {
00704         // Store the node locations
00705         std::vector<c_vector<double, SPACE_DIM> > old_node_locations;
00706         unsigned new_index = 0;
00707         for (unsigned i=0; i<this->GetNumAllNodes(); i++)
00708         {
00709             if (this->mNodes[i]->IsDeleted())
00710             {
00711                 map.SetDeleted(i);
00712             }
00713             else
00714             {
00715                 map.SetNewIndex(i, new_index);
00716                 old_node_locations.push_back(this->mNodes[i]->rGetLocation());
00717                 new_index++;
00718             }
00719         }
00720 
00721         // Remove current data
00722         Clear();
00723 
00724         // Construct the nodes and boundary nodes
00725         for (unsigned node_index=0; node_index<old_node_locations.size(); node_index++)
00726         {
00727             Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(node_index, old_node_locations[node_index], false);
00728             this->mNodes.push_back(p_node);
00729 
00730             // As we're in 1D, the boundary nodes are simply at either end of the mesh
00731             if ( node_index==0 || node_index==old_node_locations.size()-1 )
00732             {
00733                 this->mBoundaryNodes.push_back(p_node);
00734             }
00735         }
00736 
00737         // Create a map between node indices and node locations
00738         std::map<double, unsigned> location_index_map;
00739         for (unsigned i=0; i<this->mNodes.size(); i++)
00740         {
00741             location_index_map[this->mNodes[i]->rGetLocation()[0]] = this->mNodes[i]->GetIndex();
00742         }
00743 
00744         // Use this map to generate a vector of node indices that are ordered spatially
00745         std::vector<unsigned> node_indices_ordered_spatially;
00746         for (std::map<double, unsigned>::iterator iter = location_index_map.begin();
00747              iter != location_index_map.end();
00748              ++iter)
00749         {
00750             node_indices_ordered_spatially.push_back(iter->second);
00751         }
00752 
00753         // Construct the elements
00754         this->mElements.reserve(old_node_locations.size()-1);
00755         for (unsigned element_index=0; element_index<old_node_locations.size()-1; element_index++)
00756         {
00757             std::vector<Node<SPACE_DIM>*> nodes;
00758             for (unsigned j=0; j<2; j++)
00759             {
00760                 unsigned global_node_index = node_indices_ordered_spatially[element_index + j];
00761                 assert(global_node_index < this->mNodes.size());
00762                 nodes.push_back(this->mNodes[global_node_index]);
00763             }
00764             this->mElements.push_back(new Element<ELEMENT_DIM, SPACE_DIM>(element_index, nodes));
00765         }
00766 
00767         // Construct the two boundary elements - as we're in 1D, these are simply at either end of the mesh
00768         std::vector<Node<SPACE_DIM>*> nodes;
00769         nodes.push_back(this->mNodes[0]);
00770         this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(0, nodes));
00771 
00772         nodes.clear();
00773         nodes.push_back(this->mNodes[old_node_locations.size()-1]);
00774         this->mBoundaryElements.push_back(new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(1, nodes));
00775 
00776         this->RefreshJacobianCachedData();
00777     }
00778     else if (SPACE_DIM==2)  // In 2D, remesh using triangle via library calls
00779     {
00780         struct triangulateio mesher_input, mesher_output;
00781         this->InitialiseTriangulateIo(mesher_input);
00782         this->InitialiseTriangulateIo(mesher_output);
00783 
00784         this->ExportToMesher(map, mesher_input);
00785 
00786         // Library call
00787         triangulate((char*)"Qze", &mesher_input, &mesher_output, NULL);
00788 
00789         this->ImportFromMesher(mesher_output, mesher_output.numberoftriangles, mesher_output.trianglelist, mesher_output.numberofedges, mesher_output.edgelist, mesher_output.edgemarkerlist);
00790 
00791         //Tidy up triangle
00792         this->FreeTriangulateIo(mesher_input);
00793         this->FreeTriangulateIo(mesher_output);
00794     }
00795     else // in 3D, remesh using tetgen
00796     {
00797 
00798         class tetgen::tetgenio mesher_input, mesher_output;
00799 
00800         this->ExportToMesher(map, mesher_input);
00801 
00802         // Library call
00803         tetgen::tetrahedralize((char*)"Qz", &mesher_input, &mesher_output);
00804 
00805         this->ImportFromMesher(mesher_output, mesher_output.numberoftetrahedra, mesher_output.tetrahedronlist, mesher_output.numberoftrifaces, mesher_output.trifacelist, NULL);
00806     }
00807 }
00808 
00809 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00810 void MutableMesh<ELEMENT_DIM, SPACE_DIM>::ReMesh()
00811 {
00812     NodeMap map(GetNumNodes());
00813     ReMesh(map);
00814 }
00815 
00816 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00817 std::vector<c_vector<unsigned, 5> > MutableMesh<ELEMENT_DIM, SPACE_DIM>::SplitLongEdges(double cutoffLength)
00818 {
00819     assert(ELEMENT_DIM == 2);
00820     assert(SPACE_DIM == 3);
00821 
00822     std::vector<c_vector<unsigned, 5> > history;
00823 
00824 
00825     bool long_edge_exists = true;
00826 
00827     while(long_edge_exists)
00828     {
00829         std::set<std::pair<unsigned, unsigned> > long_edges;
00830 
00831         // Loop over elements to check for Long edges
00832         for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator elem_iter = this->GetElementIteratorBegin();
00833              elem_iter != this->GetElementIteratorEnd();
00834              ++elem_iter)
00835         {
00836             unsigned num_nodes = ELEMENT_DIM+1;
00837 
00838             // Loop over element vertices
00839             for (unsigned local_index=0; local_index<num_nodes; local_index++)
00840             {
00841                 // Find locations of current node (node a) and anticlockwise node (node b)
00842                 Node<SPACE_DIM>* p_node_a = elem_iter->GetNode(local_index);
00843                 unsigned local_index_plus_one = (local_index+1)%num_nodes; 
00844                 Node<SPACE_DIM>* p_node_b = elem_iter->GetNode(local_index_plus_one);
00845 
00846                 // Find distance between nodes
00847                 double distance_between_nodes = this->GetDistanceBetweenNodes(p_node_a->GetIndex(), p_node_b->GetIndex());
00848 
00849                 if (distance_between_nodes > cutoffLength)
00850                 {
00851                     if (p_node_a->GetIndex() < p_node_b->GetIndex())
00852                     {
00853                         std::pair<unsigned, unsigned> long_edge(p_node_a->GetIndex(),p_node_b->GetIndex());
00854                         long_edges.insert(long_edge);
00855                     }
00856                     else
00857                     {
00858                         std::pair<unsigned, unsigned> long_edge(p_node_b->GetIndex(),p_node_a->GetIndex());
00859                         long_edges.insert(long_edge);
00860                     }
00861                 }
00862             }
00863         }
00864 
00865         if (long_edges.size() > 0) //Split the edges in decreasing order.
00866         {
00867             while (long_edges.size() > 0)
00868             {
00869                 double longest_edge = 0.0;
00870                 std::set<std::pair<unsigned, unsigned> >::iterator longest_edge_iter;
00871 
00872                 //Find the longest edge in the set and split it
00873                 for (std::set<std::pair<unsigned, unsigned> >::iterator edge_iter = long_edges.begin();
00874                          edge_iter != long_edges.end();
00875                          ++edge_iter)
00876                 {
00877                     unsigned node_a_global_index = edge_iter->first;
00878                     unsigned node_b_global_index = edge_iter->second;
00879 
00880                     double distance_between_nodes = this->GetDistanceBetweenNodes(node_a_global_index, node_b_global_index);
00881 
00882                     if (distance_between_nodes > longest_edge)
00883                     {
00884                         longest_edge = distance_between_nodes;
00885                         longest_edge_iter = edge_iter;
00886                     }
00887                 }
00888                 assert(longest_edge >0);
00889 
00890                 c_vector<unsigned, 3> new_node_index = SplitEdge(this->GetNode(longest_edge_iter->first), this->GetNode(longest_edge_iter->second));
00891 
00892                 c_vector<unsigned, 5> node_set;
00893                 node_set(0) = new_node_index[0];
00894                 node_set(1) = longest_edge_iter->first;
00895                 node_set(2) = longest_edge_iter->second;
00896                 node_set(3) = new_node_index[1];
00897                 node_set(4) = new_node_index[2];
00898                 history.push_back(node_set);
00899 
00900                 // Delete pair from set
00901                 long_edges.erase(*longest_edge_iter);
00902             }
00903         }
00904         else
00905         {
00906             long_edge_exists = false;
00907         }
00908     }
00909 
00910     return history;
00911 }
00912 
00913 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00914 c_vector<unsigned, 3> MutableMesh<ELEMENT_DIM, SPACE_DIM>::SplitEdge(Node<SPACE_DIM>* pNodeA, Node<SPACE_DIM>* pNodeB)
00915 {
00916     c_vector<unsigned, 3> new_node_index_vector;
00917 
00918     std::set<unsigned> elements_of_node_a = pNodeA->rGetContainingElementIndices();
00919     std::set<unsigned> elements_of_node_b = pNodeB->rGetContainingElementIndices();
00920 
00921     std::set<unsigned> intersection_elements;
00922     std::set_intersection(elements_of_node_a.begin(), elements_of_node_a.end(),
00923                           elements_of_node_b.begin(), elements_of_node_b.end(),
00924                           std::inserter(intersection_elements, intersection_elements.begin()));
00925 
00926     // Create the new node
00927     c_vector<double, SPACE_DIM> new_node_location = pNodeA->rGetLocation() + 0.5*this->GetVectorFromAtoB(pNodeA->rGetLocation(), pNodeB->rGetLocation());
00928 
00929     bool is_boundary_node =  intersection_elements.size() == 1; // If only in one element then its on the boundary
00930 
00931     Node<SPACE_DIM>* p_new_node = new Node<SPACE_DIM>(0u,new_node_location,is_boundary_node); // Index is rewritten once added to mesh
00932 
00933     unsigned new_node_index = this->AddNode(p_new_node);
00934 
00935     new_node_index_vector[0] = new_node_index;
00936 
00937     unsigned counter = 1;
00938 
00939     for (std::set<unsigned>::const_iterator it = intersection_elements.begin(); it != intersection_elements.end(); ++it)
00940     {
00941         unsigned elementIndex = *it;
00942 
00943         Element<ELEMENT_DIM,SPACE_DIM>* p_original_element = this->GetElement(elementIndex);
00944 
00945         // First, make a copy of the current element and assign an unused index
00946         Element<ELEMENT_DIM,SPACE_DIM>* p_new_element = new Element<ELEMENT_DIM,SPACE_DIM>(*p_original_element, UINT_MAX);
00947 
00948         // Second, add the new element to the set of existing elements. This method will assign a proper index to the element.
00949         AddElement(p_new_element);
00950 
00951         // Third, update node a in the element with the new one
00952         p_new_element->ReplaceNode(pNodeA, this->mNodes[new_node_index]);
00953 
00954         // Last, update node b in the original element with the new one
00955         p_original_element->ReplaceNode(pNodeB, this->mNodes[new_node_index]);
00956 
00957         //Add node in both of these elements to new_node_index_vector (this enables us to add a new spring in the MeshBasedCellPopulation
00958         unsigned other_node_index = UNSIGNED_UNSET;
00959 
00960         if ( (p_original_element->GetNodeGlobalIndex(0) != new_node_index) &&
00961              (p_original_element->GetNodeGlobalIndex(0) != pNodeA->GetIndex() ) )
00962         {
00963             other_node_index = p_original_element->GetNodeGlobalIndex(0);
00964         }
00965         else if ( (p_original_element->GetNodeGlobalIndex(1) != new_node_index) &&
00966                   (p_original_element->GetNodeGlobalIndex(1) != pNodeA->GetIndex() ) )
00967         {
00968             other_node_index = p_original_element->GetNodeGlobalIndex(1);
00969         }
00970         else if ( (p_original_element->GetNodeGlobalIndex(2) != new_node_index) &&
00971                   (p_original_element->GetNodeGlobalIndex(2) != pNodeA->GetIndex() ) )
00972         {
00973             other_node_index = p_original_element->GetNodeGlobalIndex(2);
00974         }
00975         else
00976         {
00977             NEVER_REACHED;
00978         }
00979         new_node_index_vector[counter] = other_node_index;
00980         counter++;
00981     }
00982 
00983     assert(counter<4);
00984     assert(counter>1);// need to be in at least one element
00985 
00986     if (counter == 2) // only one new element
00987     {
00988         new_node_index_vector[2] = UNSIGNED_UNSET;
00989     }
00990 
00991     return new_node_index_vector;
00992 }
00993 
00994 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
00995 bool MutableMesh<ELEMENT_DIM, SPACE_DIM>::CheckIsVoronoi(Element<ELEMENT_DIM, SPACE_DIM>* pElement, double maxPenetration)
00996 {
00997     assert(ELEMENT_DIM == SPACE_DIM);
00998     unsigned num_nodes = pElement->GetNumNodes();
00999     std::set<unsigned> neighbouring_elements_indices;
01000     std::set< Element<ELEMENT_DIM,SPACE_DIM> *> neighbouring_elements;
01001     std::set<unsigned> neighbouring_nodes_indices;
01002 
01003     // Form a set of neighbouring elements via the nodes
01004     for (unsigned i=0; i<num_nodes; i++)
01005     {
01006         Node<SPACE_DIM>* p_node = pElement->GetNode(i);
01007         neighbouring_elements_indices = p_node->rGetContainingElementIndices();
01008         for (std::set<unsigned>::const_iterator it = neighbouring_elements_indices.begin();
01009              it != neighbouring_elements_indices.end();
01010              ++it)
01011         {
01012             neighbouring_elements.insert(this->GetElement(*it));
01013         }
01014     }
01015     neighbouring_elements.erase(pElement);
01016 
01017     // For each neighbouring element find the supporting nodes
01018     typedef typename std::set<Element<ELEMENT_DIM,SPACE_DIM> *>::const_iterator ElementIterator;
01019 
01020     for (ElementIterator it = neighbouring_elements.begin();
01021          it != neighbouring_elements.end();
01022          ++it)
01023     {
01024         for (unsigned i=0; i<num_nodes; i++)
01025         {
01026             neighbouring_nodes_indices.insert((*it)->GetNodeGlobalIndex(i));
01027         }
01028     }
01029 
01030     // Remove the nodes that support this element
01031     for (unsigned i = 0; i < num_nodes; i++)
01032     {
01033         neighbouring_nodes_indices.erase(pElement->GetNodeGlobalIndex(i));
01034     }
01035 
01036     // Get the circumsphere information
01037     c_vector<double, SPACE_DIM+1> this_circum_centre;
01038 
01039     this_circum_centre = pElement->CalculateCircumsphere(this->mElementJacobians[pElement->GetIndex()], this->mElementInverseJacobians[pElement->GetIndex()]);
01040 
01041     // Copy the actualy circumcentre into a smaller vector
01042     c_vector<double, ELEMENT_DIM> circum_centre;
01043     for (unsigned i=0; i<ELEMENT_DIM; i++)
01044     {
01045         circum_centre[i] = this_circum_centre[i];
01046     }
01047 
01048     for (std::set<unsigned>::const_iterator it = neighbouring_nodes_indices.begin();
01049          it != neighbouring_nodes_indices.end();
01050          ++it)
01051     {
01052         c_vector<double, ELEMENT_DIM> node_location = this->GetNode(*it)->rGetLocation();
01053 
01054         // Calculate vector from circumcenter to node
01055         node_location -= circum_centre;
01056 
01057         // This is to calculate the squared distance betweeen them
01058         double squared_distance = inner_prod(node_location, node_location);
01059 
01060         // If the squared idstance is less than the elements circum-radius(squared),
01061         // then the Voronoi property is violated.
01062         if (squared_distance < this_circum_centre[ELEMENT_DIM])
01063         {
01064             // We know the node is inside the circumsphere, but we don't know how far
01065             double radius = sqrt(this_circum_centre[ELEMENT_DIM]);
01066             double distance = radius - sqrt(squared_distance);
01067 
01068             // If the node penetration is greater than supplied maximum penetration factor
01069             if (distance/radius > maxPenetration)
01070             {
01071                 return false;
01072             }
01073         }
01074     }
01075     return true;
01076 }
01077 
01078 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
01079 bool MutableMesh<ELEMENT_DIM, SPACE_DIM>::CheckIsVoronoi(double maxPenetration)
01080 {
01081     // Looping through all the elements in the mesh
01083     for (unsigned i=0; i<this->mElements.size(); i++)
01084     {
01085         // Check if the element is not deleted
01086         if (!this->mElements[i]->IsDeleted())
01087         {
01088             // Checking the Voronoi of the Element
01089             if (CheckIsVoronoi(this->mElements[i], maxPenetration) == false)
01090             {
01091                 return false;
01092             }
01093         }
01094     }
01095     return true;
01096 }
01097 
01098 
01100 // Explicit instantiation
01102 
01103 template class MutableMesh<1,1>;
01104 template class MutableMesh<1,2>;
01105 template class MutableMesh<1,3>;
01106 template class MutableMesh<2,2>;
01107 template class MutableMesh<2,3>;
01108 template class MutableMesh<3,3>;
01109 
01110 
01111 // Serialization for Boost >= 1.36
01112 #include "SerializationExportWrapperForCpp.hpp"
01113 EXPORT_TEMPLATE_CLASS_ALL_DIMS(MutableMesh)

Generated by  doxygen 1.6.2