33 #ifndef __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV__
34 #define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV__
37 #include "../convection_diffusion_base.h"
41 namespace ConvectionDiffusionPlugin{
56 template<
typename TDomain>
81 template <
typename TElem,
typename TFVGeom>
90 template <
typename TElem,
typename TFVGeom>
94 template <
typename TElem,
typename TFVGeom>
98 template <
typename TElem,
typename TFVGeom>
102 template <
typename TElem,
typename TFVGeom>
106 template <
typename TElem,
typename TFVGeom>
110 template <
typename TElem,
typename TFVGeom>
114 template <
typename TElem,
typename TFVGeom>
119 template <
typename TElem,
typename TFVGeom>
125 template <
typename TElem,
typename TFVGeom>
131 template <
typename TElem,
typename TFVGeom>
137 template <
typename TElem,
typename TFVGeom>
139 std::vector<std::vector<number> > vvvLinDef[],
143 template <
typename TElem,
typename TFVGeom>
145 std::vector<std::vector<number> > vvvLinDef[],
149 template <
typename TElem,
typename TFVGeom>
151 std::vector<std::vector<number> > vvvLinDef[],
155 template <
typename TElem,
typename TFVGeom>
161 template <
typename TElem,
typename TFVGeom>
163 std::vector<std::vector<number> > vvvLinDef[],
167 template <
typename TElem,
typename TFVGeom>
169 std::vector<std::vector<number> > vvvLinDef[],
174 static const size_t _C_ = 0;
194 template <
typename TElem,
typename TFVGeom>
204 std::vector<std::vector<number> > vvvDeriv[]);
207 template <
typename TElem,
typename TFVGeom>
221 virtual void prepare_setting(
const std::vector<LFEID>& vLfeID,
bool bNonRegularGrid);
240 template <
typename TElem,
typename TFVGeom>
void register_func();
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_fv.h:58
void register_all_funcs(const LFEID &lfeID, const int quadOrder)
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition: convection_diffusion_fv.cpp:149
LFEID m_lfeID
current shape function set
Definition: convection_diffusion_fv.h:231
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition: convection_diffusion_fv.cpp:534
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_fv.cpp:946
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_fv.cpp:746
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
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_fv.cpp:571
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition: convection_diffusion_fv.cpp:105
ConvectionDiffusionBase< TDomain > base_type
Base class type.
Definition: convection_diffusion_fv.h:61
static const size_t _C_
abbreviation for the local solution
Definition: convection_diffusion_fv.h:174
void register_func()
Definition: convection_diffusion_fv.cpp:1175
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_fv.cpp:809
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_fv.cpp:299
static const int dim
World dimension.
Definition: convection_diffusion_fv.h:68
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_fv.cpp:683
void set_quad_order(size_t order)
sets the quad order
Definition: convection_diffusion_fv.cpp:93
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition: convection_diffusion_fv.cpp:56
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_fv.cpp:194
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_fv.cpp:717
ConvectionDiffusionFV< TDomain > this_type
Own type.
Definition: convection_diffusion_fv.h:64
virtual bool use_hanging() const
returns if hanging nodes are needed
Definition: convection_diffusion_fv.cpp:86
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 velocity
Definition: convection_diffusion_fv.cpp:608
int m_quadOrder
Definition: convection_diffusion_fv.h:235
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_fv.cpp:482
void fsh_elem_loop()
finishes the loop over all elements
Definition: convection_diffusion_fv.cpp:143
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_fv.cpp:838
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_fv.cpp:640
bool m_bQuadOrderUserDef
current integration order
Definition: convection_diffusion_fv.h:234
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_fv.cpp:775
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_fv.cpp:340
ConvectionDiffusionFV(const char *functions, const char *subsets)
Constructor.
Definition: convection_diffusion_fv.cpp:47