Plugins
convection_diffusion_fvcr.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3  * Author: Christian Wehner, 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_FVCR__
34 #define __H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FVCR__
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  ConvectionDiffusionFVCR(const char* functions, const char* subsets);
73 
75 
81 
82  private:
84 
89  template <typename TElem, typename TFVGeom>
90  void prep_elem_loop(const ReferenceObjectID roid, const int si);
91 
93 
98  template <typename TElem, typename TFVGeom>
99  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
100 
102  template <typename TElem, typename TFVGeom>
103  void fsh_elem_loop();
104 
106  template <typename TElem, typename TFVGeom>
107  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
108 
110  template <typename TElem, typename TFVGeom>
111  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
112 
114  template <typename TElem, typename TFVGeom>
115  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
116 
118  template <typename TElem, typename TFVGeom>
119  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
120 
122  template <typename TElem, typename TFVGeom>
123  void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
124 
125  protected:
127  template <typename TElem, typename TFVGeom>
128  void lin_def_velocity(const LocalVector& u,
129  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
130  const size_t nip);
131 
133  template <typename TElem, typename TFVGeom>
134  void lin_def_diffusion(const LocalVector& u,
135  std::vector<std::vector<MathMatrix<dim,dim> > > vvvLinDef[],
136  const size_t nip);
137 
139  template <typename TElem, typename TFVGeom>
140  void lin_def_reaction(const LocalVector& u,
141  std::vector<std::vector<number> > vvvLinDef[],
142  const size_t nip);
143 
145  template <typename TElem, typename TFVGeom>
146  void lin_def_reaction_rate(const LocalVector& u,
147  std::vector<std::vector<number> > vvvLinDef[],
148  const size_t nip);
149 
151  template <typename TElem, typename TFVGeom>
152  void lin_def_source(const LocalVector& u,
153  std::vector<std::vector<number> > vvvLinDef[],
154  const size_t nip);
155 
157  template <typename TElem, typename TFVGeom>
159  std::vector<std::vector<MathVector<dim> > > vvvLinDef[],
160  const size_t nip);
161 
163  template <typename TElem, typename TFVGeom>
164  void lin_def_mass_scale(const LocalVector& u,
165  std::vector<std::vector<number> > vvvLinDef[],
166  const size_t nip);
167 
169  template <typename TElem, typename TFVGeom>
170  void lin_def_mass(const LocalVector& u,
171  std::vector<std::vector<number> > vvvLinDef[],
172  const size_t nip);
173 
174  private:
176  static const size_t _C_ = 0;
177 
180  using base_type::m_imFlux;
181  using base_type::m_imSource;
189  using base_type::m_imMass;
190 
191  using base_type::m_exGrad;
192  using base_type::m_exValue;
193 
194  protected:
197 
201 
203  template <typename TElem, typename TFVGeom>
204  void ex_value(number vValue[],
205  const MathVector<dim> vGlobIP[],
206  number time, int si,
207  const LocalVector& u,
208  GridObject* elem,
209  const MathVector<dim> vCornerCoords[],
210  const MathVector<TFVGeom::dim> vLocIP[],
211  const size_t nip,
212  bool bDeriv,
213  std::vector<std::vector<number> > vvvDeriv[]);
214 
216  template <typename TElem, typename TFVGeom>
217  void ex_grad(MathVector<dim> vValue[],
218  const MathVector<dim> vGlobIP[],
219  number time, int si,
220  const LocalVector& u,
221  GridObject* elem,
222  const MathVector<dim> vCornerCoords[],
223  const MathVector<TFVGeom::dim> vLocIP[],
224  const size_t nip,
225  bool bDeriv,
226  std::vector<std::vector<MathVector<dim> > > vvvDeriv[]);
227 
228  public:
230  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
231 
233  virtual bool use_hanging() const;
234 
235  protected:
238 
241  void register_all_funcs(bool bHang);
242  template <typename TElem, typename TFVGeom> void register_func();
244 };
245 
246 // end group convection_diffusion
248 
249 } // end ConvectionDiffusionPlugin
250 } // end namespace ug
251 
252 
253 #endif /*__H__UG__LIB_DISC__CONVECTION_DIFFUSION__CONVECTION_DIFFUSION_FVCR__*/
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_fvcr.h:58
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_fvcr.cpp:514
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_fvcr.cpp:343
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_fvcr.cpp:727
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_fvcr.cpp:316
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for assembling
Definition: convection_diffusion_fvcr.cpp:146
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_fvcr.cpp:677
static const size_t _C_
abbreviation for the local solution
Definition: convection_diffusion_fvcr.h:176
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the loop over all elements
Definition: convection_diffusion_fvcr.cpp: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 velocity
Definition: convection_diffusion_fvcr.cpp:552
IConvectionShapes< dim > conv_shape_type
returns the updated convection shapes
Definition: convection_diffusion_fvcr.h:199
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition: convection_diffusion_fvcr.cpp:487
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_fvcr.cpp:627
virtual bool use_hanging() const
returns if hanging nodes are needed
Definition: convection_diffusion_fvcr.cpp:78
void fsh_elem_loop()
finishes the loop over all elements
Definition: convection_diffusion_fvcr.cpp:140
const IConvectionShapes< dim > & get_updated_conv_shapes(const FVGeometryBase &geo)
Definition: convection_diffusion_fvcr.cpp:905
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_fvcr.cpp:814
void set_upwind(SmartPtr< IConvectionShapes< dim > > shapes)
set the upwind method
Definition: convection_diffusion_fvcr.cpp:899
SmartPtr< IConvectionShapes< dim > > m_spConvShape
method to compute the upwind shapes
Definition: convection_diffusion_fvcr.h:196
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_fvcr.cpp:451
ConvectionDiffusionFVCR(const char *functions, const char *subsets)
Constructor.
Definition: convection_diffusion_fvcr.cpp:49
ConvectionDiffusionFVCR< TDomain > this_type
Own type.
Definition: convection_diffusion_fvcr.h:64
ConvectionDiffusionBase< TDomain > base_type
Base class type.
Definition: convection_diffusion_fvcr.h:61
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_fvcr.cpp:602
static const int dim
World dimension.
Definition: convection_diffusion_fvcr.h:68
bool m_bNonRegularGrid
current regular grid flag
Definition: convection_diffusion_fvcr.h:237
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_fvcr.cpp:203
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_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_fvcr.cpp:702
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition: convection_diffusion_fvcr.cpp:59
void register_func()
Definition: convection_diffusion_fvcr.cpp:987
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_fvcr.cpp:652
number time() const
double number
ReferenceObjectID