Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
approximation_space.h
Go to the documentation of this file.
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
41
42namespace ug{
43
45
62{
63 public:
66
69
70 public:
74
78 const AlgebraType& algebraType);
79
80 protected:
83 SmartPtr<grid_type> spMG, const AlgebraType& algebraType);
84
85 public:
88
90 void clear() {m_spDoFDistributionInfo->clear();}
91
92 public:
94
98 void add(const std::vector<std::string>& vName, LFEID id){
99 m_spDoFDistributionInfo->add(vName, id);
100 }
101
103 void add(const std::vector<std::string>& vName, const char* type, int order);
104
106 void add(const std::vector<std::string>& vName, const char* type);
107
109 void add(const char* name, const char* type, int order);
110
112 void add(const char* name, const char* type);
113
114 public:
116
121 void add(const std::vector<std::string>& vName, LFEID id,
122 const std::vector<std::string>& vSubset){
123 m_spDoFDistributionInfo->add(vName, id, vSubset);
124 }
125
127 void add(const std::vector<std::string>& vName, const char* type, int order,
128 const std::vector<std::string>& vSubsets);
129
131
137 void add(const char* name, const char* type, int order, const char* subsets);
138
140 void add(const std::vector<std::string>& vName, const char* type,
141 const std::vector<std::string>& vSubsets);
142
144
149 void add(const char* name, const char* type, const char* subsets);
150
151 public:
154
156 size_t num_levels() const {return m_spMGSH->num_levels();}
157
160
162 bool grouped() const {return m_bGrouped;}
163
165 bool might_contain_ghosts(int lvl) const;
166
168 bool might_contain_ghosts() const;
169
172 SmartPtr<DoFDistribution> dof_distribution(const GridLevel& gl, bool bCreate = true);
173 SmartPtr<DoFDistribution> dd(const GridLevel& gl, bool bCreate = true);
175
178 ConstSmartPtr<DoFDistribution> dof_distribution(const GridLevel& gl, bool bCreate = true) const;
179 ConstSmartPtr<DoFDistribution> dd(const GridLevel& gl, bool bCreate = true) const;
181
183 std::vector<SmartPtr<DoFDistribution> > dof_distributions() const;
184
190
192 void print_statistic(std::string flags) const;
193
195 void print_statistic() const;
196
198 void print_layout_statistic() const;
199
200
202 void init_levels();
203
205 void init_surfaces();
206
208 void init_top_surface();
209
211 const RevisionCounter& revision() const {return m_RevCnt;}
212
213 protected:
215 void create_dof_distribution(const GridLevel& gl);
216
219
222
223 protected:
225 void reinit();
226
231
234
239
242
243 protected:
246
249
252
255
258
261
264
265 protected:
267 std::vector<SmartPtr<DoFDistribution> > m_vDD;
268
274};
275
277template <typename TDomain>
279{
280 public:
282 typedef TDomain domain_type;
283
285 static const int dim = domain_type::dim;
286
288 typedef typename domain_type::subset_handler_type subset_handler_type;
289
290 public:
293
296
299
302
303 int get_dim() const { return dim; }
304
305 protected:
308};
309
310
311} // end namespace ug
312
313#endif /* __H__UG__LIB_DISC__FUNCTION_SPACE__APPROXIMATION_SPACE__ */
location name
Definition checkpoint_util.lua:128
Definition smart_pointer.h:296
Definition smart_pointer.h:108
class describing the type of an algebra
Definition algebra_type.h:49
base class for approximation spaces without type of algebra or dof distribution
Definition approximation_space.h:279
domain_type::subset_handler_type subset_handler_type
Subset Handler type.
Definition approximation_space.h:288
SmartPtr< TDomain > m_spDomain
Domain, where solution lives.
Definition approximation_space.h:307
int get_dim() const
Definition approximation_space.h:303
ConstSmartPtr< TDomain > domain() const
Return the domain.
Definition approximation_space.h:298
SmartPtr< TDomain > domain()
Return the domain.
Definition approximation_space.h:301
static const int dim
World Dimension.
Definition approximation_space.h:285
TDomain domain_type
Domain type.
Definition approximation_space.h:282
Definition dof_distribution_info.h:152
Definition grid_level.h:42
A message sent along with "GridRefinement" messages.
Definition lib_grid_messages.h:91
Definition lib_grid_messages.h:166
describes the ansatz spaces on a domain
Definition approximation_space.h:62
void init_surfaces()
initializes all surface dof distributions
Definition approximation_space.cpp:295
bool m_bAdaptionIsActive
Definition approximation_space.h:230
ConstSmartPtr< MGSubsetHandler > subset_handler() const
get underlying subset handler
Definition approximation_space.h:153
bool grouped() const
returns if dofs are grouped
Definition approximation_space.h:162
void register_at_adaption_msg_hub()
registers at message hub for grid adaption
Definition approximation_space.cpp:421
SmartPtr< SurfaceView > m_spSurfaceView
Surface View.
Definition approximation_space.h:251
~IApproximationSpace()
Destructor.
Definition approximation_space.cpp:111
void grid_changed_callback(const GridMessage_Adaption &msg)
Definition approximation_space.cpp:435
SmartPtr< DoFIndexStorage > m_spDoFIndexStrgForLevelNoGhost
Definition approximation_space.h:271
size_t num_levels() const
returns the number of level
Definition approximation_space.h:156
void surface_view_required()
creates surface SurfaceView if needed
Definition approximation_space.cpp:351
void print_statistic() const
prints statistic about DoF Distribution
Definition approximation_space.cpp:753
void print_layout_statistic() const
prints statistic on layouts
Definition approximation_space.cpp:837
AlgebraType m_algebraType
suitable algebra type for the index distribution pattern
Definition approximation_space.h:254
MessageHub::SPCallbackId m_spGridAdaptionCallbackID
message hub id
Definition approximation_space.h:228
bool might_contain_ghosts() const
returns if ghosts might be present on any level
Definition approximation_space.cpp:154
ConstSmartPtr< SurfaceView > surface_view() const
returns the approximation space
Definition approximation_space.h:159
std::vector< SmartPtr< DoFDistribution > > m_vDD
MG Level DoF Distribution.
Definition approximation_space.h:267
void grid_distribution_callback(const GridMessage_Distribution &msg)
called during parallel redistribution
Definition approximation_space.cpp:478
const RevisionCounter & revision() const
returns the current revision
Definition approximation_space.h:211
RevisionCounter m_RevCnt
revision counter
Definition approximation_space.h:263
SmartPtr< MultiGrid > m_spMG
multigrid, where elements are stored
Definition approximation_space.h:245
MGSubsetHandler subset_handler_type
Type of Subset Handler.
Definition approximation_space.h:65
SmartPtr< MGSubsetHandler > m_spMGSH
subsethandler, where elements are stored
Definition approximation_space.h:248
void init_levels()
initializes all level dof distributions
Definition approximation_space.cpp:286
ConstSmartPtr< DoFDistributionInfo > ddinfo() const
Definition approximation_space.h:187
void add(const std::vector< std::string > &vName, LFEID id)
add single solutions of LocalShapeFunctionSetID to the entire domain
Definition approximation_space.h:98
std::vector< SmartPtr< DoFDistribution > > dof_distributions() const
returns all currently created dof distributions
Definition approximation_space.cpp:280
SmartPtr< DoFIndexStorage > m_spDoFIndexStrgForLevelWithGhost
Definition approximation_space.h:272
void add(const std::vector< std::string > &vName, LFEID id, const std::vector< std::string > &vSubset)
add single solutions of LocalShapeFunctionSetID to selected subsets
Definition approximation_space.h:121
SmartPtr< DoFDistributionInfo > m_spDoFDistributionInfo
DofDistributionInfo.
Definition approximation_space.h:260
void clear()
clears functions
Definition approximation_space.h:90
bool m_bGrouped
flag if DoFs should be grouped
Definition approximation_space.h:257
void create_dof_distribution(const GridLevel &gl)
creates a dof distribution
Definition approximation_space.cpp:318
ConstSmartPtr< DoFDistributionInfo > dof_distribution_info() const
Definition approximation_space.h:188
void init_top_surface()
initializes all top surface dof distributions
Definition approximation_space.cpp:304
MultiGrid grid_type
Type of Subset Handler.
Definition approximation_space.h:68
SmartPtr< DoFDistribution > dd(const GridLevel &gl, bool bCreate=true)
Definition approximation_space.cpp:262
void reinit()
reinits all data after grid adaption
Definition approximation_space.cpp:405
void dof_distribution_info_required()
create dof distribution info
Definition approximation_space.cpp:358
MessageHub::SPCallbackId m_spGridDistributionCallbackID
Definition approximation_space.h:229
SmartPtr< DoFDistribution > dof_distribution(const GridLevel &gl, bool bCreate=true)
Definition approximation_space.cpp:247
Identifier for Local Finite Elements.
Definition local_finite_element_id.h:98
Definition multi_grid.h:72
Handles subsets on a per level basis.
Definition subset_handler_multi_grid.h:60
Class used to identify a state of adaption of a grid, approx-space, ...
Definition revision_counter.h:56
virtual void init()
the ug namespace