ug4
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/function.hpp>
46 #include <vector>
47 #include <string>
48 #include <utility> // std::pair
49 
50 // other ug4 modules
51 #include "common/common.h"
52 
53 // library intern headers
59 
60 
61 
62 namespace ug
63 {
64 
66 {
67  public:
69  static const size_t _IGNORE_ = -1;
70 };
71 
72 
75 {
76  // vector of fluxFctValues
77  std::vector<number> flux;
78  std::vector<size_t> from;
79  std::vector<size_t> to;
80 };
81 
85 {
86  // vector of fluxFctDerivValues
87  // rows: fluxIndex, cols: dependency;
88  // in the pair: first: with respect to which local function index, second: value
89  std::vector<std::vector<std::pair<size_t,number> > > fluxDeriv;
90  std::vector<size_t> from;
91  std::vector<size_t> to;
92 };
93 
94 
95 template <typename TImpl, typename TDomain>
97 : public IElemDisc<TDomain>
98 {
99  public:
102 
103  public:
105  typedef TDomain domain_type;
106 
107  private:
110 
113 
114  public:
116  static const int dim = base_type::dim;
117 
120 
123 
124 
125  public:
126 
128  FV1InnerBoundaryElemDisc(const char* functions = "", const char* subsets = "")
129  : IElemDisc<domain_type>(functions, subsets), m_bNonRegularGrid(false), m_bCurrElemIsHSlave(false), m_si(0)
130  {
132  }
133 
135  FV1InnerBoundaryElemDisc(const std::vector<std::string>& functions, const std::vector<std::string>& subsets)
136  : IElemDisc<domain_type>(functions, subsets), m_bNonRegularGrid(false), m_bCurrElemIsHSlave(false), m_si(0)
137  {
139  }
140 
143 
145  //void set_fluxFunction(UserData<number, dim>& fluxFct) {m_fluxFct.set_data(fluxFct);}
146 
148  void set_flux_scale(number scale);
150 #ifdef UG_FOR_LUA
151  void set_flux_scale(const char* luaScaleFctName);
152 #endif
153 
154  public: // inherited from IElemDisc
156  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid) override;
157 
159  virtual bool use_hanging() const override;
160 
161  protected:
162  void previous_solution_required(bool b);
163 
164  template <typename TAlgebra>
165  void prep_timestep(number future_time, number time, VectorProxyBase* upb);
166  private:
168 
173  template<typename TElem, typename TFVGeom>
174  void prep_elem_loop(const ReferenceObjectID roid, const int si);
175 
177 
182  template<typename TElem, typename TFVGeom>
183  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
184 
186  template<typename TElem, typename TFVGeom>
187  void fsh_elem_loop();
188 
190  template<typename TElem, typename TFVGeom>
191  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
192 
194  template<typename TElem, typename TFVGeom>
195  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
196 
198  template<typename TElem, typename TFVGeom>
199  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
200 
202  template<typename TElem, typename TFVGeom>
203  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
204 
206  template<typename TElem, typename TFVGeom>
207  void add_rhs_elem(LocalVector& rhs, GridObject* elem, const MathVector<dim> vCornerCoords[]);
208 
210  template <typename TElem, typename TFVGeom>
211  void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si);
212 
214  template <typename TElem, typename TFVGeom>
215  void prep_err_est_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
216 
218  template <typename TElem, typename TFVGeom>
219  void compute_err_est_A_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
220 
222  template <typename TElem, typename TFVGeom>
223  void fsh_err_est_elem_loop();
224 
225  private:
226  DataImport<number, dim> m_fluxScale; // data import for scaling of fluxes
227 
228  private:
229  void register_all_fv1_funcs();
230 
231  template <typename TElem, typename TFVGeom>
232  void register_fv1_func();
233 
235  {
236  // friend class this_type;
239  template< typename TElem > void operator()(TElem&)
240  {
243  else
245  }
246  };
247 
248  template <typename List>
250  {
252  {
253  static const bool isEmpty = boost::mpl::empty<List>::value;
254  (typename boost::mpl::if_c<isEmpty, RegEnd, RegNext>::type(p));
255  }
256  struct RegEnd
257  {
259  };
260  struct RegNext
261  {
263  {
264  typedef typename boost::mpl::front<List>::type AlgebraType;
265  typedef typename boost::mpl::pop_front<List>::type NextList;
266 
268  p->set_prep_timestep_fct(aid, &this_type::template prep_timestep<AlgebraType>);
270  }
271  };
272  };
273 
274 
276  struct ShapeValues
277  {
278  public:
279  void resize(size_t nSip, size_t _nSh)
280  {
281  nSh = _nSh;
282  sideVals.resize(nSip);
283  for (size_t i = 0; i < nSip; i++) sideVals[i].resize(nSh);
284  }
285  number& shapeAtSideIP(size_t sh, size_t ip)
286  {
287  UG_ASSERT(ip < sideVals.size(), "Requested data for IP " << ip << ", but only " << sideVals.size() << " IPs present.");
288  UG_ASSERT(sh < sideVals[ip].size(), "Requested data for shape fct " << sh << ", but only " << sideVals[ip].size() << " shape fcts present.");
289  return sideVals.at(ip).at(sh);
290  }
291  number* shapesAtSideIP(size_t ip) {return &sideVals[ip][0];}
292  size_t num_sh() {return nSh;}
293  private:
294  size_t nSh;
295  std::vector<std::vector<number> > sideVals;
297 
298  private:
301 
305 
306  // temporary vectors (avoid re-allocating heap memory over and over again)
307  std::vector<LocalVector::value_type> m_uAtCorner;
308  std::vector<LocalVector::value_type> m_uOldAtCorner;
309 
310  int m_si;
311 };
312 
313 }
314 
315 #endif /*__H__UG__LIB_DISC__SPACIAL_DISCRETIZATION__ELEM_DISC__NEUMANN_BOUNDARY__FV1__INNER_BOUNDARY__*/
parameterString p
Definition: smart_pointer.h:108
class describing the type of an algebra
Definition: algebra_type.h:49
Type based UserData.
Definition: user_data.h:501
Definition: inner_boundary.h:98
virtual ~FV1InnerBoundaryElemDisc()
Destructor.
Definition: inner_boundary.h:142
FV1InnerBoundaryElemDisc(const std::vector< std::string > &functions, const std::vector< std::string > &subsets)
Constructor with functions.
Definition: inner_boundary.h:135
MultipleSideAndElemErrEstData< domain_type > err_est_type
error estimator type
Definition: inner_boundary.h:122
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:310
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:105
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:226
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:304
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:112
LocalVector m_locUOld
Definition: inner_boundary.h:303
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:100
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:116
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:101
IElemDisc< domain_type > base_type
Base class type.
Definition: inner_boundary.h:109
std::vector< LocalVector::value_type > m_uOldAtCorner
Definition: inner_boundary.h:308
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:128
base_type::position_type position_type
Position type.
Definition: inner_boundary.h:119
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:299
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:307
bool m_bCurrElemIsHSlave
Definition: inner_boundary.h:300
bool m_bPrevSolRequired
Definition: inner_boundary.h:302
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
TDomain::position_type position_type
Position type.
Definition: elem_disc_interface.h:479
TDomain domain_type
Domain type.
Definition: elem_disc_interface.h:476
Definition: elem_disc_interface.h:800
static const int dim
Definition: elem_disc_interface.h:803
Definition: local_algebra.h:422
Definition: local_algebra.h:198
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:238
void operator()(TElem &)
Definition: inner_boundary.h:239
RegisterAssemblingFuncs(this_type *pThis)
Definition: inner_boundary.h:237
RegEnd(this_type *)
Definition: inner_boundary.h:258
RegNext(this_type *p)
Definition: inner_boundary.h:262
RegisterPrepTimestep(this_type *p)
Definition: inner_boundary.h:251
struct holding values of shape functions in IPs
Definition: inner_boundary.h:277
void resize(size_t nSip, size_t _nSh)
Definition: inner_boundary.h:279
number & shapeAtSideIP(size_t sh, size_t ip)
Definition: inner_boundary.h:285
number * shapesAtSideIP(size_t ip)
Definition: inner_boundary.h:291
size_t num_sh()
Definition: inner_boundary.h:292
std::vector< std::vector< number > > sideVals
Definition: inner_boundary.h:295
size_t nSh
Definition: inner_boundary.h:294
Definition: inner_boundary.h:66
static const size_t _IGNORE_
index value for which a flux is ignored
Definition: inner_boundary.h:69
struct that holds information about the flux densities and from where to where the flux occurs
Definition: inner_boundary.h:75
std::vector< size_t > from
Definition: inner_boundary.h:78
std::vector< number > flux
Definition: inner_boundary.h:77
std::vector< size_t > to
Definition: inner_boundary.h:79
Definition: inner_boundary.h:85
std::vector< size_t > from
Definition: inner_boundary.h:90
std::vector< std::vector< std::pair< size_t, number > > > fluxDeriv
Definition: inner_boundary.h:89
std::vector< size_t > to
Definition: inner_boundary.h:91
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