Chaste  Release::3.4
FineCoarseMeshPair.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "FineCoarseMeshPair.hpp"
37 
38 template<unsigned DIM>
40  : mrFineMesh(rFineMesh),
41  mrCoarseMesh(rCoarseMesh),
42  mpFineMeshBoxCollection(NULL),
43  mpCoarseMeshBoxCollection(NULL)
44 {
46 }
47 
48 template<unsigned DIM>
50 {
51  return mrFineMesh;
52 }
53 
54 template<unsigned DIM>
56 {
57  return mrCoarseMesh;
58 }
59 
60 template<unsigned DIM>
62 {
63  DeleteFineBoxCollection();
64  DeleteCoarseBoxCollection();
65 }
66 
67 template<unsigned DIM>
69 {
70  if (mpFineMeshBoxCollection != NULL)
71  {
72  delete mpFineMeshBoxCollection;
73  mpFineMeshBoxCollection = NULL;
74  }
75 }
76 
77 template<unsigned DIM>
79 {
80  if (mpCoarseMeshBoxCollection != NULL)
81  {
82  delete mpCoarseMeshBoxCollection;
83  mpCoarseMeshBoxCollection = NULL;
84  }
85 }
86 
88 // Setting up boxes methods
90 
91 template<unsigned DIM>
93 {
94  SetUpBoxes(mrFineMesh, boxWidth, mpFineMeshBoxCollection);
95 }
96 
97 template<unsigned DIM>
99 {
100  SetUpBoxes(mrCoarseMesh, boxWidth, mpCoarseMeshBoxCollection);
101 }
102 
103 template<unsigned DIM>
105  double boxWidth,
106  BoxCollection<DIM>*& rpBoxCollection)
107 {
108  if (rpBoxCollection)
109  {
110  delete rpBoxCollection;
111  rpBoxCollection = NULL;
112  }
113 
114  // Compute min and max values for the fine mesh nodes
115  ChasteCuboid<DIM> bounding_box = rMesh.CalculateBoundingBox();
116 
117  // Set up the boxes using a domain that is slightly larger than the fine mesh
118  c_vector<double,2*DIM> extended_min_and_max;
119  for (unsigned i=0; i<DIM; i++)
120  {
121  double width = bounding_box.GetWidth(i);
122 
123  // Subtract from the minima
124  extended_min_and_max(2*i) = bounding_box.rGetLowerCorner()[i] - 0.05*width;
125 
126  // Add to the maxima
127  extended_min_and_max(2*i+1) = bounding_box.rGetUpperCorner()[i] + 0.05*width;
128  }
129 
130  if (boxWidth < 0)
131  {
132  /*
133  * Use default value = max(max_edge_length, w20), where w20 is the width
134  * corresponding to 20 boxes in the x-direction.
135  *
136  * BoxCollection creates an extra box so divide by 19 not 20. Add a little
137  * bit on to ensure minor numerical fluctuations don't change the answer.
138  */
139  boxWidth = (extended_min_and_max(1) - extended_min_and_max(0))/19.000000001;
140 
141  // Determine the maximum edge length
142  c_vector<double, 2> min_max_edge_length = rMesh.CalculateMinMaxEdgeLengths();
143 
144  if (boxWidth < min_max_edge_length[1])
145  {
146  boxWidth = 1.1*min_max_edge_length[1];
147  }
148  }
149 
150  rpBoxCollection = new BoxCollection<DIM>(boxWidth, extended_min_and_max);
151  rpBoxCollection->SetupAllLocalBoxes();
152 
153  // For each element, if ANY of its nodes are physically in a box, put that element in that box
154  for (unsigned i=0; i<rMesh.GetNumElements(); i++)
155  {
156  Element<DIM,DIM>* p_element = rMesh.GetElement(i);
157 
158  std::set<unsigned> box_indices_each_node_this_elem;
159  for (unsigned j=0; j<DIM+1; j++) // num vertices per element
160  {
161  Node<DIM>* p_node = p_element->GetNode(j);
162  unsigned box_index = rpBoxCollection->CalculateContainingBox(p_node);
163  box_indices_each_node_this_elem.insert(box_index);
164  }
165 
166  for (std::set<unsigned>::iterator iter = box_indices_each_node_this_elem.begin();
167  iter != box_indices_each_node_this_elem.end();
168  ++iter)
169  {
170  rpBoxCollection->rGetBox( *iter ).AddElement(p_element);
171  }
172  }
173 }
174 
176 // ComputeFineElementsAndWeightsForCoarseQuadPoints()
177 // and
178 // ComputeFineElementsAndWeightsForCoarseNodes()
179 // and
180 // common method
182 
183 template<unsigned DIM>
185  bool safeMode)
186 {
187  if (mpFineMeshBoxCollection == NULL)
188  {
189  EXCEPTION("Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseQuadPoints()");
190  }
191 
192  // Get the quad point (physical) positions
193  QuadraturePointsGroup<DIM> quad_point_posns(mrCoarseMesh, rQuadRule);
194 
195  // Resize the elements and weights vector.
196  mFineMeshElementsAndWeights.resize(quad_point_posns.Size());
197 
198  #define COVERAGE_IGNORE
199  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
200  {
201  std::cout << "\nComputing fine elements and weights for coarse quad points\n";
202  }
203  #undef COVERAGE_IGNORE
204 
205 
206  ResetStatisticsVariables();
207  for (unsigned i=0; i<quad_point_posns.Size(); i++)
208  {
209  #define COVERAGE_IGNORE
210  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
211  {
212  std::cout << "\t" << i << " of " << quad_point_posns.Size() << std::flush;
213  }
214  #undef COVERAGE_IGNORE
215 
216  // Get the box this point is in
217  unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( quad_point_posns.rGet(i) );
218 
219  // A chaste point version of the c-vector is needed for the GetContainingElement call.
220  ChastePoint<DIM> point(quad_point_posns.rGet(i));
221 
222  ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
223  }
224 
225  if (mStatisticsCounters[1] > 0)
226  {
227  WARNING(mStatisticsCounters[1] << " of " << quad_point_posns.Size() << " coarse-mesh quadrature points were outside the fine mesh");
228  }
229 }
230 
231 template<unsigned DIM>
233 {
234  if (mpFineMeshBoxCollection==NULL)
235  {
236  EXCEPTION("Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseNodes()");
237  }
238 
239  // Resize the elements and weights vector.
240  mFineMeshElementsAndWeights.resize(mrCoarseMesh.GetNumNodes());
241 
242  #define COVERAGE_IGNORE
243  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
244  {
245  std::cout << "\nComputing fine elements and weights for coarse nodes\n";
246  }
247  #undef COVERAGE_IGNORE
248 
249 
250  ResetStatisticsVariables();
251  for (unsigned i=0; i<mrCoarseMesh.GetNumNodes(); i++)
252  {
253  #define COVERAGE_IGNORE
254  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
255  {
256  std::cout << "\t" << i << " of " << mrCoarseMesh.GetNumNodes() << std::flush;
257  }
258  #undef COVERAGE_IGNORE
259 
260  Node<DIM>* p_node = mrCoarseMesh.GetNode(i);
261 
262  // Get the box this point is in
263  unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( p_node->rGetModifiableLocation() );
264 
265  // A chaste point version of the c-vector is needed for the GetContainingElement call
266  ChastePoint<DIM> point(p_node->rGetLocation());
267 
268  ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
269  }
270 }
271 
278 template<unsigned DIM>
280  bool safeMode,
281  unsigned boxForThisPoint,
282  unsigned index)
283 {
284  std::set<unsigned> test_element_indices;
285 
286  /*
287  * The elements to try (initially) are those contained in the box the point is in.
288  *
289  * Note: it is possible the point to be in an element not 'in' this box, as it is
290  * possible for all element nodes to be in different boxes.
291  */
292  CollectElementsInContainingBox(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
293 
294  unsigned elem_index;
295  c_vector<double,DIM+1> weight;
296 
297  try
298  {
299  //std::cout << "\n" << "# test elements initially " << test_element_indices.size() << "\n";
300  // Try these elements only, initially
301  elem_index = mrFineMesh.GetContainingElementIndex(rPoint,
302  false,
303  test_element_indices,
304  true /* quit if not in test_elements */);
305  weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
306 
307  mStatisticsCounters[0]++;
308  }
309  catch(Exception&) //not_in_box
310  {
311  // Now try all elements, trying the elements contained in the boxes local to this element first
312  test_element_indices.clear();
313 
314  CollectElementsInLocalBoxes(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
315 
316  try
317  {
318  elem_index = mrFineMesh.GetContainingElementIndex(rPoint, false, test_element_indices, true);
319  weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
320  mStatisticsCounters[0]++;
321  }
322  catch(Exception&) //not_in_local_boxes
323  {
324  if (safeMode)
325  {
326  // Try the remaining elements
327  try
328  {
329  elem_index = mrFineMesh.GetContainingElementIndex(rPoint, false);
330  weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
331  mStatisticsCounters[0]++;
332 
333  }
334  catch (Exception&) // not_in_mesh
335  {
336  // The point is not in ANY element, so store the nearest element and corresponding weights
337  elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
338  weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
339 
340  mNotInMesh.push_back(index);
341  mNotInMeshNearestElementWeights.push_back(weight);
342  mStatisticsCounters[1]++;
343  }
344  }
345  else
346  {
347  assert(test_element_indices.size() > 0); // boxes probably too small if this fails
348 
349  /*
350  * Immediately assume it isn't in the rest of the mesh - this should be the
351  * case assuming the box width was chosen suitably. Store the nearest element
352  * and corresponding weights.
353  */
354  elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
355  weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
356 
357  mNotInMesh.push_back(index);
358  mNotInMeshNearestElementWeights.push_back(weight);
359  mStatisticsCounters[1]++;
360  }
361  }
362  }
363 
364  mFineMeshElementsAndWeights[index].ElementNum = elem_index;
365  mFineMeshElementsAndWeights[index].Weights = weight;
366 }
367 
369 // ComputeCoarseElementsForFineNodes
370 // and
371 // ComputeCoarseElementsForFineElementCentroids
372 // and
373 // common method
375 
376 template<unsigned DIM>
378 {
379  if (mpCoarseMeshBoxCollection==NULL)
380  {
381  EXCEPTION("Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineNodes()");
382  }
383 
384  #define COVERAGE_IGNORE
385  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
386  {
387  std::cout << "\nComputing coarse elements for fine nodes\n";
388  }
389  #undef COVERAGE_IGNORE
390 
391  mCoarseElementsForFineNodes.resize(mrFineMesh.GetNumNodes());
392 
393  ResetStatisticsVariables();
394  for (unsigned i=0; i<mCoarseElementsForFineNodes.size(); i++)
395  {
396  #define COVERAGE_IGNORE
397  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
398  {
399  std::cout << "\t" << i << " of " << mCoarseElementsForFineNodes.size() << std::flush;
400  }
401  #undef COVERAGE_IGNORE
402 
403  ChastePoint<DIM> point = mrFineMesh.GetNode(i)->GetPoint();
404 
405  // Get the box this point is in
406  unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox(mrFineMesh.GetNode(i)->rGetModifiableLocation());
407 
408  mCoarseElementsForFineNodes[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
409  }
410 }
411 
412 template<unsigned DIM>
414 {
415  if (mpCoarseMeshBoxCollection==NULL)
416  {
417  EXCEPTION("Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineElementCentroids()");
418  }
419 
420  #define COVERAGE_IGNORE
421  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
422  {
423  std::cout << "\nComputing coarse elements for fine element centroids\n";
424  }
425  #undef COVERAGE_IGNORE
426 
427  mCoarseElementsForFineElementCentroids.resize(mrFineMesh.GetNumElements());
428 
429  ResetStatisticsVariables();
430  for (unsigned i=0; i<mrFineMesh.GetNumElements(); i++)
431  {
432  #define COVERAGE_IGNORE
433  if(CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
434  {
435  std::cout << "\t" << i << " of " << mrFineMesh.GetNumElements() << std::flush;
436  }
437  #undef COVERAGE_IGNORE
438 
439  c_vector<double,DIM> point_cvec = mrFineMesh.GetElement(i)->CalculateCentroid();
440  ChastePoint<DIM> point(point_cvec);
441 
442  // Get the box this point is in
443  unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox( point_cvec );
444 
445  mCoarseElementsForFineElementCentroids[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
446  }
447 }
448 
449 template<unsigned DIM>
451  bool safeMode,
452  unsigned boxForThisPoint)
453 {
460  std::set<unsigned> test_element_indices;
461  CollectElementsInContainingBox(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
462 
463  unsigned elem_index;
464 
465  try
466  {
467  elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint,
468  false,
469  test_element_indices,
470  true /* quit if not in test_elements */);
471 
472  mStatisticsCounters[0]++;
473  }
474  catch(Exception&) // not_in_box
475  {
476  // Now try all elements, trying the elements contained in the boxes local to this element first
477  test_element_indices.clear();
478  CollectElementsInLocalBoxes(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
479 
480  try
481  {
482  elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint, false, test_element_indices, true);
483  mStatisticsCounters[0]++;
484  }
485  catch(Exception&) // not_in_local_boxes
486  {
487  if (safeMode)
488  {
489  // Try the remaining elements
490  try
491  {
492  elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint, false);
493 
494  mStatisticsCounters[0]++;
495  }
496  catch (Exception&) // not_in_mesh
497  {
498  // The point is not in ANY element, so store the nearest element and corresponding weights
499  elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
500  mStatisticsCounters[1]++;
501  }
502  }
503  else
504  {
505  assert(test_element_indices.size() > 0); // boxes probably too small if this fails
506 
507  /*
508  * Immediately assume it isn't in the rest of the mesh - this should be the
509  * case assuming the box width was chosen suitably. Store the nearest element
510  * and corresponding weights.
511  */
512  elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
513  mStatisticsCounters[1]++;
514  }
515  }
516  }
517 
518  return elem_index;
519 }
520 
522 // Helper methods for code
524 
525 template<unsigned DIM>
527  unsigned boxIndex,
528  std::set<unsigned>& rElementIndices)
529 {
530  for (typename std::set<Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->rGetBox(boxIndex).rGetElementsContained().begin();
531  elem_iter != rpBoxCollection->rGetBox(boxIndex).rGetElementsContained().end();
532  ++elem_iter)
533  {
534  rElementIndices.insert((*elem_iter)->GetIndex());
535  }
536 }
537 
538 template<unsigned DIM>
540  unsigned boxIndex,
541  std::set<unsigned>& rElementIndices)
542 {
543  std::set<unsigned> local_boxes = rpBoxCollection->GetLocalBoxes(boxIndex);
544  for (std::set<unsigned>::iterator local_box_iter = local_boxes.begin();
545  local_box_iter != local_boxes.end();
546  ++local_box_iter)
547  {
548  for (typename std::set<Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->rGetBox(*local_box_iter).rGetElementsContained().begin();
549  elem_iter != rpBoxCollection->rGetBox(*local_box_iter).rGetElementsContained().end();
550  ++elem_iter)
551  {
552  rElementIndices.insert((*elem_iter)->GetIndex());
553  }
554  }
555 }
556 
558 // Statistics related methods
560 
561 template<unsigned DIM>
563 {
564  mNotInMesh.clear();
565  mNotInMeshNearestElementWeights.clear();
566  mStatisticsCounters.resize(2, 0u);
567 }
568 
569 template<unsigned DIM>
571 {
572  std::cout << "\nFineCoarseMeshPair statistics for the last-called method:\n";
573 
574 // std::cout << "\tNum points for which containing element was found, using box containing that point = " << mStatisticsCounters[0] << "\n";
575 // std::cout << "\tNum points for which containing element was in local box = " << mStatisticsCounters[1] << "\n";
576 // std::cout << "\tNum points for which containing element was in an element in a non-local box = " << mStatisticsCounters[2] << "\n";
577 // std::cout << "\tNum points for which no containing element was found = " << mStatisticsCounters[3] << "\n";
578 
579  std::cout << "\tNum points for which containing element was found: " << mStatisticsCounters[0] << "\n";
580  std::cout << "\tNum points for which no containing element was found = " << mStatisticsCounters[1] << "\n";
581 
582  if (mNotInMesh.size() > 0)
583  {
584  std::cout << "\tIndices and weights for points (nodes/quad points) for which no containing element was found:\n";
585  for (unsigned i=0; i<mNotInMesh.size(); i++)
586  {
587  std::cout << "\t\t" << mNotInMesh[i] << ", " << mNotInMeshNearestElementWeights[i] << "\n";
588  }
589  }
590 }
591 
593 // Explicit instantiation
595 
596 template class FineCoarseMeshPair<1>;
597 template class FineCoarseMeshPair<2>;
598 template class FineCoarseMeshPair<3>;
void ComputeFineElementAndWeightForGivenPoint(ChastePoint< DIM > &rPoint, bool safeMode, unsigned boxForThisPoint, unsigned index)
Definition: Node.hpp:58
FineCoarseMeshPair(AbstractTetrahedralMesh< DIM, DIM > &rFineMesh, AbstractTetrahedralMesh< DIM, DIM > &rCoarseMesh)
void ComputeFineElementsAndWeightsForCoarseNodes(bool safeMode)
virtual unsigned GetNumElements() const
#define EXCEPTION(message)
Definition: Exception.hpp:143
Element< ELEMENT_DIM, SPACE_DIM > * GetElement(unsigned index) const
const AbstractTetrahedralMesh< DIM, DIM > & GetCoarseMesh() const
unsigned ComputeCoarseElementForGivenPoint(ChastePoint< DIM > &rPoint, bool safeMode, unsigned boxForThisPoint)
std::set< unsigned > GetLocalBoxes(unsigned boxIndex)
bool OptionExists(const std::string &rOption)
void ComputeFineElementsAndWeightsForCoarseQuadPoints(GaussianQuadratureRule< DIM > &rQuadRule, bool safeMode)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
void CollectElementsInContainingBox(BoxCollection< DIM > *&rpBoxCollection, unsigned boxIndex, std::set< unsigned > &rElementIndices)
const ChastePoint< SPACE_DIM > & rGetUpperCorner() const
double GetWidth(unsigned rDimension) const
unsigned CalculateContainingBox(Node< DIM > *pNode)
void SetUpBoxes(AbstractTetrahedralMesh< DIM, DIM > &rMesh, double boxWidth, BoxCollection< DIM > *&rpBoxCollection)
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition: Node.cpp:140
void SetupAllLocalBoxes()
c_vector< double, DIM > & rGet(unsigned elementIndex, unsigned quadIndex)
void SetUpBoxesOnCoarseMesh(double boxWidth=-1)
Box< DIM > & rGetBox(unsigned boxIndex)
static CommandLineArguments * Instance()
ChasteCuboid< SPACE_DIM > CalculateBoundingBox(const std::vector< Node< SPACE_DIM > * > &rNodes) const
void SetUpBoxesOnFineMesh(double boxWidth=-1)
void ComputeCoarseElementsForFineNodes(bool safeMode)
virtual c_vector< double, 2 > CalculateMinMaxEdgeLengths()
c_vector< double, SPACE_DIM > & rGetModifiableLocation()
Definition: Node.cpp:152
const AbstractTetrahedralMesh< DIM, DIM > & GetFineMesh() const
void CollectElementsInLocalBoxes(BoxCollection< DIM > *&rpBoxCollection, unsigned boxIndex, std::set< unsigned > &rElementIndices)
void ComputeCoarseElementsForFineElementCentroids(bool safeMode)
const ChastePoint< SPACE_DIM > & rGetLowerCorner() const