Plugins
navier_stokes_fe.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2014: 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__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FE__NAVIER_STOKES_FE__
34 #define __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FE__NAVIER_STOKES_FE__
35 
36 // other ug4 modules
37 #include "common/common.h"
38 #include "lib_grid/lg_base.h"
39 
40 // library intern headers
44 
45 #include "../incompressible_navier_stokes_base.h"
46 
47 namespace ug{
48 namespace NavierStokes{
49 
52 
54 template< typename TDomain>
56  : public IncompressibleNavierStokesBase<TDomain>
57 {
58  protected:
61 
64 
65  public:
67  static const int dim = base_type::dim;
68 
69  public:
72  NavierStokesFE(const char* functions, const char* subsets);
73  NavierStokesFE(const std::vector<std::string>& vFct, const std::vector<std::string>& vSubset);
75 
77  void set_quad_order(size_t order);
78 
81 
84 
87 
90 
93 
95  void set_stabilization(number alpha) {m_stabParam = alpha;}
96 
97  public:
99  void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
100 
102  virtual std::string disc_type() const {return "fe";};
103 
104  protected:
107 
111 
115 
118 
119  void init();
120 
121  private:
124 
127 
130 
132  static const size_t _P_ = dim;
133 
136  using base_type::m_bStokes;
137  using base_type::m_bLaplace;
138 
139  public:
140  template<typename TElem, typename VGeom, typename PGeom>
141  void prep_elem_loop(const ReferenceObjectID roid, const int si);
142 
143  template<typename TElem, typename VGeom, typename PGeom>
144  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
145 
147  template<typename TElem, typename VGeom, typename PGeom>
148  void fsh_elem_loop();
149 
151  template<typename TElem, typename VGeom, typename PGeom>
152  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
153 
155  template<typename TElem, typename VGeom, typename PGeom>
156  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
157 
159  template<typename TElem, typename VGeom, typename PGeom>
160  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
161 
163  template<typename TElem, typename VGeom, typename PGeom>
164  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
165 
167  template<typename TElem, typename VGeom, typename PGeom>
168  void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
169 
170  // FVHO Assemblings
171  void register_all_funcs(const LFEID& vLfeID, const LFEID& pLfeID, const int quadOrder);
172  template<typename TElem, typename VGeom, typename PGeom> void register_func();
173 
174 };
175 
177 
178 } // namespace NavierStokes
179 } // end namespace ug
180 
181 #endif /*__H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FE__NAVIER_STOKES_FE__*/
function NavierStokes(fcts, subsets, discType)
SmartPtr< CplUserData< number, dim > > user_data()
Finite Volume Element Discretization for the incompressible Navier-Stokes Equation.
Definition: incompressible_navier_stokes_base.h:146
number m_bFullNewtonFactor
factor for exact jacobian, (1 for exact jacobian, 0 for fix point)
Definition: navier_stokes_base.h:207
bool m_bLaplace
flag if using only laplace term
Definition: incompressible_navier_stokes_base.h:256
bool m_bStokes
flag if solving the Stokes equation
Definition: incompressible_navier_stokes_base.h:253
bool m_bPecletBlend
flag if using Peclet Blending
Definition: incompressible_navier_stokes_base.h:244
static const int dim
World dimension.
Definition: incompressible_navier_stokes_base.h:156
Finite Element Discretization for the incompressible Navier-Stokes Equation.
Definition: navier_stokes_fe.h:57
void init()
Definition: navier_stokes_fe.cpp:63
DataImport< number, dim > m_imKinViscosity
Data import for kinematic viscosity.
Definition: navier_stokes_fe.h:126
static const size_t _P_
abbreviation for pressure
Definition: navier_stokes_fe.h:132
GridObject * m_pElem
current element
Definition: navier_stokes_fe.h:117
virtual std::string disc_type() const
returns string identifying disc type
Definition: navier_stokes_fe.h:102
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: navier_stokes_fe.cpp:482
void prep_elem_loop(const ReferenceObjectID roid, const int si)
Definition: navier_stokes_fe.cpp:150
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: navier_stokes_fe.cpp:206
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: navier_stokes_fe.cpp:507
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
Definition: navier_stokes_fe.cpp:184
void set_kinematic_viscosity(SmartPtr< CplUserData< number, dim > > user)
sets the kinematic viscosity
Definition: navier_stokes_fe.cpp:124
bool m_bQuadOrderUserDef
quadrature order
Definition: navier_stokes_fe.h:109
void fsh_elem_loop()
finishes the loop over all elements
Definition: navier_stokes_fe.cpp:177
void set_stabilization(number alpha)
sets the stabilization parameter
Definition: navier_stokes_fe.h:95
IncompressibleNavierStokesBase< TDomain > base_type
Base class type.
Definition: navier_stokes_fe.h:60
void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition: navier_stokes_fe.cpp:98
void set_density(SmartPtr< CplUserData< number, dim > > user)
sets the density
Definition: navier_stokes_fe.cpp:131
LFEID m_vLFEID
current shape function set
Definition: navier_stokes_fe.h:113
SmartPtr< CplUserData< number, dim > > kinematic_viscosity()
returns kinematic viscosity
Definition: navier_stokes_fe.h:83
int m_quadOrder
Definition: navier_stokes_fe.h:110
NavierStokesFE< TDomain > this_type
own type
Definition: navier_stokes_fe.h:63
LFEID m_pLFEID
Definition: navier_stokes_fe.h:114
void register_all_funcs(const LFEID &vLfeID, const LFEID &pLfeID, const int quadOrder)
void set_source(SmartPtr< CplUserData< MathVector< dim >, dim > > user)
sets the source function
Definition: navier_stokes_fe.cpp:138
DataImport< number, dim > m_imDensity
Data import for density.
Definition: navier_stokes_fe.h:129
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: navier_stokes_fe.cpp:350
DataImport< MathVector< dim >, dim > m_imSource
Data import for source.
Definition: navier_stokes_fe.h:123
NavierStokesFE(const char *functions, const char *subsets)
Definition: navier_stokes_fe.cpp:46
void set_quad_order(size_t order)
sets the quad order
Definition: navier_stokes_fe.cpp:89
number m_stabParam
stabilization parameter
Definition: navier_stokes_fe.h:102
SmartPtr< CplUserData< number, dim > > density()
returns density
Definition: navier_stokes_fe.h:89
static const int dim
World dimension.
Definition: navier_stokes_fe.h:67
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
assembles the local right hand side
Definition: navier_stokes_fe.cpp:536
void register_func()
Definition: navier_stokes_fe.cpp:590
double number
ReferenceObjectID