33 #ifndef __H__UG__LIB_DISC__REFERENCE_ELEMENT__REFERENCE_MAPPING_PROVIDER__
34 #define __H__UG__LIB_DISC__REFERENCE_ELEMENT__REFERENCE_MAPPING_PROVIDER__
51 template <
int TDim,
int TWorldDim = TDim>
59 static const int dim = TDim;
86 const size_t maxIter = 1000,
87 const number tol = 1e-10)
const = 0;
92 const size_t maxIter = 1000,
93 const number tol = 1e-10)
const = 0;
98 const size_t maxIter = 1000,
99 const number tol = 1e-10)
const = 0;
144 std::vector<number>& vDet,
194 template <
int TDim,
int TWorldDim>
198 UG_STATIC_ASSERT(TWorldDim <= 3, only_implemented_for_ref_dim_smaller_equal_3);
200 UG_ASSERT(roid >= 0,
"Roid specified incorrectly.");
206 template <
int TDim,
int TWorldDim>
209 m_vvvMapping[TDim][TWorldDim][roid] =
reinterpret_cast<void*
>(&map);
224 template <
int TDim,
int TWorldDim>
229 UG_THROW(
"ReferenceMappingProvider: ReferenceMapping not found for "
230 <<roid<<
" from R^"<<TDim<<
" to R^"<<TWorldDim);
247 template <
int TDim,
int TWorldDim>
268 template <
int TDim,
int TWorldDim>
virtual base class for reference mappings
Definition: reference_mapping_provider.h:53
virtual void jacobian_transposed(MathMatrix< dim, worldDim > *vJT, const MathVector< dim > *vLocPos, size_t n) const =0
returns transposed of jacobian for n local positions
virtual void jacobian_transposed_inverse(MathMatrix< worldDim, dim > *vJTInv, number *vDet, const MathVector< dim > *vLocPos, size_t n) const =0
returns transposed of the inverse of the jacobian for n local positions
virtual void jacobian_transposed(MathMatrix< dim, worldDim > &JT, const MathVector< dim > &locPos) const =0
returns transposed of jacobian
virtual void jacobian(std::vector< MathMatrix< worldDim, dim > > &vJ, const std::vector< MathVector< dim > > &vLocPos) const =0
returns jacobian for a vector of local positions
virtual void local_to_global(std::vector< MathVector< worldDim > > &vGlobPos, const std::vector< MathVector< dim > > &vLocPos) const =0
map local coordinate to global coordinate for a vector of local positions
virtual void jacobian_transposed_inverse(std::vector< MathMatrix< worldDim, dim > > &vJTInv, const std::vector< MathVector< dim > > &vLocPos) const =0
returns transposed of the inverse of the jacobian for a vector of positions
virtual void global_to_local(std::vector< MathVector< dim > > &vLocPos, const std::vector< MathVector< worldDim > > &vGlobPos, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate for a vector of local positions
static const int dim
reference dimension (domain space dimension)
Definition: reference_mapping_provider.h:59
virtual number sqrt_gram_det(const MathVector< dim > &locPos) const =0
returns the determinate of the jacobian
static const int worldDim
world dimension (range space dimension)
Definition: reference_mapping_provider.h:56
virtual void sqrt_gram_det(number *vDet, const MathVector< dim > *vLocPos, size_t n) const =0
returns the determinate of the jacobian for n local positions
virtual void global_to_local(MathVector< dim > *vLocPos, const MathVector< worldDim > *vGlobPos, size_t n, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate for n local positions
virtual void sqrt_gram_det(std::vector< number > &vDet, const std::vector< MathVector< dim > > &vLocPos) const =0
returns the determinate of the jacobian for a vector of local positions
virtual ~DimReferenceMapping()
virtual destructor
Definition: reference_mapping_provider.h:158
virtual number jacobian_transposed_inverse(MathMatrix< worldDim, dim > &JTInv, const MathVector< dim > &locPos) const =0
returns transposed of the inverse of the jacobian and returns sqrt of gram determinante
virtual void jacobian_transposed_inverse(MathMatrix< worldDim, dim > *vJTInv, const MathVector< dim > *vLocPos, size_t n) const =0
returns transposed of the inverse of the jacobian for n local positions
virtual void jacobian_transposed(std::vector< MathMatrix< dim, worldDim > > &vJT, const std::vector< MathVector< dim > > &vLocPos) const =0
returns transposed of jacobian for a vector of positions
virtual void local_to_global(MathVector< worldDim > &globPos, const MathVector< dim > &locPos) const =0
map local coordinate to global coordinate
virtual void local_to_global(MathVector< worldDim > *vGlobPos, const MathVector< dim > *vLocPos, size_t n) const =0
map n local coordinate to global coordinate
virtual void update(const MathVector< worldDim > *vCornerCoord)=0
refresh mapping for new set of corners
virtual void update(const std::vector< MathVector< worldDim > > &vCornerCoord)=0
refresh mapping for new set of corners
virtual void jacobian_transposed_inverse(std::vector< MathMatrix< worldDim, dim > > &vJTInv, std::vector< number > &vDet, const std::vector< MathVector< dim > > &vLocPos) const =0
returns transposed of the inverse of the jacobian for a vector of positions
virtual void global_to_local(MathVector< dim > &locPos, const MathVector< worldDim > &globPos, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate
virtual bool is_linear() const =0
returns if mapping is affine
virtual void jacobian(MathMatrix< worldDim, dim > &J, const MathVector< dim > &locPos) const =0
returns jacobian
virtual void jacobian(MathMatrix< worldDim, dim > *vJ, const MathVector< dim > *vLocPos, size_t n) const =0
returns jacobian for n local positions
class to provide reference mappings
Definition: reference_mapping_provider.h:168
ReferenceMappingProvider & operator=(const ReferenceMappingProvider &)
void * m_vvvMapping[4][4][NUM_REFERENCE_OBJECTS]
Definition: reference_mapping_provider.h:185
void set_mapping(ReferenceObjectID roid, DimReferenceMapping< TDim, TWorldDim > &map)
Definition: reference_mapping_provider.h:207
DimReferenceMapping< TDim, TWorldDim > * get_mapping(ReferenceObjectID roid)
Definition: reference_mapping_provider.h:195
ReferenceMappingProvider(const ReferenceMappingProvider &)
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid, const MathVector< TWorldDim > *vCornerCoord)
returns a reference to a DimReferenceMapping with updated element corners
Definition: reference_mapping_provider.h:269
ReferenceMappingProvider()
Definition: reference_mapping_provider.cpp:222
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid, const std::vector< MathVector< TWorldDim > > &vCornerCoord)
returns a reference to a DimReferenceMapping with updated element corners
Definition: reference_mapping_provider.h:248
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid)
returns a reference to a DimReferenceMapping
Definition: reference_mapping_provider.h:225
~ReferenceMappingProvider()
Definition: reference_mapping_provider.h:178
static ReferenceMappingProvider & inst()
Definition: reference_mapping_provider.h:181
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_THROW(msg)
Definition: error.h:57
#define UG_STATIC_ASSERT(expr, msg)
Checks an expression at compile-time and raises a compile-error if the expression equals 0.
Definition: static_assert.h:63
double number
Definition: types.h:124
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74
@ NUM_REFERENCE_OBJECTS
Definition: grid_base_objects.h:85