Line data Source code
1 : /*
2 : * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Andreas Vogel
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 : #ifndef __H__UG__LIB_DISC__FUNCTION_SPACE__APPROXIMATION_SPACE__
34 : #define __H__UG__LIB_DISC__FUNCTION_SPACE__APPROXIMATION_SPACE__
35 :
36 : #include "lib_disc/common/revision_counter.h"
37 : #include "lib_disc/dof_manager/dof_distribution_info.h"
38 : #include "lib_disc/dof_manager/dof_distribution.h"
39 : #include "lib_grid/tools/surface_view.h"
40 : #include "lib_algebra/algebra_type.h"
41 :
42 : namespace ug{
43 :
44 : /// describes the ansatz spaces on a domain
45 : /**
46 : * This class provides grid function spaces on a domain.
47 : *
48 : * The Domain defines a partition of the Grid/Multigrid in terms of subsets.
49 : * The user can add discrete functions on this subsets or unions of them.
50 : *
51 : * Once finalized, this function pattern is fixed. Internally DoF indices are
52 : * created. Using this Approximation Space the user can create GridFunctions
53 : * of the following types:
54 : *
55 : * - surface grid function = grid function representing the space on the surface grid
56 : * - level grid function = grid function representing the space on a level grid
57 : * (NOTE: For a fully refined Multigrid a level grid covers the whole domain.
58 : * However for a locally/adaptively refined MultiGrid the level grid
59 : * solution is only living on a part of the domain)
60 : */
61 : class IApproximationSpace : public DoFDistributionInfoProvider
62 : {
63 : public:
64 : /// Type of Subset Handler
65 : typedef MGSubsetHandler subset_handler_type;
66 :
67 : /// Type of Subset Handler
68 : typedef MultiGrid grid_type;
69 :
70 : public:
71 : /// Constructor
72 : IApproximationSpace(SmartPtr<subset_handler_type> spMGSH,
73 : SmartPtr<grid_type> spMG);
74 :
75 : /// Constructor setting the grouping flag
76 : IApproximationSpace(SmartPtr<subset_handler_type> spMGSH,
77 : SmartPtr<grid_type>,
78 : const AlgebraType& algebraType);
79 :
80 : protected:
81 : /// initializing
82 : void init(SmartPtr<subset_handler_type> spMGSH,
83 : SmartPtr<grid_type> spMG, const AlgebraType& algebraType);
84 :
85 : public:
86 : /// Destructor
87 : ~IApproximationSpace();
88 :
89 : /// clears functions
90 0 : void clear() {m_spDoFDistributionInfo->clear();}
91 :
92 : public:
93 : /// add single solutions of LocalShapeFunctionSetID to the entire domain
94 : /**
95 : * \param[in] name name(s) of single solution (comma separated)
96 : * \param[in] id Shape Function set id
97 : */
98 : void add(const std::vector<std::string>& vName, LFEID id){
99 0 : m_spDoFDistributionInfo->add(vName, id);
100 : }
101 :
102 : /// adds function using string to indicate finite element type
103 : void add(const std::vector<std::string>& vName, const char* type, int order);
104 :
105 : /// adds function using string to indicate finite element type
106 : void add(const std::vector<std::string>& vName, const char* type);
107 :
108 : /// adds function using string to indicate finite element type
109 : void add(const char* name, const char* type, int order);
110 :
111 : /// adds function using string to indicate finite element type
112 : void add(const char* name, const char* type);
113 :
114 : public:
115 : /// add single solutions of LocalShapeFunctionSetID to selected subsets
116 : /**
117 : * \param[in] name name(s) of single solution (comma separated)
118 : * \param[in] id Shape Function set id
119 : * \param[in] subsets Subsets separated by ','
120 : */
121 : void add(const std::vector<std::string>& vName, LFEID id,
122 : const std::vector<std::string>& vSubset){
123 0 : m_spDoFDistributionInfo->add(vName, id, vSubset);
124 0 : }
125 :
126 : /// adds function using string to indicate finite element type
127 : void add(const std::vector<std::string>& vName, const char* type, int order,
128 : const std::vector<std::string>& vSubsets);
129 :
130 : /// adds function using string to indicate finite element type
131 : /**
132 : * \param[in] name name(s) of single solution (comma separated)
133 : * \param[in] type type of local finite element space
134 : * \param[in] order order of local finite element space
135 : * \param[in] subsets Subsets separated by ','
136 : */
137 : void add(const char* name, const char* type, int order, const char* subsets);
138 :
139 : /// adds function using string to indicate finite element type
140 : void add(const std::vector<std::string>& vName, const char* type,
141 : const std::vector<std::string>& vSubsets);
142 :
143 : /// adds function using string to indicate finite element type
144 : /**
145 : * \param[in] name name(s) of single solution (comma separated)
146 : * \param[in] type type of local finite element space
147 : * \param[in] subsets Subsets separated by ','
148 : */
149 : void add(const char* name, const char* type, const char* subsets);
150 :
151 : public:
152 : /// get underlying subset handler
153 : ConstSmartPtr<MGSubsetHandler> subset_handler() const {return m_spMGSH;}
154 :
155 : /// returns the number of level
156 0 : size_t num_levels() const {return m_spMGSH->num_levels();}
157 :
158 : /// returns the approximation space
159 0 : ConstSmartPtr<SurfaceView> surface_view() const {return m_spSurfaceView;}
160 :
161 : /// returns if dofs are grouped
162 : bool grouped() const {return m_bGrouped;}
163 :
164 : /// returns if ghosts might be present on a level
165 : bool might_contain_ghosts(int lvl) const;
166 :
167 : /// returns if ghosts might be present on any level
168 : bool might_contain_ghosts() const;
169 :
170 : /// returns dof distribution for a grid level
171 : /// \{
172 : SmartPtr<DoFDistribution> dof_distribution(const GridLevel& gl, bool bCreate = true);
173 : SmartPtr<DoFDistribution> dd(const GridLevel& gl, bool bCreate = true);
174 : /// \}
175 :
176 : /// returns dof distribution for a grid level
177 : /// \{
178 : ConstSmartPtr<DoFDistribution> dof_distribution(const GridLevel& gl, bool bCreate = true) const;
179 : ConstSmartPtr<DoFDistribution> dd(const GridLevel& gl, bool bCreate = true) const;
180 : /// \}
181 :
182 : /// returns all currently created dof distributions
183 : std::vector<SmartPtr<DoFDistribution> > dof_distributions() const;
184 :
185 : /// returns dof distribution info
186 : /// \{
187 : ConstSmartPtr<DoFDistributionInfo> ddinfo() const {return m_spDoFDistributionInfo;}
188 : ConstSmartPtr<DoFDistributionInfo> dof_distribution_info() const {return ddinfo();}
189 : /// \}
190 :
191 : /// prints statistic about DoF Distribution
192 : void print_statistic(std::string flags) const;
193 :
194 : /// prints statistic about DoF Distribution
195 : void print_statistic() const;
196 :
197 : /// prints statistic on layouts
198 : void print_layout_statistic() const;
199 :
200 :
201 : /// initializes all level dof distributions
202 : void init_levels();
203 :
204 : /// initializes all surface dof distributions
205 : void init_surfaces();
206 :
207 : /// initializes all top surface dof distributions
208 : void init_top_surface();
209 :
210 : /// returns the current revision
211 0 : const RevisionCounter& revision() const {return m_RevCnt;}
212 :
213 : protected:
214 : /// creates a dof distribution
215 : void create_dof_distribution(const GridLevel& gl);
216 :
217 : /// creates surface SurfaceView if needed
218 : void surface_view_required();
219 :
220 : /// create dof distribution info
221 : void dof_distribution_info_required();
222 :
223 : protected:
224 : /// reinits all data after grid adaption
225 : void reinit();
226 :
227 : /// message hub id
228 : MessageHub::SPCallbackId m_spGridAdaptionCallbackID;
229 : MessageHub::SPCallbackId m_spGridDistributionCallbackID;
230 : bool m_bAdaptionIsActive;
231 :
232 : /// registers at message hub for grid adaption
233 : void register_at_adaption_msg_hub();
234 :
235 : /** this callback is called by the message hub, when a grid change has been
236 : * performed. It will call all necessary actions in order to keep the grid
237 : * correct for computations.*/
238 : void grid_changed_callback(const GridMessage_Adaption& msg);
239 :
240 : /// called during parallel redistribution
241 : void grid_distribution_callback(const GridMessage_Distribution& msg);
242 :
243 : protected:
244 : /// multigrid, where elements are stored
245 : SmartPtr<MultiGrid> m_spMG;
246 :
247 : /// subsethandler, where elements are stored
248 : SmartPtr<MGSubsetHandler> m_spMGSH;
249 :
250 : /// Surface View
251 : SmartPtr<SurfaceView> m_spSurfaceView;
252 :
253 : /// suitable algebra type for the index distribution pattern
254 : AlgebraType m_algebraType;
255 :
256 : /// flag if DoFs should be grouped
257 : bool m_bGrouped;
258 :
259 : /// DofDistributionInfo
260 : SmartPtr<DoFDistributionInfo> m_spDoFDistributionInfo;
261 :
262 : /// revision counter
263 : RevisionCounter m_RevCnt;
264 :
265 : protected:
266 : /// MG Level DoF Distribution
267 : std::vector<SmartPtr<DoFDistribution> > m_vDD;
268 :
269 : /// Index Storage for Level (ghost / noghost)
270 : /// \{
271 : SmartPtr<DoFIndexStorage> m_spDoFIndexStrgForLevelNoGhost;
272 : SmartPtr<DoFIndexStorage> m_spDoFIndexStrgForLevelWithGhost;
273 : /// \}
274 : };
275 :
276 : /// base class for approximation spaces without type of algebra or dof distribution
277 : template <typename TDomain>
278 0 : class ApproximationSpace : public IApproximationSpace
279 : {
280 : public:
281 : /// Domain type
282 : typedef TDomain domain_type;
283 :
284 : /// World Dimension
285 : static const int dim = domain_type::dim;
286 :
287 : /// Subset Handler type
288 : typedef typename domain_type::subset_handler_type subset_handler_type;
289 :
290 : public:
291 : /// constructor
292 : ApproximationSpace(SmartPtr<TDomain> domain);
293 :
294 : /// constructor passing requested algebra type
295 : ApproximationSpace(SmartPtr<TDomain> domain, const AlgebraType& algebraType);
296 :
297 : /// Return the domain
298 0 : ConstSmartPtr<TDomain> domain() const {return m_spDomain;}
299 :
300 : /// Return the domain
301 0 : SmartPtr<TDomain> domain() {return m_spDomain;}
302 :
303 0 : int get_dim() const { return dim; }
304 :
305 : protected:
306 : /// Domain, where solution lives
307 : SmartPtr<TDomain> m_spDomain;
308 : };
309 :
310 :
311 : } // end namespace ug
312 :
313 : #endif /* __H__UG__LIB_DISC__FUNCTION_SPACE__APPROXIMATION_SPACE__ */
|