33 #ifndef __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__DENSITY_DRIVEN_FLOW__FV1__CONV_SHAPE__
34 #define __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__DENSITY_DRIVEN_FLOW__FV1__CONV_SHAPE__
36 #include <boost/mpl/range_c.hpp>
37 #include <boost/mpl/for_each.hpp>
62 static const int dim = TDim;
82 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
87 template <
typename TFVGeom>
88 bool update(
const TFVGeom* geo,
102 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1Geometry<TElem, dim> >();
103 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1CondensedGeometry<TElem, dim> >();
104 m_pThis->template register_func_for_elem_fvgeom< TElem, HFV1Geometry<TElem, dim> >();
109 template <
typename TElem,
typename TFVGeom>
113 base_type::template register_update_func<TFVGeom, TFunc>(&this_type::template update<TFVGeom>);
125 template <
int refDim>
130 base_type::template register_update_func<TGeom, TFunc>(&this_type::template update<TGeom>);
135 template <
typename TFVGeom>
138 update(
const TFVGeom* geo,
144 UG_ASSERT(Velocity != NULL,
"Null pointer");
147 const size_t numSH = geo->num_sh();
150 for(
size_t ip = 0; ip < geo->num_scvf(); ++ip)
153 const typename TFVGeom::SCVF& scvf = geo->scvf(ip);
156 const number flux =
VecDot(scvf.normal(), Velocity[ip]);
159 for(
size_t sh = 0; sh < scvf.num_sh(); sh++)
160 conv_shape(ip, sh) = flux * scvf.shape(sh);
167 for(
size_t sh = scvf.num_sh(); sh < numSH; sh++)
168 conv_shape(ip, sh) = 0.0;
172 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
173 VecScale(D_vel(ip, sh), scvf.normal(), scvf.shape(sh));
176 for(
size_t sh = scvf.num_sh(); sh < numSH; sh++)
177 VecSet(D_vel(ip, sh), 0.0);
203 static const int dim = TDim;
223 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
228 template <
typename TFVGeom>
229 bool update(
const TFVGeom* geo,
243 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1Geometry<TElem, dim> >();
244 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1CondensedGeometry<TElem, dim> >();
245 m_pThis->template register_func_for_elem_fvgeom< TElem, HFV1Geometry<TElem, dim> >();
250 template <
typename TElem,
typename TFVGeom>
254 base_type::template register_update_func<TFVGeom, TFunc>(&this_type::template update<TFVGeom>);
266 template <
int refDim>
271 base_type::template register_update_func<TGeom, TFunc>(&this_type::template update<TGeom>);
276 template <
typename TFVGeom>
278 update(
const TFVGeom* geo,
284 UG_ASSERT(Velocity != NULL,
"Null pointer");
287 const size_t numSH = geo->num_sh();
290 for(
size_t ip = 0; ip < geo->num_scvf(); ++ip)
293 const typename TFVGeom::SCVF& scvf = geo->scvf(ip);
296 const number flux =
VecDot(scvf.normal(), Velocity[ip]);
298 size_t from = scvf.from();
299 size_t to = scvf.to();
302 if((from >= scvf.num_sh()) || (to >= scvf.num_sh())){
305 for(
size_t sh = 0; sh < scvf.num_sh(); sh++)
306 conv_shape(ip, sh) = flux * scvf.shape(sh);
308 UG_ASSERT(scvf.num_sh() == numSH,
"sh's have to match!");
312 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
313 VecScale(D_vel(ip, sh), scvf.normal(), scvf.shape(sh));
321 size_t up = (flux >= 0) ? from : to;
324 for(
size_t sh = 0; sh < scvf.num_sh(); ++sh) conv_shape(ip, sh) = 0.0;
326 conv_shape(ip, up) = flux;
331 for(
size_t sh = 0; sh < numSH; ++sh)
VecSet(D_vel(ip, sh), 0.0);
332 D_vel(ip, up) = scvf.normal();
358 static const int dim = TDim;
378 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
392 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
400 template <
typename TFVGeom>
401 bool update(
const TFVGeom* geo,
419 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1Geometry<TElem, dim> >();
420 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1CondensedGeometry<TElem, dim> >();
421 m_pThis->template register_func_for_elem_fvgeom< TElem, HFV1Geometry<TElem, dim> >();
426 template <
typename TElem,
typename TFVGeom>
430 base_type::template register_update_func<TFVGeom, TFunc>(&this_type::template update<TFVGeom>);
442 template <
int refDim>
447 base_type::template register_update_func<TGeom, TFunc>(&this_type::template update<TGeom>);
452 template <
typename TFVGeom>
454 update(
const TFVGeom* geo,
460 UG_ASSERT(Velocity != NULL,
"Null pointer");
463 const size_t numSH = geo->num_sh();
466 for(
size_t ip = 0; ip < geo->num_scvf(); ++ip)
469 const typename TFVGeom::SCVF& scvf = geo->scvf(ip);
472 const number flux =
VecDot(scvf.normal(), Velocity[ip]);
474 size_t from = scvf.from();
475 size_t to = scvf.to();
478 if((from >= scvf.num_sh()) || (to >= scvf.num_sh())){
481 for(
size_t sh = 0; sh < scvf.num_sh(); sh++)
482 conv_shape(ip, sh) = flux * scvf.shape(sh);
484 UG_ASSERT(scvf.num_sh() == numSH,
"sh's have to match!");
488 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
489 VecScale(D_vel(ip, sh), scvf.normal(), scvf.shape(sh));
496 size_t up = (flux >= 0) ? from : to;
499 const number noUpFlux = (1.-m_weight)*flux;
500 for(
size_t sh = 0; sh < scvf.num_sh(); ++sh)
501 conv_shape(ip, sh) = noUpFlux * scvf.shape(sh);
504 conv_shape(ip, up) += m_weight * flux;
510 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
511 VecScale(D_vel(ip, sh), scvf.normal(),
512 (1.-m_weight)*scvf.shape(sh));
514 for(
size_t sh = scvf.num_sh(); sh < numSH; sh++)
515 VecSet(D_vel(ip, sh), 0.0);
545 static const int dim = TDim;
561 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
566 template <
typename TFVGeom>
567 bool update(
const TFVGeom* geo,
581 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1Geometry<TElem, dim> >();
582 m_pThis->template register_func_for_elem_fvgeom< TElem, FV1CondensedGeometry<TElem, dim> >();
587 template <
typename TElem,
typename TFVGeom>
591 base_type::template register_update_func<TFVGeom, TFunc>(&this_type::template update<TFVGeom>);
603 template <
int refDim>
608 base_type::template register_update_func<TGeom, TFunc>(&this_type::template update<TGeom>);
613 template <
typename TFVGeom>
616 update(
const TFVGeom* geo,
622 UG_ASSERT(Velocity != NULL,
"Null pointer");
627 const number vol = ElementSize<dim>(geo->roid(), geo->corners());
630 for(
size_t i = 0; i < geo->num_scvf(); ++i)
633 const typename TFVGeom::SCVF& scvf = geo->scvf(i);
636 const size_t from = scvf.from();
637 const size_t to = scvf.to();
659 lambda = - GradDiffGrad * vol;
666 for(
size_t sh = 0; sh < scvf.num_sh(); ++sh)
667 conv_shape(i, sh) = 0.0;
669 for(
size_t sh = 0; sh < scvf.num_sh(); ++sh)
670 VecSet(D_vel(i, sh), 0.0);
675 if(from >= scvf.num_sh() || to >= scvf.num_sh())
676 UG_THROW(
"PartialUpwind: Currently not implemented for hanging nodes.")
682 if(lambda <= 0 || DiffDisp == NULL)
685 const size_t up = (flux >= 0) ? scvf.from() : scvf.to();
688 conv_shape(i, up) = flux;
694 D_vel(i, up) = scvf.normal();
697 set_non_zero_deriv_diffusion_flag(
false);
708 if (2 * lambda > fabs(flux))
710 conv_shape(i, from) = flux / 2.0;
711 conv_shape(i, to) = flux / 2.0;
715 set_non_zero_deriv_diffusion_flag(
false);
717 VecScale(D_vel(i,from), scvf.normal(), 1.0/2.0);
718 VecScale(D_vel(i, to), scvf.normal(), 1.0/2.0);
729 set_non_zero_deriv_diffusion_flag(
true);
732 conv_shape(i, from) = flux - lambda;
733 conv_shape(i, to) = lambda;
736 D_vel(i,from) = scvf.normal();
740 conv_shape(i, from) = - lambda;
741 conv_shape(i, to) = flux + lambda;
744 D_vel(i,to) = scvf.normal();
749 for(
size_t sh = 0; sh < scvf.num_sh(); ++sh)
750 MatSet(conv_shape_diffusion(i, sh), 0.0);
752 for (
size_t k = 0; k < (size_t)
dim; k++)
753 for (
size_t l = 0; l < (size_t)
dim; l++)
755 conv_shape_diffusion(i, from)(k,l) = gradFrom[k]*gradTo[l]*vol;
756 conv_shape_diffusion(i, to)(k,l) = - gradFrom[k]*gradTo[l]*vol;
781 static const int dim = TDim;
801 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
806 template <
typename TFVGeom>
807 bool update(
const TFVGeom *geo,
818 template <
typename TElem>
821 m_pThis->template register_func_for_elem_fvgeom<TElem, FV1Geometry<TElem, dim>>();
828 template <
typename TElem,
typename TFVGeom>
832 base_type::template register_update_func<TFVGeom, TFunc>(&this_type::template update<TFVGeom>);
840 template <
typename TRefDim>
845 template <
int refDim>
855 template <
typename TRefElem,
int TWorldDim>
868 if (!ElementSideRayIntersection<TRefElem, TWorldDim>(side, globalIntersection, localIntersection,
869 IP, IPVel,
false , vCornerCoords))
870 UG_THROW(
"GetNodeNextToCut: ElementSideRayIntersection failed");
877 number min = std::numeric_limits<number>::max();
880 for (
size_t i = 0; i < rRefElem.num(
dim - 1, side, 0); ++i)
883 const size_t co = rRefElem.id(
dim - 1, side, 0, i);
898 template <
typename TFVGeom>
900 update(
const TFVGeom *geo,
906 UG_ASSERT(Velocity != NULL,
"Null pointer");
909 const size_t numSH = geo->num_sh();
912 for (
size_t ip = 0; ip < geo->num_scvf(); ++ip)
915 const typename TFVGeom::SCVF &scvf = geo->scvf(ip);
918 const number flux =
VecDot(scvf.normal(), Velocity[ip]);
924 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
925 conv_shape(ip, sh) = flux * scvf.shape(sh);
926 for (
size_t sh = scvf.num_sh(); sh < numSH; sh++)
927 conv_shape(ip, sh) = 0.0;
930 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
931 VecScale(D_vel(ip, sh), scvf.normal(), scvf.shape(sh));
933 for (
size_t sh = scvf.num_sh(); sh < numSH; sh++)
934 VecSet(D_vel(ip, sh), 0.0);
940 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
941 conv_shape(ip, sh) = 0.0;
949 for (
size_t sh = scvf.num_sh(); sh < numSH; sh++)
950 conv_shape(ip, sh) = 0.0;
956 GetNodeNextToCut<typename TFVGeom::ref_elem_type, dim>(up, scvf.global_ip(), Velocity[ip], vCornerCoords);
960 conv_shape(ip, up) = flux;
965 for (
size_t sh = 0; sh < numSH; ++sh)
966 VecSet(D_vel(ip, sh), 0.0);
967 D_vel(ip, up) = scvf.normal();
993 static const int dim = TDim;
1013 boost::mpl::for_each<typename domain_traits<dim>::AllElemList>(
RegisterElemFunc(
this));
1014 boost::mpl::for_each<boost::mpl::range_c<int, 1, dim + 1>>(
RegisterRefDimFunc(
this));
1018 template <
typename TFVGeom>
1019 bool update(
const TFVGeom *geo,
1030 template <
typename TElem>
1033 m_pThis->template register_func_for_elem_fvgeom<TElem, FV1Geometry<TElem, dim>>();
1040 template <
typename TElem,
typename TFVGeom>
1044 base_type::template register_update_func<TFVGeom, TFunc>(&this_type::template update<TFVGeom>);
1052 template <
typename TRefDim>
1057 template <
int refDim>
1067 template <
typename TFVGeom>
1069 update(
const TFVGeom *geo,
1075 UG_ASSERT(Velocity != NULL,
"Null pointer");
1078 const size_t numSH = geo->num_sh();
1081 for (
size_t ip = 0; ip < geo->num_scvf(); ++ip)
1084 const typename TFVGeom::SCVF &scvf = geo->scvf(ip);
1087 const number flux =
VecDot(scvf.normal(), Velocity[ip]);
1093 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
1094 conv_shape(ip, sh) = flux * scvf.shape(sh);
1095 for (
size_t sh = scvf.num_sh(); sh < numSH; sh++)
1096 conv_shape(ip, sh) = 0.0;
1099 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
1100 VecScale(D_vel(ip, sh), scvf.normal(), scvf.shape(sh));
1102 for (
size_t sh = scvf.num_sh(); sh < numSH; sh++)
1103 VecSet(D_vel(ip, sh), 0.0);
1109 for (
size_t sh = 0; sh < scvf.num_sh(); sh++)
1110 conv_shape(ip, sh) = 0.0;
1118 for (
size_t sh = scvf.num_sh(); sh < numSH; sh++)
1119 conv_shape(ip, sh) = 0.0;
1123 typedef typename TFVGeom::ref_elem_type TRefElem;
1130 ElementSideRayIntersection<TRefElem, dim>(side, globalIntersection, localIntersection,
1131 scvf.global_ip(), Velocity[ip],
false , vCornerCoords);
1133 UG_CATCH_THROW(
"GetLinearProfileSkewedUpwindShapes: Cannot find cut side.");
1144 std::vector<number> vShape;
1146 std::vector<MathVector<TRefElem::dim>> vShapeGrad;
1148 TrialSpace.
shapes(vShape, localIntersection);
1149 TrialSpace.
grads(vShapeGrad, localIntersection);
1151 size_t num_corners_inters_side = rRefElem.num(
dim - 1, side, 0);
1154 for (
size_t j = 0; j < num_corners_inters_side; ++j)
1157 const size_t co = rRefElem.id(
dim-1, side, 0, j);
1160 conv_shape(ip, co) = flux * vShape[co];
Definition: conv_shape.h:194
IConvectionShapes< TDim > base_type
Base class.
Definition: conv_shape.h:197
void register_func_for_elem_fvgeom()
registers the update function for an element type and a FV geometry
Definition: conv_shape.h:251
void register_func_for_refDim()
registers the update function for a reference dimension
Definition: conv_shape.h:267
static const int dim
Dimension.
Definition: conv_shape.h:203
bool update(const TFVGeom *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
update of values for FV1Geometry
Definition: conv_shape.h:278
ConvectionShapesFullUpwind()
constructor
Definition: conv_shape.h:216
ConvectionShapesFullUpwind< TDim > this_type
This class.
Definition: conv_shape.h:200
Definition: conv_shape.h:984
ConvectionShapesLinearProfileSkewedUpwind< TDim > this_type
This class.
Definition: conv_shape.h:990
void register_func_for_refDim()
registers the update function for a reference dimension
Definition: conv_shape.h:1058
ConvectionShapesLinearProfileSkewedUpwind()
constructor
Definition: conv_shape.h:1006
void register_func_for_elem_fvgeom()
registers the update function for an element type and a FV geometry
Definition: conv_shape.h:1041
static const int dim
Dimension.
Definition: conv_shape.h:993
IConvectionShapes< TDim > base_type
Base class.
Definition: conv_shape.h:987
bool update(const TFVGeom *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
update of values for FV geometry
Definition: conv_shape.h:1069
Definition: conv_shape.h:53
static const int dim
Dimension.
Definition: conv_shape.h:62
void register_func_for_elem_fvgeom()
registers the update function for an element type and a FV geometry
Definition: conv_shape.h:110
bool update(const TFVGeom *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
update of values for FV geometry
Definition: conv_shape.h:138
IConvectionShapes< TDim > base_type
Base class.
Definition: conv_shape.h:56
ConvectionShapesNoUpwind()
constructor
Definition: conv_shape.h:75
void register_func_for_refDim()
registers the update function for a reference dimension
Definition: conv_shape.h:126
ConvectionShapesNoUpwind< TDim > this_type
This class.
Definition: conv_shape.h:59
Definition: conv_shape.h:536
ConvectionShapesPartialUpwind< TDim > this_type
This class.
Definition: conv_shape.h:542
IConvectionShapes< TDim > base_type
Base class.
Definition: conv_shape.h:539
void register_func_for_elem_fvgeom()
registers the update function for an element type and a FV geometry
Definition: conv_shape.h:588
void register_func_for_refDim()
registers the update function for a reference dimension
Definition: conv_shape.h:604
ConvectionShapesPartialUpwind()
constructor
Definition: conv_shape.h:558
bool update(const TFVGeom *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
update of values for FV1Geometry
Definition: conv_shape.h:616
static const int dim
Dimension.
Definition: conv_shape.h:545
Definition: conv_shape.h:772
static const int dim
Dimension.
Definition: conv_shape.h:781
ConvectionShapesSkewedUpwind< TDim > this_type
This class.
Definition: conv_shape.h:778
ConvectionShapesSkewedUpwind()
constructor
Definition: conv_shape.h:794
void register_func_for_refDim()
registers the update function for a reference dimension
Definition: conv_shape.h:846
bool update(const TFVGeom *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
update of values for FV geometry
Definition: conv_shape.h:900
void register_func_for_elem_fvgeom()
registers the update function for an element type and a FV geometry
Definition: conv_shape.h:829
IConvectionShapes< TDim > base_type
Base class.
Definition: conv_shape.h:775
Definition: conv_shape.h:349
void register_func_for_elem_fvgeom()
registers the update function for an element type and a FV geometry
Definition: conv_shape.h:427
void register_func_for_refDim()
registers the update function for a reference dimension
Definition: conv_shape.h:443
number m_weight
Definition: conv_shape.h:408
ConvectionShapesWeightedUpwind< TDim > this_type
This class.
Definition: conv_shape.h:355
ConvectionShapesWeightedUpwind(number weight)
constructor
Definition: conv_shape.h:383
IConvectionShapes< TDim > base_type
Base class.
Definition: conv_shape.h:352
bool update(const TFVGeom *geo, const MathVector< dim > *Velocity, const MathMatrix< dim, dim > *DiffDisp, bool computeDeriv)
update of values for FV1Geometry
Definition: conv_shape.h:454
void set_weight(number weight)
set weighting between full upwind (1.0) and no upwind (0.0)
Definition: conv_shape.h:397
ConvectionShapesWeightedUpwind()
constructor
Definition: conv_shape.h:371
static const int dim
Dimension.
Definition: conv_shape.h:358
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition: fv1_geom.h:665
Interface class for upwind methods.
Definition: conv_shape_interface.h:100
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 register_update_func(TAssFunc func)
register a update function for a Geometry
Definition: conv_shape_interface.h:249
void set_non_zero_deriv_diffusion_flag(bool flag)
sets the shape ip flag
Definition: conv_shape_interface.h:174
bool non_zero_deriv_diffusion() const
returns if upwind shape w.r.t. ip vel is non-zero
Definition: conv_shape_interface.h:149
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
Identifier for Local Finite Elements.
Definition: local_finite_element_id.h:98
@ LAGRANGE
Definition: local_finite_element_id.h:104
virtual base class for local shape function sets
Definition: local_shape_function_set.h:70
virtual void grads(grad_type *vGrad, const MathVector< dim > &x) const =0
returns all gradients evaluated at a point
virtual void shapes(shape_type *vShape, const MathVector< dim > &x) const =0
returns all shape functions evaluated at a point
static TClass & get()
returns a singleton based on the identifier
Definition: provider.h:54
void MatSet(matrix_t &mInOut, typename matrix_t::value_type s)
Set each matrix entry to a scalar (componentwise)
Definition: math_matrix_functions_common_impl.hpp:778
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_CATCH_THROW(msg)
Definition: error.h:64
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
void MatVecMult(vector_t_out &vOut, const matrix_t &m, const vector_t_in &v)
Matrix - Vector Multiplication.
Definition: math_matrix_vector_functions_common_impl.hpp:49
void VecScaleAppend(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1)
Scales a Vector and adds it to a second vector.
Definition: math_vector_functions_common_impl.hpp:126
void VecSet(vector_t &vInOut, typename vector_t::value_type s)
Set each vector component to scalar (componentwise)
Definition: math_vector_functions_common_impl.hpp:539
vector_t::value_type VecDistanceSq(const vector_t &v1, const vector_t &v2)
returns the squared distance of two vector_ts.
Definition: math_vector_functions_common_impl.hpp:351
void VecScale(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
scales a MathVector<N>
Definition: math_vector_functions_common_impl.hpp:252
vector_t::value_type VecTwoNorm(const vector_t &v)
Definition: math_vector_functions_common_impl.hpp:577
vector_t::value_type VecDot(const vector_t &v1, const vector_t &v2)
returns the dot-product of two vector_ts
Definition: math_vector_functions_common_impl.hpp:385
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74
void GetNodeNextToCut(size_t &coOut, const MathVector< TWorldDim > &IP, const MathVector< TWorldDim > &IPVel, const MathVector< TWorldDim > *vCornerCoords)
computes the closest node to a elem side ray intersection
Definition: conv_shape.h:856
functor for registering the shapes for the element-templated FV geometries
Definition: conv_shape.h:238
void operator()(TElem &)
Definition: conv_shape.h:241
this_type * m_pThis
Definition: conv_shape.h:239
RegisterElemFunc(this_type *pThis)
Definition: conv_shape.h:240
functor for registering the shapes for the reference-dimension-templated FV geometries
Definition: conv_shape.h:259
RegisterRefDimFunc(this_type *pThis)
Definition: conv_shape.h:261
this_type * m_pThis
Definition: conv_shape.h:260
void operator()(TRefDim &)
Definition: conv_shape.h:262
functor for registering the shapes for the element-templated FV geometries
Definition: conv_shape.h:1027
void operator()(TElem &)
Definition: conv_shape.h:1031
RegisterElemFunc(this_type *pThis)
Definition: conv_shape.h:1029
this_type * m_pThis
Definition: conv_shape.h:1028
functor for registering the shapes for the reference-dimension-templated FV geometries
Definition: conv_shape.h:1049
RegisterRefDimFunc(this_type *pThis)
Definition: conv_shape.h:1051
void operator()(TRefDim &)
Definition: conv_shape.h:1053
this_type * m_pThis
Definition: conv_shape.h:1050
functor for registering the shapes for the element-templated FV geometries
Definition: conv_shape.h:97
this_type * m_pThis
Definition: conv_shape.h:98
RegisterElemFunc(this_type *pThis)
Definition: conv_shape.h:99
void operator()(TElem &)
Definition: conv_shape.h:100
functor for registering the shapes for the reference-dimension-templated FV geometries
Definition: conv_shape.h:118
void operator()(TRefDim &)
Definition: conv_shape.h:121
RegisterRefDimFunc(this_type *pThis)
Definition: conv_shape.h:120
this_type * m_pThis
Definition: conv_shape.h:119
functor for registering the shapes for the element-templated FV geometries
Definition: conv_shape.h:576
void operator()(TElem &)
Definition: conv_shape.h:579
this_type * m_pThis
Definition: conv_shape.h:577
RegisterElemFunc(this_type *pThis)
Definition: conv_shape.h:578
functor for registering the shapes for the reference-dimension-templated FV geometries
Definition: conv_shape.h:596
this_type * m_pThis
Definition: conv_shape.h:597
void operator()(TRefDim &)
Definition: conv_shape.h:599
RegisterRefDimFunc(this_type *pThis)
Definition: conv_shape.h:598
functor for registering the shapes for the element-templated FV geometries
Definition: conv_shape.h:815
this_type * m_pThis
Definition: conv_shape.h:816
void operator()(TElem &)
Definition: conv_shape.h:819
RegisterElemFunc(this_type *pThis)
Definition: conv_shape.h:817
functor for registering the shapes for the reference-dimension-templated FV geometries
Definition: conv_shape.h:837
this_type * m_pThis
Definition: conv_shape.h:838
RegisterRefDimFunc(this_type *pThis)
Definition: conv_shape.h:839
void operator()(TRefDim &)
Definition: conv_shape.h:841
functor for registering the shapes for the element-templated FV geometries
Definition: conv_shape.h:414
RegisterElemFunc(this_type *pThis)
Definition: conv_shape.h:416
void operator()(TElem &)
Definition: conv_shape.h:417
this_type * m_pThis
Definition: conv_shape.h:415
functor for registering the shapes for the reference-dimension-templated FV geometries
Definition: conv_shape.h:435
void operator()(TRefDim &)
Definition: conv_shape.h:438
RegisterRefDimFunc(this_type *pThis)
Definition: conv_shape.h:437
this_type * m_pThis
Definition: conv_shape.h:436