Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
compressible_navier_stokes_fv1.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2014: G-CSC, Goethe University Frankfurt
3 * Author: Raphael Prohl
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/*
34 * main parts are copied from the discretization of the incompressible Navier-Stokes Equations
35 * of Andreas Vogel and Christian Wehner
36 */
37
38#ifndef __H__UG__PLUGINS__NAVIER_STOKES__COMPRESSIBLE__COMPRESSIBLE_NAVIER_STOKES_FV1_H__
39#define __H__UG__PLUGINS__NAVIER_STOKES__COMPRESSIBLE__COMPRESSIBLE_NAVIER_STOKES_FV1_H__
40
41// other ug4 modules
42#include "common/common.h"
43#include "lib_grid/lg_base.h"
44
45// library intern headers
49
50#include "../compressible_navier_stokes_base.h"
51#include "../../upwind_interface.h"
52
53namespace ug{
54namespace NavierStokes{
55
58
60
129template< typename TDomain>
131 : public CompressibleNavierStokesBase<TDomain>
132{
133 protected:
136
139
140 public:
142 static const int dim = base_type::dim;
143
144 public:
147 CompressibleNavierStokesFV1(const char* functions, const char* subsets);
148 CompressibleNavierStokesFV1(const std::vector<std::string>& vFct, const std::vector<std::string>& vSubset);
150
153
156
159
162
165
168
171 {m_spConvUpwind = spUpwind;}
172
174 void set_upwind(const std::string& name)
175 {m_spConvUpwind = CreateNavierStokesUpwind<dim>(name);}
176
177 public:
179 void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
180
182 virtual std::string disc_type() const {return "fv1";};
183
184 public:
186
190 template<typename TElem, typename TFVGeom>
191 void prep_timestep_elem(const number time, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
192
194
203 template <typename TElem, typename TFVGeom>
204 void prep_elem_loop(const ReferenceObjectID roid, const int si);
205
207
218 template <typename TElem, typename TFVGeom>
219 void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
220
222 template <typename TElem, typename TFVGeom>
223 void fsh_elem_loop();
224
226
289 template <typename TElem, typename TFVGeom>
290 void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
291
293
359 template <typename TElem, typename TFVGeom>
360 void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
361
363
385 template <typename TElem, typename TFVGeom>
386 void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
387
389
413 template <typename TElem, typename TFVGeom>
414 void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
415
417
436 template <typename TElem, typename TFVGeom>
437 void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
438
440
471 template <typename TFVGeom>
472 inline number mach_number_blending(MathVector<dim>& UpwindVel, const TFVGeom& geo, size_t ip,
473 const MathVector<dim>& StdVel, number kinVisco);
474
475 protected:
479
482
485
488
490 static const size_t _P_ = dim;
492 static const size_t _Rho_ = dim + 1;
493
496
497 void init();
498
504
507
508 private:
511 void register_all_funcs(bool bHang);
512 template <typename TElem, typename TFVGeom> void register_func();
514};
515
517
518} // namespace NavierStokes
519} // end namespace ug
520
521#endif /*__H__UG__PLUGINS__NAVIER_STOKES__COMPRESSIBLE__COMPRESSIBLE_NAVIER_STOKES_FV1_H__*/
function NavierStokes(fcts, subsets, discType)
SmartPtr< CplUserData< TData, dim > > user_data()
number time() const
Finite Volume Element Discretization for the compressible Navier-Stokes Equation.
Definition compressible_navier_stokes_base.h:126
number m_bFullNewtonFactor
factor for exact jacobian, (1 for exact jacobian, 0 for fix point)
Definition navier_stokes_base.h:207
static const int dim
World dimension.
Definition compressible_navier_stokes_base.h:136
bool m_bMachNrBlend
flag if using Mach-number Blending
Definition compressible_navier_stokes_base.h:187
Finite Volume Element Discretization for the compressible Navier-Stokes Equation.
Definition compressible_navier_stokes_fv1.h:132
virtual std::string disc_type() const
returns string identifying disc type
Definition compressible_navier_stokes_fv1.h:182
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
adds the source part to the local defect
Definition compressible_navier_stokes_fv1.cpp:837
void set_upwind(const std::string &name)
sets the upwind based on a string identifier
Definition compressible_navier_stokes_fv1.h:174
void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition compressible_navier_stokes_fv1.cpp:93
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the element loop
Definition compressible_navier_stokes_fv1.cpp:177
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
adds the mass part to the local defect
Definition compressible_navier_stokes_fv1.cpp:794
SmartPtr< CplUserData< number, dim > > kinematic_viscosity()
returns kinematic viscosity
Definition compressible_navier_stokes_fv1.h:155
CompressibleNavierStokesBase< TDomain > base_type
Base class type.
Definition compressible_navier_stokes_fv1.h:135
void register_func()
Definition compressible_navier_stokes_fv1.cpp:955
DataImport< MathVector< dim >, dim > sourceSCVF()
returns scvf source
Definition compressible_navier_stokes_fv1.h:167
static const size_t _P_
abbreviation for pressure
Definition compressible_navier_stokes_fv1.h:490
static const int dim
World dimension.
Definition compressible_navier_stokes_fv1.h:142
SmartPtr< INavierStokesUpwind< dim > > m_spConvUpwind
Upwinding for velocity in convective term of momentum equation.
Definition compressible_navier_stokes_fv1.h:487
number m_maxPressure
Definition compressible_navier_stokes_fv1.h:501
static const size_t _Rho_
abbreviation for density
Definition compressible_navier_stokes_fv1.h:492
CompressibleNavierStokesFV1< TDomain > this_type
own type
Definition compressible_navier_stokes_fv1.h:138
void set_upwind(SmartPtr< INavierStokesUpwind< dim > > spUpwind)
sets an upwinding for the convective term of momentum equation
Definition compressible_navier_stokes_fv1.h:170
void add_jac_A_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
adds the stiffness part to the local jacobian
Definition compressible_navier_stokes_fv1.cpp:279
DataImport< MathVector< dim >, dim > m_imSourceSCVF
Definition compressible_navier_stokes_fv1.h:478
void set_adiabatic_index(SmartPtr< CplUserData< number, dim > > user)
sets the adiabatic index
Definition compressible_navier_stokes_fv1.cpp:120
SmartPtr< CplUserData< number, dim > > adiabatic_index()
returns adiabatic index
Definition compressible_navier_stokes_fv1.h:161
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for evaluation
Definition compressible_navier_stokes_fv1.cpp:240
void add_def_A_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
adds the stiffness part to the local defect
Definition compressible_navier_stokes_fv1.cpp:552
number m_maxDensity
Definition compressible_navier_stokes_fv1.h:502
number m_refMachNrSq
Definition compressible_navier_stokes_fv1.h:505
void prep_timestep_elem(const number time, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
prepares the element loop
Definition compressible_navier_stokes_fv1.cpp:142
void set_source(SmartPtr< CplUserData< MathVector< dim >, dim > > user)
sets the source function
Definition compressible_navier_stokes_fv1.cpp:127
number m_maxVel
reference quantities for getting dimensionless equations
Definition compressible_navier_stokes_fv1.h:500
DataImport< number, dim > m_imAdiabaticIndex
Data import for adiabatic index.
Definition compressible_navier_stokes_fv1.h:484
void fsh_elem_loop()
finishes the element loop
Definition compressible_navier_stokes_fv1.cpp:233
void set_kinematic_viscosity(SmartPtr< CplUserData< number, dim > > user)
sets the kinematic viscosity
Definition compressible_navier_stokes_fv1.cpp:113
DataImport< MathVector< dim >, dim > m_imSourceSCV
Data import for source.
Definition compressible_navier_stokes_fv1.h:477
void init()
Definition compressible_navier_stokes_fv1.cpp:68
DataImport< number, dim > m_imKinViscosity
Data import for kinematic viscosity.
Definition compressible_navier_stokes_fv1.h:481
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
adds the mass part to the local jacobian
Definition compressible_navier_stokes_fv1.cpp:760
number mach_number_blending(MathVector< dim > &UpwindVel, const TFVGeom &geo, size_t ip, const MathVector< dim > &StdVel, number kinVisco)
computes the Mach-number blended Upwind velocity
Definition compressible_navier_stokes_fv1.cpp:871
number m_numeratorOfRefReynoldsNr
Definition compressible_navier_stokes_fv1.h:506
number m_refLength
Definition compressible_navier_stokes_fv1.h:503
Definition upwind_interface.h:57
double number
ReferenceObjectID