Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
navier_stokes_fvcr.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2010-2017: G-CSC, Goethe University Frankfurt
3 * Author: Christian Wehner
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__FVCR__NAVIER_STOKES_FVCR__
34#define __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FVCR__NAVIER_STOKES_FVCR__
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#include "../../upwind_interface.h"
47
48namespace ug{
49namespace NavierStokes{
50
53
54template< typename TDomain>
56 : public IncompressibleNavierStokesBase<TDomain>
57{
58 private:
61
64
65 public:
67 static const int dim = base_type::dim;
68
69 public:
72 NavierStokesFVCR(const char* functions, const char* subsets);
73 NavierStokesFVCR(const std::vector<std::string>& vFct, const std::vector<std::string>& vSubset);
75
78
81
84
87
90
91 void set_defect_upwind(bool defectUpwind) { m_bDefectUpwind = defectUpwind;}
93
97
99 void set_upwind(const std::string& name)
100 {m_spConvUpwind = CreateNavierStokesUpwind<dim>(name);}
101
102 public:
104 void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
105
107 virtual bool use_hanging() const;
108
110 virtual std::string disc_type() const {return "fvcr";};
111
112 protected:
113 void init();
114
115 public:
117 template <typename TFVGeom>
118 inline number peclet_blend(MathVector<dim>& UpwindVel, const TFVGeom& geo, size_t ip,
119 const MathVector<dim>& StdVel, number kinVisco);
120
121 protected:
124
127
130
134
137
140
143
145 static const size_t _P_ = dim;
146
152
153 public:
154 template <typename TElem, typename TFVGeom>
155 void prep_elem_loop(const ReferenceObjectID roid, const int si);
156
157 template <typename TElem, typename TFVGeom>
158 void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
159
160 template <typename TElem, typename TFVGeom>
161 void fsh_elem_loop();
162
163 template <typename TElem, typename TFVGeom>
164 void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
165
166 template <typename TElem, typename TFVGeom>
167 void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
168
169 template <typename TElem, typename TFVGeom>
170 void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
171
172 template <typename TElem, typename TFVGeom>
173 void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
174
175 template <typename TElem, typename TFVGeom>
176 void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
177
178 private:
179 template <typename TElem, typename TFVGeom>
180 void register_func();
181 void register_all_funcs(bool bHang);
182};
183
185
186} // namespace NavierStokes
187} // end namespace ug
188
189#endif /*__H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FVCR__NAVIER_STOKES_FVCR__*/
function NavierStokes(fcts, subsets, discType)
SmartPtr< CplUserData< TData, dim > > user_data()
Definition upwind_interface.h:57
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
number m_gradDivFactor
factor for div grad stabilization
Definition incompressible_navier_stokes_base.h:250
static const int dim
World dimension.
Definition incompressible_navier_stokes_base.h:156
Definition navier_stokes_fvcr.h:57
virtual std::string disc_type() const
returns string identifying disc type
Definition navier_stokes_fvcr.h:110
NavierStokesFVCR< TDomain > this_type
own type
Definition navier_stokes_fvcr.h:63
SmartPtr< CplUserData< number, dim > > density()
returns density
Definition navier_stokes_fvcr.h:86
void init()
Definition navier_stokes_fvcr.cpp:64
void add_jac_A_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition navier_stokes_fvcr.cpp:271
void register_func()
Definition navier_stokes_fvcr.cpp:829
void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition navier_stokes_fvcr.cpp:92
bool m_bDefectUpwind
flag if using upwind in defect computation
Definition navier_stokes_fvcr.h:123
void prep_elem_loop(const ReferenceObjectID roid, const int si)
Definition navier_stokes_fvcr.cpp:148
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition navier_stokes_fvcr.cpp:735
void set_kinematic_viscosity(SmartPtr< CplUserData< number, dim > > user)
sets the kinematic viscosity
Definition navier_stokes_fvcr.cpp:120
void set_upwind(const std::string &name)
sets the upwind based on a string identifier
Definition navier_stokes_fvcr.h:99
void set_defect_upwind(bool defectUpwind)
Definition navier_stokes_fvcr.h:91
DataImport< MathMatrix< dim, dim >, dim > m_imCentralGradient
Data import for central gradient.
Definition navier_stokes_fvcr.h:136
bool get_defect_upwind()
Definition navier_stokes_fvcr.h:92
DataImport< MathVector< dim >, dim > m_imSource
Data import for source.
Definition navier_stokes_fvcr.h:126
SmartPtr< CplUserData< number, dim > > kinematic_viscosity()
returns kinematic viscosity
Definition navier_stokes_fvcr.h:80
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition navier_stokes_fvcr.cpp:675
virtual bool use_hanging() const
returns if hanging nodes are needed
Definition navier_stokes_fvcr.cpp:113
void set_source(SmartPtr< CplUserData< MathVector< dim >, dim > > user)
sets the source function
Definition navier_stokes_fvcr.cpp:135
number peclet_blend(MathVector< dim > &UpwindVel, const TFVGeom &geo, size_t ip, const MathVector< dim > &StdVel, number kinVisco)
mixed upwind for Crouzeix-Raviart elements
Definition navier_stokes_fvcr.cpp:249
void fsh_elem_loop()
Definition navier_stokes_fvcr.cpp:201
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
Definition navier_stokes_fvcr.cpp:208
DataImport< number, dim > m_imKinViscosity
Data import for kinematic viscosity.
Definition navier_stokes_fvcr.h:129
void set_density(SmartPtr< CplUserData< number, dim > > user)
sets the density
Definition navier_stokes_fvcr.cpp:127
void add_def_A_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition navier_stokes_fvcr.cpp:509
static const size_t _P_
abbreviation for pressure
Definition navier_stokes_fvcr.h:145
void set_upwind(SmartPtr< INavierStokesUpwind< dim > > spUpwind)
sets an upwinding for the convective term of momentum equation
Definition navier_stokes_fvcr.h:95
IncompressibleNavierStokesBase< TDomain > base_type
Base class type.
Definition navier_stokes_fvcr.h:60
static const int dim
World dimension.
Definition navier_stokes_fvcr.h:67
SmartPtr< INavierStokesUpwind< dim > > m_spConvUpwind
Upwinding for velocity in convective term of momentum equation.
Definition navier_stokes_fvcr.h:142
DataImport< number, dim > m_imDensitySCV
Definition navier_stokes_fvcr.h:133
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition navier_stokes_fvcr.cpp:705
DataImport< number, dim > m_imDensitySCVF
Data import for density.
Definition navier_stokes_fvcr.h:132
DataImport< MathVector< dim >, dim > m_imPressureGradient
Data import for central gradient.
Definition navier_stokes_fvcr.h:139
double number
ReferenceObjectID