Line data Source code
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 0 : unsigned long MarkOutOfRangeElems
47 : (
48 : SmartPtr<IRefiner> refiner,
49 : ConstSmartPtr<TGridFunction> u,
50 : size_t cmp,
51 : number lowerBnd,
52 : number upperBnd
53 : )
54 : {
55 : typedef typename TGridFunction::template traits<TBaseElem>::const_iterator elem_it;
56 : typedef typename domain_traits<TGridFunction::domain_type::dim>::element_type elem_type;
57 : typedef typename Grid::traits<elem_type>::secure_container elem_list;
58 :
59 0 : 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 0 : 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 0 : for (size_t i = 0; i < nInd; ++i)
75 : {
76 : const number& val = DoFRef(*u, vDI[i]);
77 0 : 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 0 : for (size_t e = 0; e < elSz; ++e)
83 : {
84 0 : if (refiner->get_mark(el[e]) == RM_NONE)
85 : {
86 0 : refiner->mark(el[e], RM_FULL);
87 0 : ++nMarked;
88 : }
89 : }
90 : }
91 : }
92 : }
93 :
94 0 : return nMarked;
95 : }
96 :
97 : template <typename TGridFunction>
98 0 : void MarkOutOfRangeElems
99 : (
100 : SmartPtr<IRefiner> refiner,
101 : ConstSmartPtr<TGridFunction> u,
102 : size_t cmp,
103 : number lowerBnd,
104 : number upperBnd
105 : )
106 : {
107 : unsigned long nMarked = 0;
108 0 : if (u->max_fct_dofs(cmp, 0))
109 0 : nMarked += MarkOutOfRangeElems<TGridFunction, Vertex>(refiner, u, cmp, lowerBnd, upperBnd);
110 0 : if (u->max_fct_dofs(cmp, 1))
111 0 : nMarked += MarkOutOfRangeElems<TGridFunction, Edge>(refiner, u, cmp, lowerBnd, upperBnd);
112 0 : if (u->max_fct_dofs(cmp, 2))
113 0 : nMarked += MarkOutOfRangeElems<TGridFunction, Face>(refiner, u, cmp, lowerBnd, upperBnd);
114 0 : if (u->max_fct_dofs(cmp, 3))
115 0 : nMarked += MarkOutOfRangeElems<TGridFunction, Volume>(refiner, u, cmp, lowerBnd, upperBnd);
116 :
117 : #ifdef UG_PARALLEL
118 : if (pcl::NumProcs() > 1)
119 : {
120 : pcl::ProcessCommunicator pc;
121 : nMarked = pc.allreduce(nMarked, PCL_RO_SUM);
122 : }
123 : #endif
124 :
125 :
126 0 : if (nMarked)
127 : UG_LOGN(" +++ Marked for refinement: " << nMarked << " elements.");
128 0 : }
129 :
130 : } // namespace ug
|