Plugins
symmetric_boundary_fv1.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2014: 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__FV1__BND__SYMMETRIC___
34 #define __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FV1__BND__SYMMETRIC___
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 "../stabilization.h"
47 
48 
49 namespace ug{
50 namespace NavierStokes{
51 
54 template< typename TDomain>
56  : public IElemDisc<TDomain>
57 {
58  private:
61 
64 
65  public:
68 
70  static const int dim = base_type::dim;
71 
74 
75  public:
78 
80  void add(const char* subsets);
81 
83 
89  {m_imKinViscosity.set_data(data);}
90 
92 
98  {m_imDensity.set_data(data);}
99 
100  public:
102  virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
103 
104  public:
106  virtual bool requests_local_time_series() {return true;}
107 
109  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
110  void prep_elem_loop(const ReferenceObjectID roid, const int si);
111 
113  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
114  void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
115 
117  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
118  void fsh_elem_loop();
119 
121  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
122  void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
123 
125  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
126  void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
127 
128  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
129  void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
130  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
131  void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
132  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
133  void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
134 
135  private:
138 
140  std::vector<std::string> m_vScheduledBndSubSets; // names
141  std::vector<int> m_vBndSubSetIndex; // indices
142 
143  void extract_scheduled_data(); // convert m_vScheduledBndSubSets -> m_vBndSubSetIndex
144 
150  std::vector<MathVector<dim> > m_vLocIP;
151  std::vector<MathVector<dim> > m_vGloIP;
152 
155 
157  static const size_t _P_ = dim;
158 
159  private:
160  void register_all_funcs(bool bHang);
161 
162  template <template <class Elem, int WorldDim> class TFVGeom>
163  struct RegisterFV1 {
164  RegisterFV1(this_type* pThis) : m_pThis(pThis){}
166  template< typename TElem > void operator()(TElem&)
167  {m_pThis->register_func<TElem, TFVGeom>();}
168  };
169 
170  template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
171  void register_func();
172 };
173 
174 } // namespace NavierStokes
175 } // end namespace ug
176 
177 #endif /* __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FV1__BND__SYMMETRIC___ */
function NavierStokes(fcts, subsets, discType)
void set_data(SmartPtr< CplUserData< number, dim > > spData)
TDomain::position_type position_type
static const int dim
Finite Volume Element Discretization for the incompressible Navier-Stokes Equation.
Definition: incompressible_navier_stokes_base.h:146
Definition: symmetric_boundary_fv1.h:57
base_type::domain_type domain_type
Domain type.
Definition: symmetric_boundary_fv1.h:67
virtual bool requests_local_time_series()
returns if local time series is needed
Definition: symmetric_boundary_fv1.h:106
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition: symmetric_boundary_fv1.cpp:46
void register_func()
Definition: symmetric_boundary_fv1.cpp:452
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition: symmetric_boundary_fv1.h:131
DataImport< MathVector< dim >, dim > m_imSourceSCVF
Data import for source.
Definition: symmetric_boundary_fv1.h:154
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition: symmetric_boundary_fv1.h:133
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the element loop
Definition: symmetric_boundary_fv1.cpp:130
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for evaluation
Definition: symmetric_boundary_fv1.cpp:185
NavierStokesSymBCFV1(SmartPtr< IncompressibleNavierStokesBase< TDomain > > spMaster)
Constructor (setting default values)
Definition: symmetric_boundary_fv1.cpp:404
DataImport< number, dim > m_imKinViscosity
Data import for kinematic viscosity.
Definition: symmetric_boundary_fv1.h:146
static const size_t _P_
abbreviation for pressure
Definition: symmetric_boundary_fv1.h:157
std::vector< MathVector< dim > > m_vGloIP
Definition: symmetric_boundary_fv1.h:151
void register_all_funcs(bool bHang)
Definition: symmetric_boundary_fv1.cpp:436
NavierStokesSymBCFV1< TDomain > this_type
own type
Definition: symmetric_boundary_fv1.h:63
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: symmetric_boundary_fv1.cpp:329
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: symmetric_boundary_fv1.cpp:227
void extract_scheduled_data()
Definition: symmetric_boundary_fv1.cpp:67
std::vector< std::string > m_vScheduledBndSubSets
The boundary subsets:
Definition: symmetric_boundary_fv1.h:140
std::vector< int > m_vBndSubSetIndex
Definition: symmetric_boundary_fv1.h:141
void fsh_elem_loop()
finishes the element loop
Definition: symmetric_boundary_fv1.cpp:167
base_type::position_type position_type
Position type.
Definition: symmetric_boundary_fv1.h:73
SmartPtr< IncompressibleNavierStokesBase< TDomain > > m_spMaster
The master discretization:
Definition: symmetric_boundary_fv1.h:137
static const int dim
World dimension.
Definition: symmetric_boundary_fv1.h:70
DataImport< number, dim > m_imDensity
Data import for density.
Definition: symmetric_boundary_fv1.h:148
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition: symmetric_boundary_fv1.h:129
void add(const char *subsets)
adds a boundary segment
Definition: symmetric_boundary_fv1.cpp:116
void set_kinematic_viscosity(SmartPtr< CplUserData< number, dim > > data)
sets the kinematic viscosity
Definition: symmetric_boundary_fv1.h:88
std::vector< MathVector< dim > > m_vLocIP
Boundary integration points of the viscosity and the density.
Definition: symmetric_boundary_fv1.h:150
void set_density(SmartPtr< CplUserData< number, dim > > data)
sets the density
Definition: symmetric_boundary_fv1.h:97
IElemDisc< TDomain > base_type
Base class type.
Definition: symmetric_boundary_fv1.h:60
ReferenceObjectID
Definition: symmetric_boundary_fv1.h:163
void operator()(TElem &)
Definition: symmetric_boundary_fv1.h:166
RegisterFV1(this_type *pThis)
Definition: symmetric_boundary_fv1.h:164
this_type * m_pThis
Definition: symmetric_boundary_fv1.h:165