ug4
marking_utils_impl.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2011-2018: G-CSC, Goethe University Frankfurt
3  * Author: Markus Breit
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 
34 #include "lib_disc/common/multi_index.h" // for DoFIndex
35 #include "lib_disc/domain_traits.h" // for domain_traits
36 #include "lib_grid/refinement/refiner_interface.h" // for IRefiner
37 
38 
39 namespace ug {
40 
41 
42 // //////////////////////////////////////////////////////////////////////////////
43 // Check values are within bounds
44 // //////////////////////////////////////////////////////////////////////////////
45 template <typename TGridFunction, typename TBaseElem>
46 unsigned long MarkOutOfRangeElems
47 (
48  SmartPtr<IRefiner> refiner,
50  size_t cmp,
51  number lowerBnd,
52  number upperBnd
53 )
54 {
55  typedef typename TGridFunction::template traits<TBaseElem>::const_iterator elem_it;
57  typedef typename Grid::traits<elem_type>::secure_container elem_list;
58 
59  Grid& grid = *refiner->grid();
60 
61  unsigned long nMarked = 0;
62  elem_list el;
63 
64  // loop the elements in the subset
65  std::vector<DoFIndex> vDI;
66  elem_it it = u->template begin<TBaseElem>();
67  elem_it itEnd = u->template end<TBaseElem>();
68  for (; it != itEnd; ++it)
69  {
70  TBaseElem* elem = *it;
71 
72  // loop indices at this element
73  const size_t nInd = u->inner_dof_indices(elem, cmp, vDI, true);
74  for (size_t i = 0; i < nInd; ++i)
75  {
76  const number& val = DoFRef(*u, vDI[i]);
77  if (val < lowerBnd || val > upperBnd)
78  {
79  // mark neighbors for refinement
80  grid.associated_elements(el, elem);
81  const size_t elSz = el.size();
82  for (size_t e = 0; e < elSz; ++e)
83  {
84  if (refiner->get_mark(el[e]) == RM_NONE)
85  {
86  refiner->mark(el[e], RM_FULL);
87  ++nMarked;
88  }
89  }
90  }
91  }
92  }
93 
94  return nMarked;
95 }
96 
97 template <typename TGridFunction>
99 (
100  SmartPtr<IRefiner> refiner,
102  size_t cmp,
103  number lowerBnd,
104  number upperBnd
105 )
106 {
107  unsigned long nMarked = 0;
108  if (u->max_fct_dofs(cmp, 0))
109  nMarked += MarkOutOfRangeElems<TGridFunction, Vertex>(refiner, u, cmp, lowerBnd, upperBnd);
110  if (u->max_fct_dofs(cmp, 1))
111  nMarked += MarkOutOfRangeElems<TGridFunction, Edge>(refiner, u, cmp, lowerBnd, upperBnd);
112  if (u->max_fct_dofs(cmp, 2))
113  nMarked += MarkOutOfRangeElems<TGridFunction, Face>(refiner, u, cmp, lowerBnd, upperBnd);
114  if (u->max_fct_dofs(cmp, 3))
115  nMarked += MarkOutOfRangeElems<TGridFunction, Volume>(refiner, u, cmp, lowerBnd, upperBnd);
116 
117 #ifdef UG_PARALLEL
118  if (pcl::NumProcs() > 1)
119  {
121  nMarked = pc.allreduce(nMarked, PCL_RO_SUM);
122  }
123 #endif
124 
125 
126  if (nMarked)
127  UG_LOGN(" +++ Marked for refinement: " << nMarked << " elements.");
128 }
129 
130 } // namespace ug
Definition: smart_pointer.h:296
Definition: smart_pointer.h:108
Definition: pcl_process_communicator.h:70
void allreduce(const void *sendBuf, void *recBuf, int count, DataType type, ReduceOperation op) const
performs MPI_Allreduce on the processes of the communicator.
Definition: pcl_process_communicator.cpp:318
Manages the elements of a grid and their interconnection.
Definition: grid.h:132
void associated_elements(traits< Vertex >::secure_container &elemsOut, TElem *e)
Puts all elements of type TAss which are contained in 'e' or which contain 'e' into elemsOut.
Definition: grid_impl.hpp:466
Container which holds an array of pointers.
Definition: pointer_const_array.h:84
@ RM_NONE
no refinement is performed
Definition: refiner_interface.h:49
@ RM_FULL
Fully refines an element and all associated sides and edges.
Definition: refiner_interface.h:54
#define PCL_RO_SUM
Definition: pcl_methods.h:63
int NumProcs()
returns the number of processes
Definition: pcl_base.cpp:91
base_type::TBaseElem TBaseElem
#define UG_LOGN(msg)
Definition: log.h:369
double number
Definition: types.h:124
the ug namespace
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition: multi_index.h:276
void MarkOutOfRangeElems(SmartPtr< IRefiner > refiner, ConstSmartPtr< TGridFunction > u, size_t cmp, number lowerBnd, number upperBnd)
Definition: marking_utils_impl.h:99
Definition: domain_traits.h:53