33 #ifndef __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV1__
34 #define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV1__
40 #include "../convection_diffusion_base.h"
41 #include "../convection_diffusion_sss.h"
44 namespace ConvectionDiffusionPlugin{
59 template<
typename TDomain>
110 template <
typename TElem,
typename TFVGeom>
119 template <
typename TElem,
typename TFVGeom>
123 template <
typename TElem,
typename TFVGeom>
127 template <
typename TElem,
typename TFVGeom>
131 template <
typename TElem,
typename TFVGeom>
135 template <
typename TElem,
typename TFVGeom>
139 template <
typename TElem,
typename TFVGeom>
143 template <
typename TElem,
typename TFVGeom>
147 template <
typename TElem,
typename TFVGeom>
151 template <
typename TElem,
typename TFVGeom>
155 template <
typename TElem,
typename TFVGeom>
159 template <
typename TElem,
typename TFVGeom>
163 template <
typename TElem,
typename TFVGeom>
167 template <
typename TElem,
typename TFVGeom>
171 template <
typename TElem,
typename TFVGeom>
177 template<
typename TElem,
typename TFVGeom>
189 template<
typename TElem,
typename TFVGeom>
202 template <
typename TElem,
typename TFVGeom>
208 template <
typename TElem,
typename TFVGeom>
214 template <
typename TElem,
typename TFVGeom>
220 template <
typename TElem,
typename TFVGeom>
222 std::vector<std::vector<number> > vvvLinDef[],
226 template <
typename TElem,
typename TFVGeom>
228 std::vector<std::vector<number> > vvvLinDef[],
232 template <
typename TElem,
typename TFVGeom>
234 std::vector<std::vector<number> > vvvLinDef[],
238 template <
typename TElem,
typename TFVGeom>
244 template <
typename TElem,
typename TFVGeom>
246 std::vector<std::vector<number> > vvvLinDef[],
250 template <
typename TElem,
typename TFVGeom>
252 std::vector<std::vector<number> > vvvLinDef[],
257 static const size_t _C_ = 0;
287 template <
typename TElem,
typename TFVGeom>
297 std::vector<std::vector<number> > vvvDeriv[]);
300 template <
typename TElem,
typename TFVGeom>
314 virtual void prepare_setting(
const std::vector<LFEID>& vLfeID,
bool bNonRegularGrid);
330 template <
typename TElem,
typename TFVGeom>
void register_func();
345 void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
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
SmartPtr< CDSingularSourcesAndSinks< dim > > sss_manager()
get singular sources and sinks
Definition: convection_diffusion_fv1.h:92
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
ConvectionDiffusionFV1(const char *functions, const char *subsets)
Constructor.
Definition: convection_diffusion_fv1.cpp:51
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
ConvectionDiffusionBase< TDomain > base_type
Base class type.
Definition: convection_diffusion_fv1.h:64
Definition: convection_diffusion_fv1.h:332
void operator()(TElem &)
Definition: convection_diffusion_fv1.h:335
bool m_bHang
Definition: convection_diffusion_fv1.h:333
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
void resize(std::size_t nEip, std::size_t nSip, std::size_t _nSh)
Definition: convection_diffusion_fv1.h:345
number * shapesAtElemIP(std::size_t ip)
Definition: convection_diffusion_fv1.h:355
std::size_t nSh
Definition: convection_diffusion_fv1.h:359
number * shapesAtSideIP(std::size_t ip)
Definition: convection_diffusion_fv1.h:356
number & shapeAtElemIP(std::size_t sh, std::size_t ip)
Definition: convection_diffusion_fv1.h:353
std::vector< std::vector< number > > sideVals
Definition: convection_diffusion_fv1.h:361
number & shapeAtSideIP(std::size_t sh, std::size_t ip)
Definition: convection_diffusion_fv1.h:354