12#ifndef MFEM_TEMPLATE_EVALUATOR
13#define MFEM_TEMPLATE_EVALUATOR
32template <
class FE,
class IR,
bool TP,
typename real_t>
36template <
class FE,
class IR,
typename real_t>
40 static const int DOF = FE::dofs;
41 static const int NIP = IR::qpts;
42 static const int DIM = FE::dim;
53 fe.CalcShapes(IR::GetIntRule(), B.
data, G.
data);
63 template <
typename dof_layout_t,
typename dof_data_t,
64 typename qpt_layout_t,
typename qpt_data_t>
65 inline MFEM_ALWAYS_INLINE
66 void Calc(
const dof_layout_t &dof_layout,
const dof_data_t &dof_data,
67 const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data)
const
69 MFEM_STATIC_ASSERT(dof_layout_t::rank == 2 &&
70 dof_layout_t::dim_1 == DOF,
71 "invalid dof_layout_t.");
72 MFEM_STATIC_ASSERT(qpt_layout_t::rank == 2 &&
73 qpt_layout_t::dim_1 == NIP,
74 "invalid qpt_layout_t.");
75 MFEM_STATIC_ASSERT(dof_layout_t::dim_2 == qpt_layout_t::dim_2,
76 "incompatible dof- and qpt- layouts.");
80 qpt_layout, qpt_data);
86 typename qpt_layout_t,
typename qpt_data_t,
87 typename dof_layout_t,
typename dof_data_t>
88 inline MFEM_ALWAYS_INLINE
89 void CalcT(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
90 const dof_layout_t &dof_layout, dof_data_t &dof_data)
const
92 MFEM_STATIC_ASSERT(dof_layout_t::rank == 2 &&
93 dof_layout_t::dim_1 == DOF,
94 "invalid dof_layout_t.");
95 MFEM_STATIC_ASSERT(qpt_layout_t::rank == 2 &&
96 qpt_layout_t::dim_1 == NIP,
97 "invalid qpt_layout_t.");
98 MFEM_STATIC_ASSERT(dof_layout_t::dim_2 == qpt_layout_t::dim_2,
99 "incompatible dof- and qpt- layouts.");
102 qpt_layout, qpt_data,
103 dof_layout, dof_data);
108 template <
typename dof_layout_t,
typename dof_data_t,
109 typename grad_layout_t,
typename grad_data_t>
110 inline MFEM_ALWAYS_INLINE
112 const dof_data_t &dof_data,
113 const grad_layout_t &grad_layout,
114 grad_data_t &grad_data)
const
116 MFEM_STATIC_ASSERT(dof_layout_t::rank == 2 &&
117 dof_layout_t::dim_1 == DOF,
118 "invalid dof_layout_t.");
119 MFEM_STATIC_ASSERT(grad_layout_t::rank == 3 &&
120 grad_layout_t::dim_1 == NIP &&
121 grad_layout_t::dim_2 ==
DIM,
122 "invalid grad_layout_t.");
123 MFEM_STATIC_ASSERT(dof_layout_t::dim_2 == grad_layout_t::dim_3,
124 "incompatible dof- and grad- layouts.");
127 dof_layout, dof_data,
128 grad_layout.merge_12(), grad_data);
134 typename grad_layout_t,
typename grad_data_t,
135 typename dof_layout_t,
typename dof_data_t>
136 inline MFEM_ALWAYS_INLINE
138 const grad_data_t &grad_data,
139 const dof_layout_t &dof_layout,
140 dof_data_t &dof_data)
const
142 MFEM_STATIC_ASSERT(dof_layout_t::rank == 2 &&
143 dof_layout_t::dim_1 == DOF,
144 "invalid dof_layout_t.");
145 MFEM_STATIC_ASSERT(grad_layout_t::rank == 3 &&
146 grad_layout_t::dim_1 == NIP &&
147 grad_layout_t::dim_2 ==
DIM,
148 "invalid grad_layout_t.");
149 MFEM_STATIC_ASSERT(dof_layout_t::dim_2 == grad_layout_t::dim_3,
150 "incompatible dof- and grad- layouts.");
153 grad_layout.merge_12(), grad_data,
154 dof_layout, dof_data);
160 template <
typename qpt_layout_t,
typename qpt_data_t,
161 typename M_layout_t,
typename M_data_t>
162 inline MFEM_ALWAYS_INLINE
163 void Assemble(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
164 const M_layout_t &M_layout, M_data_t &M_data)
const
171 qpt_layout.template split_1<1,NIP>(), qpt_data,
172 M_layout.template split_1<DOF,1>(), M_data);
176 qpt_layout.template split_1<1,NIP>(), qpt_data,
177 M_layout.template split_1<DOF,1>(), M_data);
184 template <
typename qpt_layout_t,
typename qpt_data_t,
185 typename D_layout_t,
typename D_data_t>
186 inline MFEM_ALWAYS_INLINE
188 const qpt_data_t &qpt_data,
189 const D_layout_t &D_layout,
190 D_data_t &D_data)
const
192 const int NC = qpt_layout_t::dim_4;
193 typedef typename qpt_data_t::data_type entry_type;
195 for (
int k = 0; k < NC; k++)
199 for (
int j = 0; j < NIP; j++)
213template <
int Dim,
int DOF,
int NIP,
typename real_t>
217template <
int DOF,
int NIP,
typename real_t>
221 static const int TDOF = DOF;
231 template <
typename dof_layout_t,
typename dof_data_t,
232 typename qpt_layout_t,
typename qpt_data_t>
233 inline MFEM_ALWAYS_INLINE
234 void Calc(
const dof_layout_t &dof_layout,
const dof_data_t &dof_data,
235 const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data)
const
238 dof_layout, dof_data,
239 qpt_layout, qpt_data);
245 typename qpt_layout_t,
typename qpt_data_t,
246 typename dof_layout_t,
typename dof_data_t>
247 inline MFEM_ALWAYS_INLINE
248 void CalcT(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
249 const dof_layout_t &dof_layout, dof_data_t &dof_data)
const
252 qpt_layout, qpt_data,
253 dof_layout, dof_data);
258 template <
typename dof_layout_t,
typename dof_data_t,
259 typename grad_layout_t,
typename grad_data_t>
260 inline MFEM_ALWAYS_INLINE
262 const dof_data_t &dof_data,
263 const grad_layout_t &grad_layout,
264 grad_data_t &grad_data)
const
268 dof_layout, dof_data,
269 grad_layout.merge_12(), grad_data);
275 typename grad_layout_t,
typename grad_data_t,
276 typename dof_layout_t,
typename dof_data_t>
277 inline MFEM_ALWAYS_INLINE
279 const grad_data_t &grad_data,
280 const dof_layout_t &dof_layout,
281 dof_data_t &dof_data)
const
286 grad_layout.merge_12(), grad_data,
287 dof_layout, dof_data);
292 template <
typename qpt_layout_t,
typename qpt_data_t,
293 typename M_layout_t,
typename M_data_t>
294 inline MFEM_ALWAYS_INLINE
295 void Assemble(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
296 const M_layout_t &M_layout, M_data_t &M_data)
const
303 qpt_layout.template split_1<1,NIP>(), qpt_data,
304 M_layout.template split_1<DOF,1>(), M_data);
308 qpt_layout.template split_1<1,NIP>(), qpt_data,
309 M_layout.template split_1<DOF,1>(), M_data);
316 template <
typename qpt_layout_t,
typename qpt_data_t,
317 typename D_layout_t,
typename D_data_t>
318 inline MFEM_ALWAYS_INLINE
320 const qpt_data_t &qpt_data,
321 const D_layout_t &D_layout,
322 D_data_t &D_data)
const
329 qpt_layout.merge_12().merge_23().template split_1<1,NIP>(), qpt_data,
330 D_layout.template split_1<DOF,1>(), D_data);
334 qpt_layout.merge_12().merge_23().template split_1<1,NIP>(), qpt_data,
335 D_layout.template split_1<DOF,1>(), D_data);
341template <
int DOF,
int NIP,
typename real_t>
349 static const int TDOF = DOF*DOF;
350 static const int TNIP = NIP*NIP;
354 template <
bool Dx,
bool Dy,
355 typename dof_layout_t,
typename dof_data_t,
356 typename qpt_layout_t,
typename qpt_data_t>
357 inline MFEM_ALWAYS_INLINE
358 void Calc(
const dof_layout_t &dof_layout,
const dof_data_t &dof_data,
359 const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data)
const
361 const int NC = dof_layout_t::dim_2;
362 typedef typename qpt_data_t::data_type entry_type;
368 dof_layout.
template split_1<DOF,DOF>(), dof_data,
373 qpt_layout.template split_1<NIP,NIP>(), qpt_data);
378 template <
typename dof_layout_t,
typename dof_data_t,
379 typename qpt_layout_t,
typename qpt_data_t>
380 inline MFEM_ALWAYS_INLINE
381 void Calc(
const dof_layout_t &dof_layout,
const dof_data_t &dof_data,
382 const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data)
const
384 Calc<false,false>(dof_layout, dof_data, qpt_layout, qpt_data);
387 template <
bool Dx,
bool Dy,
bool Add,
388 typename qpt_layout_t,
typename qpt_data_t,
389 typename dof_layout_t,
typename dof_data_t>
390 inline MFEM_ALWAYS_INLINE
391 void CalcT(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
392 const dof_layout_t &dof_layout, dof_data_t &dof_data)
const
394 const int NC = dof_layout_t::dim_2;
395 typedef typename qpt_data_t::data_type entry_type;
401 qpt_layout.template split_1<NIP,NIP>(), qpt_data,
406 dof_layout.template split_1<DOF,DOF>(), dof_data);
412 typename qpt_layout_t,
typename qpt_data_t,
413 typename dof_layout_t,
typename dof_data_t>
414 inline MFEM_ALWAYS_INLINE
415 void CalcT(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
416 const dof_layout_t &dof_layout, dof_data_t &dof_data)
const
418 CalcT<false,false,Add>(qpt_layout, qpt_data, dof_layout, dof_data);
423 template <
typename dof_layout_t,
typename dof_data_t,
424 typename grad_layout_t,
typename grad_data_t>
425 inline MFEM_ALWAYS_INLINE
427 const dof_data_t &dof_data,
428 const grad_layout_t &grad_layout,
429 grad_data_t &grad_data)
const
431 Calc<true,false>(dof_layout, dof_data,
432 grad_layout.ind2(0), grad_data);
433 Calc<false,true>(dof_layout, dof_data,
434 grad_layout.ind2(1), grad_data);
441 typename grad_layout_t,
typename grad_data_t,
442 typename dof_layout_t,
typename dof_data_t>
443 inline MFEM_ALWAYS_INLINE
445 const grad_data_t &grad_data,
446 const dof_layout_t &dof_layout,
447 dof_data_t &dof_data)
const
449 CalcT<true,false, Add>(grad_layout.ind2(0), grad_data,
450 dof_layout, dof_data);
451 CalcT<false,true,true>(grad_layout.ind2(1), grad_data,
452 dof_layout, dof_data);
457 template <
typename qpt_layout_t,
typename qpt_data_t,
458 typename M_layout_t,
typename M_data_t>
459 inline MFEM_ALWAYS_INLINE
460 void Assemble(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
461 const M_layout_t &M_layout, M_data_t &M_data)
const
463 const int NC = qpt_layout_t::dim_2;
464 typedef typename qpt_data_t::data_type entry_type;
473 qpt_layout.template split_1<NIP,NIP>(), qpt_data,
479 M_layout.merge_23().template split_12<DOF,DOF,DOF,DOF*NC>(), M_data);
485 qpt_layout.template split_1<NIP,NIP>(), qpt_data,
491 M_layout.merge_23().template split_12<DOF,DOF,DOF,DOF*NC>(), M_data);
496 for (
int k = 0; k < NC; k++)
500 qpt_layout.ind2(k).template split_1<NIP,NIP>().
501 template split_1<1,NIP>(), qpt_data,
502 B_1d.
layout, B_1d, F3.
layout.template split_1<1,NIP>(), F3);
513 M_layout.ind3(k).template split_1<DOF,DOF>(),
519 template <
int D1,
int D2,
bool Add,
520 typename qpt_layout_t,
typename qpt_data_t,
521 typename D_layout_t,
typename D_data_t>
522 inline MFEM_ALWAYS_INLINE
524 const qpt_data_t &qpt_data,
525 const D_layout_t &D_layout,
526 D_data_t &D_data)
const
528 const int NC = qpt_layout_t::dim_2;
529 typedef typename qpt_data_t::data_type entry_type;
536 Bt_1d.
layout, D1 == 0 ? Bt_1d : Gt_1d,
537 B_1d.
layout, D2 == 0 ? B_1d : G_1d,
538 qpt_layout.template split_1<NIP,NIP>(), qpt_data,
542 Bt_1d.
layout, D1 == 1 ? Bt_1d : Gt_1d,
543 B_1d.
layout, D2 == 1 ? B_1d : G_1d,
545 D_layout.merge_23().template split_12<DOF,DOF,DOF,DOF*NC>(), D_data);
551 template <
typename qpt_layout_t,
typename qpt_data_t,
552 typename D_layout_t,
typename D_data_t>
553 inline MFEM_ALWAYS_INLINE
555 const qpt_data_t &qpt_data,
556 const D_layout_t &D_layout,
557 D_data_t &D_data)
const
560 Assemble<0,0,false>(qpt_layout.ind23(0,0), qpt_data, D_layout, D_data);
561 Assemble<1,0,true >(qpt_layout.ind23(1,0), qpt_data, D_layout, D_data);
562 Assemble<0,1,true >(qpt_layout.ind23(0,1), qpt_data, D_layout, D_data);
563 Assemble<1,1,true >(qpt_layout.ind23(1,1), qpt_data, D_layout, D_data);
565 const int NC = qpt_layout_t::dim_4;
570 for (
int k = 0; k < NC; k++)
572 for (
int d1 = 0; d1 < 2; d1++)
578 template split_1<NIP,NIP>().
579 template split_1<1,NIP>(), qpt_data,
581 F3.
layout.template split_1<1,NIP>(), F3);
588 template split_1<NIP,NIP>().
589 template split_1<1,NIP>(), qpt_data,
591 F3.
layout.template split_1<1,NIP>(), F3);
604 D_layout.ind3(k).template split_1<DOF,DOF>(),
611 D_layout.ind3(k).template split_1<DOF,DOF>(),
621template <
int DOF,
int NIP,
typename real_t>
629 static const int TDOF = DOF*DOF*DOF;
630 static const int TNIP = NIP*NIP*NIP;
634 template <
bool Dx,
bool Dy,
bool Dz,
635 typename dof_layout_t,
typename dof_data_t,
636 typename qpt_layout_t,
typename qpt_data_t>
637 inline MFEM_ALWAYS_INLINE
638 void Calc(
const dof_layout_t &dof_layout,
const dof_data_t &dof_data,
639 const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data)
const
641 const int NC = dof_layout_t::dim_2;
642 typedef typename qpt_data_t::data_type entry_type;
648 dof_layout.template split_1<DOF,DOF*DOF>(), dof_data,
657 qpt_layout.template split_1<NIP*NIP,NIP>(), qpt_data);
662 template <
typename dof_layout_t,
typename dof_data_t,
663 typename qpt_layout_t,
typename qpt_data_t>
664 inline MFEM_ALWAYS_INLINE
665 void Calc(
const dof_layout_t &dof_layout,
const dof_data_t &dof_data,
666 const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data)
const
668 Calc<false,false,false>(dof_layout, dof_data, qpt_layout, qpt_data);
671 template <
bool Dx,
bool Dy,
bool Dz,
bool Add,
672 typename qpt_layout_t,
typename qpt_data_t,
673 typename dof_layout_t,
typename dof_data_t>
674 inline MFEM_ALWAYS_INLINE
675 void CalcT(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
676 const dof_layout_t &dof_layout, dof_data_t &dof_data)
const
678 const int NC = dof_layout_t::dim_2;
679 typedef typename qpt_data_t::data_type entry_type;
685 qpt_layout.template split_1<NIP*NIP,NIP>(), qpt_data,
694 dof_layout.template split_1<DOF,DOF*DOF>(), dof_data);
700 typename qpt_layout_t,
typename qpt_data_t,
701 typename dof_layout_t,
typename dof_data_t>
702 inline MFEM_ALWAYS_INLINE
703 void CalcT(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
704 const dof_layout_t &dof_layout, dof_data_t &dof_data)
const
706 CalcT<false,false,false,Add>(qpt_layout, qpt_data, dof_layout, dof_data);
711 template <
typename dof_layout_t,
typename dof_data_t,
712 typename grad_layout_t,
typename grad_data_t>
713 inline MFEM_ALWAYS_INLINE
715 const dof_data_t &dof_data,
716 const grad_layout_t &grad_layout,
717 grad_data_t &grad_data)
const
719 Calc<true,false,false>(dof_layout, dof_data,
720 grad_layout.ind2(0), grad_data);
721 Calc<false,true,false>(dof_layout, dof_data,
722 grad_layout.ind2(1), grad_data);
723 Calc<false,false,true>(dof_layout, dof_data,
724 grad_layout.ind2(2), grad_data);
733 typename grad_layout_t,
typename grad_data_t,
734 typename dof_layout_t,
typename dof_data_t>
735 inline MFEM_ALWAYS_INLINE
737 const grad_data_t &grad_data,
738 const dof_layout_t &dof_layout,
739 dof_data_t &dof_data)
const
741 CalcT<true,false,false, Add>(grad_layout.ind2(0), grad_data,
742 dof_layout, dof_data);
743 CalcT<false,true,false,true>(grad_layout.ind2(1), grad_data,
744 dof_layout, dof_data);
745 CalcT<false,false,true,true>(grad_layout.ind2(2), grad_data,
746 dof_layout, dof_data);
751 template <
typename qpt_layout_t,
typename qpt_data_t,
752 typename M_layout_t,
typename M_data_t>
753 inline MFEM_ALWAYS_INLINE
754 void Assemble(
const qpt_layout_t &qpt_layout,
const qpt_data_t &qpt_data,
755 const M_layout_t &M_layout, M_data_t &M_data)
const
757 const int NC = qpt_layout_t::dim_2;
758 typedef typename qpt_data_t::data_type entry_type;
768 qpt_layout.template split_1<NIP*NIP,NIP>(), qpt_data,
779 M_layout.merge_23().template split_12<DOF,DOF*DOF,DOF,DOF*DOF*NC>(),
785 qpt_layout.template split_1<NIP*NIP,NIP>(), qpt_data,
796 M_layout.merge_23().template split_12<DOF,DOF*DOF,DOF,DOF*DOF*NC>(),
801 template <
int D1,
int D2,
bool Add,
802 typename qpt_layout_t,
typename qpt_data_t,
803 typename D_layout_t,
typename D_data_t>
804 inline MFEM_ALWAYS_INLINE
806 const qpt_data_t &qpt_data,
807 const D_layout_t &D_layout,
808 D_data_t &D_data)
const
810 const int NC = qpt_layout_t::dim_2;
811 typedef typename qpt_data_t::data_type entry_type;
819 Bt_1d.
layout, D1 != 2 ? Bt_1d : Gt_1d,
820 B_1d.
layout, D2 != 2 ? B_1d : G_1d,
821 qpt_layout.template split_1<NIP*NIP,NIP>(), qpt_data,
825 Bt_1d.
layout, D1 != 1 ? Bt_1d : Gt_1d,
826 B_1d.
layout, D2 != 1 ? B_1d : G_1d,
831 Bt_1d.
layout, D1 != 0 ? Bt_1d : Gt_1d,
832 B_1d.
layout, D2 != 0 ? B_1d : G_1d,
834 D_layout.merge_23().template split_12<DOF,DOF*DOF,DOF,DOF*DOF*NC>(),
839 template <
typename qpt_layout_t,
typename qpt_data_t,
840 typename D_layout_t,
typename D_data_t>
841 inline MFEM_ALWAYS_INLINE
843 const qpt_layout_t &qpt_layout,
844 const qpt_data_t &qpt_data,
845 const D_layout_t &D_layout,
846 D_data_t &D_data)
const
848 const int NC = qpt_layout_t::dim_2;
856 Bt_1d.
layout, D1 != 2 ? Bt_1d : Gt_1d,
857 B_1d.
layout, D2 != 2 ? B_1d : G_1d,
858 qpt_layout.template split_1<NIP*NIP,NIP>(), qpt_data,
862 Bt_1d.
layout, D1 != 1 ? Bt_1d : Gt_1d,
863 B_1d.
layout, D2 != 1 ? B_1d : G_1d,
868 Bt_1d.
layout, D1 != 0 ? Bt_1d : Gt_1d,
869 B_1d.
layout, D2 != 0 ? B_1d : G_1d,
871 D_layout.merge_23().template split_12<DOF,DOF*DOF,DOF,DOF*DOF*NC>(),
879 template <
typename qpt_layout_t,
typename qpt_data_t,
880 typename D_layout_t,
typename D_data_t>
881 inline MFEM_ALWAYS_INLINE
883 const qpt_data_t &qpt_data,
884 const D_layout_t &D_layout,
885 D_data_t &D_data)
const
889 Assemble<0,0,false>(qpt_layout.ind23(0,0), qpt_data, D_layout, D_data);
890 Assemble<1,0,true >(qpt_layout.ind23(1,0), qpt_data, D_layout, D_data);
891 Assemble<2,0,true >(qpt_layout.ind23(2,0), qpt_data, D_layout, D_data);
892 Assemble<0,1,true >(qpt_layout.ind23(0,1), qpt_data, D_layout, D_data);
893 Assemble<1,1,true >(qpt_layout.ind23(1,1), qpt_data, D_layout, D_data);
894 Assemble<2,1,true >(qpt_layout.ind23(2,1), qpt_data, D_layout, D_data);
895 Assemble<0,2,true >(qpt_layout.ind23(0,2), qpt_data, D_layout, D_data);
896 Assemble<1,2,true >(qpt_layout.ind23(1,2), qpt_data, D_layout, D_data);
897 Assemble<2,2,true >(qpt_layout.ind23(2,2), qpt_data, D_layout, D_data);
900 for (
int d2 = 0; d2 < 3; d2++)
902 for (
int d1 = 0; d1 < 3; d1++)
904 Assemble(d1, d2, qpt_layout.ind23(d1,d2), qpt_data,
913template <
class FE,
class IR,
typename real_t>
920 using base_class::B_1d;
921 using base_class::Bt_1d;
922 using base_class::G_1d;
923 using base_class::Gt_1d;
928 fe.Calc1DShapes(IR::Get1DIntRule(), B_1d.data, G_1d.data);
930 B_1d.layout.transpose_12(), B_1d);
932 G_1d.layout.transpose_12(), G_1d);
939template <
class FE,
class IR,
typename real_t>
945 static const int dim = FE::dim;
946 static const int qpts = IR::qpts;
947 static const bool tensor_prod = FE::tensor_prod && IR::tensor_prod;
952 using base_class::Calc;
953 using base_class::CalcT;
954 using base_class::CalcGrad;
955 using base_class::CalcGradT;
966template <
typename FESpace_t,
typename VecLayout_t,
typename IR,
967 typename complex_t,
typename real_t>
979 inline MFEM_ALWAYS_INLINE
988 inline MFEM_ALWAYS_INLINE
995template <
typename FESpace_t,
typename VecLayout_t,
typename IR,
1010 static const int dofs = FE_type::dofs;
1011 static const int dim = FE_type::dim;
1013 static const int vdim = VecLayout_t::vec_dim;
1028 inline MFEM_ALWAYS_INLINE
1031 const complex_t *global_data_in, complex_t *global_data_out)
1038 inline MFEM_ALWAYS_INLINE
1040 const complex_t *global_data_in, complex_t *global_data_out)
1047 inline MFEM_ALWAYS_INLINE
1049 const complex_t *global_data_in, complex_t *global_data_out)
1061 inline MFEM_ALWAYS_INLINE
1068 template <
typename val_layout_t,
typename val_data_t>
1069 inline MFEM_ALWAYS_INLINE
1070 void GetValues(
int el,
const val_layout_t &l, val_data_t &vals)
1072 const int ne = val_layout_t::dim_3;
1080 template <
typename grad_layout_t,
typename grad_data_t>
1081 inline MFEM_ALWAYS_INLINE
1084 const int ne = grad_layout_t::dim_4;
1089 l.merge_34(), grad);
1094 template <
typename DataType>
1095 inline MFEM_ALWAYS_INLINE
1102 template <
typename DataType>
1103 inline MFEM_ALWAYS_INLINE
1110 template <
bool Add,
typename DataType>
1111 inline MFEM_ALWAYS_INLINE
1119 template <
bool Add,
typename DataType>
1120 inline MFEM_ALWAYS_INLINE
1127#ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
1128 template <
typename DataType>
1129 inline MFEM_ALWAYS_INLINE
1136 template <
bool Add,
typename DataType>
1137 inline MFEM_ALWAYS_INLINE
1139 typename DataType::vcomplex_t *loc_dofs)
1162 template<
int IOData,
typename impl_traits_t>
struct AData;
1164 template <
typename it_t>
struct AData<0,it_t>
1169 template <
typename it_t>
struct AData<1,it_t>
1171 static const int ne = it_t::batch_size;
1173#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1182 template <
typename it_t>
struct AData<2,it_t>
1184 static const int ne = it_t::batch_size;
1186#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1195 template <
typename it_t>
struct AData<3,it_t>
1197 static const int ne = it_t::batch_size;
1199#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1211 template <
int IData,
int OData,
typename it_t>
1222 template <
int Ops,
bool dummy>
struct Action;
1231 template <
typename vec_layout_t,
typename AData_t>
1232 static inline MFEM_ALWAYS_INLINE
1235#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1236 typename AData_t::val_dofs_t &val_dofs = D.val_dofs;
1238 typename AData_t::val_dofs_t val_dofs;
1240 T.
fespace.VectorExtract(l, T.
data_in, val_dofs.layout, val_dofs);
1241 T.
shapeEval.Calc(val_dofs.layout.merge_23(), val_dofs,
1242 D.val_qpts.layout.merge_23(), D.val_qpts);
1245 template <
bool Add,
typename vec_layout_t,
typename AData_t>
1246 static inline MFEM_ALWAYS_INLINE
1250#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1251 typename AData_t::val_dofs_t &val_dofs = D.val_dofs;
1253 typename AData_t::val_dofs_t val_dofs;
1256 D.val_qpts.layout.merge_23(), D.val_qpts,
1257 val_dofs.layout.merge_23(), val_dofs);
1258 T.
fespace.template VectorAssemble<Op>(
1259 val_dofs.layout, val_dofs, l, T.
data_out);
1262#ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
1263 template <
typename AData_t>
1264 static inline MFEM_ALWAYS_INLINE
1266 const typename AData_t::vcomplex_t *loc_dofs,
1269 T.
shapeEval.Calc(AData_t::val_dofs_t::layout.merge_23(), loc_dofs,
1270 D.val_qpts.layout.merge_23(), D.val_qpts);
1273 template <
bool Add,
typename AData_t>
1274 static inline MFEM_ALWAYS_INLINE
1276 typename AData_t::vcomplex_t *loc_dofs)
1279 D.val_qpts.layout.merge_23(), D.val_qpts,
1280 AData_t::val_dofs_t::layout.merge_23(), loc_dofs);
1287 template <
typename vec_layout_t,
typename AData_t>
1288 static inline MFEM_ALWAYS_INLINE
1291#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1292 typename AData_t::val_dofs_t &val_dofs = D.val_dofs;
1294 typename AData_t::val_dofs_t val_dofs;
1296 T.
fespace.VectorExtract(l, T.
data_in, val_dofs.layout, val_dofs);
1297 T.
shapeEval.CalcGrad(val_dofs.layout.merge_23(), val_dofs,
1298 D.grad_qpts.layout.merge_34(), D.grad_qpts);
1301 template <
bool Add,
typename vec_layout_t,
typename AData_t>
1302 static inline MFEM_ALWAYS_INLINE
1306#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1307 typename AData_t::val_dofs_t &val_dofs = D.val_dofs;
1309 typename AData_t::val_dofs_t val_dofs;
1312 D.grad_qpts.layout.merge_34(), D.grad_qpts,
1313 val_dofs.layout.merge_23(), val_dofs);
1314 T.
fespace.template VectorAssemble<Op>(
1315 val_dofs.layout, val_dofs, l, T.
data_out);
1318#ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
1319 template <
typename AData_t>
1320 static inline MFEM_ALWAYS_INLINE
1322 const typename AData_t::vcomplex_t *loc_dofs,
1325 T.
shapeEval.CalcGrad(AData_t::val_dofs_t::layout.merge_23(), loc_dofs,
1326 D.grad_qpts.layout.merge_34(), D.grad_qpts);
1329 template <
bool Add,
typename AData_t>
1330 static inline MFEM_ALWAYS_INLINE
1332 typename AData_t::vcomplex_t *loc_dofs)
1335 D.grad_qpts.layout.merge_34(), D.grad_qpts,
1336 AData_t::val_dofs_t::layout.merge_23(), loc_dofs);
1343 template <
typename vec_layout_t,
typename AData_t>
1344 static inline MFEM_ALWAYS_INLINE
1347#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1348 typename AData_t::val_dofs_t &val_dofs = D.val_dofs;
1350 typename AData_t::val_dofs_t val_dofs;
1352 T.
fespace.VectorExtract(l, T.
data_in, val_dofs.layout, val_dofs);
1353 T.
shapeEval.Calc(val_dofs.layout.merge_23(), val_dofs,
1354 D.val_qpts.layout.merge_23(), D.val_qpts);
1355 T.
shapeEval.CalcGrad(val_dofs.layout.merge_23(), val_dofs,
1356 D.grad_qpts.layout.merge_34(), D.grad_qpts);
1359 template <
bool Add,
typename vec_layout_t,
typename AData_t>
1360 static inline MFEM_ALWAYS_INLINE
1364#ifdef MFEM_TEMPLATE_FIELD_EVAL_DATA_HAS_DOFS
1365 typename AData_t::val_dofs_t &val_dofs = D.val_dofs;
1367 typename AData_t::val_dofs_t val_dofs;
1370 D.val_qpts.layout.merge_23(), D.val_qpts,
1371 val_dofs.layout.merge_23(), val_dofs);
1373 D.grad_qpts.layout.merge_34(), D.grad_qpts,
1374 val_dofs.layout.merge_23(), val_dofs);
1375 T.
fespace.template VectorAssemble<Op>(
1376 val_dofs.layout, val_dofs, l, T.
data_out);
1379#ifdef MFEM_TEMPLATE_ENABLE_SERIALIZE
1380 template <
typename AData_t>
1381 static inline MFEM_ALWAYS_INLINE
1383 const typename AData_t::vcomplex_t *loc_dofs,
1386 T.
shapeEval.Calc(AData_t::val_dofs_t::layout.merge_23(), loc_dofs,
1387 D.val_qpts.layout.merge_23(), D.val_qpts);
1388 T.
shapeEval.CalcGrad(AData_t::val_dofs_t::layout.merge_23(), loc_dofs,
1389 D.grad_qpts.layout.merge_34(), D.grad_qpts);
1392 template <
bool Add,
typename AData_t>
1393 static inline MFEM_ALWAYS_INLINE
1395 typename AData_t::vcomplex_t *loc_dofs)
1398 D.val_qpts.layout.merge_23(), D.val_qpts,
1399 AData_t::val_dofs_t::layout.merge_23(), loc_dofs);
1401 D.grad_qpts.layout.merge_34(), D.grad_qpts,
1402 AData_t::val_dofs_t::layout.merge_23(), loc_dofs);
1416 template <
typename qpt_layout_t,
typename qpt_data_t,
1417 typename M_layout_t,
typename M_data_t>
1418 static inline MFEM_ALWAYS_INLINE
1419 void Compute(
const qpt_layout_t &
a,
const qpt_data_t &A,
1422 ev.Assemble(
a.template split_1<qpts,1>(), A,
1423 m.template split_2<dofs,1>(), M);
1439 template <
typename qpt_layout_t,
typename qpt_data_t,
1440 typename M_layout_t,
typename M_data_t>
1441 static inline MFEM_ALWAYS_INLINE
1442 void Compute(
const qpt_layout_t &
a,
const qpt_data_t &A,
1445 ev.AssembleGradGrad(
a.template split_3<dim,1>(), A,
1446 m.template split_2<dofs,1>(), M);
1450 template <
typename kernel_t,
typename impl_traits_t>
struct Spec
Field evaluators – values of a given global FE grid function This is roughly speaking a templated ver...
MFEM_ALWAYS_INLINE FieldEvaluator_base(const FESpace_t &tfes, const ShapeEval_type &shape_eval, const VecLayout_t &vec_layout)
With this constructor, fespace is a shallow copy.
FESpace_t::FE_type FE_type
MFEM_ALWAYS_INLINE FieldEvaluator_base(const FE_type &fe, const FiniteElementSpace &fes)
This constructor creates new fespace, not a shallow copy.
ShapeEvaluator< FE_type, IR, real_t > ShapeEval_type
complex_t - dof/qpt data type, real_t - ShapeEvaluator (FE basis) data type
MFEM_ALWAYS_INLINE void GetValues(int el, const val_layout_t &l, val_data_t &vals)
val_layout_t is (qpts x vdim x NE)
MFEM_ALWAYS_INLINE VecLayout_type & VecLayout()
InOutData
Enumeration for the data type used by the Eval() and Assemble() methods. The types can be obtained by...
MFEM_ALWAYS_INLINE void EvalSerialized(const typename DataType::vcomplex_t *loc_dofs, DataType &F)
FESpace_t::FE_type FE_type
MFEM_ALWAYS_INLINE void Eval(int el, DataType &F)
MFEM_ALWAYS_INLINE FieldEvaluator(const FESpace_t &tfes, const ShapeEval_type &shape_eval, const VecLayout_type &vec_layout, const complex_t *global_data_in, complex_t *global_data_out)
With this constructor, fespace is a shallow copy of tfes.
FieldEvaluator< FESpace_t, VecLayout_t, IR, complex_t, real_t > T_type
VecLayout_t VecLayout_type
FieldEvaluator_base< FESpace_t, VecLayout_t, IR, complex_t, real_t > base_class
MFEM_ALWAYS_INLINE FieldEvaluator(const FiniteElementSpace &fes, const complex_t *global_data_in, complex_t *global_data_out)
This constructor creates a new fespace, not a shallow copy.
MFEM_ALWAYS_INLINE ShapeEval_type & ShapeEval()
MFEM_ALWAYS_INLINE void Assemble(int el, DataType &F)
const complex_t * data_in
MFEM_ALWAYS_INLINE void AssembleSerialized(const DataType &F, typename DataType::vcomplex_t *loc_dofs)
MFEM_ALWAYS_INLINE FESpace_type & FESpace()
MFEM_ALWAYS_INLINE void Eval(DataType &F)
MFEM_ALWAYS_INLINE void Assemble(DataType &F)
MFEM_ALWAYS_INLINE FieldEvaluator(const FieldEvaluator &f, const complex_t *global_data_in, complex_t *global_data_out)
With this constructor, fespace is a shallow copy of f.fespace.
ShapeEvaluator< FE_type, IR, real_t > ShapeEval_type
MFEM_ALWAYS_INLINE void SetElement(int el)
MFEM_ALWAYS_INLINE void GetGradients(int el, const grad_layout_t &l, grad_data_t &grad)
grad_layout_t is (qpts x dim x vdim x NE)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
MFEM_ALWAYS_INLINE void CalcGrad(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const grad_layout_t &grad_layout, grad_data_t &grad_data) const
Multi-component gradient evaluation from DOFs to quadrature points. dof_layout is (DOF x NumComp) and...
TMatrix< NIP, DOF, real_t, true > B
MFEM_ALWAYS_INLINE void CalcGradT(const grad_layout_t &grad_layout, const grad_data_t &grad_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component gradient evaluation transpose from quadrature points to DOFs. grad_layout is (NIP x D...
TTensor3< NIP, DIM, DOF, real_t, true > G
ShapeEvaluator_base(const FE &fe)
MFEM_ALWAYS_INLINE void CalcT(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component shape evaluation transpose from quadrature points to DOFs. qpt_layout is (NIP x NumCo...
MFEM_ALWAYS_INLINE void Calc(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data) const
Multi-component shape evaluation from DOFs to quadrature points. dof_layout is (DOF x NumComp) and qp...
TTensor3< DOF, NIP, DIM, real_t > Gt
MFEM_ALWAYS_INLINE void AssembleGradGrad(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
Multi-component assemble of grad-grad element matrices. qpt_layout is (NIP x DIM x DIM x NumComp),...
TMatrix< DOF, NIP, real_t, true > Bt
MFEM_ALWAYS_INLINE void Assemble(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const M_layout_t &M_layout, M_data_t &M_data) const
Multi-component assemble. qpt_layout is (NIP x NumComp), M_layout is (DOF x DOF x NumComp)
TProductShapeEvaluator< FE::dim, FE::dofs_1d, IR::qpts_1d, real_t > base_class
ShapeEvaluator_base(const FE &fe)
Shape evaluators – values of basis functions on the reference element.
General ShapeEvaluator for any scalar FE type (L2 or H1)
ShapeEvaluator_base< FE, IR, tensor_prod, real_t > base_class
static const bool tensor_prod
ShapeEvaluator(const FE &fe)
TMatrix< DOF, NIP, real_t, true > Bt_1d
MFEM_ALWAYS_INLINE void CalcGradT(const grad_layout_t &grad_layout, const grad_data_t &grad_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component gradient evaluation transpose from quadrature points to DOFs. grad_layout is (NIP x D...
MFEM_ALWAYS_INLINE void CalcT(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component shape evaluation transpose from quadrature points to DOFs. qpt_layout is (NIP x NumCo...
TMatrix< NIP, DOF, real_t, true > B_1d
MFEM_ALWAYS_INLINE void Calc(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data) const
Multi-component shape evaluation from DOFs to quadrature points. dof_layout is (DOF x NumComp) and qp...
MFEM_ALWAYS_INLINE void Assemble(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const M_layout_t &M_layout, M_data_t &M_data) const
Multi-component assemble. qpt_layout is (NIP x NumComp), M_layout is (DOF x DOF x NumComp)
MFEM_ALWAYS_INLINE void CalcGrad(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const grad_layout_t &grad_layout, grad_data_t &grad_data) const
Multi-component gradient evaluation from DOFs to quadrature points. dof_layout is (DOF x NumComp) and...
MFEM_ALWAYS_INLINE void AssembleGradGrad(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
Multi-component assemble of grad-grad element matrices. qpt_layout is (NIP x DIM x DIM x NumComp),...
MFEM_ALWAYS_INLINE void Assemble(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
TMatrix< DOF, NIP, real_t, true > Bt_1d
TMatrix< NIP, DOF, real_t, true > B_1d
MFEM_ALWAYS_INLINE void CalcT(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component shape evaluation transpose from quadrature points to DOFs. qpt_layout is (TNIP x NumC...
MFEM_ALWAYS_INLINE void CalcT(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
MFEM_ALWAYS_INLINE void AssembleGradGrad(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
Multi-component assemble of grad-grad element matrices. qpt_layout is (TNIP x DIM x DIM x NumComp),...
MFEM_ALWAYS_INLINE void Calc(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data) const
MFEM_ALWAYS_INLINE void CalcGrad(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const grad_layout_t &grad_layout, grad_data_t &grad_data) const
Multi-component gradient evaluation from DOFs to quadrature points. dof_layout is (TDOF x NumComp) an...
MFEM_ALWAYS_INLINE void Assemble(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const M_layout_t &M_layout, M_data_t &M_data) const
Multi-component assemble. qpt_layout is (TNIP x NumComp), M_layout is (TDOF x TDOF x NumComp)
MFEM_ALWAYS_INLINE void Calc(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data) const
Multi-component shape evaluation from DOFs to quadrature points. dof_layout is (TDOF x NumComp) and q...
MFEM_ALWAYS_INLINE void CalcGradT(const grad_layout_t &grad_layout, const grad_data_t &grad_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component gradient evaluation transpose from quadrature points to DOFs. grad_layout is (TNIP x ...
MFEM_ALWAYS_INLINE void Assemble(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
TMatrix< DOF, NIP, real_t, true > Bt_1d
MFEM_ALWAYS_INLINE void CalcGrad(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const grad_layout_t &grad_layout, grad_data_t &grad_data) const
Multi-component gradient evaluation from DOFs to quadrature points. dof_layout is (TDOF x NumComp) an...
MFEM_ALWAYS_INLINE void CalcT(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
MFEM_ALWAYS_INLINE void Assemble(int D1, int D2, const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
MFEM_ALWAYS_INLINE void AssembleGradGrad(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const D_layout_t &D_layout, D_data_t &D_data) const
Multi-component assemble of grad-grad element matrices. qpt_layout is (TNIP x DIM x DIM x NumComp),...
MFEM_ALWAYS_INLINE void Assemble(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const M_layout_t &M_layout, M_data_t &M_data) const
Multi-component assemble. qpt_layout is (TNIP x NumComp), M_layout is (TDOF x TDOF x NumComp)
MFEM_ALWAYS_INLINE void Calc(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data) const
MFEM_ALWAYS_INLINE void Calc(const dof_layout_t &dof_layout, const dof_data_t &dof_data, const qpt_layout_t &qpt_layout, qpt_data_t &qpt_data) const
Multi-component shape evaluation from DOFs to quadrature points. dof_layout is (TDOF x NumComp) and q...
MFEM_ALWAYS_INLINE void CalcGradT(const grad_layout_t &grad_layout, const grad_data_t &grad_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component gradient evaluation transpose from quadrature points to DOFs. grad_layout is (TNIP x ...
TMatrix< NIP, DOF, real_t, true > B_1d
MFEM_ALWAYS_INLINE void CalcT(const qpt_layout_t &qpt_layout, const qpt_data_t &qpt_data, const dof_layout_t &dof_layout, dof_data_t &dof_data) const
Multi-component shape evaluation transpose from quadrature points to DOFs. qpt_layout is (TNIP x NumC...
MFEM_ALWAYS_INLINE void Mult_1_2(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)
MFEM_ALWAYS_INLINE void Mult_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)
MFEM_ALWAYS_INLINE void Mult_2_1(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)
MFEM_ALWAYS_INLINE void TensorProduct(const A_layout_t &a, const A_data_t &A, const B_layout_t &b, const B_data_t &B, const C_layout_t &c, C_data_t &C)
void TAssign(const A_layout_t &A_layout, A_data_t &A_data, const scalar_t value)
MFEM_ALWAYS_INLINE void TensorAssemble(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)
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
TTensor3< dofs, vdim, ne, vcomplex_t, true > val_dofs_t
TTensor3< qpts, vdim, ne, vcomplex_t > val_qpts
it_t::vcomplex_t vcomplex_t
TTensor3< dofs, vdim, ne, vcomplex_t > val_dofs_t
TTensor3< dofs, vdim, ne, vcomplex_t, true > val_dofs_t
it_t::vcomplex_t vcomplex_t
TTensor3< dofs, vdim, ne, vcomplex_t > val_dofs_t
TTensor4< qpts, dim, vdim, ne, vcomplex_t > grad_qpts
it_t::vcomplex_t vcomplex_t
TTensor3< dofs, vdim, ne, vcomplex_t > val_dofs_t
TTensor3< dofs, vdim, ne, vcomplex_t, true > val_dofs_t
TTensor3< qpts, vdim, ne, vcomplex_t, true > val_qpts
TTensor4< qpts, dim, vdim, ne, vcomplex_t > grad_qpts
Auxiliary templated struct AData, used by the Eval() and Assemble() methods.
static MFEM_ALWAYS_INLINE void Eval(const vec_layout_t &l, T_type &T, AData_t &D)
static MFEM_ALWAYS_INLINE void AssembleSerialized(T_type &T, const AData_t &D, typename AData_t::vcomplex_t *loc_dofs)
static MFEM_ALWAYS_INLINE void Assemble(const vec_layout_t &l, T_type &T, AData_t &D)
static MFEM_ALWAYS_INLINE void EvalSerialized(T_type &T, const typename AData_t::vcomplex_t *loc_dofs, AData_t &D)
static MFEM_ALWAYS_INLINE void AssembleSerialized(T_type &T, const AData_t &D, typename AData_t::vcomplex_t *loc_dofs)
static MFEM_ALWAYS_INLINE void EvalSerialized(T_type &T, const typename AData_t::vcomplex_t *loc_dofs, AData_t &D)
static MFEM_ALWAYS_INLINE void Assemble(const vec_layout_t &l, T_type &T, AData_t &D)
static MFEM_ALWAYS_INLINE void Eval(const vec_layout_t &l, T_type &T, AData_t &D)
static MFEM_ALWAYS_INLINE void Eval(const vec_layout_t &l, T_type &T, AData_t &D)
static MFEM_ALWAYS_INLINE void AssembleSerialized(T_type &T, const AData_t &D, typename AData_t::vcomplex_t *loc_dofs)
static MFEM_ALWAYS_INLINE void Assemble(const vec_layout_t &l, T_type &T, AData_t &D)
static MFEM_ALWAYS_INLINE void EvalSerialized(T_type &T, const typename AData_t::vcomplex_t *loc_dofs, AData_t &D)
This struct implements the input (Eval, EvalSerialized) and output (Assemble, AssembleSerialized) ope...
This struct is similar to struct AData, adding separate static data members for the input (InData) an...
TElementMatrix< InData, OutData, impl_traits_t > ElementMatrix
BData< InData, OutData, impl_traits_t > DataType
static MFEM_ALWAYS_INLINE void Compute(const qpt_layout_t &a, const qpt_data_t &A, const M_layout_t &m, M_data_t &M, ShapeEval_type &ev)
static MFEM_ALWAYS_INLINE void Compute(const qpt_layout_t &a, const qpt_data_t &A, const M_layout_t &m, M_data_t &M, ShapeEval_type &ev)
Assemble element mass matrix.
This struct implements element matrix computation for some combinations of input (InOps) and output (...
static StridedLayout2D< N2, S2, N1, S1 > transpose_12()
static StridedLayout2D< N1, S1, N2 *N3, S2 > merge_23()
static StridedLayout2D< N1 *N2, S1, N3, S3 > merge_12()
static OffsetStridedLayout2D< N2, S2, N3, S3 > ind1(int i1)
static OffsetStridedLayout2D< N2, S2, N3, S3 > ind14(int i1, int i4)
static StridedLayout3D< N1, S1, N2, S2, N3 *N4, S3 > merge_34()
static StridedLayout3D< N1 *N2, S1, N3, S3, N4, S4 > merge_12()
static const layout_type layout
data_t data[aligned_size >0?aligned_size:1]
data_t data[aligned_size >0?aligned_size:1]
static const layout_type layout
static const layout_type layout