27template<
int T_D1D = 0,
int T_Q1D = 0>
28static void EAHdivAssemble2D(
const int NE,
29 const Array<real_t> &Bo_,
30 const Array<real_t> &Bc_,
32 const Vector &pa_data,
38 const int D1D = T_D1D ? T_D1D : d1d;
39 const int Q1D = T_Q1D ? T_Q1D : q1d;
42 const int NDOF = 2*(D1D-1)*D1D;
43 const auto Bo =
Reshape(Bo_.Read(), Q1D, D1D-1);
44 const auto Bc =
Reshape(Bc_.Read(), Q1D, D1D);
45 const auto D =
Reshape(pa_data.Read(), Q1D, Q1D, coeff_dim, NE);
46 const bool symmetric = (coeff_dim == 3);
47 auto M =
Reshape(
add ? ea_data.ReadWrite() : ea_data.
Write(), NDOF, NDOF, NE);
50 constexpr int MD1 = T_D1D ? T_D1D : DofQuadLimits::HDIV_MAX_D1D;
51 constexpr int MQ1 = T_Q1D ? T_Q1D : DofQuadLimits::HDIV_MAX_Q1D;
55 for (
int d = 0; d < D1D; d++)
57 for (
int q = 0; q < Q1D; q++)
59 if (d < D1D - 1) { r_Bo[q][d] = Bo(q,d); }
64 MFEM_SHARED
real_t s_D[4][MQ1][MQ1];
65 MFEM_FOREACH_THREAD(idx_q, x, Q1D*Q1D)
67 const int qx = idx_q % Q1D;
68 const int qy = idx_q / Q1D;
71 const real_t val = D(qx, qy, 0, e);
72 for (
int i = 0; i < 4; ++i) { s_D[i][qx][qy] = val; }
76 s_D[0][qx][qy] = D(qx, qy, 0, e);
77 s_D[1][qx][qy] = D(qx, qy, 1, e);
78 s_D[2][qx][qy] = (symmetric) ? s_D[1][qx][qy] : D(qx, qy, 2, e);
79 s_D[3][qx][qy] = (symmetric) ? D(qx, qy, 2, e) : D(qx, qy, 3, e);
84 MFEM_FOREACH_THREAD(idx_i, x, NDOF)
86 const int ic = idx_i / D1D / (D1D-1);
87 const int idx_ii = idx_i % (D1D * (D1D-1));
88 const int ix = (ic == 0) ? idx_ii%D1D : idx_ii%(D1D-1);
89 const int iy = (ic == 0) ? idx_ii/D1D : idx_ii/(D1D-1);
91 const real_t (&Bi1)[MQ1][MD1] = (ic == 0) ? r_Bc : r_Bo;
92 const real_t (&Bi2)[MQ1][MD1] = (ic == 0) ? r_Bo : r_Bc;
94 for (
int idx_j = 0; idx_j < NDOF; ++idx_j)
96 const int jc = idx_j / (D1D*(D1D-1));
97 const int idx_jj = idx_j % (D1D * (D1D-1));
98 const int jx = (jc == 0) ? idx_jj%D1D : idx_jj%(D1D-1);
99 const int jy = (jc == 0) ? idx_jj/D1D : idx_jj/(D1D-1);
101 const real_t (&Bj1)[MQ1][MD1] = (jc == 0) ? r_Bc : r_Bo;
102 const real_t (&Bj2)[MQ1][MD1] = (jc == 0) ? r_Bo : r_Bc;
105 for (
int qx = 0; qx < Q1D; ++qx)
107 for (
int qy = 0; qy < Q1D; ++qy)
109 const double coeff = s_D[ic + jc*2][qx][qy];
110 val += coeff*Bi1[qx][ix]*Bi2[qy][iy]*Bj1[qx][jx]*Bj2[qy][jy];
115 M(idx_i, idx_j, e) += val;
119 M(idx_i, idx_j, e) = val;
134template<
int T_D1D = 0,
int T_Q1D = 0>
135static void EAHdivAssemble3D(
const int NE,
136 const Array<real_t> &Bo_,
137 const Array<real_t> &Bc_,
139 const Vector &pa_data,
145 const int D1D = T_D1D ? T_D1D : d1d;
146 const int Q1D = T_Q1D ? T_Q1D : q1d;
149 const int NDOF_C = (D1D-1)*(D1D-1)*D1D;
150 const int NDOF = 3*NDOF_C;
151 const auto Bo =
Reshape(Bo_.Read(), Q1D, D1D-1);
152 const auto Bc =
Reshape(Bc_.Read(), Q1D, D1D);
153 const auto D =
Reshape(pa_data.Read(), Q1D, Q1D, Q1D, coeff_dim, NE);
154 const bool symmetric = (coeff_dim == 6);
155 auto M =
Reshape(
add ? ea_data.ReadWrite() : ea_data.
Write(), NDOF, NDOF, NE);
158 constexpr int MD1 = T_D1D ? T_D1D : DofQuadLimits::HDIV_MAX_D1D;
159 constexpr int MQ1 = T_Q1D ? T_Q1D : DofQuadLimits::HDIV_MAX_Q1D;
163 for (
int d = 0; d < D1D; d++)
165 for (
int q = 0; q < Q1D; q++)
167 if (d < D1D - 1) { r_Bo[q][d] = Bo(q,d); }
168 r_Bc[q][d] = Bc(q,d);
172 MFEM_SHARED
real_t s_D[9][MQ1][MQ1][MQ1];
173 MFEM_FOREACH_THREAD(idx_q, x, Q1D*Q1D*Q1D)
175 const int qx = idx_q % Q1D;
176 const int qy = (idx_q / Q1D) % Q1D;
177 const int qz = (idx_q / Q1D) / Q1D;
180 const real_t val = D(qx,qy,qz,0,e);
181 for (
int i = 0; i < 9; ++i) { s_D[i][qx][qy][qz] = val; }
185 s_D[0][qx][qy][qz] = D(qx,qy,qz,0,e);
186 s_D[1][qx][qy][qz] = D(qx,qy,qz,1,e);
187 s_D[2][qx][qy][qz] = D(qx,qy,qz,2,e);
188 s_D[3][qx][qy][qz] = symmetric ? s_D[1][qx][qy][qz] : D(qx,qy,qz,3,e);
189 s_D[4][qx][qy][qz] = symmetric ? D(qx,qy,qz,3,e) : D(qx,qy,qz,4,e);
190 s_D[5][qx][qy][qz] = symmetric ? D(qx,qy,qz,4,e) : D(qx,qy,qz,5,e);
191 s_D[6][qx][qy][qz] = symmetric ? s_D[2][qx][qy][qz] : D(qx,qy,qz,6,e);
192 s_D[7][qx][qy][qz] = symmetric ? s_D[5][qx][qy][qz] : D(qx,qy,qz,7,e);
193 s_D[8][qx][qy][qz] = symmetric ? D(qx,qy,qz,5,e) : D(qx,qy,qz,8,e);
198 MFEM_FOREACH_THREAD(idx_i, x, NDOF)
200 const int ic = idx_i / NDOF_C;
201 const int idx_ii = idx_i % NDOF_C;
203 const int nx_i = (ic == 0) ? D1D : D1D-1;
204 const int ny_i = (ic == 1) ? D1D : D1D-1;
206 const int ix = idx_ii % nx_i;
207 const int iy = (idx_ii / nx_i) % ny_i;
208 const int iz = (idx_ii / nx_i) / ny_i;
210 const real_t (&Bi1)[MQ1][MD1] = (ic == 0) ? r_Bc : r_Bo;
211 const real_t (&Bi2)[MQ1][MD1] = (ic == 1) ? r_Bc : r_Bo;
212 const real_t (&Bi3)[MQ1][MD1] = (ic == 2) ? r_Bc : r_Bo;
214 for (
int idx_j = 0; idx_j < NDOF; ++idx_j)
216 const int jc = idx_j / NDOF_C;
217 const int idx_jj = idx_j % NDOF_C;
219 const int nx_j = (jc == 0) ? D1D : D1D-1;
220 const int ny_j = (jc == 1) ? D1D : D1D-1;
222 const int jx = idx_jj % nx_j;
223 const int jy = (idx_jj / nx_j) % ny_j;
224 const int jz = (idx_jj / nx_j) / ny_j;
226 const real_t (&Bj1)[MQ1][MD1] = (jc == 0) ? r_Bc : r_Bo;
227 const real_t (&Bj2)[MQ1][MD1] = (jc == 1) ? r_Bc : r_Bo;
228 const real_t (&Bj3)[MQ1][MD1] = (jc == 2) ? r_Bc : r_Bo;
231 for (
int qx = 0; qx < Q1D; ++qx)
233 for (
int qy = 0; qy < Q1D; ++qy)
235 for (
int qz = 0; qz < Q1D; ++qz)
237 const double coeff = s_D[ic + jc*3][qx][qy][qz];
238 val += coeff*Bi1[qx][ix]*Bi2[qy][iy]*Bi3[qz][iz]*
239 Bj1[qx][jx]*Bj2[qy][jy]*Bj3[qz][jz];
245 M(idx_i, idx_j, e) += val;
249 M(idx_i, idx_j, e) = val;
265 MFEM_ABORT(
"Unsupported kernel.");
274 auto kernel = EAHdivAssemble2D<0,0>;
277 case 0x22: kernel = EAHdivAssemble2D<2,2>;
break;
278 case 0x33: kernel = EAHdivAssemble2D<3,3>;
break;
279 case 0x44: kernel = EAHdivAssemble2D<4,4>;
break;
280 case 0x55: kernel = EAHdivAssemble2D<5,5>;
break;
287 auto kernel = EAHdivAssemble3D<0,0>;
290 case 0x23: kernel = EAHdivAssemble3D<2,3>;
break;
291 case 0x34: kernel = EAHdivAssemble3D<3,4>;
break;
292 case 0x45: kernel = EAHdivAssemble3D<4,5>;
break;
293 case 0x56: kernel = EAHdivAssemble3D<5,6>;
break;
297 MFEM_ABORT(
"Unknown kernel.");
311 auto kernel = EAHdivAssemble2D<0,0>;
312 switch ((dofs1D << 4 ) | quad1D)
314 case 0x22: kernel = EAHdivAssemble2D<2,2>;
break;
315 case 0x33: kernel = EAHdivAssemble2D<3,3>;
break;
316 case 0x44: kernel = EAHdivAssemble2D<4,4>;
break;
317 case 0x55: kernel = EAHdivAssemble2D<5,5>;
break;
319 return kernel(ne,Bo,Gc,1,pa_data,ea_data,
add,dofs1D,quad1D);
323 auto kernel = EAHdivAssemble3D<0,0>;
324 switch ((dofs1D << 4 ) | quad1D)
326 case 0x23: kernel = EAHdivAssemble3D<2,3>;
break;
327 case 0x34: kernel = EAHdivAssemble3D<3,4>;
break;
328 case 0x45: kernel = EAHdivAssemble3D<4,5>;
break;
329 case 0x56: kernel = EAHdivAssemble3D<5,6>;
break;
331 return kernel(ne,Bo,Gc,1,pa_data,ea_data,
add,dofs1D,quad1D);
333 MFEM_ABORT(
"Unknown kernel.");
void AssembleEA(const FiniteElementSpace &fes, Vector &emat, const bool add) override
Method defining element assembly.
void AssemblePA(const FiniteElementSpace &fes) override
Method defining partial assembly.
Array< real_t > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Array< real_t > B
Basis functions evaluated at quadrature points.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
@ DIV
Implements CalcDivShape methods.
const DofToQuad * mapsO
Not owned. DOF-to-quad map, open.
void AssembleEA(const FiniteElementSpace &fes, Vector &emat, const bool add) override
Method defining element assembly.
void AssemblePA(const FiniteElementSpace &fes) override
Method defining partial assembly.
bool symmetric
False if using a nonsymmetric matrix coefficient.
const DofToQuad * mapsC
Not owned. DOF-to-quad map, closed.
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
void add(const Vector &v1, const Vector &v2, Vector &v)
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(int N, int X, int Y, lambda &&body)
static const DeviceDofQuadLimits & Get()
Return a const reference to the DeviceDofQuadLimits singleton.