13#include <kernel/shape.hpp>
14#include <kernel/util/tiny_algebra.hpp>
16#include <kernel/util/math.hpp>
18#include <kernel/util/cuda_math.cuh>
46 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
47 typename HessType_,
typename HessInvType_,
typename Shape_,
int img_dim>
51 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
52 typename HessType_,
typename HessInvType_,
int img_dim>
53 struct EvalHelper<DataType_, DomPointType_, ImgPointType_, JacMatType_, JacMatInvType_, JacDetType_, HessType_, HessInvType_, Shape::Simplex<1>, img_dim>
75 static constexpr int image_dim = img_dim;
92 template<
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
95 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
97 const VertexType& v0 = vertex_set[index_set(cell_index, 0)];
98 const VertexType& v1 = vertex_set[index_set(cell_index, 1)];
101 for(
int i(0); i < image_dim; ++i)
104 coeffs[i][1] =
DataType(v1[i] - v0[i]);
123 for(
int i(0); i < image_dim; ++i)
143 for(
int i(0); i < img_dim; ++i)
178 #ifndef __CUDA_ARCH__
179 for(
int i(0); i < image_dim; ++i)
183 for(
int i(0); i < image_dim; ++i)
184 v += CudaMath::cuda_sqr(coeffs[i][1]);
185 return CudaMath::cuda_sqrt(v);
206 return volume(coeffs);
210 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
211 typename HessType_,
typename HessInvType_,
int img_dim>
212 struct EvalHelper<DataType_, DomPointType_, ImgPointType_, JacMatType_, JacMatInvType_, JacDetType_, HessType_, HessInvType_, Shape::Simplex<2>, img_dim>
234 static constexpr int image_dim = img_dim;
252 template<
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
255 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
257 const VertexType& v0 = vertex_set[index_set(cell_index, 0)];
258 const VertexType& v1 = vertex_set[index_set(cell_index, 1)];
259 const VertexType& v2 = vertex_set[index_set(cell_index, 2)];
262 for(
int i(0); i < image_dim; ++i)
265 coeffs[i][1] =
DataType(v1[i] - v0[i]);
266 coeffs[i][2] =
DataType(v2[i] - v0[i]);
285 for(
int i(0); i < image_dim; ++i)
305 for(
int i(0); i < image_dim; ++i)
341 for(
int i(0); i < image_dim; ++i)
378 for(
int i(0); i < image_dim; ++i)
386 ref_ray.set_mat_vec_mult(
jac_inv, ray);
397 #ifndef __CUDA_ARCH__
405 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
406 typename HessType_,
typename HessInvType_,
int img_dim>
407 struct EvalHelper<DataType_, DomPointType_, ImgPointType_, JacMatType_, JacMatInvType_, JacDetType_, HessType_, HessInvType_, Shape::Simplex<3>, img_dim>
429 static constexpr int image_dim = img_dim;
446 template<
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
450 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
451 const VertexType& v0 = vertex_set[index_set(cell_index, 0)];
452 const VertexType& v1 = vertex_set[index_set(cell_index, 1)];
453 const VertexType& v2 = vertex_set[index_set(cell_index, 2)];
454 const VertexType& v3 = vertex_set[index_set(cell_index, 3)];
457 for(
int i(0); i < image_dim; ++i)
460 coeffs[i][1] =
DataType(v1[i] - v0[i]);
461 coeffs[i][2] =
DataType(v2[i] - v0[i]);
462 coeffs[i][3] =
DataType(v3[i] - v0[i]);
481 for(
int i(0); i < image_dim; ++i)
503 for(
int i(0); i < image_dim; ++i)
540 for(
int i(0); i < image_dim; ++i)
573 for(
int i(0); i < image_dim; ++i)
582 ref_ray.set_mat_vec_mult(
jac_inv, ray);
593 #ifndef __CUDA_ARCH__
597 return DataType(2) * CudaMath::cuda_rsqrt(CudaMath::sqr(
DataType(2) * ref_ray[0] + ref_ray[1] + ref_ray[2]) +
598 CudaMath::cuda_sqr(ref_ray[1] - ref_ray[2]) +
DataType(2) * CudaMath::cuda_sqr(ref_ray[1] + ref_ray[2]));
603 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
604 typename HessType_,
typename HessInvType_,
int img_dim>
605 struct EvalHelper<DataType_, DomPointType_, ImgPointType_, JacMatType_, JacMatInvType_, JacDetType_, HessType_, HessInvType_, Shape::Hypercube<1>, img_dim>
627 static constexpr int image_dim = img_dim;
644 template<
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
647 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
649 const VertexType& v0 = vertex_set[index_set(cell_index, 0)];
650 const VertexType& v1 = vertex_set[index_set(cell_index, 1)];
653 for(
int i(0); i < img_dim; ++i)
675 for(
int i(0); i < image_dim; ++i)
695 for(
int i(0); i < img_dim; ++i)
730 #ifndef __CUDA_ARCH__
731 for(
int i(0); i < image_dim; ++i)
735 for(
int i(0); i < image_dim; ++i)
736 v += CudaMath::cuda_sqr(coeffs[i][1]);
737 return DataType(2) * CudaMath::cuda_sqrt(v);
758 return volume(coeffs);
762 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
763 typename HessType_,
typename HessInvType_,
int img_dim>
764 struct EvalHelper<DataType_, DomPointType_, ImgPointType_, JacMatType_, JacMatInvType_, JacDetType_, HessType_, HessInvType_, Shape::Hypercube<2>, img_dim>
786 static constexpr int image_dim = img_dim;
803 template<
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
806 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
808 const VertexType& v0 = vertex_set[index_set(cell_index, 0)];
809 const VertexType& v1 = vertex_set[index_set(cell_index, 1)];
810 const VertexType& v2 = vertex_set[index_set(cell_index, 2)];
811 const VertexType& v3 = vertex_set[index_set(cell_index, 3)];
814 for(
int i(0); i < image_dim; ++i)
835 template<
typename ThreadGroup_,
typename VertexSetType_,
typename IndexSetType_>
836 CUDA_DEVICE
static void __forceinline__ _group_gather_vertex_set(
const ThreadGroup_& tg,
typename VertexSetType_::VertexType* __restrict__ vn,
const VertexSetType_& vertex_set,
const IndexSetType_& index_set,
const Index cell_index)
839 for(
int i = tg.thread_rank(); i < num_verts; i += tg.num_threads())
840 vn[i] = vertex_set[index_set(cell_index, i)];
843 template<
typename ThreadGroup_,
typename VertexType_>
844 CUDA_DEVICE
static void __forceinline__ _group_set_coeffs(
const ThreadGroup_& tg, DataType* __restrict__ coeffs,
const VertexType_* __restrict__ vn)
850 for(
int i = tg.thread_rank(); i < image_dim; i += tg.num_threads())
852 coeffs[4*i+0] = DataType(0.25) * DataType( vn[0][i] + vn[1][i] + vn[2][i] + vn[3][i]);
853 coeffs[4*i+1] = DataType(0.25) * DataType(-vn[0][i] + vn[1][i] - vn[2][i] + vn[3][i]);
854 coeffs[4*i+2] = DataType(0.25) * DataType(-vn[0][i] - vn[1][i] + vn[2][i] + vn[3][i]);
855 coeffs[4*i+3] = DataType(0.25) * DataType( vn[0][i] - vn[1][i] - vn[2][i] + vn[3][i]);
878 template<
typename ThreadGroup_,
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
879 CUDA_DEVICE
static void inline grouped_set_coefficients(
const ThreadGroup_& tg, DataType* coeffs,
const VertexSetType_& vertex_set,
const IndexSetType_& index_set,
const IT_ cell_index)
881 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
883 __shared__ VertexType vn[8];
885 _group_gather_vertex_set(tg, vn, vertex_set, index_set, cell_index);
889 _group_set_coeffs(tg, coeffs, vn);
909 for(
int i(0); i < image_dim; ++i)
930 for(
int i(0); i < image_dim; ++i)
957 template<
typename ThreadGroup_>
958 CUDA_DEVICE
static inline void grouped_calc_jac_mat(
const ThreadGroup_& tg, JacobianMatrixType&
jac_mat,
const DomainPointType&
dom_point,
const DataType* coeffs)
960 for(
int idx = tg.thread_rank(); idx < 2*image_dim; idx += tg.num_threads())
962 __builtin_assume(idx >= 0);
963 const int i = idx/image_dim;
964 const int curr = idx%image_dim;
966 jac_mat(i,curr) = coeffs[i*num_verts + curr+1] +
dom_point[1-curr] * coeffs[i*num_verts + 3];
987 for(
int i(0); i < image_dim; ++i)
1019 for(
int i(0); i < image_dim; ++i)
1050 for(
int i(0); i < image_dim; ++i)
1058 ref_ray.set_mat_vec_mult(
jac_inv, ray);
1061 return DataType(2) / ref_ray.norm_euclid();
1065 template<
typename DataType_,
typename DomPointType_,
typename ImgPointType_,
typename JacMatType_,
typename JacMatInvType_,
typename JacDetType_,
1066 typename HessType_,
typename HessInvType_,
int img_dim>
1067 struct EvalHelper<DataType_, DomPointType_, ImgPointType_, JacMatType_, JacMatInvType_, JacDetType_, HessType_, HessInvType_, Shape::Hypercube<3>, img_dim>
1089 static constexpr int image_dim = img_dim;
1106 template<
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
1109 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
1111 const VertexType& v0 = vertex_set[index_set(cell_index, 0)];
1112 const VertexType& v1 = vertex_set[index_set(cell_index, 1)];
1113 const VertexType& v2 = vertex_set[index_set(cell_index, 2)];
1114 const VertexType& v3 = vertex_set[index_set(cell_index, 3)];
1115 const VertexType& v4 = vertex_set[index_set(cell_index, 4)];
1116 const VertexType& v5 = vertex_set[index_set(cell_index, 5)];
1117 const VertexType& v6 = vertex_set[index_set(cell_index, 6)];
1118 const VertexType& v7 = vertex_set[index_set(cell_index, 7)];
1123 for(
int i(0); i < image_dim; ++i)
1125 coeffs[i][0] =
DataType(0.125) *
DataType( + v0[i] + v1[i] + v2[i] + v3[i] + v4[i] + v5[i] + v6[i] + v7[i]);
1126 coeffs[i][1] =
DataType(0.125) *
DataType( - v0[i] + v1[i] - v2[i] + v3[i] - v4[i] + v5[i] - v6[i] + v7[i]);
1127 coeffs[i][2] =
DataType(0.125) *
DataType( - v0[i] - v1[i] + v2[i] + v3[i] - v4[i] - v5[i] + v6[i] + v7[i]);
1128 coeffs[i][3] =
DataType(0.125) *
DataType( - v0[i] - v1[i] - v2[i] - v3[i] + v4[i] + v5[i] + v6[i] + v7[i]);
1129 coeffs[i][4] =
DataType(0.125) *
DataType( + v0[i] - v1[i] - v2[i] + v3[i] + v4[i] - v5[i] - v6[i] + v7[i]);
1130 coeffs[i][5] =
DataType(0.125) *
DataType( + v0[i] - v1[i] + v2[i] - v3[i] - v4[i] + v5[i] - v6[i] + v7[i]);
1131 coeffs[i][6] =
DataType(0.125) *
DataType( + v0[i] + v1[i] - v2[i] - v3[i] - v4[i] - v5[i] + v6[i] + v7[i]);
1132 coeffs[i][7] =
DataType(0.125) *
DataType( - v0[i] + v1[i] + v2[i] - v3[i] + v4[i] - v5[i] - v6[i] + v7[i]);
1138 template<
typename ThreadGroup_,
typename VertexSetType_,
typename IndexSetType_>
1139 CUDA_DEVICE
static void __forceinline__ _group_gather_vertex_set(
const ThreadGroup_& tg,
typename VertexSetType_::VertexType* __restrict__ vn,
const VertexSetType_& vertex_set,
const IndexSetType_& index_set,
const Index cell_index)
1142 for(
int i = tg.thread_rank(); i < num_verts; i += tg.num_threads())
1143 vn[i] = vertex_set[index_set(cell_index, i)];
1146 template<
typename ThreadGroup_,
typename VertexType_>
1147 CUDA_DEVICE
static void __forceinline__ _group_set_coeffs(
const ThreadGroup_& tg, DataType* __restrict__ coeffs,
const VertexType_* __restrict__ vn)
1153 for(
int i = tg.thread_rank(); i < image_dim; i += tg.num_threads())
1155 coeffs[8*i+0] = DataType(0.125) * DataType( + vn[0][i] + vn[1][i] + vn[2][i] + vn[3][i] + vn[4][i] + vn[5][i] + vn[6][i] + vn[7][i]);
1156 coeffs[8*i+1] = DataType(0.125) * DataType( - vn[0][i] + vn[1][i] - vn[2][i] + vn[3][i] - vn[4][i] + vn[5][i] - vn[6][i] + vn[7][i]);
1157 coeffs[8*i+2] = DataType(0.125) * DataType( - vn[0][i] - vn[1][i] + vn[2][i] + vn[3][i] - vn[4][i] - vn[5][i] + vn[6][i] + vn[7][i]);
1158 coeffs[8*i+3] = DataType(0.125) * DataType( - vn[0][i] - vn[1][i] - vn[2][i] - vn[3][i] + vn[4][i] + vn[5][i] + vn[6][i] + vn[7][i]);
1159 coeffs[8*i+4] = DataType(0.125) * DataType( + vn[0][i] - vn[1][i] - vn[2][i] + vn[3][i] + vn[4][i] - vn[5][i] - vn[6][i] + vn[7][i]);
1160 coeffs[8*i+5] = DataType(0.125) * DataType( + vn[0][i] - vn[1][i] + vn[2][i] - vn[3][i] - vn[4][i] + vn[5][i] - vn[6][i] + vn[7][i]);
1161 coeffs[8*i+6] = DataType(0.125) * DataType( + vn[0][i] + vn[1][i] - vn[2][i] - vn[3][i] - vn[4][i] - vn[5][i] + vn[6][i] + vn[7][i]);
1162 coeffs[8*i+7] = DataType(0.125) * DataType( - vn[0][i] + vn[1][i] + vn[2][i] - vn[3][i] + vn[4][i] - vn[5][i] - vn[6][i] + vn[7][i]);
1185 template<
typename ThreadGroup_,
typename VertexSetType_,
typename IndexSetType_,
typename IT_>
1186 CUDA_DEVICE
static void inline grouped_set_coefficients(
const ThreadGroup_& tg, DataType* coeffs,
const VertexSetType_& vertex_set,
const IndexSetType_& index_set,
const IT_ cell_index)
1188 typedef std::remove_const_t<std::remove_reference_t<
decltype(vertex_set[0])>> VertexType;
1190 __shared__ VertexType vn[8];
1192 _group_gather_vertex_set(tg, vn, vertex_set, index_set, cell_index);
1196 _group_set_coeffs(tg, coeffs, vn);
1216 for(
int i(0); i < image_dim; ++i)
1240 for(
int i(0); i < image_dim; ++i)
1268 template<
typename ThreadGroup_>
1269 CUDA_DEVICE
static inline void grouped_calc_jac_mat(
const ThreadGroup_& tg, JacobianMatrixType&
jac_mat,
const DomainPointType&
dom_point,
const DataType* coeffs)
1271 for(
int idx = tg.thread_rank(); idx < 3*image_dim; idx += tg.num_threads())
1273 __builtin_assume(idx >= 0);
1274 const int i = idx/image_dim;
1275 const int curr = idx%image_dim;
1277 jac_mat(i,curr) = coeffs[i*num_verts + curr+1] +
dom_point[1-(curr+1)/2] * coeffs[i* num_verts + 4+curr/2] +
dom_point[2-curr/2] * (coeffs[i*num_verts + 5+(curr+1)/2] +
dom_point[1-(curr+1)/2] * coeffs[i*num_verts + 7]);
1296 for(
int i(0); i < image_dim; ++i)
1323 const DataType cx =
DataType(FEAT_F128C(0.57735026918962576450914878050195745564760175127));
1330 for(
int i(0); i < 8; ++i)
1333 for(
int j(0); j < 3; ++j)
1334 cub_pt[j] =
DataType((((i >> j) & 1) << 1) - 1) * cx;
1337 calc_jac_mat(
jac_mat, cub_pt, coeffs);
1365 for(
int i(0); i < image_dim; ++i)
1374 ref_ray.set_mat_vec_mult(
jac_inv, ray);
1377 return DataType(2) / ref_ray.norm_euclid();
Tiny Matrix class template.
T_ sqrt(T_ x)
Returns the square-root of a value.
T_ sqr(T_ x)
Returns the square of a value.
std::uint64_t Index
Index data type.
@ img_point
specifies whether the trafo should supply image point coordinates
@ dom_point
specifies whether the trafo should supply domain point coordinates
@ hess_ten
specifies whether the trafo should supply hessian tensors
@ jac_inv
specifies whether the trafo should supply inverse jacobian matrices
@ jac_mat
specifies whether the trafo should supply jacobian matrices
Face traits tag struct template.
Hypercube shape tag struct template.
Simplex shape tag struct template.
HessType_ HessianTensorType
hessian tensor type
static CUDA_HOST_DEVICE void set_coefficients(Tiny::Matrix< DataType, image_dim, num_verts > &coeffs, const VertexSetType_ &vertex_set, const IndexSetType_ &index_set, IT_ cell_index)
Sets the coefficient vector based on the underlying standard transformation.
static CUDA_HOST_DEVICE void map_point(ImagePointType &img_point, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Maps a point from the reference cell to the selected cell.
static CUDA_HOST_DEVICE void calc_jac_mat(JacobianMatrixType &jac_mat, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Calculates the jacobian matrix for a given point.
Shape::Hypercube< 2 > ShapeType
the shape type
DomPointType_ DomainPointType
domain point type
ImgPointType_ ImagePointType
image point type
static CUDA_HOST_DEVICE void calc_hess_ten(HessianTensorType &hess_ten, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes the hessian tensor for a given domain point.
static CUDA_HOST_DEVICE DataType volume(const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the volume of the current cell.
JacDetType_ JacobianDeterminantType
jacobian determinant type
JacMatType_ JacobianMatrixType
jacobian matrix type
DataType_ DataType
evaluation data type
static CUDA_HOST_DEVICE DataType width_directed(const ImagePointType &ray, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the directed mesh width.
JacMatInvType_ JacobianInverseType
jacobian inverse matrix type
HessInvType_ HessianInverseType
hessian inverse tensor type
HessType_ HessianTensorType
hessian tensor type
static CUDA_HOST_DEVICE void calc_hess_ten(HessianTensorType &hess_ten, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes the hessian tensor for a given domain point.
DomPointType_ DomainPointType
domain point type
DataType_ DataType
evaluation data type
static CUDA_HOST_DEVICE DataType width_directed(const ImagePointType &ray, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the directed mesh width.
JacMatInvType_ JacobianInverseType
jacobian inverse matrix type
Shape::Simplex< 1 > ShapeType
the shape type
static CUDA_HOST_DEVICE void set_coefficients(Tiny::Matrix< DataType, image_dim, num_verts > &coeffs, const VertexSetType_ &vertex_set, const IndexSetType_ &index_set, IT_ cell_index)
Sets the coefficient vector based on the underlying standard transformation.
static CUDA_HOST_DEVICE void map_point(ImagePointType &img_point, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Maps a point from the reference cell to the selected cell.
static CUDA_HOST_DEVICE DataType volume(const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the volume of the current cell.
ImgPointType_ ImagePointType
image point type
JacMatType_ JacobianMatrixType
jacobian matrix type
HessInvType_ HessianInverseType
hessian inverse tensor type
JacDetType_ JacobianDeterminantType
jacobian determinant type
static CUDA_HOST_DEVICE void calc_jac_mat(JacobianMatrixType &jac_mat, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Calculates the jacobian matrix for a given point.
static CUDA_HOST_DEVICE void set_coefficients(Tiny::Matrix< DataType, image_dim, num_verts > &coeffs, const VertexSetType_ &vertex_set, const IndexSetType_ &index_set, IT_ cell_index)
Sets the coefficient vector based on the underlying standard transformation.
ImgPointType_ ImagePointType
image point type
JacMatInvType_ JacobianInverseType
jacobian inverse matrix type
HessInvType_ HessianInverseType
hessian inverse tensor type
static CUDA_HOST_DEVICE void calc_hess_ten(HessianTensorType &hess_ten, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes the hessian tensor for a given domain point.
DomPointType_ DomainPointType
domain point type
JacMatType_ JacobianMatrixType
jacobian matrix type
static CUDA_HOST_DEVICE DataType volume(const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the volume of the current cell.
HessType_ HessianTensorType
hessian tensor type
static CUDA_HOST_DEVICE void map_point(ImagePointType &img_point, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Maps a point from the reference cell to the selected cell.
Shape::Simplex< 2 > ShapeType
the shape type
static CUDA_HOST_DEVICE DataType width_directed(const ImagePointType &ray, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the directed mesh width.
DataType_ DataType
evaluation data type
static CUDA_HOST_DEVICE void calc_jac_mat(JacobianMatrixType &jac_mat, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Calculates the jacobian matrix for a given point.
JacDetType_ JacobianDeterminantType
jacobian determinant type
DataType_ DataType
evaluation data type
static CUDA_HOST_DEVICE void set_coefficients(Tiny::Matrix< DataType, image_dim, num_verts > &coeffs, const VertexSetType_ &vertex_set, const IndexSetType_ &index_set, IT_ cell_index)
Sets the coefficient vector based on the underlying standard transformation.
JacMatType_ JacobianMatrixType
jacobian matrix type
HessType_ HessianTensorType
hessian tensor type
static CUDA_HOST_DEVICE void map_point(ImagePointType &img_point, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Maps a point from the reference cell to the selected cell.
HessInvType_ HessianInverseType
hessian inverse tensor type
Shape::Simplex< 3 > ShapeType
the shape type
static CUDA_HOST_DEVICE void calc_hess_ten(HessianTensorType &hess_ten, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes the hessian tensor for a given domain point.
ImgPointType_ ImagePointType
image point type
static CUDA_HOST_DEVICE DataType width_directed(const ImagePointType &ray, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the directed mesh width.
static CUDA_HOST_DEVICE DataType volume(const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the volume of the current cell.
DomPointType_ DomainPointType
domain point type
JacDetType_ JacobianDeterminantType
jacobian determinant type
static CUDA_HOST_DEVICE void calc_jac_mat(JacobianMatrixType &jac_mat, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Calculates the jacobian matrix for a given point.
JacMatInvType_ JacobianInverseType
jacobian inverse matrix type
HessType_ HessianTensorType
hessian tensor type
static CUDA_HOST_DEVICE void calc_hess_ten(HessianTensorType &hess_ten, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes the hessian tensor for a given domain point.
static CUDA_HOST_DEVICE void set_coefficients(Tiny::Matrix< DataType, img_dim, num_verts > &coeffs, const VertexSetType_ &vertex_set, const IndexSetType_ &index_set, IT_ cell_index)
Sets the coefficient vector based on the underlying standard transformation.
JacDetType_ JacobianDeterminantType
jacobian determinant type
JacMatType_ JacobianMatrixType
jacobian matrix type
DataType_ DataType
evaluation data type
Shape::Hypercube< 3 > ShapeType
the shape type
static CUDA_HOST_DEVICE DataType volume(const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the volume of the current cell.
ImgPointType_ ImagePointType
image point type
static CUDA_HOST_DEVICE void map_point(ImagePointType &img_point, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Maps a point from the reference cell to the selected cell.
HessInvType_ HessianInverseType
hessian inverse tensor type
static CUDA_HOST_DEVICE void calc_jac_mat(JacobianMatrixType &jac_mat, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Calculates the jacobian matrix for a given point.
JacMatInvType_ JacobianInverseType
jacobian inverse matrix type
DomPointType_ DomainPointType
domain point type
static CUDA_HOST_DEVICE DataType width_directed(const ImagePointType &ray, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the directed mesh width.
JacMatType_ JacobianMatrixType
jacobian matrix type
static CUDA_HOST_DEVICE void map_point(ImagePointType &img_point, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Maps a point from the reference cell to the selected cell.
HessInvType_ HessianInverseType
hessian inverse tensor type
HessType_ HessianTensorType
hessian tensor type
static CUDA_HOST_DEVICE void set_coefficients(Tiny::Matrix< DataType, image_dim, num_verts > &coeffs, const VertexSetType_ &vertex_set, const IndexSetType_ &index_set, IT_ cell_index)
Sets the coefficient vector based on the underlying standard transformation.
DataType_ DataType
evaluation data type
static CUDA_HOST_DEVICE DataType width_directed(const ImagePointType &ray, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the directed mesh width.
DomPointType_ DomainPointType
domain point type
static CUDA_HOST_DEVICE DataType volume(const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes and returns the volume of the current cell.
Shape::Hypercube< 1 > ShapeType
the shape type
JacDetType_ JacobianDeterminantType
jacobian determinant type
ImgPointType_ ImagePointType
image point type
static CUDA_HOST_DEVICE void calc_jac_mat(JacobianMatrixType &jac_mat, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Calculates the jacobian matrix for a given point.
JacMatInvType_ JacobianInverseType
jacobian inverse matrix type
static CUDA_HOST_DEVICE void calc_hess_ten(HessianTensorType &hess_ten, const DomainPointType &dom_point, const Tiny::Matrix< DataType, image_dim, num_verts > &coeffs)
Computes the hessian tensor for a given domain point.
Evalautator helper class for standard trafo.