Plugins
convection_diffusion_fv.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_FV__
34 #define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV__
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 
56 template< typename TDomain>
58 {
59  private:
62 
65 
66  public:
68  static const int dim = base_type::dim;
69 
70  public:
72  ConvectionDiffusionFV(const char* functions, const char* subsets);
73 
74  private:
76 
81  template <typename TElem, typename TFVGeom>
82  void prep_elem_loop(const ReferenceObjectID roid, const int si);
83 
85 
90  template <typename TElem, typename TFVGeom>
91  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
92 
94  template <typename TElem, typename TFVGeom>
95  void fsh_elem_loop();
96 
98  template <typename TElem, typename TFVGeom>
99  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
100 
102  template <typename TElem, typename TFVGeom>
103  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
104 
106  template <typename TElem, typename TFVGeom>
107  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
108 
110  template <typename TElem, typename TFVGeom>
111  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
112 
114  template <typename TElem, typename TFVGeom>
115  void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
116 
117  protected:
119  template <typename TElem, typename TFVGeom>
120  void lin_def_velocity(const LocalVector& u,
121  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
122  const size_t nip);
123 
125  template <typename TElem, typename TFVGeom>
126  void lin_def_diffusion(const LocalVector& u,
127  std::vector<std::vector<MathMatrix<dim,dim> > > vvvLinDef[],
128  const size_t nip);
129 
131  template <typename TElem, typename TFVGeom>
132  void lin_def_flux(const LocalVector& u,
133  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
134  const size_t nip);
135 
137  template <typename TElem, typename TFVGeom>
138  void lin_def_reaction(const LocalVector& u,
139  std::vector<std::vector<number> > vvvLinDef[],
140  const size_t nip);
141 
143  template <typename TElem, typename TFVGeom>
144  void lin_def_reaction_rate(const LocalVector& u,
145  std::vector<std::vector<number> > vvvLinDef[],
146  const size_t nip);
147 
149  template <typename TElem, typename TFVGeom>
150  void lin_def_source(const LocalVector& u,
151  std::vector<std::vector<number> > vvvLinDef[],
152  const size_t nip);
153 
155  template <typename TElem, typename TFVGeom>
157  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
158  const size_t nip);
159 
161  template <typename TElem, typename TFVGeom>
162  void lin_def_mass_scale(const LocalVector& u,
163  std::vector<std::vector<number> > vvvLinDef[],
164  const size_t nip);
165 
167  template <typename TElem, typename TFVGeom>
168  void lin_def_mass(const LocalVector& u,
169  std::vector<std::vector<number> > vvvLinDef[],
170  const size_t nip);
171 
172  private:
174  static const size_t _C_ = 0;
175 
178  using base_type::m_imFlux;
179  using base_type::m_imSource;
187  using base_type::m_imMass;
188 
189  using base_type::m_exGrad;
190  using base_type::m_exValue;
191 
192  protected:
194  template <typename TElem, typename TFVGeom>
195  void ex_value(number vValue[],
196  const MathVector<dim> vGlobIP[],
197  number time, int si,
198  const LocalVector& u,
199  GridObject* elem,
200  const MathVector<dim> vCornerCoords[],
201  const MathVector<TFVGeom::dim> vLocIP[],
202  const size_t nip,
203  bool bDeriv,
204  std::vector<std::vector<number> > vvvDeriv[]);
205 
207  template <typename TElem, typename TFVGeom>
208  void ex_grad(MathVector<dim> vValue[],
209  const MathVector<dim> vGlobIP[],
210  number time, int si,
211  const LocalVector& u,
212  GridObject* elem,
213  const MathVector<dim> vCornerCoords[],
214  const MathVector<TFVGeom::dim> vLocIP[],
215  const size_t nip,
216  bool bDeriv,
217  std::vector<std::vector<MathVector<dim> > > vvvDeriv[]);
218 
219  public:
221  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
222 
224  virtual bool use_hanging() const;
225 
227  void set_quad_order(size_t order);
228 
229  protected:
232 
236 
239  void register_all_funcs(const LFEID& lfeID, const int quadOrder);
240  template <typename TElem, typename TFVGeom> void register_func();
242 };
243 
244 // end group convection_diffusion
246 
247 } // end ConvectionDiffusionPlugin
248 } // end namespace ug
249 
250 
251 #endif /*__H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FV__*/
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
number time() const
double number
ReferenceObjectID