Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
symmetric_boundary_fvcr.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2012-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__FVCR__BND__SYMMETRIC__
34#define __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FVCR__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
47namespace ug{
48namespace NavierStokes{
49
52template< typename TDomain>
54 : public IElemDisc<TDomain>
55{
56 private:
59
62
63 public:
66
68 static const int dim = base_type::dim;
69
72
73 public:
76
78 void add(const char* subsets);
79
81
88
90
97
98 public:
100 virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
101
102 public:
104 virtual bool requests_local_time_series() {return true;}
105
107 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
108 void prep_elem_loop(const ReferenceObjectID roid, const int si);
109
111 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
112 void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
113
115 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
116 void fsh_elem_loop();
117
119 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
120 void add_jac_A_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
121
123 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
124 void add_def_A_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
125
126 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
127 void add_jac_M_elem(LocalMatrix& J, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
128 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
129 void add_def_M_elem(LocalVector& d, const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
130 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
131 void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]) {}
132
133 private:
136
138 std::vector<std::string> m_vScheduledBndSubSets; // names
139 std::vector<int> m_vBndSubSetIndex; // indices
140
141 void extract_scheduled_data(); // convert m_vScheduledBndSubSets -> m_vBndSubSetIndex
142
148 std::vector<MathVector<dim> > m_vLocIP;
149 std::vector<MathVector<dim> > m_vGloIP;
150
152 static const size_t _P_ = dim;
153
154 private:
155 void register_all_funcs(bool bHang);
156
157 template <typename TElem, template <class Elem, int WorldDim> class TFVGeom>
158 void register_func();
159};
160
161} // namespace NavierStokes
162} // end namespace ug
163
164#endif /* __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FVCR__BND__SYMMETRIC__ */
function NavierStokes(fcts, subsets, discType)
void set_data(SmartPtr< CplUserData< TData, dim > > spData)
TDomain::position_type position_type
static const int dim
Definition symmetric_boundary_fvcr.h:55
base_type::domain_type domain_type
Domain type.
Definition symmetric_boundary_fvcr.h:65
std::vector< int > m_vBndSubSetIndex
Definition symmetric_boundary_fvcr.h:139
std::vector< MathVector< dim > > m_vLocIP
Boundary integration points of the viscosity and the density.
Definition symmetric_boundary_fvcr.h:148
void set_density(SmartPtr< CplUserData< number, dim > > data)
sets the density
Definition symmetric_boundary_fvcr.h:95
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_fvcr.cpp:223
virtual bool requests_local_time_series()
returns if local time series is needed
Definition symmetric_boundary_fvcr.h:104
std::vector< MathVector< dim > > m_vGloIP
Definition symmetric_boundary_fvcr.h:149
void extract_scheduled_data()
Definition symmetric_boundary_fvcr.cpp:63
static const size_t _P_
abbreviation for pressure
Definition symmetric_boundary_fvcr.h:152
std::vector< std::string > m_vScheduledBndSubSets
The boundary subsets:
Definition symmetric_boundary_fvcr.h:138
CRNavierStokesSymBC< TDomain > this_type
own type
Definition symmetric_boundary_fvcr.h:61
DataImport< number, dim > m_imDensity
Data import for density.
Definition symmetric_boundary_fvcr.h:146
void add_jac_M_elem(LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition symmetric_boundary_fvcr.h:127
base_type::position_type position_type
Position type.
Definition symmetric_boundary_fvcr.h:71
IElemDisc< TDomain > base_type
Base class type.
Definition symmetric_boundary_fvcr.h:58
DataImport< number, dim > m_imKinViscosity
Data import for kinematic viscosity.
Definition symmetric_boundary_fvcr.h:144
void fsh_elem_loop()
finishes the element loop
Definition symmetric_boundary_fvcr.cpp:163
void prep_elem_loop(const ReferenceObjectID roid, const int si)
prepares the element loop
Definition symmetric_boundary_fvcr.cpp:126
void add_def_M_elem(LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition symmetric_boundary_fvcr.h:129
static const int dim
World dimension.
Definition symmetric_boundary_fvcr.h:68
SmartPtr< IncompressibleNavierStokesBase< TDomain > > m_spMaster
The master discretization:
Definition symmetric_boundary_fvcr.h:135
void set_kinematic_viscosity(SmartPtr< CplUserData< number, dim > > data)
sets the kinematic viscosity
Definition symmetric_boundary_fvcr.h:86
void add(const char *subsets)
adds a boundary segment
Definition symmetric_boundary_fvcr.cpp:112
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_fvcr.cpp:295
void register_func()
Definition symmetric_boundary_fvcr.cpp:458
void add_rhs_elem(LocalVector &d, GridObject *elem, const MathVector< dim > vCornerCoords[])
Definition symmetric_boundary_fvcr.h:131
void prep_elem(const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
prepares the element for evaluation
Definition symmetric_boundary_fvcr.cpp:181
virtual void prepare_setting(const std::vector< LFEID > &vLfeID, bool bNonRegularGrid)
type of trial space for each function used
Definition symmetric_boundary_fvcr.cpp:46
Finite Volume Element Discretization for the incompressible Navier-Stokes Equation.
Definition incompressible_navier_stokes_base.h:146
ReferenceObjectID