Plugins
navier_stokes_fv1.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__FV1__NAVIER_STOKES_FV1__
34 #define __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FV1__NAVIER_STOKES_FV1__
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 #include "stabilization.h"
48 
49 namespace ug{
50 namespace NavierStokes{
51 
54 
56 
145 template< typename TDomain>
147  : public IncompressibleNavierStokesBase<TDomain>
148 {
149  protected:
152 
155 
156  public:
158  static const int dim = base_type::dim;
159 
160  public:
163  NavierStokesFV1(const char* functions, const char* subsets);
164  NavierStokesFV1(const std::vector<std::string>& vFct, const std::vector<std::string>& vSubset);
166 
169 
172 
175 
178 
181 
184 
187  {m_spStab = spStab;}
188 
190  void set_stabilization(const std::string& name)
191  {m_spStab = CreateNavierStokesStabilization<dim>(name);
192  if(m_spConvUpwind.valid()) m_spStab->set_upwind(m_spConvUpwind);}
193 
195  void set_stabilization(const std::string& name, const std::string& diffLength)
196  {SmartPtr<INavierStokesSRFV1Stabilization<dim> > spStab = CreateNavierStokesStabilization<dim>(name);
197  spStab->set_diffusion_length(diffLength);
198  m_spStab = spStab;
199  if(m_spConvUpwind.valid()) m_spStab->set_upwind(m_spConvUpwind);}
200 
203 
206  {m_spConvStab = spStab; m_spConvUpwind = SPNULL;}
207 
210  {m_spConvStab = SPNULL; m_spConvUpwind = spUpwind;}
211 
213  void set_upwind(const std::string& name)
214  {m_spConvStab = SPNULL; m_spConvUpwind = CreateNavierStokesUpwind<dim>(name);
215  if(m_spStab.valid() && m_spStab->upwind().invalid()) m_spStab->set_upwind(m_spConvUpwind);}
216 
217  void set_pac_upwind(bool bPac)
218  {
219  if (bPac==true){
220  if (m_spConvUpwind.valid()) m_spStab->set_upwind(m_spConvUpwind);
221  else UG_THROW("Upwind must be specified previously.\n");
222  if (m_spStab.valid()) set_upwind(m_spStab);
223  else UG_THROW("Stabilization must be specified previously.\n");
224  }
225  }
226 
227  protected:
228 
231 
232  public:
234  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
235 
237  virtual std::string disc_type() const {return "fv1";};
238 
239  public:
241 
250  template <typename TElem, typename TFVGeom>
251  void prep_elem_loop(const ReferenceObjectID roid, const int si);
252 
254 
265  template <typename TElem, typename TFVGeom>
266  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
267 
269  template <typename TElem, typename TFVGeom>
270  void fsh_elem_loop();
271 
273 
336  template <typename TElem, typename TFVGeom>
337  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
338 
340 
406  template <typename TElem, typename TFVGeom>
407  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
408 
410 
432  template <typename TElem, typename TFVGeom>
433  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
434 
436 
460  template <typename TElem, typename TFVGeom>
461  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
462 
464 
483  template <typename TElem, typename TFVGeom>
484  void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
485 
487 
518  template <typename TFVGeom>
519  inline number peclet_blend(MathVector<dim>& UpwindVel, const TFVGeom& geo, size_t ip,
520  const MathVector<dim>& StdVel, number kinVisco);
521 
523  template <typename TElem, typename TFVGeom>
524  void ex_nodal_velocity(MathVector<dim> vValue[],
525  const MathVector<dim> vGlobIP[],
526  number time, int si,
527  const LocalVector& u,
528  GridObject* elem,
529  const MathVector<dim> vCornerCoords[],
530  const MathVector<TFVGeom::dim> vLocIP[],
531  const size_t nip,
532  bool bDeriv,
533  std::vector<std::vector<MathVector<dim> > > vvvDeriv[]);
534 
536  template <typename TElem, typename TFVGeom>
538  const MathVector<dim> vGlobIP[],
539  number time, int si,
540  const LocalVector& u,
541  GridObject* elem,
542  const MathVector<dim> vCornerCoords[],
543  const MathVector<TFVGeom::dim> vLocIP[],
544  const size_t nip,
545  bool bDeriv,
546  std::vector<std::vector<MathMatrix<dim,dim> > > vvvDeriv[]);
547 
548  protected:
552 
555 
559 
562 
566 
569 
571  static const size_t _P_ = dim;
572 
575  using base_type::m_bStokes;
576  using base_type::m_bLaplace;
577 
578  virtual void init();
579 
580  protected:
583  virtual void register_all_funcs(bool bHang);
584  template <typename TElem, typename TFVGeom> void register_func();
586 };
587 
589 
590 } // namespace NavierStokes
591 } // end namespace ug
592 
593 #endif /*__H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FV1__NAVIER_STOKES_FV1__*/
function NavierStokes(fcts, subsets, discType)
SmartPtr< CplUserData< number, dim > > user_data()
number time() const
Definition: upwind_interface.h:57
Finite Volume Element Discretization for the incompressible Navier-Stokes Equation.
Definition: incompressible_navier_stokes_base.h:146
SmartPtr< DataExport< MathMatrix< dim, dim >, dim > > m_exVelocityGrad
Export for the velocity gradient.
Definition: incompressible_navier_stokes_base.h:262
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
SmartPtr< DataExport< MathVector< dim >, dim > > m_exVelocity
Export for the velocity.
Definition: incompressible_navier_stokes_base.h:259
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 Volume Element Discretization for the incompressible Navier-Stokes Equation.
Definition: navier_stokes_fv1.h:148
void ex_nodal_velocity(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[])
export value of the velocity
Definition: navier_stokes_fv1.cpp:898
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the element loop
Definition: navier_stokes_fv1.cpp:139
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for evaluation
Definition: navier_stokes_fv1.cpp:211
void set_upwind(const std::string &name)
sets the upwind based on a string identifier
Definition: navier_stokes_fv1.h:213
SmartPtr< CplUserData< number, dim > > kinematic_viscosity()
returns kinematic viscosity
Definition: navier_stokes_fv1.h:171
virtual void init()
Definition: navier_stokes_fv1.cpp:63
void fsh_elem_loop()
finishes the element loop
Definition: navier_stokes_fv1.cpp:204
number peclet_blend(MathVector< dim > &UpwindVel, const TFVGeom &geo, size_t ip, const MathVector< dim > &StdVel, number kinVisco)
computes the pecled blended Upwind veloctity
Definition: navier_stokes_fv1.cpp:876
DataImport< MathVector< dim >, dim > m_imSourceSCVF
Definition: navier_stokes_fv1.h:551
DataImport< MathVector< dim >, dim > m_imSourceSCV
Data import for source.
Definition: navier_stokes_fv1.h:550
IncompressibleNavierStokesBase< TDomain > base_type
Base class type.
Definition: navier_stokes_fv1.h:151
DataImport< MathVector< dim >, dim > sourceSCVF()
returns scvf source
Definition: navier_stokes_fv1.h:183
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition: navier_stokes_fv1.cpp:91
void set_stabilization(const std::string &name)
sets stabilization based on string identifier
Definition: navier_stokes_fv1.h:190
DataImport< number, dim > m_imKinViscosity
Data import for kinematic viscosity.
Definition: navier_stokes_fv1.h:554
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: navier_stokes_fv1.cpp:600
void register_func()
Definition: navier_stokes_fv1.cpp:1164
NavierStokesFV1< TDomain > this_type
own type
Definition: navier_stokes_fv1.h:154
void set_source(SmartPtr< CplUserData< MathVector< dim >, dim > > user)
sets the source function
Definition: navier_stokes_fv1.cpp:126
void set_upwind(SmartPtr< INavierStokesFV1Stabilization< dim > > spStab)
sets a stabilization for upwinding (Physical Advection Correction)
Definition: navier_stokes_fv1.h:205
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
adds the source part to the local defect
Definition: navier_stokes_fv1.cpp:844
SmartPtr< INavierStokesFV1Stabilization< dim > > m_spConvStab
Definition: navier_stokes_fv1.h:565
SmartPtr< CplUserData< number, dim > > density()
returns density
Definition: navier_stokes_fv1.h:177
virtual void register_all_funcs(bool bHang)
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: navier_stokes_fv1.cpp:253
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: navier_stokes_fv1.cpp:814
DataImport< number, dim > m_imDensitySCV
Definition: navier_stokes_fv1.h:558
NavierStokesFV1(const char *functions, const char *subsets)
Definition: navier_stokes_fv1.cpp:46
void set_stabilization(SmartPtr< INavierStokesFV1Stabilization< dim > > spStab)
sets the stabilization used to compute the stabilized velocities
Definition: navier_stokes_fv1.h:186
static const int dim
World dimension.
Definition: navier_stokes_fv1.h:158
void set_stabilization(const std::string &name, const std::string &diffLength)
sets stabilization and diff length method based on string identifier (only for Schneider-Raw stabiliz...
Definition: navier_stokes_fv1.h:195
SmartPtr< INavierStokesUpwind< dim > > m_spConvUpwind
Upwinding for velocity in convective term of momentum equation.
Definition: navier_stokes_fv1.h:568
void set_pac_upwind(bool bPac)
Definition: navier_stokes_fv1.h:217
virtual std::string disc_type() const
returns string identifying disc type
Definition: navier_stokes_fv1.h:237
static const size_t _P_
abbreviation for pressure
Definition: navier_stokes_fv1.h:571
void ex_velocity_grad(MathMatrix< dim, 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< MathMatrix< dim, dim > > > vvvDeriv[])
computes the value of the gradient of the velocity
Definition: navier_stokes_fv1.cpp:997
void set_density(SmartPtr< CplUserData< number, dim > > user)
sets the density
Definition: navier_stokes_fv1.cpp:118
SmartPtr< INavierStokesFV1Stabilization< dim > > stabilization()
returns stabilization
Definition: navier_stokes_fv1.h:202
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: navier_stokes_fv1.cpp:784
void set_kinematic_viscosity(SmartPtr< CplUserData< number, dim > > user)
sets the kinematic viscosity
Definition: navier_stokes_fv1.cpp:111
DataImport< number, dim > m_imDensitySCVF
Data import for density.
Definition: navier_stokes_fv1.h:557
void set_upwind(SmartPtr< INavierStokesUpwind< dim > > spUpwind)
sets an upwinding for the convective term of momentum equation
Definition: navier_stokes_fv1.h:209
SmartPtr< INavierStokesFV1Stabilization< dim > > m_spStab
Stabilization for velocity in continuity equation.
Definition: navier_stokes_fv1.h:561
const NullSmartPtr SPNULL
#define UG_THROW(msg)
double number
ReferenceObjectID