12#ifndef MFEM_TEMPLATE_BILININTEG
13#define MFEM_TEMPLATE_BILININTEG
25template <
typename coeff_t,
template<
int,
int,
typename>
class kernel_t>
31 template <
int SDim,
int Dim,
typename complex_t>
32 struct kernel {
typedef kernel_t<SDim,Dim,complex_t>
type; };
41template <
int SDim,
int Dim,
typename complex_t>
69 template <
typename IR,
typename coeff_t,
typename impl_traits_t>
82 template <
typename T_result_t,
typename Q_t,
typename q_t,
84 static inline MFEM_ALWAYS_INLINE
85 void Action(
const int k,
const T_result_t &F,
86 const Q_t &Q,
const q_t &q, S_data_t &R)
88 typedef typename T_result_t::Jt_type::data_type real_t_;
89 const int M = S_data_t::eval_type::qpts;
90 const int NC = S_data_t::eval_type::vdim;
91 MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
92 "incompatible dimensions");
93 MFEM_FLOPS_ADD(M*(1+NC));
94 for (
int i = 0; i < M; i++)
98 for (
int j = 0; j < NC; j++)
100 R.val_qpts(i,j,k) *= wi;
115 template <
typename T_result_t,
typename Q_t,
typename q_t,
int qpts>
116 static inline MFEM_ALWAYS_INLINE
120 typedef typename T_result_t::Jt_type::data_type real_t_;
122 const int M = T_result_t::Jt_type::layout_type::dim_1;
123 MFEM_STATIC_ASSERT(qpts == M,
"incompatible dimensions");
125 for (
int i = 0; i < M; i++)
127 A[i] = Q.get(q,i,k) *
TDet<real_t_>(F.Jt.layout.ind14(i,k), F.Jt);
137 template <
int qpts,
typename S_data_t>
138 static inline MFEM_ALWAYS_INLINE
141 const int M = S_data_t::eval_type::qpts;
142 const int NC = S_data_t::eval_type::vdim;
143 MFEM_STATIC_ASSERT(qpts == M,
"incompatible dimensions");
144 MFEM_FLOPS_ADD(M*NC);
145 for (
int i = 0; i < M; i++)
147 for (
int j = 0; j < NC; j++)
149 R.val_qpts(i,j,k) *= A[i];
159template <
int SDim,
int Dim,
typename complex_t>
163template <
typename complex_t>
169 static const bool uses_Jacobians =
true;
173 static const bool in_values =
false;
174 static const bool in_gradients =
true;
175 static const bool out_values =
false;
176 static const bool out_gradients =
true;
191 template <
typename IR,
typename coeff_t,
typename impl_traits_t>
204 template <
typename T_result_t,
typename Q_t,
typename q_t,
206 static inline MFEM_ALWAYS_INLINE
207 void Action(
const int k,
const T_result_t &F,
208 const Q_t &Q,
const q_t &q, S_data_t &R)
210 const int M = S_data_t::eval_type::qpts;
211 const int NC = S_data_t::eval_type::vdim;
212 MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
213 "incompatible dimensions");
214 MFEM_FLOPS_ADD(M*(1+NC));
215 for (
int i = 0; i < M; i++)
217 const complex_t wi = Q.get(q,i,k) / F.Jt(i,0,0,k);
218 for (
int j = 0; j < NC; j++)
220 R.grad_qpts(i,0,j,k) *= wi;
239 template <
typename T_result_t,
typename Q_t,
typename q_t,
typename asm_type>
240 static inline MFEM_ALWAYS_INLINE
242 const Q_t &Q,
const q_t &q, asm_type &A)
244 const int M = T_result_t::Jt_type::layout_type::dim_1;
245 MFEM_STATIC_ASSERT(asm_type::layout_type::dim_1 == M,
246 "incompatible dimensions");
248 for (
int i = 0; i < M; i++)
251 A[i] = Q.get(q,i,k) / F.Jt(i,0,0,k);
261 template <
int qpts,
typename S_data_t>
262 static inline MFEM_ALWAYS_INLINE
266 const int M = S_data_t::eval_type::qpts;
267 const int NC = S_data_t::eval_type::vdim;
268 MFEM_STATIC_ASSERT(qpts == M,
"incompatible dimensions");
269 MFEM_FLOPS_ADD(M*NC);
270 for (
int i = 0; i < M; i++)
272 for (
int j = 0; j < NC; j++)
274 R.grad_qpts(i,0,j,k) *= A(i,0);
281template <
typename complex_t>
287 static const bool uses_Jacobians =
true;
291 static const bool in_values =
false;
292 static const bool in_gradients =
true;
293 static const bool out_values =
false;
294 static const bool out_gradients =
true;
309 template <
typename IR,
typename coeff_t,
typename impl_traits_t>
323 template <
typename T_result_t,
typename Q_t,
typename q_t,
325 static inline MFEM_ALWAYS_INLINE
326 void Action(
const int k,
const T_result_t &F,
327 const Q_t &Q,
const q_t &q, S_data_t &R)
329 const int M = S_data_t::eval_type::qpts;
330 const int NC = S_data_t::eval_type::vdim;
331 MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
332 "incompatible dimensions");
333 MFEM_FLOPS_ADD(M*(4+NC*14));
334 for (
int i = 0; i < M; i++)
336 typedef typename T_result_t::Jt_type::data_type real_t_;
337 const real_t_ J11 = F.Jt(i,0,0,k);
338 const real_t_ J12 = F.Jt(i,1,0,k);
339 const real_t_ J21 = F.Jt(i,0,1,k);
340 const real_t_ J22 = F.Jt(i,1,1,k);
341 const complex_t w_det_J = Q.get(q,i,k) / (J11 * J22 - J21 * J12);
342 for (
int j = 0; j < NC; j++)
344 const complex_t x1 = R.grad_qpts(i,0,j,k);
345 const complex_t x2 = R.grad_qpts(i,1,j,k);
347 const complex_t z1 = J22 * x1 - J21 * x2;
348 const complex_t z2 = J11 * x2 - J12 * x1;
349 R.grad_qpts(i,0,j,k) = w_det_J * (J22 * z1 - J12 * z2);
350 R.grad_qpts(i,1,j,k) = w_det_J * (J11 * z2 - J21 * z1);
367 template <
typename T_result_t,
typename Q_t,
typename q_t,
typename asm_type>
368 static inline MFEM_ALWAYS_INLINE
370 const Q_t &Q,
const q_t &q, asm_type &A)
372 typedef typename T_result_t::Jt_type::data_type real_t_;
373 const int M = T_result_t::Jt_type::layout_type::dim_1;
374 MFEM_STATIC_ASSERT(asm_type::layout_type::dim_1 == M,
375 "incompatible dimensions");
376 MFEM_FLOPS_ADD(16*M);
377 const bool Symm = (asm_type::layout_type::rank == 2);
378 for (
int i = 0; i < M; i++)
380 const real_t_ J11 = F.Jt(i,0,0,k);
381 const real_t_ J12 = F.Jt(i,1,0,k);
382 const real_t_ J21 = F.Jt(i,0,1,k);
383 const real_t_ J22 = F.Jt(i,1,1,k);
384 const complex_t w_det_J = Q.get(q,i,k) / (J11 * J22 - J21 * J12);
385 internal::MatrixOps<2,2>::Symm<Symm>::Set(
387 + w_det_J * (J12*J12 + J22*J22),
388 - w_det_J * (J11*J12 + J21*J22),
389 + w_det_J * (J11*J11 + J21*J21)
401 template <
int qpts,
typename S_data_t>
402 static inline MFEM_ALWAYS_INLINE
406 const int M = S_data_t::eval_type::qpts;
407 const int NC = S_data_t::eval_type::vdim;
408 MFEM_STATIC_ASSERT(qpts == M,
"incompatible dimensions");
409 MFEM_FLOPS_ADD(6*M*NC);
410 for (
int i = 0; i < M; i++)
412 const complex_t A11 = A(i,0);
413 const complex_t A21 = A(i,1);
414 const complex_t A22 = A(i,2);
415 for (
int j = 0; j < NC; j++)
417 const complex_t x1 = R.grad_qpts(i,0,j,k);
418 const complex_t x2 = R.grad_qpts(i,1,j,k);
419 R.grad_qpts(i,0,j,k) = A11 * x1 + A21 * x2;
420 R.grad_qpts(i,1,j,k) = A21 * x1 + A22 * x2;
427template <
typename complex_t>
433 static const bool uses_Jacobians =
true;
437 static const bool in_values =
false;
438 static const bool in_gradients =
true;
439 static const bool out_values =
false;
440 static const bool out_gradients =
true;
455 template <
typename IR,
typename coeff_t,
typename impl_traits_t>
468 template <
typename T_result_t,
typename Q_t,
typename q_t,
470 static inline MFEM_ALWAYS_INLINE
471 void Action(
const int k,
const T_result_t &F,
472 const Q_t &Q,
const q_t &q, S_data_t &R)
474 const int M = S_data_t::eval_type::qpts;
475 const int NC = S_data_t::eval_type::vdim;
476 MFEM_STATIC_ASSERT(T_result_t::Jt_type::layout_type::dim_1 == M,
477 "incompatible dimensions");
479 for (
int i = 0; i < M; i++)
481 typedef typename T_result_t::Jt_type::data_type real_t_;
483 const complex_t w_det_J =
489 R.grad_qpts.
layout.ind14(i,k), R.grad_qpts,
494 R.grad_qpts.
layout.ind14(i,k), R.grad_qpts);
510 template <
typename T_result_t,
typename Q_t,
typename q_t,
typename asm_type>
511 static inline MFEM_ALWAYS_INLINE
513 const Q_t &Q,
const q_t &q, asm_type &A)
515 typedef typename T_result_t::Jt_type::data_type real_t_;
516 const int M = T_result_t::Jt_type::layout_type::dim_1;
517 MFEM_STATIC_ASSERT(asm_type::layout_type::dim_1 == M,
518 "incompatible dimensions");
519 MFEM_FLOPS_ADD(37*M);
520 const bool Symm = (asm_type::layout_type::rank == 2);
521 for (
int i = 0; i < M; i++)
528 internal::MatrixOps<3,3>::Symm<Symm>::Set(
530 u*(B(0,0)*B(0,0)+B(0,1)*B(0,1)+B(0,2)*B(0,2)),
531 u*(B(0,0)*B(1,0)+B(0,1)*B(1,1)+B(0,2)*B(1,2)),
532 u*(B(0,0)*B(2,0)+B(0,1)*B(2,1)+B(0,2)*B(2,2)),
533 u*(B(1,0)*B(1,0)+B(1,1)*B(1,1)+B(1,2)*B(1,2)),
534 u*(B(1,0)*B(2,0)+B(1,1)*B(2,1)+B(1,2)*B(2,2)),
535 u*(B(2,0)*B(2,0)+B(2,1)*B(2,1)+B(2,2)*B(2,2))
546 template <
int qpts,
typename S_data_t>
547 static inline MFEM_ALWAYS_INLINE
551 const int M = S_data_t::eval_type::qpts;
552 const int NC = S_data_t::eval_type::vdim;
553 MFEM_STATIC_ASSERT(qpts == M,
"incompatible dimensions");
554 MFEM_FLOPS_ADD(15*M*NC);
555 for (
int i = 0; i < M; i++)
557 const complex_t A11 = A(i,0);
558 const complex_t A21 = A(i,1);
559 const complex_t A31 = A(i,2);
560 const complex_t A22 = A(i,3);
561 const complex_t A32 = A(i,4);
562 const complex_t A33 = A(i,5);
563 for (
int j = 0; j < NC; j++)
565 const complex_t x1 = R.grad_qpts(i,0,j,k);
566 const complex_t x2 = R.grad_qpts(i,1,j,k);
567 const complex_t x3 = R.grad_qpts(i,2,j,k);
568 R.grad_qpts(i,0,j,k) = A11*x1 + A21*x2 + A31*x3;
569 R.grad_qpts(i,1,j,k) = A21*x1 + A22*x2 + A32*x3;
570 R.grad_qpts(i,2,j,k) = A31*x1 + A32*x2 + A33*x3;
The Integrator class combines a kernel and a coefficient.
TIntegrator(const coefficient_type &c)
MFEM_ALWAYS_INLINE void sMult_AB(const A_layout_t &A_layout, const A_data_t &A_data, const B_layout_t &B_layout, const B_data_t &B_data, const C_layout_t &C_layout, C_data_t &C_data)
real_t u(const Vector &xvec)
scalar_t TDet(const layout_t &a, const data_t &A)
scalar_t TAdjDet(const A_layout_t &a, const A_data_t &A, const B_layout_t &b, B_data_t &B)
static StridedLayout2D< N2, S2, N1, S1 > transpose_12()
IntRuleCoefficient< IR, coeff_t, impl_traits_t >::Type Type
Partially assembled data type for one element with the given number of quadrature points....
TTensor3< qpts, 1, 1, complex_t > type
Partially assembled data type for one element with the given number of quadrature points....
TMatrix< qpts, 1, complex_t > type
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TMatrix< qpts, 1, complex_t > &A, S_data_t &R)
Method for partially assembled action.
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
Method used for un-assembled (matrix free) action.
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, asm_type &A)
Method defining partial assembly. The pointwise Dim x Dim matrices are stored as symmetric (when asm_...
IntRuleCoefficient< IR, coeff_t, impl_traits_t >::Type Type
Partially assembled data type for one element with the given number of quadrature points....
TTensor3< qpts, 2, 2, complex_t > type
Partially assembled data type for one element with the given number of quadrature points....
TMatrix< qpts, 3, complex_t > type
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
Method used for un-assembled (matrix free) action.
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TMatrix< qpts, 3, complex_t > &A, S_data_t &R)
Method for partially assembled action.
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, asm_type &A)
Method defining partial assembly. The pointwise Dim x Dim matrices are stored as symmetric (when asm_...
IntRuleCoefficient< IR, coeff_t, impl_traits_t >::Type Type
Partially assembled data type for one element with the given number of quadrature points....
TTensor3< qpts, 3, 3, complex_t > type
Partially assembled data type for one element with the given number of quadrature points....
TMatrix< qpts, 6, complex_t > type
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TMatrix< qpts, 6, complex_t > &A, S_data_t &R)
Method for partially assembled action. A [M x Dim*(Dim+1)/2] - partially assembled Dim x Dim symmetri...
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
Method used for un-assembled (matrix free) action. grad_qpts = (w/det(J)) adj(J) adj(J)^t grad_qpts J...
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, asm_type &A)
Method defining partial assembly. The pointwise Dim x Dim matrices are stored as symmetric (when asm_...
kernel_t< SDim, Dim, complex_t > type
IntRuleCoefficient< IR, coeff_t, impl_traits_t >::Type Type
Partially assembled data type for one element with the given number of quadrature points....
TVector< qpts, complex_t > type
Partially assembled data type for one element with the given number of quadrature points....
TVector< qpts, complex_t > type
static MFEM_ALWAYS_INLINE void Assemble(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, TVector< qpts, complex_t > &A)
Method defining partial assembly. Result in A is the quadrature-point dependent part of element matri...
static MFEM_ALWAYS_INLINE void MultAssembled(const int k, const TVector< qpts, complex_t > &A, S_data_t &R)
Method for partially assembled action.
static const bool in_values
static const bool out_gradients
static const bool out_values
static const bool in_gradients
static MFEM_ALWAYS_INLINE void Action(const int k, const T_result_t &F, const Q_t &Q, const q_t &q, S_data_t &R)
Method used for un-assembled (matrix free) action.
static const bool uses_Jacobians
Needed for the TElementTransformation::Result class.
static const layout_type layout
void Scale(const data_t scale)