Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
convection_diffusion_fv1.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_FV1__
34#define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV1__
35
36// ug4 headers
38
39// plugin's internal headers
40#include "../convection_diffusion_base.h"
41#include "../convection_diffusion_sss.h"
42
43namespace ug{
44namespace ConvectionDiffusionPlugin{
45
46// \ingroup lib_disc_elem_disc
49
51
59template< typename TDomain>
61{
62 private:
65
68
71
72 public:
74 static const int dim = base_type::dim;
75
76 public:
78 ConvectionDiffusionFV1(const char* functions, const char* subsets);
79
81
87
90
93
95 void set_condensed_FV(bool condensed) {m_bCondensedFV = condensed;}
96
99
100 private:
102 virtual void prep_assemble_loop();
103
105
110 template <typename TElem, typename TFVGeom>
111 void prep_elem_loop(const ReferenceObjectID roid, const int si);
112
114
119 template <typename TElem, typename TFVGeom>
120 void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
121
123 template <typename TElem, typename TFVGeom>
124 void fsh_elem_loop();
125
127 template <typename TElem, typename TFVGeom>
128 void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
129
131 template <typename TElem, typename TFVGeom>
132 void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
133
135 template <typename TElem, typename TFVGeom>
136 void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
137
139 template <typename TElem, typename TFVGeom>
140 void add_def_A_expl_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
141
143 template <typename TElem, typename TFVGeom>
144 void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
145
147 template <typename TElem, typename TFVGeom>
148 void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
149
151 template <typename TElem, typename TFVGeom>
152 void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si);
153
155 template <typename TElem, typename TFVGeom>
156 void prep_err_est_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
157
159 template <typename TElem, typename TFVGeom>
160 void compute_err_est_A_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
161
163 template <typename TElem, typename TFVGeom>
164 void compute_err_est_M_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
165
167 template <typename TElem, typename TFVGeom>
168 void compute_err_est_rhs_elem(GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
169
171 template <typename TElem, typename TFVGeom>
173
174 private:
175
177 template<typename TElem, typename TFVGeom>
179 (
180 LocalMatrix& J,
181 const LocalVector& u,
182 GridObject* elem,
183 const TFVGeom& geo,
184 size_t i,
185 number flux
186 );
187
189 template<typename TElem, typename TFVGeom>
191 (
192 LocalVector& d,
193 const LocalVector& u,
194 GridObject* elem,
195 const TFVGeom& geo,
196 size_t i,
197 number flux
198 );
199
200 protected:
202 template <typename TElem, typename TFVGeom>
203 void lin_def_velocity(const LocalVector& u,
204 std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
205 const size_t nip);
206
208 template <typename TElem, typename TFVGeom>
209 void lin_def_diffusion(const LocalVector& u,
210 std::vector<std::vector<MathMatrix<dim,dim> > > vvvLinDef[],
211 const size_t nip);
212
214 template <typename TElem, typename TFVGeom>
215 void lin_def_flux(const LocalVector& u,
216 std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
217 const size_t nip);
218
220 template <typename TElem, typename TFVGeom>
221 void lin_def_reaction(const LocalVector& u,
222 std::vector<std::vector<number> > vvvLinDef[],
223 const size_t nip);
224
226 template <typename TElem, typename TFVGeom>
228 std::vector<std::vector<number> > vvvLinDef[],
229 const size_t nip);
230
232 template <typename TElem, typename TFVGeom>
233 void lin_def_source(const LocalVector& u,
234 std::vector<std::vector<number> > vvvLinDef[],
235 const size_t nip);
236
238 template <typename TElem, typename TFVGeom>
240 std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
241 const size_t nip);
242
244 template <typename TElem, typename TFVGeom>
245 void lin_def_mass_scale(const LocalVector& u,
246 std::vector<std::vector<number> > vvvLinDef[],
247 const size_t nip);
248
250 template <typename TElem, typename TFVGeom>
251 void lin_def_mass(const LocalVector& u,
252 std::vector<std::vector<number> > vvvLinDef[],
253 const size_t nip);
254
255 private:
257 static const size_t _C_ = 0;
258
261
274
277
278 protected:
281
284 const IConvectionShapes<dim>& get_updated_conv_shapes(const FVGeometryBase& geo, bool compute_deriv);
285
287 template <typename TElem, typename TFVGeom>
288 void ex_value(number vValue[],
289 const MathVector<dim> vGlobIP[],
290 number time, int si,
291 const LocalVector& u,
292 GridObject* elem,
293 const MathVector<dim> vCornerCoords[],
294 const MathVector<TFVGeom::dim> vLocIP[],
295 const size_t nip,
296 bool bDeriv,
297 std::vector<std::vector<number> > vvvDeriv[]);
298
300 template <typename TElem, typename TFVGeom>
301 void ex_grad(MathVector<dim> vValue[],
302 const MathVector<dim> vGlobIP[],
303 number time, int si,
304 const LocalVector& u,
305 GridObject* elem,
306 const MathVector<dim> vCornerCoords[],
307 const MathVector<TFVGeom::dim> vLocIP[],
308 const size_t nip,
309 bool bDeriv,
310 std::vector<std::vector<MathVector<dim> > > vvvDeriv[]);
311
312 public:
314 virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
315
317 virtual bool use_hanging() const;
318
319 protected:
322
325
328 void register_all_funcs(bool bHang);
329 template <typename TElem> void register_func_for_(bool bHang);
330 template <typename TElem, typename TFVGeom> void register_func();
332 {
334 RegisterLocalDiscr(this_type * pThis, bool bHang) : m_pThis(pThis), m_bHang(bHang) {}
335 template< typename TElem > void operator() (TElem &)
337 };
339
340 private:
343 {
344 public:
345 void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
346 {
347 nSh = _nSh;
348 elemVals.resize(nEip);
349 sideVals.resize(nSip);
350 for (std::size_t i = 0; i < nEip; i++) elemVals[i].resize(nSh);
351 for (std::size_t i = 0; i < nSip; i++) sideVals[i].resize(nSh);
352 }
353 number& shapeAtElemIP(std::size_t sh, std::size_t ip) {return elemVals[ip][sh];}
354 number& shapeAtSideIP(std::size_t sh, std::size_t ip) {return sideVals[ip][sh];}
355 number* shapesAtElemIP(std::size_t ip) {return &elemVals[ip][0];}
356 number* shapesAtSideIP(std::size_t ip) {return &sideVals[ip][0];}
357 std::size_t num_sh() {return nSh;}
358 private:
359 std::size_t nSh;
360 std::vector<std::vector<number> > elemVals;
361 std::vector<std::vector<number> > sideVals;
363};
364
365// end group convection_diffusion
367
368} // end ConvectionDiffusionPlugin
369} // end namespace ug
370
371
372#endif /*__H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV1__*/
Definition convection_diffusion_sss.h:194
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
FV Discretization for the Convection-Diffusion Equation.
Definition convection_diffusion_fv1.h:61
SmartPtr< IConvectionShapes< dim > > m_spConvShape
method to compute the upwind shapes
Definition convection_diffusion_fv1.h:280
void register_all_funcs(bool bHang)
Definition convection_diffusion_fv1.cpp:1725
void add_sss_def_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const TFVGeom &geo, size_t i, number flux)
adds contributions of a singular source or sink to the defect
Definition convection_diffusion_fv1.cpp:430
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 diffusion
Definition convection_diffusion_fv1.cpp:1239
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_fv1.cpp:452
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition convection_diffusion_fv1.cpp:700
bool m_bNonRegularGrid
current regular grid flag
Definition convection_diffusion_fv1.h:321
SmartPtr< CDSingularSourcesAndSinks< dim > > m_sss_mngr
singular sources and sinks manager
Definition convection_diffusion_fv1.h:260
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_fv1.cpp:402
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_fv1.cpp:1317
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_fv1.cpp:1288
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition convection_diffusion_fv1.cpp:101
void set_condensed_FV(bool condensed)
activates/deactivates the condensed scvf ip's for the FV scheme
Definition convection_diffusion_fv1.h:95
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_fv1.cpp:1342
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_fv1.cpp:252
void ex_value(number vValue[], const MathVector< dim > vGlobIP[], number time, int si, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< TFVGeom::dim > vLocIP[], const size_t nip, bool bDeriv, std::vector< std::vector< number > > vvvDeriv[])
computes the concentration
Definition convection_diffusion_fv1.cpp:1471
void fsh_elem_loop()
finishes the loop over all elements
Definition convection_diffusion_fv1.cpp:142
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_fv1.cpp:1446
void add_sss_jac_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const TFVGeom &geo, size_t i, number flux)
adds contributions of a singular source or sink to the matrix
Definition convection_diffusion_fv1.cpp:233
const IConvectionShapes< dim > & get_updated_conv_shapes(const FVGeometryBase &geo, bool compute_deriv)
Definition convection_diffusion_fv1.cpp:1692
virtual bool use_hanging() const
returns if hanging nodes are needed
Definition convection_diffusion_fv1.cpp:80
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition convection_diffusion_fv1.cpp:61
void register_func()
Definition convection_diffusion_fv1.cpp:1762
struct ug::ConvectionDiffusionPlugin::ConvectionDiffusionFV1::ShapeValues m_shapeValues
void register_func_for_(bool bHang)
Definition convection_diffusion_fv1.cpp:1741
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_fv1.cpp:818
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_fv1.cpp:1421
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_fv1.cpp:1071
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_fv1.cpp:664
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_fv1.cpp:1393
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_fv1.cpp:742
SideAndElemErrEstData< TDomain > err_est_type
error estimator type
Definition convection_diffusion_fv1.h:70
virtual void prep_assemble_loop()
prepares assembling
Definition convection_diffusion_fv1.cpp:91
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< TFVGeom::dim > vLocIP[], const size_t nip, bool bDeriv, std::vector< std::vector< MathVector< dim > > > vvvDeriv[])
computes the gradient of the concentration
Definition convection_diffusion_fv1.cpp:1585
IConvectionShapes< dim > conv_shape_type
returns the updated convection shapes
Definition convection_diffusion_fv1.h:283
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_fv1.cpp:1367
static const size_t _C_
abbreviation for the local solution
Definition convection_diffusion_fv1.h:257
static const int dim
World dimension.
Definition convection_diffusion_fv1.h:74
bool m_bCondensedFV
if to use the 'condensed' FV scvf ip's
Definition convection_diffusion_fv1.h:324
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_fv1.cpp:1188
void set_upwind(SmartPtr< IConvectionShapes< dim > > shapes)
set the upwind method
Definition convection_diffusion_fv1.cpp:1686
bool condensed_FV()
returns the 'condensed scvf ip' flag
Definition convection_diffusion_fv1.h:98
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_fv1.cpp:909
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_fv1.cpp:1201
void add_def_A_expl_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the stiffness part of the local defect explicit reaction, reaction_rate and source
Definition convection_diffusion_fv1.cpp:605
ConvectionDiffusionFV1< TDomain > this_type
Own type.
Definition convection_diffusion_fv1.h:67
void set_sss_manager(SmartPtr< CDSingularSourcesAndSinks< dim > > sss_mngr)
set singular sources and sinks
Definition convection_diffusion_fv1.h:89
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition convection_diffusion_fv1.cpp: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_fv1.cpp:1122
SmartPtr< CDSingularSourcesAndSinks< dim > > sss_manager()
get singular sources and sinks
Definition convection_diffusion_fv1.h:92
ConvectionDiffusionBase< TDomain > base_type
Base class type.
Definition convection_diffusion_fv1.h:64
number time() const
double number
ReferenceObjectID
void operator()(TElem &)
Definition convection_diffusion_fv1.h:335
RegisterLocalDiscr(this_type *pThis, bool bHang)
Definition convection_diffusion_fv1.h:334
this_type * m_pThis
Definition convection_diffusion_fv1.h:333
struct holding values of shape functions in IPs
Definition convection_diffusion_fv1.h:343
std::size_t num_sh()
Definition convection_diffusion_fv1.h:357
std::vector< std::vector< number > > elemVals
Definition convection_diffusion_fv1.h:360
number * shapesAtSideIP(std::size_t ip)
Definition convection_diffusion_fv1.h:356
void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
Definition convection_diffusion_fv1.h:345
std::size_t nSh
Definition convection_diffusion_fv1.h:359
number * shapesAtElemIP(std::size_t ip)
Definition convection_diffusion_fv1.h:355
number & shapeAtSideIP(std::size_t sh, std::size_t ip)
Definition convection_diffusion_fv1.h:354
std::vector< std::vector< number > > sideVals
Definition convection_diffusion_fv1.h:361
number & shapeAtElemIP(std::size_t sh, std::size_t ip)
Definition convection_diffusion_fv1.h:353