ug4
Loading...
Searching...
No Matches
inner_boundary.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2011-2015: 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 * Finite Volume Element Discretization for an inner BndCond that depends on the unknowns (on the bnd)
35 *
36 * This class implements the IElemDisc interface to provide element local
37 * assemblings for the unknown-dependent Neumann-flux over an inner boundary.
38 * The equation of this flux and its derivative should be given
39 * in a concretization of this class.
40 */
41
42#ifndef __H__UG__LIB_DISC__SPACIAL_DISCRETIZATION__ELEM_DISC__NEUMANN_BOUNDARY__FV1__INNER_BOUNDARY__
43#define __H__UG__LIB_DISC__SPACIAL_DISCRETIZATION__ELEM_DISC__NEUMANN_BOUNDARY__FV1__INNER_BOUNDARY__
44
45#include <boost/mpl/empty.hpp>
46#include <boost/mpl/if.hpp>
47#include <boost/mpl/front.hpp>
48#include <boost/mpl/pop_front.hpp>
49
50#include <vector>
51#include <string>
52#include <utility> // std::pair
53
54// other ug4 modules
55#include "common/common.h"
56
57// library intern headers
63
64
65
66namespace ug
67{
68
70{
71 public:
73 static const size_t _IGNORE_ = -1;
74};
75
76
79{
80 // vector of fluxFctValues
81 std::vector<number> flux;
82 std::vector<size_t> from;
83 std::vector<size_t> to;
84};
85
89{
90 // vector of fluxFctDerivValues
91 // rows: fluxIndex, cols: dependency;
92 // in the pair: first: with respect to which local function index, second: value
93 std::vector<std::vector<std::pair<size_t,number> > > fluxDeriv;
94 std::vector<size_t> from;
95 std::vector<size_t> to;
96};
97
98
99template <typename TImpl, typename TDomain>
101: public IElemDisc<TDomain>
102{
103 public:
106
107 public:
109 typedef TDomain domain_type;
110
111 private:
114
117
118 public:
120 static const int dim = base_type::dim;
121
123 typedef typename base_type::position_type position_type;
124
127
128
129 public:
130
132 FV1InnerBoundaryElemDisc(const char* functions = "", const char* subsets = "")
133 : IElemDisc<domain_type>(functions, subsets), m_bNonRegularGrid(false), m_bCurrElemIsHSlave(false), m_si(0)
134 {
136 }
137
139 FV1InnerBoundaryElemDisc(const std::vector<std::string>& functions, const std::vector<std::string>& subsets)
140 : IElemDisc<domain_type>(functions, subsets), m_bNonRegularGrid(false), m_bCurrElemIsHSlave(false), m_si(0)
141 {
143 }
144
147
149 //void set_fluxFunction(UserData<number, dim>& fluxFct) {m_fluxFct.set_data(fluxFct);}
150
152 void set_flux_scale(number scale);
154#ifdef UG_FOR_LUA
155 void set_flux_scale(const char* luaScaleFctName);
156#endif
157
158 public: // inherited from IElemDisc
160 virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid) override;
161
163 virtual bool use_hanging() const override;
164
165 protected:
166 void previous_solution_required(bool b);
167
168 template <typename TAlgebra>
169 void prep_timestep(number future_time, number time, VectorProxyBase* upb);
170 private:
172
177 template<typename TElem, typename TFVGeom>
178 void prep_elem_loop(const ReferenceObjectID roid, const int si);
179
181
186 template<typename TElem, typename TFVGeom>
187 void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
188
190 template<typename TElem, typename TFVGeom>
191 void fsh_elem_loop();
192
194 template<typename TElem, typename TFVGeom>
195 void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
196
198 template<typename TElem, typename TFVGeom>
199 void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
200
202 template<typename TElem, typename TFVGeom>
203 void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
204
206 template<typename TElem, typename TFVGeom>
207 void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
208
210 template<typename TElem, typename TFVGeom>
211 void add_rhs_elem(LocalVector& rhs, GridObject* elem, const MathVector<dim> vCornerCoords[]);
212
214 template <typename TElem, typename TFVGeom>
215 void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si);
216
218 template <typename TElem, typename TFVGeom>
219 void prep_err_est_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
220
222 template <typename TElem, typename TFVGeom>
223 void compute_err_est_A_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
224
226 template <typename TElem, typename TFVGeom>
228
229 private:
230 DataImport<number, dim> m_fluxScale; // data import for scaling of fluxes
231
232 private:
234
235 template <typename TElem, typename TFVGeom>
236 void register_fv1_func();
237
239 {
240 // friend class this_type;
243 template< typename TElem > void operator()(TElem&)
244 {
245 if (m_pThis->m_bNonRegularGrid)
246 m_pThis->register_fv1_func<TElem, HFV1ManifoldGeometry<TElem, dim> >();
247 else
248 m_pThis->register_fv1_func<TElem, FV1ManifoldGeometry<TElem, dim> >();
249 }
250 };
251
252 template <typename List>
254 {
256 {
257 static const bool isEmpty = boost::mpl::empty<List>::value;
258 (typename boost::mpl::if_c<isEmpty, RegEnd, RegNext>::type(p));
259 }
260 struct RegEnd
261 {
263 };
264 struct RegNext
265 {
267 {
268 typedef typename boost::mpl::front<List>::type AlgebraType;
269 typedef typename boost::mpl::pop_front<List>::type NextList;
270
272 p->set_prep_timestep_fct(aid, &this_type::template prep_timestep<AlgebraType>);
274 }
275 };
276 };
277
278
281 {
282 public:
283 void resize(size_t nSip, size_t _nSh)
284 {
285 nSh = _nSh;
286 sideVals.resize(nSip);
287 for (size_t i = 0; i < nSip; i++) sideVals[i].resize(nSh);
288 }
289 number& shapeAtSideIP(size_t sh, size_t ip)
290 {
291 UG_ASSERT(ip < sideVals.size(), "Requested data for IP " << ip << ", but only " << sideVals.size() << " IPs present.");
292 UG_ASSERT(sh < sideVals[ip].size(), "Requested data for shape fct " << sh << ", but only " << sideVals[ip].size() << " shape fcts present.");
293 return sideVals.at(ip).at(sh);
294 }
295 number* shapesAtSideIP(size_t ip) {return &sideVals[ip][0];}
296 size_t num_sh() {return nSh;}
297 private:
298 size_t nSh;
299 std::vector<std::vector<number> > sideVals;
301
302 private:
305
309
310 // temporary vectors (avoid re-allocating heap memory over and over again)
311 std::vector<LocalVector::value_type> m_uAtCorner;
312 std::vector<LocalVector::value_type> m_uOldAtCorner;
313
314 int m_si;
315};
316
317}
318
319#endif /*__H__UG__LIB_DISC__SPACIAL_DISCRETIZATION__ELEM_DISC__NEUMANN_BOUNDARY__FV1__INNER_BOUNDARY__*/
parameterString p
Definition smart_pointer.h:107
class describing the type of an algebra
Definition algebra_type.h:49
Type based UserData.
Definition user_data.h:502
Data import.
Definition data_import.h:180
Definition inner_boundary.h:102
virtual ~FV1InnerBoundaryElemDisc()
Destructor.
Definition inner_boundary.h:146
FV1InnerBoundaryElemDisc(const std::vector< std::string > &functions, const std::vector< std::string > &subsets)
Constructor with functions.
Definition inner_boundary.h:139
MultipleSideAndElemErrEstData< domain_type > err_est_type
error estimator type
Definition inner_boundary.h:126
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition inner_boundary_impl.h:124
int m_si
Definition inner_boundary.h:314
void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements of one type for the computation of the error estimator
Definition inner_boundary_impl.h:367
TDomain domain_type
Domain type.
Definition inner_boundary.h:109
void register_all_fv1_funcs()
Definition inner_boundary_impl.h:581
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid) override
type of trial space for each function used
Definition inner_boundary_impl.h:74
DataImport< number, dim > m_fluxScale
Definition inner_boundary.h:230
void register_fv1_func()
Definition inner_boundary_impl.h:594
void prep_timestep(number future_time, number time, VectorProxyBase *upb)
prepare the time step
Definition inner_boundary_impl.h:104
SmartPtr< VectorProxyBase > m_spOldSolutionProxy
Definition inner_boundary.h:308
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the mass part of the local defect
Definition inner_boundary_impl.h:350
FV1InnerBoundaryElemDisc< TImpl, TDomain > this_type
own type
Definition inner_boundary.h:116
LocalVector m_locUOld
Definition inner_boundary.h:307
void add_def_A_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the stiffness part of the local defect
Definition inner_boundary_impl.h:281
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local mass matrix using a finite volume scheme
Definition inner_boundary_impl.h:274
InnerBoundaryFluxCond FluxCond
Definition inner_boundary.h:104
void previous_solution_required(bool b)
Definition inner_boundary_impl.h:95
virtual bool use_hanging() const override
returns if hanging nodes are used
Definition inner_boundary_impl.h:88
static const int dim
World dimension.
Definition inner_boundary.h:120
void add_rhs_elem(LocalVector &rhs, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition inner_boundary_impl.h:356
void fsh_elem_loop()
finishes the loop over all elements
Definition inner_boundary_impl.h:144
InnerBoundaryFluxDerivCond FluxDerivCond
Definition inner_boundary.h:105
IElemDisc< domain_type > base_type
Base class type.
Definition inner_boundary.h:113
std::vector< LocalVector::value_type > m_uOldAtCorner
Definition inner_boundary.h:312
void add_jac_A_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local stiffness matrix using a finite volume scheme
Definition inner_boundary_impl.h:205
void set_flux_scale(number scale)
Setting the flux function.
Definition inner_boundary_impl.h:47
struct ug::FV1InnerBoundaryElemDisc::ShapeValues m_shapeValues
FV1InnerBoundaryElemDisc(const char *functions="", const char *subsets="")
Constructor with c-strings.
Definition inner_boundary.h:132
base_type::position_type position_type
Position type.
Definition inner_boundary.h:123
void compute_err_est_A_elem(const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
computes the error estimator contribution for one element
Definition inner_boundary_impl.h:491
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition inner_boundary_impl.h:151
bool m_bNonRegularGrid
Definition inner_boundary.h:303
void prep_err_est_elem(const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
prepares the element for assembling the error estimator
Definition inner_boundary_impl.h:430
std::vector< LocalVector::value_type > m_uAtCorner
Definition inner_boundary.h:311
bool m_bCurrElemIsHSlave
Definition inner_boundary.h:304
bool m_bPrevSolRequired
Definition inner_boundary.h:306
void fsh_err_est_elem_loop()
postprocesses the loop over all elements of one type in the computation of the error estimator
Definition inner_boundary_impl.h:565
Definition fv1_geom.h:1194
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
Definition hfv1_geom.h:976
number time() const
returns currently set timepoint
Definition elem_disc_interface.h:677
Definition elem_disc_interface.h:800
Definition local_algebra.h:422
Definition local_algebra.h:198
a mathematical Vector with N entries.
Definition math_vector.h:97
Definition err_est_data.h:506
#define UG_ASSERT(expr, msg)
Definition assert.h:70
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
this_type * m_pThis
Definition inner_boundary.h:242
void operator()(TElem &)
Definition inner_boundary.h:243
RegisterAssemblingFuncs(this_type *pThis)
Definition inner_boundary.h:241
RegEnd(this_type *)
Definition inner_boundary.h:262
RegNext(this_type *p)
Definition inner_boundary.h:266
RegisterPrepTimestep(this_type *p)
Definition inner_boundary.h:255
struct holding values of shape functions in IPs
Definition inner_boundary.h:281
number & shapeAtSideIP(size_t sh, size_t ip)
Definition inner_boundary.h:289
number * shapesAtSideIP(size_t ip)
Definition inner_boundary.h:295
void resize(size_t nSip, size_t _nSh)
Definition inner_boundary.h:283
size_t num_sh()
Definition inner_boundary.h:296
std::vector< std::vector< number > > sideVals
Definition inner_boundary.h:299
size_t nSh
Definition inner_boundary.h:298
Definition inner_boundary.h:70
static const size_t _IGNORE_
index value for which a flux is ignored
Definition inner_boundary.h:73
struct that holds information about the flux densities and from where to where the flux occurs
Definition inner_boundary.h:79
std::vector< size_t > from
Definition inner_boundary.h:82
std::vector< number > flux
Definition inner_boundary.h:81
std::vector< size_t > to
Definition inner_boundary.h:83
Definition inner_boundary.h:89
std::vector< size_t > from
Definition inner_boundary.h:94
std::vector< std::vector< std::pair< size_t, number > > > fluxDeriv
Definition inner_boundary.h:93
std::vector< size_t > to
Definition inner_boundary.h:95
Proxy struct for generic passing of any vector type.
Definition elem_disc_interface.h:70
size_t id()
Definition util_algebra_dependent.h:130
static AlgebraTypeIDProvider & instance()
Definition util_algebra_dependent.h:107