27namespace quadrature_interpolator
30template<QVectorLayout Q_LAYOUT>
31static void Values1D(
const int NE,
39 const auto b =
Reshape(b_, q1d, d1d);
40 const auto x =
Reshape(x_, d1d, vdim, NE);
47 for (
int c = 0; c < vdim; c++)
49 for (
int q = 0; q < q1d; q++)
52 for (
int d = 0; d < d1d; d++)
54 u +=
b(q, d) * x(d, c, e);
65 int T_VDIM = 0,
int T_D1D = 0,
int T_Q1D = 0,
66 int T_NBZ = 1,
int MAX_D1D = 0,
int MAX_Q1D = 0>
67static void Values2D(
const int NE,
75 static constexpr int NBZ = T_NBZ ? T_NBZ : 1;
77 const int D1D = T_D1D ? T_D1D : d1d;
78 const int Q1D = T_Q1D ? T_Q1D : q1d;
79 const int VDIM = T_VDIM ? T_VDIM : vdim;
81 const auto b =
Reshape(b_, Q1D, D1D);
82 const auto x =
Reshape(x_, D1D, D1D, VDIM, NE);
84 Reshape(y_, Q1D, Q1D, VDIM, NE):
85 Reshape(y_, VDIM, Q1D, Q1D, NE);
89 const int D1D = T_D1D ? T_D1D : d1d;
90 const int Q1D = T_Q1D ? T_Q1D : q1d;
91 const int VDIM = T_VDIM ? T_VDIM : vdim;
92 constexpr int MQ1 = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
93 constexpr int MD1 = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
94 constexpr int MDQ = (MQ1 > MD1) ? MQ1 : MD1;
95 const int tidz = MFEM_THREAD_ID(z);
97 MFEM_SHARED
real_t sB[MQ1*MD1];
98 MFEM_SHARED
real_t sm0[NBZ][MDQ*MDQ];
99 MFEM_SHARED
real_t sm1[NBZ][MDQ*MDQ];
101 kernels::internal::LoadB<MD1,MQ1>(D1D,Q1D,
b,sB);
108 for (
int c = 0; c < VDIM; c++)
110 kernels::internal::LoadX(e,D1D,c,x,DD);
111 kernels::internal::EvalX(D1D,Q1D,B,DD,DQ);
112 kernels::internal::EvalY(D1D,Q1D,B,DQ,QQ);
113 MFEM_FOREACH_THREAD(qy,y,Q1D)
115 MFEM_FOREACH_THREAD(qx,x,Q1D)
129 int T_VDIM = 0,
int T_D1D = 0,
int T_Q1D = 0>
130static void Values3D(
const int NE,
138 const int D1D = T_D1D ? T_D1D : d1d;
139 const int Q1D = T_Q1D ? T_Q1D : q1d;
140 const int VDIM = T_VDIM ? T_VDIM : vdim;
142 const auto b =
Reshape(b_, Q1D, D1D);
143 const auto x =
Reshape(x_, D1D, D1D, D1D, VDIM, NE);
145 Reshape(y_, Q1D, Q1D, Q1D, VDIM, NE):
146 Reshape(y_, VDIM, Q1D, Q1D, Q1D, NE);
150 const int D1D = T_D1D ? T_D1D : d1d;
151 const int Q1D = T_Q1D ? T_Q1D : q1d;
152 const int VDIM = T_VDIM ? T_VDIM : vdim;
153 constexpr int MQ1 = T_Q1D ? T_Q1D : DofQuadLimits::MAX_INTERP_1D;
154 constexpr int MD1 = T_D1D ? T_D1D : DofQuadLimits::MAX_INTERP_1D;
155 constexpr int MDQ = (MQ1 > MD1) ? MQ1 : MD1;
157 MFEM_SHARED
real_t sB[MQ1*MD1];
158 MFEM_SHARED
real_t sm0[MDQ*MDQ*MDQ];
159 MFEM_SHARED
real_t sm1[MDQ*MDQ*MDQ];
161 kernels::internal::LoadB<MD1,MQ1>(D1D,Q1D,
b,sB);
169 for (
int c = 0; c < VDIM; c++)
171 kernels::internal::LoadX(e,D1D,c,x,DDD);
172 kernels::internal::EvalX(D1D,Q1D,B,DDD,DDQ);
173 kernels::internal::EvalY(D1D,Q1D,B,DDQ,DQQ);
174 kernels::internal::EvalZ(D1D,Q1D,B,DQQ,QQQ);
175 MFEM_FOREACH_THREAD(qz,z,Q1D)
177 MFEM_FOREACH_THREAD(qy,y,Q1D)
179 MFEM_FOREACH_THREAD(qx,x,Q1D)
181 const real_t u = QQQ(qz,qy,qx);
real_t u(const Vector &xvec)
DeviceTensor< 2, real_t > DeviceMatrix
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
void forall_2D_batch(int N, int X, int Y, int BZ, lambda &&body)
DeviceTensor< 2, const real_t > ConstDeviceMatrix
void forall_3D(int N, int X, int Y, int Z, lambda &&body)
QVectorLayout
Type describing possible layouts for Q-vectors.
@ byNODES
NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
@ byVDIM
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
DeviceTensor< 3, real_t > DeviceCube
void forall(int N, lambda &&body)