Plugins
convection_diffusion_stab_fe.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2018: G-CSC, Goethe University Frankfurt
3  * Author: Arne Naegel
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__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_STAB_FE__
34 #define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_STAB_FE__
35 
36 // library intern headers
37 #include "../convection_diffusion_base.h"
39 
40 namespace ug{
41 namespace ConvectionDiffusionPlugin{
42 
43 // \ingroup lib_disc_elem_disc
46 
48 
54 template< typename TDomain>
55 class ConvectionDiffusionStabFE : public IElemDisc<TDomain>
56 {
57  private:
60 
63 
64  public:
66  static const int dim = base_type::dim;
67 
68  public:
70  ConvectionDiffusionStabFE(const char* functions, const char* subsets);
71  ConvectionDiffusionStabFE(const char* functions, const char* subsets, number stabM);
72  ConvectionDiffusionStabFE(const char* functions, const char* subsets, number stabM, number stabA);
73 
76 
78  void set_quad_order(size_t order);
79 
80  private:
82 
87  template <typename TElem, typename TFEGeom>
88  void prep_elem_loop(const ReferenceObjectID roid, const int si);
89 
91 
96  template <typename TElem, typename TFEGeom>
97  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
98 
100  template <typename TElem, typename TFEGeom>
101  void fsh_elem_loop();
102 
103  template <typename TElem, typename TFEGeom>
104  void add_jac_X_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], double stab);
105 
107  template <typename TElem, typename TFEGeom>
108  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
109 
111  template <typename TElem, typename TFEGeom>
112  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]){};
113 
115  template <typename TElem, typename TFEGeom>
116  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
117 
119  template <typename TElem, typename TFEGeom>
120  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
121 
123  template <typename TElem, typename TFEGeom>
124  void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
125 
127  template <typename TElem, typename TFEGeom>
128  void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si) {}
129 
131  template <typename TElem, typename TFEGeom>
132  void prep_err_est_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
133 
135  template <typename TElem, typename TFEGeom>
136  void compute_err_est_A_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale) {}
137 
139  template <typename TElem, typename TFEGeom>
140  void compute_err_est_M_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale) {}
141 
143  template <typename TElem, typename TFEGeom>
144  void compute_err_est_rhs_elem(GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale) {}
145 
147  template <typename TElem, typename TFEGeom>
149 
150  private:
152  static const size_t _C_ = 0;
153 
154 
155 
156  public:
158  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
159 
161  virtual bool use_hanging() const;
162 
163  protected:
167 
170 
173  void register_all_funcs(const LFEID& lfeid, const int quadOrder);
174  template <typename TElem, typename TFEGeom> void register_func();
176 
177  private:
179  /* struct ShapeValues
180  {
181  public:
182  void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
183  {
184  nSh = _nSh;
185  elemVals.resize(nEip);
186  sideVals.resize(nSip);
187  for (std::size_t i = 0; i < nEip; i++) elemVals[i].resize(nSh);
188  for (std::size_t i = 0; i < nSip; i++) sideVals[i].resize(nSh);
189  }
190  number& shapeAtElemIP(std::size_t sh, std::size_t ip) {return elemVals[ip][sh];}
191  number& shapeAtSideIP(std::size_t sh, std::size_t ip) {return sideVals[ip][sh];}
192  number* shapesAtElemIP(std::size_t ip) {return &elemVals[ip][0];}
193  number* shapesAtSideIP(std::size_t ip) {return &sideVals[ip][0];}
194  std::size_t num_sh() {return nSh;}
195  private:
196  std::size_t nSh;
197  std::vector<std::vector<number> > elemVals;
198  std::vector<std::vector<number> > sideVals;
199  } m_shapeValues;*/
200 
201  double m_stabParamM;
202  double m_stabParamA;
203 };
204 
205 // end group convection_diffusion
207 
208 } // end ConvectionDiffusionPlugin
209 } // end namespace ug
210 
211 
212 #endif /*__H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FE__*/
Discretization for the Convection-Diffusion Equation.
Definition: convection_diffusion_stab_fe.h:56
LFEID m_lfeID
current shape function set
Definition: convection_diffusion_stab_fe.h:169
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition: convection_diffusion_stab_fe.cpp:90
void add_jac_X_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], double stab)
Definition: convection_diffusion_stab_fe.cpp:213
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition: convection_diffusion_stab_fe.cpp:149
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: convection_diffusion_stab_fe.h:112
void set_quad_order(size_t order)
sets the quad order
Definition: convection_diffusion_stab_fe.cpp:82
void prep_err_est_elem(const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
prepares the element for assembling the error estimator
Definition: convection_diffusion_stab_fe.h:132
static const int dim
World dimension.
Definition: convection_diffusion_stab_fe.h:66
int m_quadOrder
Definition: convection_diffusion_stab_fe.h:166
bool m_bQuadOrderUserDef
current integration order
Definition: convection_diffusion_stab_fe.h:165
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition: convection_diffusion_stab_fe.cpp:122
double m_stabParamM
struct holding values of shape functions in IPs
Definition: convection_diffusion_stab_fe.h:201
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: convection_diffusion_stab_fe.cpp:172
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: convection_diffusion_stab_fe.h:128
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: convection_diffusion_stab_fe.cpp:256
void register_func()
Definition: convection_diffusion_stab_fe.cpp:388
void fsh_err_est_elem_loop()
postprocesses the loop over all elements of one type in the computation of the error estimator
Definition: convection_diffusion_stab_fe.h:148
void compute_err_est_rhs_elem(GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
computes the error estimator contribution for one element
Definition: convection_diffusion_stab_fe.h:144
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition: convection_diffusion_stab_fe.h:124
double m_stabParamA
Definition: convection_diffusion_stab_fe.h:202
virtual ~ConvectionDiffusionStabFE()
Destructor.
Definition: convection_diffusion_stab_fe.h:75
virtual bool use_hanging() const
returns if hanging nodes are needed
Definition: convection_diffusion_stab_fe.cpp:110
void register_all_funcs(const LFEID &lfeid, const int quadOrder)
ConvectionDiffusionStabFE< TDomain > this_type
Own type.
Definition: convection_diffusion_stab_fe.h:62
static const size_t _C_
abbreviation for the local solution
Definition: convection_diffusion_stab_fe.h:152
void compute_err_est_M_elem(const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
computes the error estimator contribution for one element
Definition: convection_diffusion_stab_fe.h:140
void fsh_elem_loop()
finishes the loop over all elements
Definition: convection_diffusion_stab_fe.cpp:143
IElemDisc< TDomain > base_type
Base class type.
Definition: convection_diffusion_stab_fe.h:59
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: convection_diffusion_stab_fe.cpp:248
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: convection_diffusion_stab_fe.h:136
ConvectionDiffusionStabFE(const char *functions, const char *subsets)
Constructor.
Definition: convection_diffusion_stab_fe.cpp:53
static const int dim
double number
ReferenceObjectID