Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
err_est_data.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-2015: G-CSC, Goethe University Frankfurt
3 * Authors: Dmitry Logashenko, 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 * Data shared by element discretizations for a-posteriori error estimation
35 */
36
37#ifndef __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__ERR_EST_DATA__
38#define __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__ERR_EST_DATA__
39
40// extern headers
41#include <vector>
42#include <string>
43#include <limits>
44
45// intern headers
49
50#ifdef UG_PARALLEL
53#endif
54
55#include <boost/mpl/for_each.hpp>
56
57
58namespace ug{
59
61
69template <typename TDomain>
71{
72 public:
74 static const int dim = TDomain::dim;
75
77 IErrEstData () : m_consider(true), m_scale(1.0) {};
78
80 virtual ~IErrEstData () {};
81
83 virtual void alloc_err_est_data (ConstSmartPtr<SurfaceView> spSV, const GridLevel& gl) = 0;
84
86 virtual void summarize_err_est_data (SmartPtr<TDomain> spDomain) = 0;
87
89 virtual number get_elem_error_indicator(GridObject* elem, const MathVector<dim> vCornerCoords[]) = 0;
90
92 virtual void release_err_est_data () = 0;
93
95 bool consider_me() const {return m_consider;};
96
98
106 void set_consider_me(bool b) {m_consider = b;};
107
109 void set_scaling_factor(number scale) {m_scale = scale;};
110
112
122
123 private:
126};
127
129
136template <typename TDomain>
137class SideFluxErrEstData : public IErrEstData<TDomain>
138{
139public:
141 typedef TDomain domain_type;
142
144 static const int dim = TDomain::dim;
145
148
149public:
152
155
156// Functions to access data
157
159 number& operator()
160 (
161 side_type* pSide
162 )
163 {
164 return m_aaFluxJump[pSide];
165 };
166
169
170// Interface virtual functions inherited from IErrEstData
171
173 virtual void alloc_err_est_data (ConstSmartPtr<SurfaceView> spSV, const GridLevel& gl);
174
176 virtual void summarize_err_est_data (SmartPtr<TDomain> spDomain);
177
179 virtual number get_elem_error_indicator(GridObject* elem, const MathVector<dim> vCornerCoords[]) {return 0;};
180
182 virtual void release_err_est_data ();
183
184private:
187
190
193
196};
197
198
199
201
219template <typename TDomain>
220class SideAndElemErrEstData : public IErrEstData<TDomain>
221{
222public:
224 typedef TDomain domain_type;
225
227 static const int dim = TDomain::dim;
228
232
235
238
240 static const int MAX_NUM_SIDES = 8;
241
242public:
244 SideAndElemErrEstData(size_t _sideOrder, size_t _elemOrder, const char* subsets);
245 SideAndElemErrEstData(size_t _sideOrder, size_t _elemOrder,
246 std::vector<std::string> subsets = std::vector<std::string>(0));
247
250
251 // Functions to access data
252
254 size_t side_order() const {return sideOrder;}
255
257 size_t elem_order() const {return elemOrder;}
258
260 number& operator()
261 (
262 side_type* pSide,
263 size_t ip
264 );
265
267 number& operator()
268 (
269 elem_type* pElem,
270 size_t ip
271 );
272
274 template <int refDim>
276
278 template <int refDim>
280
283
286
289
291 size_t num_side_ips(const side_type* pSide);
292
294 size_t num_side_ips(const ReferenceObjectID roid);
295
297 size_t first_side_ips(const ReferenceObjectID roid, const size_t side);
298
300 size_t num_all_side_ips(const ReferenceObjectID roid);
301
303 size_t num_elem_ips(const ReferenceObjectID roid);
304
306 size_t side_ip_index(const ReferenceObjectID roid, const size_t side, const size_t ip);
307
310
311 // virtual functions inherited from IErrEstData
313 virtual void alloc_err_est_data (ConstSmartPtr<SurfaceView> spSV, const GridLevel& gl);
314
316 virtual void summarize_err_est_data (SmartPtr<TDomain> spDomain);
317
319 virtual number get_elem_error_indicator(GridObject* elem, const MathVector<dim> vCornerCoords[]);
320
322 virtual void release_err_est_data ();
323
326
327protected:
329 void init_quadrature();
330
332 template<int refDim>
334 {
335 GetQuadRules(QuadratureRule<refDim>** ppQuadRule, size_t quadOrder) :
336 m_ppQuadRule(ppQuadRule), m_quadOrder(quadOrder) {}
339 template< typename TElem > void operator()(TElem&)
340 {
342 m_ppQuadRule[roid] =
344 }
345 };
346
347private:
349 size_t sideOrder;
350 size_t elemOrder;
351
353 std::vector<std::string> m_vSs;
355
359
361 std::vector<MathVector<TDomain::dim> > m_SideIPcoords[NUM_REFERENCE_OBJECTS];
362
364 std::vector<MathVector<TDomain::dim> > m_sideGlobalIPcoords;
365 std::vector<MathVector<TDomain::dim> > m_singleSideGlobalIPcoords; // not the most elegant solution...
366 std::vector<MathVector<TDomain::dim> > m_elemGlobalIPcoords;
367
370
373
376
379
382
385
388
392};
393
394
396
421template <typename TDomain, typename TErrEstData>
422class MultipleErrEstData : public IErrEstData<TDomain>
423{
424 public:
426 static const int dim = TDomain::dim;
427
430 : IErrEstData<TDomain>(), m_spApprox(approx),
431 m_fctGrp(m_spApprox->dof_distribution_info())
432 {
433 this->set_consider_me(false);
434 }
435
438
440 virtual void add(SmartPtr<TErrEstData> spEed, const char* fct)
441 {
442 // check that fct is not already contained in fctGrp
443 size_t uid = m_spApprox->fct_id_by_name(fct);
444 if (m_fctGrp.contains(uid))
445 {
446 UG_THROW("Error estimator for function '" << fct << "' can not be added\n"
447 "as another error estimator object for the same function is already\n"
448 "contained here.")
449 }
450
451 // add to function group
452 try
453 {
454 m_fctGrp.add(fct);
455 }
456 UG_CATCH_THROW("Error estimator data object for function '"
457 << fct << "' could not be added.");
458
459 m_vEed.push_back(spEed.get());
460 }
461
463 size_t num() const {return m_vEed.size();};
464
466 TErrEstData* get(size_t uid)
467 {
468 if (!m_fctGrp.contains(uid))
469 {
470 std::string name = m_spApprox->name(uid);
471 UG_THROW("Trying to access error estimator data object "
472 "for unique function index " << uid << "(aka '" << name << "')\n"
473 "which is not present in this collection.")
474 }
475
476 return m_vEed[m_fctGrp.local_index(uid)];
477 }
478
479 // inherited from IErrEstData
481 virtual void alloc_err_est_data(ConstSmartPtr<SurfaceView> spSV, const GridLevel& gl);
482
484 virtual void summarize_err_est_data(SmartPtr<TDomain> spDomain);
485
487 virtual number get_elem_error_indicator(GridObject* elem, const MathVector<dim> vCornerCoords[]);
488
490 virtual void release_err_est_data();
491
492 protected:
493 std::vector<TErrEstData*> m_vEed;
494
497
500};
501
502
503template <typename TDomain>
505 : public MultipleErrEstData<TDomain, SideAndElemErrEstData<TDomain> >
506{
507 public:
509 static const int dim = TDomain::dim;
510
514
519
522
525 virtual void add(SmartPtr<SideAndElemErrEstData<TDomain> > spEed, const char* fct);
526
528 bool equal_side_order() const {return m_bEqSideOrder;}
529
531 bool equal_elem_order() const {return m_bEqElemOrder;}
532
533 protected:
536 void check_equal_order();
537
541
545
546 private:
549
550};
551
552} // end of namespace ug
553
554#include "err_est_data_impl.h"
555
556#endif /* __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__ERR_EST_DATA__ */
557
558/* End of File */
location name
Definition checkpoint_util.lua:128
Definition smart_pointer.h:296
Definition smart_pointer.h:108
T * get()
returns encapsulated pointer
Definition smart_pointer.h:197
base class for approximation spaces without type of algebra or dof distribution
Definition approximation_space.h:279
Definition function_group.h:52
size_t local_index(size_t uniqueID) const
return index in Function group for a function
Definition function_group.cpp:302
bool contains(size_t uniqueID) const
returns true if unique id is contained in this group
Definition function_group.cpp:279
void add(size_t fct)
adds a function by id to this group
Definition function_group.cpp:79
the generic attachment-accessor for access to grids attachment pipes.
Definition grid.h:182
Definition grid_level.h:42
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
Base class for error estimator data.
Definition err_est_data.h:71
void set_scaling_factor(number scale)
set scaling factor for final error calculation
Definition err_est_data.h:109
virtual void release_err_est_data()=0
virtual function to release data structures for the error estimator
bool m_consider
Definition err_est_data.h:124
virtual void summarize_err_est_data(SmartPtr< TDomain > spDomain)=0
virtual function called after the computation of the error estimator data in all the elements
static const int dim
world dimension
Definition err_est_data.h:74
virtual ~IErrEstData()
virtual class destructor
Definition err_est_data.h:80
bool consider_me() const
virtual function granting get access to the m_consider member
Definition err_est_data.h:95
IErrEstData()
class constructor
Definition err_est_data.h:77
void set_consider_me(bool b)
whether or not this instance is to be considered by domainDisc
Definition err_est_data.h:106
virtual void alloc_err_est_data(ConstSmartPtr< SurfaceView > spSV, const GridLevel &gl)=0
virtual function to allocate data structures for the error estimator
number scaling_factor()
get scaling factor
Definition err_est_data.h:121
number m_scale
Definition err_est_data.h:125
virtual number get_elem_error_indicator(GridObject *elem, const MathVector< dim > vCornerCoords[])=0
calculate L2 integrals
a mathematical Vector with N entries.
Definition math_vector.h:97
Error estimator data class for discretizations with more than one unknown.
Definition err_est_data.h:423
virtual number get_elem_error_indicator(GridObject *elem, const MathVector< dim > vCornerCoords[])
calculate L2 integrals
Definition err_est_data_impl.h:949
virtual void alloc_err_est_data(ConstSmartPtr< SurfaceView > spSV, const GridLevel &gl)
virtual function to allocate data structures for the error estimator
Definition err_est_data_impl.h:931
TErrEstData * get(size_t uid)
accessing the underlying error estimator data objects via function id
Definition err_est_data.h:466
virtual void release_err_est_data()
virtual function to release data structures for the error estimator
Definition err_est_data_impl.h:961
static const int dim
world dimension
Definition err_est_data.h:426
ConstSmartPtr< ApproximationSpace< TDomain > > m_spApprox
approx space
Definition err_est_data.h:496
virtual ~MultipleErrEstData()
virtual class destructor
Definition err_est_data.h:437
FunctionGroup m_fctGrp
function group (in order to map fcts to error estimator objects)
Definition err_est_data.h:499
virtual void summarize_err_est_data(SmartPtr< TDomain > spDomain)
virtual function called after the computation of the error estimator data in all the elements
Definition err_est_data_impl.h:940
size_t num() const
getting the number of underlying error estimator data objects
Definition err_est_data.h:463
virtual void add(SmartPtr< TErrEstData > spEed, const char *fct)
adding error estimator data objects
Definition err_est_data.h:440
std::vector< TErrEstData * > m_vEed
Definition err_est_data.h:493
MultipleErrEstData(ConstSmartPtr< ApproximationSpace< TDomain > > approx)
class constructor
Definition err_est_data.h:429
Definition err_est_data.h:506
void check_equal_order()
Definition err_est_data_impl.h:981
MultipleSideAndElemErrEstData(ConstSmartPtr< ApproximationSpace< TDomain > > approx)
constructor
Definition err_est_data.h:516
SideAndElemErrEstData< TDomain >::side_type side_type
type of the sides (face, edge) and the elems (volume, face)
Definition err_est_data.h:512
void check_equal_elem_order()
Definition err_est_data_impl.h:1007
void check_equal_side_order()
Definition err_est_data_impl.h:988
bool m_bEqElemOrder
Definition err_est_data.h:548
SideAndElemErrEstData< TDomain >::elem_type elem_type
Definition err_est_data.h:513
bool equal_elem_order() const
returns whether all underlying err ests have the same elem and side integration orders
Definition err_est_data.h:531
virtual void add(SmartPtr< SideAndElemErrEstData< TDomain > > spEed, const char *fct)
Definition err_est_data_impl.h:974
static const int dim
world dimension
Definition err_est_data.h:509
virtual ~MultipleSideAndElemErrEstData()
destructor
Definition err_est_data.h:521
bool m_bEqSideOrder
Definition err_est_data.h:547
bool equal_side_order() const
returns whether all underlying err ests have the same elem and side integration orders
Definition err_est_data.h:528
provides quadrature rule for a Reference Dimension
Definition quadrature.h:70
static const QuadratureRule< TDim > & get(size_t order, QuadType type=BEST)
gets quadrature rule of requested order
Definition quadrature_provider_impl.h:40
Error estimator data class storing a number vector per side and per element.
Definition err_est_data.h:221
static const int MAX_NUM_SIDES
maximal number of sides of any element
Definition err_est_data.h:240
void set_type(int type)
select L2/H1 Estimator
Definition err_est_data.h:325
QuadratureRule< dim-1 > * quadRuleSide[NUM_REFERENCE_OBJECTS]
storage for integration rules
Definition err_est_data.h:357
MultiGrid::AttachmentAccessor< side_type, attachment_type > m_aaSide
vector of side attachment accessors
Definition err_est_data.h:378
int m_type
Definition err_est_data.h:391
virtual void summarize_err_est_data(SmartPtr< TDomain > spDomain)
virtual function called after the computation of the error estimator data in all the elements
Definition err_est_data_impl.h:636
size_t side_order() const
getting the side integration order
Definition err_est_data.h:254
ConstSmartPtr< SurfaceView > & surface_view()
get the surface view
Definition err_est_data.h:309
attachment_type m_aElem
vector of attachments for elems
Definition err_est_data.h:375
SideAndElemErrEstData< TDomain > this_type
this class
Definition err_est_data.h:237
TDomain domain_type
domain type
Definition err_est_data.h:224
std::vector< MathVector< TDomain::dim > > m_singleSideGlobalIPcoords
Definition err_est_data.h:365
std::vector< MathVector< TDomain::dim > > m_elemGlobalIPcoords
Definition err_est_data.h:366
QuadratureRule< dim > * quadRuleElem[NUM_REFERENCE_OBJECTS]
Definition err_est_data.h:358
size_t m_sideIPsStartIndex[NUM_REFERENCE_OBJECTS][MAX_NUM_SIDES]
the first index for IPs of a specific side in the sideIP series for a roid
Definition err_est_data.h:369
size_t sideOrder
order of side and elem function approximations for integrating
Definition err_est_data.h:349
virtual void release_err_est_data()
virtual function to release data structures of the error estimator
Definition err_est_data_impl.h:912
size_t num_all_side_ips(const ReferenceObjectID roid)
get number of side IPs
Definition err_est_data_impl.h:442
void init_quadrature()
initialization of quadrature (to be called during construction)
Definition err_est_data_impl.h:176
std::vector< std::string > m_vSs
the subsets this error estimator will produce values for
Definition err_est_data.h:353
size_t elemOrder
Definition err_est_data.h:350
MultiGrid::AttachmentAccessor< elem_type, attachment_type > m_aaElem
vector of elem attachment accessors
Definition err_est_data.h:381
domain_traits< dim >::side_type side_type
type of the sides (face, edge) and the elems (volume, face)
Definition err_est_data.h:230
SubsetGroup m_ssg
Definition err_est_data.h:354
MathVector< TDomain::dim > * all_side_global_ips(GridObject *elem, const MathVector< dim > vCornerCoords[])
get all global side integration points
Definition err_est_data_impl.h:333
const MathVector< refDim > * elem_local_ips(const ReferenceObjectID roid)
get the local elem integration points for a specific roid
Definition err_est_data_impl.h:316
virtual void alloc_err_est_data(ConstSmartPtr< SurfaceView > spSV, const GridLevel &gl)
virtual function to allocate data structures for the error estimator
Definition err_est_data_impl.h:471
static const int dim
world dimension
Definition err_est_data.h:227
MathVector< TDomain::dim > * elem_global_ips(GridObject *elem, const MathVector< dim > vCornerCoords[])
get the global elem integration points for a specific roid
Definition err_est_data_impl.h:391
Attachment< std::vector< number > > attachment_type
attachment type
Definition err_est_data.h:234
std::vector< MathVector< TDomain::dim > > m_SideIPcoords[NUM_REFERENCE_OBJECTS]
extra storage for local side IPs (elem IPs are contained in elem quad rules)
Definition err_est_data.h:361
type
set type
Definition err_est_data.h:390
@ L2_ERROR_TYPE
Definition err_est_data.h:390
@ H1_ERROR_TYPE
Definition err_est_data.h:390
attachment_type m_aSide
vector of attachments for sides
Definition err_est_data.h:372
const MathVector< refDim > * side_local_ips(const ReferenceObjectID roid)
get the local side integration points for a specific roid
Definition err_est_data_impl.h:282
size_t side_ip_index(const ReferenceObjectID roid, const size_t side, const size_t ip)
get index of specific side IP in sideIP array returned by side_local_ips
Definition err_est_data_impl.h:458
std::vector< MathVector< TDomain::dim > > m_sideGlobalIPcoords
storage for global elem and side IPs
Definition err_est_data.h:364
domain_traits< dim >::element_type elem_type
Definition err_est_data.h:231
size_t elem_order() const
getting the elem integration order
Definition err_est_data.h:257
size_t first_side_ips(const ReferenceObjectID roid, const size_t side)
get number of first IP belonging to a specific side
Definition err_est_data_impl.h:435
virtual ~SideAndElemErrEstData()
virtual class destructor
Definition err_est_data.h:249
size_t num_side_ips(const side_type *pSide)
get number of side IPs of a specific side
Definition err_est_data_impl.h:419
MathVector< TDomain::dim > * side_global_ips(GridObject *elem, const MathVector< dim > vCornerCoords[])
get the global side integration points for a specific side roid
Definition err_est_data_impl.h:362
virtual number get_elem_error_indicator(GridObject *elem, const MathVector< dim > vCornerCoords[])
calculate L2 integrals
Definition err_est_data_impl.h:813
size_t num_elem_ips(const ReferenceObjectID roid)
get number of elem IPs
Definition err_est_data_impl.h:448
ConstSmartPtr< SurfaceView > m_spSV
Grid for the attachment.
Definition err_est_data.h:384
GridLevel m_errEstGL
Finest grid level.
Definition err_est_data.h:387
Error estimator data class storing one scalar number per side.
Definition err_est_data.h:138
ANumber m_aFluxJumpOverSide
Flux jumps for the error estimator.
Definition err_est_data.h:186
ConstSmartPtr< SurfaceView > & surface_view()
get the surface view
Definition err_est_data.h:168
ConstSmartPtr< SurfaceView > m_spSV
Grid for the attachment.
Definition err_est_data.h:192
virtual void release_err_est_data()
virtual function to release data structures of the error estimator
Definition err_est_data_impl.h:92
GridLevel m_errEstGL
Finest grid level.
Definition err_est_data.h:195
virtual void alloc_err_est_data(ConstSmartPtr< SurfaceView > spSV, const GridLevel &gl)
virtual function to allocate data structures for the error estimator
Definition err_est_data_impl.h:40
static const int dim
world dimension
Definition err_est_data.h:144
virtual void summarize_err_est_data(SmartPtr< TDomain > spDomain)
virtual function called after the computation of the error estimator data in all the elements
Definition err_est_data_impl.h:63
MultiGrid::AttachmentAccessor< side_type, ANumber > m_aaFluxJump
Attachment accessor.
Definition err_est_data.h:189
virtual number get_elem_error_indicator(GridObject *elem, const MathVector< dim > vCornerCoords[])
calculate L2 integrals
Definition err_est_data.h:179
virtual ~SideFluxErrEstData()
virtual class destructor
Definition err_est_data.h:154
SideFluxErrEstData()
constructor
Definition err_est_data.h:151
TDomain domain_type
domain type
Definition err_est_data.h:141
domain_traits< dim >::side_type side_type
type of the sides (face, edge) and the elems (volume, face)
Definition err_est_data.h:147
Group of subsets.
Definition subset_group.h:51
Definition grid_base_object_traits.h:68
#define UG_CATCH_THROW(msg)
Definition error.h:64
#define UG_THROW(msg)
Definition error.h:57
double number
Definition types.h:124
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
@ NUM_REFERENCE_OBJECTS
Definition grid_base_objects.h:85
helper struct for getting quadrature rules by use of mpl::lists
Definition err_est_data.h:334
GetQuadRules(QuadratureRule< refDim > **ppQuadRule, size_t quadOrder)
Definition err_est_data.h:335
size_t m_quadOrder
Definition err_est_data.h:338
QuadratureRule< refDim > ** m_ppQuadRule
Definition err_est_data.h:337
void operator()(TElem &)
Definition err_est_data.h:339
Definition domain_traits.h:53