12 #ifndef MFEM_TEMPLATE_COEFFICIENT 13 #define MFEM_TEMPLATE_COEFFICIENT 15 #include "../config/tconfig.hpp" 16 #include "../linalg/ttensor.hpp" 17 #include "../linalg/tlayout.hpp" 18 #include "../linalg/vector.hpp" 37 template <
typename complex_t =
double>
50 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
51 inline MFEM_ALWAYS_INLINE
52 void Eval(
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
const 54 TAssign<AssignOp::Set>(l, c,
value);
66 template <
typename Func,
typename complex_t =
double>
76 template <
int dim,
bool dummy>
struct Dim;
77 template <
bool dummy>
struct Dim<1,dummy>
79 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
80 static inline MFEM_ALWAYS_INLINE
81 void Eval(Func &
F,
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
83 const int qpts = T_result_t::x_type::layout_type::dim_1;
84 const int ne = T_result_t::x_type::layout_type::dim_3;
85 const int vs =
sizeof(T.x[0])/
sizeof(T.x[0][0]);
86 for (
int k = 0; k < ne; k++)
88 for (
int i = 0; i < qpts; i++)
90 for (
int s = 0;
s < vs;
s++)
92 c[l.ind(i,k)][
s] =
F.Eval1D(T.x(i,0,k)[
s]);
98 template <
bool dummy>
struct Dim<2,dummy>
100 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
101 static inline MFEM_ALWAYS_INLINE
102 void Eval(Func &
F,
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
104 const int qpts = T_result_t::x_type::layout_type::dim_1;
105 const int ne = T_result_t::x_type::layout_type::dim_3;
106 const int vs =
sizeof(T.x[0])/
sizeof(T.x[0][0]);
107 for (
int k = 0; k < ne; k++)
109 for (
int i = 0; i < qpts; i++)
111 for (
int s = 0;
s < vs;
s++)
113 c[l.ind(i,k)][
s] =
F.Eval2D(T.x(i,0,k)[
s], T.x(i,1,k)[
s]);
119 template <
bool dummy>
struct Dim<3,dummy>
121 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
122 static inline MFEM_ALWAYS_INLINE
123 void Eval(Func &
F,
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
125 const int qpts = T_result_t::x_type::layout_type::dim_1;
126 const int ne = T_result_t::x_type::layout_type::dim_3;
127 const int vs =
sizeof(T.x[0])/
sizeof(T.x[0][0]);
128 for (
int k = 0; k < ne; k++)
130 for (
int i = 0; i < qpts; i++)
132 for (
int s = 0;
s < vs;
s++)
135 F.Eval3D(T.x(i,0,k)[
s], T.x(i,1,k)[
s], T.x(i,2,k)[
s]);
149 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
150 inline MFEM_ALWAYS_INLINE
151 void Eval(
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
153 const int qpts = T_result_t::x_type::layout_type::dim_1;
154 const int sdim = T_result_t::x_type::layout_type::dim_2;
155 const int ne = T_result_t::x_type::layout_type::dim_3;
156 MFEM_STATIC_ASSERT(c_layout_t::rank == 2 && c_layout_t::dim_1 == qpts &&
157 c_layout_t::dim_2 == ne,
"invalid c_layout_t");
166 template <
typename complex_t =
double>
182 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
183 inline MFEM_ALWAYS_INLINE
184 void Eval(
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
186 const int ne = T_result_t::ne;
187 const int vs =
sizeof(T.attrib[0])/
sizeof(T.attrib[0][0]);
188 MFEM_STATIC_ASSERT(vs ==
sizeof(c[0])/
sizeof(c[0][0]),
"");
189 for (
int i = 0; i < ne; i++)
191 typename c_data_t::data_type ci;
192 for (
int s = 0;
s < vs;
s++)
196 TAssign<AssignOp::Set>(l.ind2(i), c, ci);
202 template <
typename FieldEval>
218 inline MFEM_ALWAYS_INLINE
225 inline MFEM_ALWAYS_INLINE
230 :
fieldEval(tfes, shapeEval, vec_layout, data, NULL)
234 inline MFEM_ALWAYS_INLINE
241 inline MFEM_ALWAYS_INLINE
243 :
fieldEval(*func.FESpace(), func.GetData(), NULL)
248 template <
typename T_result_t,
typename c_layout_t,
typename c_data_t>
249 inline MFEM_ALWAYS_INLINE
250 void Eval(
const T_result_t &T,
const c_layout_t &l, c_data_t &c)
252 const int ne = T_result_t::ne;
253 const int vdim = FieldEval::vdim;
254 const int qpts = FieldEval::qpts;
255 MFEM_STATIC_ASSERT(c_layout_t::rank == 2,
"tensor rank must be 2");
256 MFEM_STATIC_ASSERT(c_layout_t::dim_1 == qpts,
"incompatible quadrature");
257 MFEM_STATIC_ASSERT(c_layout_t::dim_2 == ne,
"");
258 MFEM_STATIC_ASSERT(vdim == 1,
"vdim != 1 is not supported");
260 fieldEval.GetValues(T.first_elem_idx, l.template split_2<1,ne>(), c);
267 template <
typename IR,
typename coeff_t,
typename impl_traits_t>
270 static const int qpts = IR::qpts;
271 static const int ne = impl_traits_t::batch_size;
275 template <
bool is_const,
bool dummy>
struct Aux;
278 template <
bool dummy>
struct Aux<true,dummy>
280 typedef struct { } result_t;
283 inline MFEM_ALWAYS_INLINE
Aux(
const IR &int_rule,
const coeff_t &c)
285 c.Eval(
true, cw.
layout, cw);
286 int_rule.ApplyWeights(cw);
289 template <
typename T_result_t>
290 inline MFEM_ALWAYS_INLINE
291 void Eval(
const T_result_t &F, result_t &res) { }
293 inline MFEM_ALWAYS_INLINE
300 template <
bool dummy>
struct Aux<false,dummy>
303 #ifdef MFEM_TEMPLATE_INTRULE_COEFF_PRECOMP 310 #ifdef MFEM_TEMPLATE_INTRULE_COEFF_PRECOMP 311 inline MFEM_ALWAYS_INLINE
Aux(
const IR &int_rule,
const coeff_t &c)
314 int_rule.template AssignWeights<AssignOp::Set>(w.
layout, w);
317 inline MFEM_ALWAYS_INLINE
Aux(
const IR &int_rule,
const coeff_t &c)
318 : int_rule(int_rule), c(c) { }
320 template <
typename T_result_t>
321 inline MFEM_ALWAYS_INLINE
324 c.Eval(F, res.
layout, res);
325 #ifdef MFEM_TEMPLATE_INTRULE_COEFF_PRECOMP 326 for (
int i = 0; i <
ne; i++)
328 TAssign<AssignOp::Mult>(res.
layout.
ind2(i), res,
332 int_rule.template AssignWeights<AssignOp::Mult>(res.
layout, res);
336 inline MFEM_ALWAYS_INLINE
348 #endif // MFEM_TEMPLATE_COEFFICIENT static const bool uses_Jacobians
TFunctionCoefficient(Func &F_)
Constructor for the case when Func has data members.
MFEM_ALWAYS_INLINE void Eval(const T_result_t &F, result_t &res)
Class for grid function - Vector with associated FE space.
Templated coefficient classes, cf. coefficient.?pp.
FieldEval::FESpace_type FESpace_type
TMatrix< qpts, ne, vcomplex_t > result_t
GridFunction coefficient class.
static const layout_type layout
static const bool uses_coordinates
MFEM_ALWAYS_INLINE void Eval(const T_result_t &T, const c_layout_t &l, c_data_t &c) const
static MFEM_ALWAYS_INLINE void Eval(Func &F, const T_result_t &T, const c_layout_t &l, c_data_t &c)
MFEM_ALWAYS_INLINE TGridFunctionCoefficient(const GridFunction &func)
FieldEval::complex_type complex_type
static const bool uses_attributes
static const bool uses_element_idxs
TFunctionCoefficient()
Constructor for the case when Func has no data members.
FieldEval::ShapeEval_type ShapeEval_type
static const bool uses_element_idxs
static const bool uses_coordinates
static const bool is_const
MFEM_ALWAYS_INLINE void Eval(const T_result_t &T, const c_layout_t &l, c_data_t &c)
impl_traits_t::vcomplex_t vcomplex_t
TMatrix< qpts, 1, complex_type > cw
MFEM_ALWAYS_INLINE Aux(const IR &int_rule, const coeff_t &c)
MFEM_ALWAYS_INLINE Aux(const IR &int_rule, const coeff_t &c)
static MFEM_ALWAYS_INLINE void Eval(Func &F, const T_result_t &T, const c_layout_t &l, c_data_t &c)
MFEM_ALWAYS_INLINE TGridFunctionCoefficient(const FESpace_type &tfes, const ShapeEval_type &shapeEval, const VecLayout_type &vec_layout, const complex_type *data)
TPiecewiseConstCoefficient(const Vector &constants)
Note: in the input array index i corresponds to mesh attribute i+1.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
static const bool is_const
TConstantCoefficient(complex_t val)
TMatrix< qpts, 1, typename IR::real_type > w
FieldEval::VecLayout_type VecLayout_type
Aux< coeff_t::is_const, true > Type
MFEM_ALWAYS_INLINE void Eval(const T_result_t &F, result_t &res)
static MFEM_ALWAYS_INLINE void Eval(Func &F, const T_result_t &T, const c_layout_t &l, c_data_t &c)
static OffsetStridedLayout1D< N1, S1 > ind2(int i2)
MFEM_ALWAYS_INLINE void Eval(const T_result_t &T, const c_layout_t &l, c_data_t &c)
coeff_t::complex_type complex_type
MFEM_ALWAYS_INLINE TGridFunctionCoefficient(const FieldEval &fE, const complex_type *data)
static const bool uses_attributes
MFEM_ALWAYS_INLINE TGridFunctionCoefficient(const FiniteElementSpace &fes, const complex_type *data)
MFEM_ALWAYS_INLINE void Eval(const T_result_t &T, const c_layout_t &l, c_data_t &c)