Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
convection_diffusion_fe.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-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__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FE__
34#define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FE__
35
36// library intern headers
37#include "../convection_diffusion_base.h"
39
40namespace ug{
41namespace ConvectionDiffusionPlugin{
42
43// \ingroup lib_disc_elem_disc
46
48
56template< typename TDomain>
58{
59 private:
62
65
68
69 public:
71 static const int dim = base_type::dim;
72
73 public:
75 ConvectionDiffusionFE(const char* functions, const char* subsets);
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
104 template <typename TElem, typename TFEGeom>
105 void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
106
108 template <typename TElem, typename TFEGeom>
109 void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
110
112 template <typename TElem, typename TFEGeom>
113 void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
114
116 template <typename TElem, typename TFEGeom>
117 void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
118
120 template <typename TElem, typename TFEGeom>
121 void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
122
124 template <typename TElem, typename TFEGeom>
125 void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si);
126
128 template <typename TElem, typename TFEGeom>
129 void prep_err_est_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
130
132 template <typename TElem, typename TFEGeom>
133 void compute_err_est_A_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
134
136 template <typename TElem, typename TFEGeom>
137 void compute_err_est_M_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
138
140 template <typename TElem, typename TFEGeom>
141 void compute_err_est_rhs_elem(GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
142
144 template <typename TElem, typename TFEGeom>
146
147 protected:
149 template <typename TElem, typename TFEGeom>
150 void lin_def_velocity(const LocalVector& u,
151 std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
152 const size_t nip);
153
155 template <typename TElem, typename TFEGeom>
156 void lin_def_diffusion(const LocalVector& u,
157 std::vector<std::vector<MathMatrix<dim,dim> > > vvvLinDef[],
158 const size_t nip);
159
161 template <typename TElem, typename TFEGeom>
162 void lin_def_flux(const LocalVector& u,
163 std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
164 const size_t nip);
165
167 template <typename TElem, typename TFEGeom>
168 void lin_def_reaction(const LocalVector& u,
169 std::vector<std::vector<number> > vvvLinDef[],
170 const size_t nip);
171
173 template <typename TElem, typename TFEGeom>
175 std::vector<std::vector<number> > vvvLinDef[],
176 const size_t nip);
177
179 template <typename TElem, typename TFEGeom>
180 void lin_def_source(const LocalVector& u,
181 std::vector<std::vector<number> > vvvLinDef[],
182 const size_t nip);
183
185 template <typename TElem, typename TFEGeom>
187 std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
188 const size_t nip);
189
191 template <typename TElem, typename TFEGeom>
192 void lin_def_mass_scale(const LocalVector& u,
193 std::vector<std::vector<number> > vvvLinDef[],
194 const size_t nip);
195
197 template <typename TElem, typename TFEGeom>
198 void lin_def_mass(const LocalVector& u,
199 std::vector<std::vector<number> > vvvLinDef[],
200 const size_t nip);
201
202 private:
204 static const size_t _C_ = 0;
205
218
221
222 protected:
224 template <typename TElem, typename TFEGeom>
225 void ex_value(number vValue[],
226 const MathVector<dim> vGlobIP[],
227 number time, int si,
228 const LocalVector& u,
229 GridObject* elem,
230 const MathVector<dim> vCornerCoords[],
231 const MathVector<TFEGeom::dim> vLocIP[],
232 const size_t nip,
233 bool bDeriv,
234 std::vector<std::vector<number> > vvvDeriv[]);
235
237 template <typename TElem, typename TFEGeom>
238 void ex_grad(MathVector<dim> vValue[],
239 const MathVector<dim> vGlobIP[],
240 number time, int si,
241 const LocalVector& u,
242 GridObject* elem,
243 const MathVector<dim> vCornerCoords[],
244 const MathVector<TFEGeom::dim> vLocIP[],
245 const size_t nip,
246 bool bDeriv,
247 std::vector<std::vector<MathVector<dim> > > vvvDeriv[]);
248
249 public:
251 virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
252
254 virtual bool use_hanging() const;
255
256 protected:
260
263
266 void register_all_funcs(const LFEID& lfeid, const int quadOrder);
267 template <typename TElem, typename TFEGeom> void register_func();
269
270 private:
273 {
274 public:
275 void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
276 {
277 nSh = _nSh;
278 elemVals.resize(nEip);
279 sideVals.resize(nSip);
280 for (std::size_t i = 0; i < nEip; i++) elemVals[i].resize(nSh);
281 for (std::size_t i = 0; i < nSip; i++) sideVals[i].resize(nSh);
282 }
283 number& shapeAtElemIP(std::size_t sh, std::size_t ip) {return elemVals[ip][sh];}
284 number& shapeAtSideIP(std::size_t sh, std::size_t ip) {return sideVals[ip][sh];}
285 number* shapesAtElemIP(std::size_t ip) {return &elemVals[ip][0];}
286 number* shapesAtSideIP(std::size_t ip) {return &sideVals[ip][0];}
287 std::size_t num_sh() {return nSh;}
288 private:
289 std::size_t nSh;
290 std::vector<std::vector<number> > elemVals;
291 std::vector<std::vector<number> > sideVals;
293
294};
295
296// end group convection_diffusion
298
299} // end ConvectionDiffusionPlugin
300} // end namespace ug
301
302
303#endif /*__H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FE__*/
Discretization for the Convection-Diffusion Equation.
Definition convection_diffusion_base.h:79
DataImport< number, dim > m_imReactionRate
Data import for the reaction term.
Definition convection_diffusion_base.h:246
DataImport< number, dim > m_imSource
Data import for the right-hand side (volume)
Definition convection_diffusion_base.h:261
DataImport< number, dim > m_imReactionExpl
Data import for the reaction term explicit.
Definition convection_diffusion_base.h:255
DataImport< number, dim > m_imSourceExpl
Data import for the source term explicit.
Definition convection_diffusion_base.h:258
DataImport< number, dim > m_imMass
Data import for the mass scale.
Definition convection_diffusion_base.h:270
DataImport< MathVector< dim >, dim > m_imVectorSource
Data import for the right-hand side (vector)
Definition convection_diffusion_base.h:264
SmartPtr< DataExport< number, dim > > m_exValue
Export for the concentration.
Definition convection_diffusion_base.h:288
DataImport< MathVector< dim >, dim > m_imFlux
Data import for the Flux.
Definition convection_diffusion_base.h:243
DataImport< MathMatrix< dim, dim >, dim > m_imDiffusion
Data import for Diffusion.
Definition convection_diffusion_base.h:237
DataImport< MathVector< dim >, dim > m_imVelocity
Data import for the Velocity field.
Definition convection_diffusion_base.h:240
DataImport< number, dim > m_imMassScale
Data import for the mass scale.
Definition convection_diffusion_base.h:267
DataImport< number, dim > m_imReaction
Data import for the reaction term.
Definition convection_diffusion_base.h:249
SmartPtr< DataExport< MathVector< dim >, dim > > m_exGrad
Export for the gradient of concentration.
Definition convection_diffusion_base.h:291
DataImport< number, dim > m_imReactionRateExpl
Data import for the reaction_rate term explicit.
Definition convection_diffusion_base.h:252
static const int dim
World dimension.
Definition convection_diffusion_base.h:86
Discretization for the Convection-Diffusion Equation.
Definition convection_diffusion_fe.h:58
static const size_t _C_
abbreviation for the local solution
Definition convection_diffusion_fe.h:204
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition convection_diffusion_fe.cpp:343
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_fe.cpp:440
virtual bool use_hanging() const
returns if hanging nodes are needed
Definition convection_diffusion_fe.cpp:85
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_fe.cpp:681
bool m_bQuadOrderUserDef
current integration order
Definition convection_diffusion_fe.h:258
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition convection_diffusion_fe.cpp:95
ConvectionDiffusionFE< TDomain > this_type
Own type.
Definition convection_diffusion_fe.h:64
void lin_def_mass_scale(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the mass scale term
Definition convection_diffusion_fe.cpp:953
void lin_def_vector_source(const LocalVector &u, std::vector< std::vector< MathVector< dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the vector source term
Definition convection_diffusion_fe.cpp:929
ConvectionDiffusionBase< TDomain > base_type
Base class type.
Definition convection_diffusion_fe.h:61
void lin_def_mass(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the mass scale term
Definition convection_diffusion_fe.cpp:984
void fsh_elem_loop()
finishes the loop over all elements
Definition convection_diffusion_fe.cpp:126
void lin_def_source(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the source term
Definition convection_diffusion_fe.cpp:906
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_fe.cpp:383
void lin_def_flux(const LocalVector &u, std::vector< std::vector< MathVector< dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the flux
Definition convection_diffusion_fe.cpp:792
SideAndElemErrEstData< TDomain > err_est_type
error estimator type
Definition convection_diffusion_fe.h:67
int m_quadOrder
Definition convection_diffusion_fe.h:259
void register_all_funcs(const LFEID &lfeid, const int quadOrder)
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_fe.cpp:239
void lin_def_velocity(const LocalVector &u, std::vector< std::vector< MathVector< dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the velocity
Definition convection_diffusion_fe.cpp:763
void set_quad_order(size_t order)
sets the quad order
Definition convection_diffusion_fe.cpp:58
void lin_def_reaction(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the reaction
Definition convection_diffusion_fe.cpp:848
static const int dim
World dimension.
Definition convection_diffusion_fe.h:71
struct ug::ConvectionDiffusionPlugin::ConvectionDiffusionFE::ShapeValues m_shapeValues
void ex_grad(MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< TFEGeom::dim > vLocIP[], const size_t nip, bool bDeriv, std::vector< std::vector< MathVector< dim > > > vvvDeriv[])
computes the gradient of the concentration
Definition convection_diffusion_fe.cpp:1091
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_fe.cpp:629
void lin_def_diffusion(const LocalVector &u, std::vector< std::vector< MathMatrix< dim, dim > > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the velocity
Definition convection_diffusion_fe.cpp:816
void ex_value(number vValue[], const MathVector< dim > vGlobIP[], number time, int si, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< TFEGeom::dim > vLocIP[], const size_t nip, bool bDeriv, std::vector< std::vector< number > > vvvDeriv[])
computes the concentration
Definition convection_diffusion_fe.cpp:1010
void lin_def_reaction_rate(const LocalVector &u, std::vector< std::vector< number > > vvvLinDef[], const size_t nip)
computes the linearized defect w.r.t to the reaction
Definition convection_diffusion_fe.cpp:874
void register_func()
Definition convection_diffusion_fe.cpp:1297
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_fe.cpp:158
LFEID m_lfeID
current shape function set
Definition convection_diffusion_fe.h:262
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_fe.cpp:301
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition convection_diffusion_fe.cpp:132
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_fe.cpp:486
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_fe.cpp:209
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_fe.cpp:749
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition convection_diffusion_fe.cpp:66
number time() const
double number
ReferenceObjectID
struct holding values of shape functions in IPs
Definition convection_diffusion_fe.h:273
std::vector< std::vector< number > > elemVals
Definition convection_diffusion_fe.h:290
std::vector< std::vector< number > > sideVals
Definition convection_diffusion_fe.h:291
std::size_t num_sh()
Definition convection_diffusion_fe.h:287
number & shapeAtElemIP(std::size_t sh, std::size_t ip)
Definition convection_diffusion_fe.h:283
void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
Definition convection_diffusion_fe.h:275
number * shapesAtSideIP(std::size_t ip)
Definition convection_diffusion_fe.h:286
number & shapeAtSideIP(std::size_t sh, std::size_t ip)
Definition convection_diffusion_fe.h:284
std::size_t nSh
Definition convection_diffusion_fe.h:289
number * shapesAtElemIP(std::size_t ip)
Definition convection_diffusion_fe.h:285