ug4
conv_shape_interface.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2011-2015: 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__LIB_DISC__SPATIAL_DISC__ELEM_DISC__DENSITY_DRIVEN_FLOW__FV1__CONV_SHAPE_INTERFACE__
34 #define __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__DENSITY_DRIVEN_FLOW__FV1__CONV_SHAPE_INTERFACE__
35 
36 #include "common/common.h"
37 #include "common/util/provider.h"
40 
41 namespace ug{
42 
44 // Interface for Convection shapes
46 
48 
98 template <int dim>
100 {
101  public:
104 
106  typedef bool (this_type::*UpdateFunc)
107  (const FVGeometryBase* geo,
108  const MathVector<dim>* Velocity,
109  const MathMatrix<dim, dim>* Diffusion,
110  bool computeDeriv);
111 
112  public:
116  {
117  m_vUpdateFunc.clear();
118  m_vConvShape.clear();
119  m_vConvShapeVel.clear();
120  m_vConvShapeDiffusion.clear();
121  }
122 
124  virtual ~IConvectionShapes() {};
125 
127  size_t num_sh() const {return m_numSh;}
128 
130  size_t num_scvf() const {return m_numScvf;}
131 
133  number operator()(size_t scvf, size_t sh) const
134  {
135  UG_ASSERT(scvf < m_vConvShape.size(), "Invalid index");
136  UG_ASSERT(sh < m_vConvShape[scvf].size(), "Invalid index");
137  return m_vConvShape[scvf][sh];
138  }
139 
141  const MathVector<dim>& D_vel(size_t scvf, size_t sh) const
142  {
143  UG_ASSERT(scvf < m_vConvShapeVel.size(), "Invalid index");
144  UG_ASSERT(sh < m_vConvShapeVel[scvf].size(), "Invalid index");
145  return m_vConvShapeVel[scvf][sh];
146  }
147 
150 
152  const MathMatrix<dim,dim>& D_diffusion(size_t scvf, size_t sh) const
153  {
154  UG_ASSERT(scvf < m_vConvShapeDiffusion.size(), "Invalid index");
155  UG_ASSERT(sh < m_vConvShapeDiffusion[scvf].size(), "Invalid index");
156  return m_vConvShapeDiffusion[scvf][sh];
157  }
158 
159  bool update(const FVGeometryBase* geo,
160  const MathVector<dim>* Velocity,
161  const MathMatrix<dim, dim>* DiffDisp,
162  bool computeDeriv)
163  {return (this->*(m_vUpdateFunc[m_id])) (geo, Velocity, DiffDisp, computeDeriv);}
164 
166  // internal handling
168 
169  protected:
171  void set_sizes(size_t numScvf, size_t numSh);
172 
175 
177  number& conv_shape(size_t scvf, size_t sh)
178  {
179  UG_ASSERT(scvf < m_vConvShape.size(), "Invalid index");
180  UG_ASSERT(sh < m_vConvShape[scvf].size(), "Invalid index");
181  return m_vConvShape[scvf][sh];
182  }
183 
185  MathVector<dim>& D_vel(size_t scvf, size_t sh)
186  {
187  UG_ASSERT(scvf < m_vConvShapeVel.size(), "Invalid index");
188  UG_ASSERT(sh < m_vConvShapeVel[scvf].size(), "Invalid index");
189  return m_vConvShapeVel[scvf][sh];
190  }
191 
193  MathMatrix<dim,dim>& conv_shape_diffusion(size_t scvf, size_t sh)
194  {
195  UG_ASSERT(scvf < m_vConvShapeDiffusion.size(), "Invalid index");
196  UG_ASSERT(sh < m_vConvShapeDiffusion[scvf].size(), "Invalid index");
197  return m_vConvShapeDiffusion[scvf][sh];
198  }
199 
200 
202  size_t m_numScvf;
203 
205  size_t m_numSh;
206 
208  std::vector<std::vector<number> > m_vConvShape;
209 
211  std::vector<std::vector<MathVector<dim> > > m_vConvShapeVel;
212 
214  std::vector<std::vector<MathMatrix<dim,dim> > > m_vConvShapeDiffusion;
215 
218 
220  // registering process
222 
223  public:
225  template <typename TFVGeom, typename TAssFunc>
226  void register_update_func(TAssFunc func);
227 
229  template <typename TFVGeom>
230  bool set_geometry_type(const TFVGeom& geo);
231 
232  protected:
234  std::vector<UpdateFunc> m_vUpdateFunc;
235 
237  int m_id;
238 };
239 
241 // Interface for Stabilization
243 
244 // register a update function for a Geometry
245 template <int dim>
246 template <typename TFVGeom, typename TAssFunc>
247 void
249 register_update_func(TAssFunc func)
250 {
251 // get unique geometry id
252  size_t id = GetUniqueFVGeomID<TFVGeom>();
253 
254 // make sure that there is enough space
255  if((size_t)id >= m_vUpdateFunc.size())
256  m_vUpdateFunc.resize(id+1, NULL);
257 
258 // set pointer
259  m_vUpdateFunc[id] = (UpdateFunc)func;
260 }
261 
262 // set the Geometry type to use for next updates
263 template <int dim>
264 template <typename TFVGeom>
265 bool
267 set_geometry_type(const TFVGeom& geo)
268 {
269 // get unique geometry id
270  size_t id = GetUniqueFVGeomID<TFVGeom>();
271 
272 // check that function exists
273  if(id >= m_vUpdateFunc.size() || m_vUpdateFunc[id] == NULL)
274  {
275  UG_LOG("ERROR in 'IConvectionShapes::set_geometry_type':"
276  " No update function registered for this Geometry.\n");
277  return false;
278  }
279 
280 // set current geometry
281  m_id = id;
282 
283 // set sizes
284  set_sizes(geo.num_scvf(), geo.num_sh());
285 
286 // we're done
287  return true;
288 }
289 
290 
291 // resize the data arrays
292 template <int dim>
293 void
295 set_sizes(size_t numScvf, size_t numSh)
296 {
297 // remember sizes
298  m_numScvf = numScvf;
299  m_numSh = numSh;
300 
301 // adjust arrays
302  m_vConvShape.resize(m_numScvf);
303  m_vConvShapeVel.resize(m_numScvf);
304  m_vConvShapeDiffusion.resize(m_numScvf);
305  for(size_t scvf = 0; scvf < m_numScvf; ++scvf)
306  {
307  m_vConvShape[scvf].resize(m_numSh);
308  m_vConvShapeVel[scvf].resize(m_numSh);
309  m_vConvShapeDiffusion[scvf].resize(m_numSh);
310  }
311 }
312 
313 
314 
315 } // end namespace ug
316 
317 #endif /* __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__DENSITY_DRIVEN_FLOW__FV1__CONV_SHAPE_INTERFACE__ */
base class for all FVGeometries
Definition: fv_geom_base.h:76
Interface class for upwind methods.
Definition: conv_shape_interface.h:100
std::vector< std::vector< MathMatrix< dim, dim > > > m_vConvShapeDiffusion
upwind shapes for corners shape functions
Definition: conv_shape_interface.h:214
std::vector< UpdateFunc > m_vUpdateFunc
Vector holding all update functions.
Definition: conv_shape_interface.h:234
IConvectionShapes< dim > this_type
Abbreviation for own type.
Definition: conv_shape_interface.h:103
bool(this_type::* UpdateFunc)(const FVGeometryBase *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *Diffusion, bool computeDeriv)
type of update function
Definition: conv_shape_interface.h:107
const MathVector< dim > & D_vel(size_t scvf, size_t sh) const
upwind shape for corner vel
Definition: conv_shape_interface.h:141
number & conv_shape(size_t scvf, size_t sh)
non-const access to ip velocity (i.e. interpolated velocity at ip)
Definition: conv_shape_interface.h:177
void set_sizes(size_t numScvf, size_t numSh)
resize the data arrays
Definition: conv_shape_interface.h:295
std::vector< std::vector< number > > m_vConvShape
upwind shapes for corners shape functions
Definition: conv_shape_interface.h:208
IConvectionShapes()
constructor
Definition: conv_shape_interface.h:114
std::vector< std::vector< MathVector< dim > > > m_vConvShapeVel
upwind shapes for corners shape functions
Definition: conv_shape_interface.h:211
void register_update_func(TAssFunc func)
register a update function for a Geometry
Definition: conv_shape_interface.h:249
bool update(const FVGeometryBase *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
Definition: conv_shape_interface.h:159
bool set_geometry_type(const TFVGeom &geo)
set the Geometry type to use for next updates
Definition: conv_shape_interface.h:267
void set_non_zero_deriv_diffusion_flag(bool flag)
sets the shape ip flag
Definition: conv_shape_interface.h:174
virtual ~IConvectionShapes()
destructor
Definition: conv_shape_interface.h:124
bool m_bNonZeroDerivDiffusion
flag if ip shapes are non-zero
Definition: conv_shape_interface.h:217
size_t num_scvf() const
returns number of sub control volume faces
Definition: conv_shape_interface.h:130
int m_id
id of current geometry type
Definition: conv_shape_interface.h:237
size_t m_numSh
number of current shape functions (usually in corners)
Definition: conv_shape_interface.h:205
MathVector< dim > & D_vel(size_t scvf, size_t sh)
non-const access to upwind shapes for corner vel
Definition: conv_shape_interface.h:185
size_t num_sh() const
returns number of shapes
Definition: conv_shape_interface.h:127
bool non_zero_deriv_diffusion() const
returns if upwind shape w.r.t. ip vel is non-zero
Definition: conv_shape_interface.h:149
number operator()(size_t scvf, size_t sh) const
shape value
Definition: conv_shape_interface.h:133
size_t m_numScvf
number of current scvf
Definition: conv_shape_interface.h:202
MathMatrix< dim, dim > & conv_shape_diffusion(size_t scvf, size_t sh)
non-const access to upwind shapes for ip vel
Definition: conv_shape_interface.h:193
const MathMatrix< dim, dim > & D_diffusion(size_t scvf, size_t sh) const
upwind shapes for ip vel
Definition: conv_shape_interface.h:152
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_LOG(msg)
Definition: log.h:367
double number
Definition: types.h:124
the ug namespace
function func(x, y, z, t, si)
typedefs for ugmath