FineCoarseMeshPair.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include "FineCoarseMeshPair.hpp"
00030 
00031 template<unsigned DIM>
00032 FineCoarseMeshPair<DIM>::FineCoarseMeshPair(TetrahedralMesh<DIM,DIM>& rFineMesh, TetrahedralMesh<DIM,DIM>& rCoarseMesh)
00033     : mrFineMesh(rFineMesh),
00034       mrCoarseMesh(rCoarseMesh)
00035 {
00036     mpFineMeshBoxCollection = NULL;
00037     mpCoarseMeshBoxCollection = NULL;
00038     ResetStatisticsVariables();
00039 }
00040 
00041 template<unsigned DIM>
00042 FineCoarseMeshPair<DIM>::~FineCoarseMeshPair()
00043 {
00044     DeleteFineBoxCollection();
00045     DeleteCoarseBoxCollection();
00046 }
00047 
00048 template<unsigned DIM>
00049 void FineCoarseMeshPair<DIM>::DeleteFineBoxCollection()
00050 {
00051     if (mpFineMeshBoxCollection != NULL)
00052     {
00053         delete mpFineMeshBoxCollection;
00054         mpFineMeshBoxCollection = NULL;
00055     }
00056 }
00057 
00058 template<unsigned DIM>
00059 void FineCoarseMeshPair<DIM>::DeleteCoarseBoxCollection()
00060 {
00061     if (mpCoarseMeshBoxCollection != NULL)
00062     {
00063         delete mpCoarseMeshBoxCollection;
00064         mpCoarseMeshBoxCollection = NULL;
00065     }
00066 }
00067 
00069 //   Setting up boxes methods
00071 
00072 template<unsigned DIM>
00073 void FineCoarseMeshPair<DIM>::SetUpBoxesOnFineMesh(double boxWidth)
00074 {
00075     SetUpBoxes(mrFineMesh, boxWidth, mpFineMeshBoxCollection);
00076 }
00077 
00078 template<unsigned DIM>
00079 void FineCoarseMeshPair<DIM>::SetUpBoxesOnCoarseMesh(double boxWidth)
00080 {
00081     SetUpBoxes(mrCoarseMesh, boxWidth, mpCoarseMeshBoxCollection);
00082 }
00083 
00084 template<unsigned DIM>
00085 void FineCoarseMeshPair<DIM>::SetUpBoxes(TetrahedralMesh<DIM,DIM>& rMesh,
00086                                          double boxWidth,
00087                                          BoxCollection<DIM>*& rpBoxCollection)
00088 {
00089     if (rpBoxCollection)
00090     {
00091         delete rpBoxCollection;
00092         rpBoxCollection = NULL;
00093     }
00094 
00095     // Compute min and max values for the fine mesh nodes
00096     ChasteCuboid<DIM> bounding_box = rMesh.CalculateBoundingBox();
00097 
00098     // Set up the boxes using a domain that is slightly larger than the fine mesh
00099     c_vector<double,2*DIM> extended_min_and_max;
00100     for (unsigned i=0; i<DIM; i++)
00101     {
00102         double width = bounding_box.GetWidth(i);
00103 
00104         // Subtract from the minima
00105         extended_min_and_max(2*i) = bounding_box.rGetLowerCorner()[i] - 0.05*width;
00106 
00107         // Add to the maxima
00108         extended_min_and_max(2*i+1) = bounding_box.rGetUpperCorner()[i] + 0.05*width;
00109     }
00110 
00111     if (boxWidth < 0)
00112     {
00113         /*
00114          * Use default value = max(max_edge_length, w20),  where w20 is the width
00115          * corresponding to 20 boxes in the x-direction.
00116          *
00117          * BoxCollection creates an extra box so divide by 19 not 20. Add a little
00118          * bit on to ensure minor numerical fluctuations don't change the answer.
00119          */
00120         boxWidth = (extended_min_and_max(1) - extended_min_and_max(0))/19.000000001;
00121 
00122         // Determine the maximum edge length
00123         double max_edge_length = -1;
00124 
00125         for (typename TetrahedralMesh<DIM,DIM>::EdgeIterator edge_iterator = rMesh.EdgesBegin();
00126              edge_iterator!=rMesh.EdgesEnd();
00127              ++edge_iterator)
00128         {
00129             c_vector<double, 3> location1 = edge_iterator.GetNodeA()->rGetLocation();
00130             c_vector<double, 3> location2 = edge_iterator.GetNodeB()->rGetLocation();
00131             double edge_length = norm_2(location1-location2);
00132 
00133             if (edge_length>max_edge_length)
00134             {
00135                 max_edge_length = edge_length;
00136             }
00137         }
00138 
00139         if (boxWidth < max_edge_length)
00140         {
00141             boxWidth = 1.1*max_edge_length;
00142         }
00143     }
00144 
00145     rpBoxCollection = new BoxCollection<DIM>(boxWidth, extended_min_and_max);
00146     rpBoxCollection->SetupAllLocalBoxes();
00147 
00148     // For each element, if ANY of its nodes are physically in a box, put that element in that box
00149     for (unsigned i=0; i<rMesh.GetNumElements(); i++)
00150     {
00151         Element<DIM,DIM>* p_element = rMesh.GetElement(i);
00152 
00153         std::set<unsigned> box_indices_each_node_this_elem;
00154         for (unsigned j=0; j<DIM+1; j++) // num vertices per element
00155         {
00156             Node<DIM>* p_node = p_element->GetNode(j);
00157             unsigned box_index = rpBoxCollection->CalculateContainingBox(p_node);
00158             box_indices_each_node_this_elem.insert(box_index);
00159         }
00160 
00161         for (std::set<unsigned>::iterator iter = box_indices_each_node_this_elem.begin();
00162             iter != box_indices_each_node_this_elem.end();
00163             ++iter)
00164         {
00165             rpBoxCollection->rGetBox( *iter ).AddElement(p_element);
00166         }
00167     }
00168 }
00169 
00171 // ComputeFineElementsAndWeightsForCoarseQuadPoints()
00172 // and
00173 // ComputeFineElementsAndWeightsForCoarseNodes()
00174 // and
00175 // common method
00177 
00178 template<unsigned DIM>
00179 void FineCoarseMeshPair<DIM>::ComputeFineElementsAndWeightsForCoarseQuadPoints(GaussianQuadratureRule<DIM>& rQuadRule,
00180                                                                                bool safeMode)
00181 {
00182     if (mpFineMeshBoxCollection == NULL)
00183     {
00184         EXCEPTION("Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseQuadPoints()");
00185     }
00186 
00187     // Get the quad point (physical) positions
00188     QuadraturePointsGroup<DIM> quad_point_posns(mrCoarseMesh, rQuadRule);
00189 
00190     // Resize the elements and weights vector.
00191     mFineMeshElementsAndWeights.resize(quad_point_posns.Size());
00192 
00193     #ifdef FINECOARSEMESHPAIR_VERBOSE
00194     std::cout << "\nComputing fine elements and weights for coarse quad points\n";
00195     #endif
00196 
00197     ResetStatisticsVariables();
00198     for (unsigned i=0; i<quad_point_posns.Size(); i++)
00199     {
00200         #ifdef FINECOARSEMESHPAIR_VERBOSE
00201         std::cout << "\r " << i << " of " << quad_point_posns.Size() << std::flush;
00202         #endif
00203 
00204         // Get the box this point is in
00205         unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( quad_point_posns.Get(i) );
00206 
00207         // A chaste point version of the c-vector is needed for the GetContainingElement call.
00208         ChastePoint<DIM> point(quad_point_posns.Get(i));
00209 
00210         ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
00211     }
00212 
00213     if (mStatisticsCounters[1] > 0)
00214     {
00215         WARNING(mStatisticsCounters[1] << " of " << quad_point_posns.Size() << " coarse-mesh quadrature points were outside the fine mesh");
00216     }
00217 }
00218 
00219 template<unsigned DIM>
00220 void FineCoarseMeshPair<DIM>::ComputeFineElementsAndWeightsForCoarseNodes(bool safeMode)
00221 {
00222     if (mpFineMeshBoxCollection==NULL)
00223     {
00224         EXCEPTION("Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseNodes()");
00225     }
00226 
00227     // Resize the elements and weights vector.
00228     mFineMeshElementsAndWeights.resize(mrCoarseMesh.GetNumNodes());
00229 
00230     #ifdef FINECOARSEMESHPAIR_VERBOSE
00231     std::cout << "\nComputing fine elements and weights for coarse nodes\n";
00232     #endif
00233 
00234     ResetStatisticsVariables();
00235     for (unsigned i=0; i<mrCoarseMesh.GetNumNodes(); i++)
00236     {
00237         #ifdef FINECOARSEMESHPAIR_VERBOSE
00238         std::cout << "\r " << i << " of " << mrCoarseMesh.GetNumNodes() << std::flush;
00239         #endif
00240 
00241         Node<DIM>* p_node = mrCoarseMesh.GetNode(i);
00242 
00243         // Get the box this point is in
00244         unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( p_node->rGetModifiableLocation() );
00245 
00246         // A chaste point version of the c-vector is needed for the GetContainingElement call
00247         ChastePoint<DIM> point(p_node->rGetLocation());
00248 
00249         ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
00250     }
00251 }
00252 
00259 template<unsigned DIM>
00260 void FineCoarseMeshPair<DIM>::ComputeFineElementAndWeightForGivenPoint(ChastePoint<DIM>& rPoint,
00261                                                                        bool safeMode,
00262                                                                        unsigned boxForThisPoint,
00263                                                                        unsigned index)
00264 {
00265     std::set<unsigned> test_element_indices;
00266 
00267     /*
00268      * The elements to try (initially) are those contained in the box the point is in.
00269      *
00270      * Note: it is possible the point to be in an element not 'in' this box, as it is
00271      * possible for all element nodes to be in different boxes.
00272      */
00273     CollectElementsInContainingBox(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
00274 
00275     unsigned elem_index;
00276     c_vector<double,DIM+1> weight;
00277 
00278     try
00279     {
00280         //std::cout << "\n" << "# test elements initially " << test_element_indices.size() << "\n";
00281         // Try these elements only, initially
00282         elem_index = mrFineMesh.GetContainingElementIndex(rPoint,
00283                                                           false,
00284                                                           test_element_indices,
00285                                                           true /* quit if not in test_elements */);
00286         weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
00287 
00288         mStatisticsCounters[0]++;
00289     }
00290     catch(Exception& e)
00291     {
00292         // Now try all elements, trying the elements contained in the boxes local to this element first
00293         std::set<unsigned> test_element_indices;
00294 
00295         CollectElementsInLocalBoxes(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
00296 
00297         try
00298         {
00299             elem_index = mrFineMesh.GetContainingElementIndex(rPoint, false, test_element_indices, true);
00300             weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
00301             mStatisticsCounters[0]++;
00302         }
00303         catch(Exception& e)
00304         {
00305             if (safeMode)
00306             {
00307                 // Try the remaining elements
00308                 try
00309                 {
00310                     elem_index = mrFineMesh.GetContainingElementIndex(rPoint, false);
00311                     weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
00312                     mStatisticsCounters[0]++;
00313 
00314                 }
00315                 catch (Exception& e)
00316                 {
00317                     // The point is not in ANY element, so store the nearest element and corresponding weights
00318                     elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
00319                     weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
00320 
00321                     mNotInMesh.push_back(index);
00322                     mNotInMeshNearestElementWeights.push_back(weight);
00323                     mStatisticsCounters[1]++;
00324                 }
00325             }
00326             else
00327             {
00328                 assert(test_element_indices.size() > 0); // boxes probably too small if this fails
00329 
00330                 /*
00331                  * Immediately assume it isn't in the rest of the mesh - this should be the
00332                  * case assuming the box width was chosen suitably. Store the nearest element
00333                  * and corresponding weights.
00334                  */
00335                 elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
00336                 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
00337 
00338                 mNotInMesh.push_back(index);
00339                 mNotInMeshNearestElementWeights.push_back(weight);
00340                 mStatisticsCounters[1]++;
00341             }
00342         }
00343     }
00344 
00345     mFineMeshElementsAndWeights[index].ElementNum = elem_index;
00346     mFineMeshElementsAndWeights[index].Weights = weight;
00347 }
00348 
00350 // ComputeCoarseElementsForFineNodes
00351 // and
00352 // ComputeCoarseElementsForFineElementCentroids
00353 // and
00354 // common method
00356 
00357 template<unsigned DIM>
00358 void FineCoarseMeshPair<DIM>::ComputeCoarseElementsForFineNodes(bool safeMode)
00359 {
00360     if (mpCoarseMeshBoxCollection==NULL)
00361     {
00362         EXCEPTION("Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineNodes()");
00363     }
00364 
00365     #ifdef FINECOARSEMESHPAIR_VERBOSE
00366     std::cout << "\nComputing coarse elements for fine nodes\n";
00367     #endif
00368 
00369     mCoarseElementsForFineNodes.resize(mrFineMesh.GetNumNodes());
00370 
00371     ResetStatisticsVariables();
00372     for (unsigned i=0; i<mCoarseElementsForFineNodes.size(); i++)
00373     {
00374         #ifdef FINECOARSEMESHPAIR_VERBOSE
00375         std::cout << "\r " << i << " of " << mCoarseElementsForFineNodes.size() << std::flush;
00376         #endif
00377 
00378         ChastePoint<DIM> point = mrFineMesh.GetNode(i)->GetPoint();
00379 
00380         // Get the box this point is in
00381         unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox(mrFineMesh.GetNode(i)->rGetModifiableLocation());
00382 
00383         mCoarseElementsForFineNodes[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
00384     }
00385 }
00386 
00387 template<unsigned DIM>
00388 void FineCoarseMeshPair<DIM>::ComputeCoarseElementsForFineElementCentroids(bool safeMode)
00389 {
00390     if (mpCoarseMeshBoxCollection==NULL)
00391     {
00392         EXCEPTION("Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineElementCentroids()");
00393     }
00394 
00395     #ifdef FINECOARSEMESHPAIR_VERBOSE
00396     std::cout << "\nComputing coarse elements for fine element centroids\n";
00397     #endif
00398 
00399     mCoarseElementsForFineElementCentroids.resize(mrFineMesh.GetNumElements());
00400 
00401     ResetStatisticsVariables();
00402     for (unsigned i=0; i<mrFineMesh.GetNumElements(); i++)
00403     {
00404         #ifdef FINECOARSEMESHPAIR_VERBOSE
00405         std::cout << "\r " << i << " of " << mrFineMesh.GetNumElements() << std::flush;
00406         #endif
00407 
00408         c_vector<double,DIM> point_cvec = mrFineMesh.GetElement(i)->CalculateCentroid();
00409         ChastePoint<DIM> point(point_cvec);
00410 
00411         // Get the box this point is in
00412         unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox( point_cvec );
00413 
00414         mCoarseElementsForFineElementCentroids[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
00415     }
00416 }
00417 
00424 template<unsigned DIM>
00425 unsigned FineCoarseMeshPair<DIM>::ComputeCoarseElementForGivenPoint(ChastePoint<DIM>& rPoint,
00426                                                                     bool safeMode,
00427                                                                     unsigned boxForThisPoint)
00428 {
00429     std::set<unsigned> test_element_indices;
00430     CollectElementsInContainingBox(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
00431 
00432     unsigned elem_index;
00433 
00434     try
00435     {
00436         elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint,
00437                                                             false,
00438                                                             test_element_indices,
00439                                                             true /* quit if not in test_elements */);
00440 
00441         mStatisticsCounters[0]++;
00442     }
00443     catch(Exception& e)
00444     {
00445         // Now try all elements, trying the elements contained in the boxes local to this element first
00446         std::set<unsigned> test_element_indices;
00447         CollectElementsInLocalBoxes(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
00448 
00449         try
00450         {
00451             elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint, false, test_element_indices, true);
00452             mStatisticsCounters[0]++;
00453         }
00454         catch(Exception& e)
00455         {
00456             if (safeMode)
00457             {
00458                 // Try the remaining elements
00459                 try
00460                 {
00461                     elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint, false);
00462 
00463                     mStatisticsCounters[0]++;
00464                 }
00465                 catch (Exception& e)
00466                 {
00467                     // The point is not in ANY element, so store the nearest element and corresponding weights
00468                     elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
00469                     mStatisticsCounters[1]++;
00470                 }
00471             }
00472             else
00473             {
00474                 assert(test_element_indices.size() > 0); // boxes probably too small if this fails
00475 
00476                 /*
00477                  * Immediately assume it isn't in the rest of the mesh - this should be the
00478                  * case assuming the box width was chosen suitably. Store the nearest element
00479                  * and corresponding weights.
00480                  */
00481                 elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
00482                 mStatisticsCounters[1]++;
00483             }
00484         }
00485     }
00486 
00487     return elem_index;
00488 }
00489 
00491 // Helper methods for code
00493 
00494 template<unsigned DIM>
00495 void FineCoarseMeshPair<DIM>::CollectElementsInContainingBox(BoxCollection<DIM>*& rpBoxCollection,
00496                                                              unsigned boxIndex,
00497                                                              std::set<unsigned>& rElementIndices)
00498 {
00499     for (typename std::set<Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->rGetBox(boxIndex).rGetElementsContained().begin();
00500          elem_iter != rpBoxCollection->rGetBox(boxIndex).rGetElementsContained().end();
00501          ++elem_iter)
00502     {
00503         rElementIndices.insert((*elem_iter)->GetIndex());
00504     }
00505 }
00506 
00507 template<unsigned DIM>
00508 void FineCoarseMeshPair<DIM>::CollectElementsInLocalBoxes(BoxCollection<DIM>*& rpBoxCollection,
00509                                                           unsigned boxIndex,
00510                                                           std::set<unsigned>& rElementIndices)
00511 {
00512     std::set<unsigned> local_boxes = rpBoxCollection->GetLocalBoxes(boxIndex);
00513     for (std::set<unsigned>::iterator local_box_iter = local_boxes.begin();
00514          local_box_iter != local_boxes.end();
00515          ++local_box_iter)
00516     {
00517         for (typename std::set<Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->rGetBox(*local_box_iter).rGetElementsContained().begin();
00518              elem_iter != rpBoxCollection->rGetBox(*local_box_iter).rGetElementsContained().end();
00519              ++elem_iter)
00520         {
00521             rElementIndices.insert((*elem_iter)->GetIndex());
00522         }
00523     }
00524 }
00525 
00527 // Statistics related methods
00529 
00530 template<unsigned DIM>
00531 void FineCoarseMeshPair<DIM>::ResetStatisticsVariables()
00532 {
00533     mNotInMesh.clear();
00534     mNotInMeshNearestElementWeights.clear();
00535     mStatisticsCounters.resize(2, 0u);
00536 }
00537 
00538 template<unsigned DIM>
00539 void FineCoarseMeshPair<DIM>::PrintStatistics()
00540 {
00541     std::cout << "\nFineCoarseMeshPair statistics for the last-called method:\n";
00542 
00543 //    std::cout << "\tNum points for which containing element was found, using box containing that point = " << mStatisticsCounters[0] << "\n";
00544 //    std::cout << "\tNum points for which containing element was in local box = " << mStatisticsCounters[1] << "\n";
00545 //    std::cout << "\tNum points for which containing element was in an element in a non-local box = " << mStatisticsCounters[2] << "\n";
00546 //    std::cout << "\tNum points for which no containing element was found = " << mStatisticsCounters[3] << "\n";
00547 
00548     std::cout << "\tNum points for which containing element was found: " << mStatisticsCounters[0] << "\n";
00549     std::cout << "\tNum points for which no containing element was found = " << mStatisticsCounters[1] << "\n";
00550 
00551     if (mNotInMesh.size() > 0)
00552     {
00553         std::cout << "\tIndices and weights for points (nodes/quad points) for which no containing element was found:\n";
00554         for (unsigned i=0; i<mNotInMesh.size(); i++)
00555         {
00556             std::cout << "\t\t" << mNotInMesh[i] << ", " << mNotInMeshNearestElementWeights[i] << "\n";
00557         }
00558     }
00559 }
00560 
00562 // Explicit instantiation
00564 
00565 template class FineCoarseMeshPair<1>;
00566 template class FineCoarseMeshPair<2>;
00567 template class FineCoarseMeshPair<3>;
Generated on Thu Dec 22 13:00:17 2011 for Chaste by  doxygen 1.6.3