MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_vecdiffusion_pa.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "../bilininteg.hpp"
15
16#include "./bilininteg_vecdiffusion_pa.hpp" // IWYU pragma: keep
17
18// #include "bilininteg_vecdiffusion_kernels.hpp"
19// #include "bilininteg_vecdiffusion_pa.hpp"
20
21namespace mfem
22{
23
29
35
38{
39 vdim = vector_dimension;
40}
41
48
50 int vector_dimension)
52{
53 Q = &q;
54 vdim = vector_dimension;
55}
56
63
70
72{
73 Mesh *mesh = fes.GetMesh();
74 const FiniteElement &el = *fes.GetTypicalFE();
75 const auto *ir = IntRule ? IntRule : &DiffusionIntegrator::GetRule(el, el);
76
77 if (DeviceCanUseCeed())
78 {
79 delete ceedOp;
80 const bool mixed =
81 mesh->GetNumGeometries(mesh->Dimension()) > 1 || fes.IsVariableOrder();
82 if (mixed) { ceedOp = new ceed::MixedPADiffusionIntegrator(*this, fes, Q); }
83 else { ceedOp = new ceed::PADiffusionIntegrator(fes, *ir, Q); }
84 return;
85 }
86
87 // If vdim is not set, set it to the space dimension
88 vdim = (vdim == -1) ? fes.GetVDim() : vdim;
89 MFEM_VERIFY(vdim == fes.GetVDim(), "vdim != fes.GetVDim()");
90
93 : pa_mt;
94
95 ne = fes.GetNE();
96 dim = mesh->Dimension();
97 sdim = mesh->SpaceDimension();
98 const int nq = ir->GetNPoints();
101 dofs1D = maps->ndof;
102 quad1D = maps->nqpt;
103 const int q1d = quad1D;
104
105 if (!(dim == 2 || dim == 3)) { MFEM_ABORT("Dimension not supported."); }
106
107 QuadratureSpace qs(*mesh, *ir);
109
110 if (Q)
111 {
112 coeff.Project(*Q);
113 }
114 else if (VQ)
115 {
116 coeff.Project(*VQ);
117 MFEM_VERIFY(VQ->GetVDim() == vdim, "VQ vdim vs. vdim error");
118 }
119 else if (MQ)
120 {
121 coeff.ProjectTranspose(*MQ);
122 MFEM_VERIFY(MQ->GetVDim() == vdim, "MQ dimension vs. vdim error");
123 MFEM_VERIFY(coeff.Size() == (vdim*vdim) * ne * nq, "MQ size error");
124 }
125 else { coeff.SetConstant(1.0); }
126
127 coeff_vdim = coeff.GetVDim();
128 const bool scalar_coeff = coeff_vdim == 1;
129 const bool vector_coeff = coeff_vdim == vdim;
130 const bool matrix_coeff = coeff_vdim == vdim * vdim;
131 MFEM_VERIFY(scalar_coeff + vector_coeff + matrix_coeff == 1, "");
132
133 const int pa_size = dim * dim;
134 pa_data.SetSize(nq * pa_size * vdim * (matrix_coeff ? dim : 1) * ne, mt);
135
136 if (dim == 2 && sdim == 3)
137 {
138 MFEM_VERIFY(scalar_coeff, "");
139 const int nc = vdim;
140 const auto W = Reshape(ir->GetWeights().Read(), q1d, q1d);
141 const auto J = Reshape(geom->J.Read(), q1d, q1d, sdim, dim, ne);
142 const auto C = Reshape(coeff.Read(), coeff_vdim, q1d, q1d, ne);
143 auto D = Reshape(pa_data.Write(), q1d, q1d, pa_size,
144 vdim * (matrix_coeff ? dim : 1), ne);
145
146 mfem::forall_2D(ne, q1d, q1d, [=] MFEM_HOST_DEVICE(int e)
147 {
148 MFEM_FOREACH_THREAD(qy, y, q1d)
149 {
150 MFEM_FOREACH_THREAD(qx, x, q1d)
151 {
152 for (int i = 0; i < nc; ++i)
153 {
154 const real_t wq = W(qx, qy);
155 const real_t J11 = J(qx, qy, 0, 0, e);
156 const real_t J21 = J(qx, qy, 1, 0, e);
157 const real_t J31 = J(qx, qy, 2, 0, e);
158 const real_t J12 = J(qx, qy, 0, 1, e);
159 const real_t J22 = J(qx, qy, 1, 1, e);
160 const real_t J32 = J(qx, qy, 2, 1, e);
161 const real_t E = J11*J11 + J21*J21 + J31*J31;
162 const real_t G = J12*J12 + J22*J22 + J32*J32;
163 const real_t F = J11*J12 + J21*J22 + J31*J32;
164 const real_t iw = 1.0 / sqrt(E*G - F*F);
165 const auto C0 = C(0, qx, qy, e);
166 const real_t alpha = wq * C0 * iw;
167 D(qx, qy, 0, i, e) = alpha * G; // 1,1
168 D(qx, qy, 1, i, e) = -alpha * F; // 1,2
169 D(qx, qy, 2, i, e) = -alpha * F; // 2,1 == 1,2
170 D(qx, qy, 3, i, e) = alpha * E; // 2,2
171 }
172 }
173 }
174 });
175 }
176 else if (dim == 2 && sdim == 2)
177 {
178 const int nc = vdim, cvdim = coeff_vdim;
179 const auto W = Reshape(ir->GetWeights().Read(), q1d, q1d);
180 const auto J = Reshape(geom->J.Read(), q1d, q1d, sdim, dim, ne);
181 const auto C = Reshape(coeff.Read(), coeff_vdim, q1d, q1d, ne);
182 auto DE = Reshape(pa_data.Write(), q1d, q1d, pa_size,
183 vdim * (matrix_coeff ? dim : 1), ne);
184 mfem::forall_2D(ne, q1d, q1d, [=] MFEM_HOST_DEVICE(int e)
185 {
186 MFEM_FOREACH_THREAD(qy, y, q1d)
187 {
188 MFEM_FOREACH_THREAD(qx, x, q1d)
189 {
190 const real_t J11 = J(qx, qy, 0, 0, e);
191 const real_t J21 = J(qx, qy, 1, 0, e);
192 const real_t J12 = J(qx, qy, 0, 1, e);
193 const real_t J22 = J(qx, qy, 1, 1, e);
194 const real_t w_detJ = W(qx, qy) / ((J11*J22)-(J21*J12));
195 const real_t D0 = w_detJ * (J12*J12 + J22*J22);
196 const real_t D1 = -w_detJ * (J12*J11 + J22*J21);
197 const real_t D2 = w_detJ * (J11*J11 + J21*J21);
198 const int map[4] = {0, 2, 1, 3};
199
200 for (int i = 0; i < (matrix_coeff ? cvdim : nc); ++i)
201 {
202 const auto k = matrix_coeff ? map[i] : (vector_coeff ? i : 0);
203 const auto Cc = C(k, qx, qy, e);
204 DE(qx, qy, 0, i, e) = D0 * Cc;
205 DE(qx, qy, 1, i, e) = D1 * Cc;
206 DE(qx, qy, 2, i, e) = D1 * Cc;
207 DE(qx, qy, 3, i, e) = D2 * Cc;
208 }
209 }
210 }
211 });
212 }
213 else if (dim == 3 && sdim == 3)
214 {
215 const int nc = vdim, cvdim = coeff_vdim;
216 const auto W = Reshape(ir->GetWeights().Read(), q1d, q1d, q1d);
217 const auto J = Reshape(geom->J.Read(), q1d, q1d, q1d, sdim, dim, ne);
218 const auto C = Reshape(coeff.Read(), coeff_vdim, q1d, q1d, q1d, ne);
219 auto DE = Reshape(pa_data.Write(), q1d, q1d, q1d, pa_size,
220 vdim * (matrix_coeff ? dim : 1), ne);
221
222 mfem::forall_3D(ne, q1d, q1d, q1d, [=] MFEM_HOST_DEVICE(int e)
223 {
224 MFEM_FOREACH_THREAD(qz, z, q1d)
225 {
226 MFEM_FOREACH_THREAD(qy, y, q1d)
227 {
228 MFEM_FOREACH_THREAD(qx, x, q1d)
229 {
230 const real_t J11 = J(qx, qy, qz, 0, 0, e);
231 const real_t J21 = J(qx, qy, qz, 1, 0, e);
232 const real_t J31 = J(qx, qy, qz, 2, 0, e);
233 const real_t J12 = J(qx, qy, qz, 0, 1, e);
234 const real_t J22 = J(qx, qy, qz, 1, 1, e);
235 const real_t J32 = J(qx, qy, qz, 2, 1, e);
236 const real_t J13 = J(qx, qy, qz, 0, 2, e);
237 const real_t J23 = J(qx, qy, qz, 1, 2, e);
238 const real_t J33 = J(qx, qy, qz, 2, 2, e);
239 const real_t detJ = J11 * (J22 * J33 - J32 * J23) -
240 J21 * (J12 * J33 - J32 * J13) +
241 J31 * (J12 * J23 - J22 * J13);
242 const real_t c_detJ = W(qx, qy, qz) / detJ;
243 // adj(J)
244 const real_t A11 = (J22 * J33) - (J23 * J32);
245 const real_t A12 = (J32 * J13) - (J12 * J33);
246 const real_t A13 = (J12 * J23) - (J22 * J13);
247 const real_t A21 = (J31 * J23) - (J21 * J33);
248 const real_t A22 = (J11 * J33) - (J13 * J31);
249 const real_t A23 = (J21 * J13) - (J11 * J23);
250 const real_t A31 = (J21 * J32) - (J31 * J22);
251 const real_t A32 = (J31 * J12) - (J11 * J32);
252 const real_t A33 = (J11 * J22) - (J12 * J21);
253 // detJ J^{-1} J^{-T} = (1/detJ) adj(J) adj(J)^T
254 const real_t D11 = c_detJ * (A11*A11 + A12*A12 + A13*A13); // 1,1
255 const real_t D21 = c_detJ * (A11*A21 + A12*A22 + A13*A23); // 2,1
256 const real_t D31 = c_detJ * (A11*A31 + A12*A32 + A13*A33); // 3,1
257 const real_t D22 = c_detJ * (A21*A21 + A22*A22 + A23*A23); // 2,2
258 const real_t D32 = c_detJ * (A21*A31 + A22*A32 + A23*A33); // 3,2
259 const real_t D33 = c_detJ * (A31*A31 + A32*A32 + A33*A33); // 3,3
260 const int map[9] = {0, 3, 6, 1, 4, 7, 2, 5, 8};
261
262 for (int i = 0; i < (matrix_coeff ? cvdim : nc); ++i)
263 {
264 const auto k = matrix_coeff ? map[i] : vector_coeff ? i : 0;
265 const auto Ck = C(k, qx, qy, qz, e);
266 DE(qx, qy, qz, 0, i, e) = D11 * Ck;
267 DE(qx, qy, qz, 1, i, e) = D21 * Ck;
268 DE(qx, qy, qz, 2, i, e) = D31 * Ck;
269 DE(qx, qy, qz, 3, i, e) = D22 * Ck;
270 DE(qx, qy, qz, 4, i, e) = D32 * Ck;
271 DE(qx, qy, qz, 5, i, e) = D33 * Ck;
272 }
273 }
274 }
275 }
276 });
277 }
278 else
279 {
280 MFEM_ABORT("Unknown VectorDiffusionIntegrator::AssemblePA kernel for"
281 << " dim:" << dim << ", vdim:" << vdim << ", sdim:" << sdim);
282 }
283}
284
285// PA Diffusion Apply kernel
287{
288 // Use CEED backend if available
289 if (DeviceCanUseCeed()) { return ceedOp->AddMult(x, y); }
290
291 // Add the VectorDiffusionAddMultPA specializations
292 static const auto vector_diffusion_kernel_specializations =
293 (
294 // 2D, SDIM = 2
295 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 2,2>::Add(),
296 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 3,3>::Add(),
297 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 4,4>::Add(),
298 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 5,5>::Add(),
299 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 6,6>::Add(),
300 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 7,7>::Add(),
301 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 8,8>::Add(),
302 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,2, 9,9>::Add(),
303 // 2D, SDIM = 3
304 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,3, 2,2>::Add(),
305 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,3, 3,3>::Add(),
306 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,3, 4,4>::Add(),
307 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<2,3, 5,5>::Add(),
308 // 3D, SDIM = 3
309 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 2,2>::Add(),
310 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 2,3>::Add(),
311 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 3,4>::Add(),
312 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 4,5>::Add(),
313 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 4,6>::Add(),
314 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 5,6>::Add(),
315 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 5,8>::Add(),
316 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 6,7>::Add(),
317 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 7,8>::Add(),
318 VectorDiffusionIntegrator::ApplyPAKernels::Specialization<3,3, 8,9>::Add(),
319 true);
320 MFEM_CONTRACT_VAR(vector_diffusion_kernel_specializations);
321
322 ApplyPAKernels::Run(dim, sdim, dofs1D, quad1D,
323 ne, coeff_vdim, maps->B, maps->G, pa_data, x, y,
324 sdim, dofs1D, quad1D);
325
326}
327
328template<int T_D1D = 0, int T_Q1D = 0>
329static void PAVectorDiffusionDiagonal2D(const int NE,
330 const Array<real_t> &b,
331 const Array<real_t> &g,
332 const Vector &d,
333 Vector &y,
334 const int d1d = 0,
335 const int q1d = 0)
336{
337 const int D1D = T_D1D ? T_D1D : d1d;
338 const int Q1D = T_Q1D ? T_Q1D : q1d;
339 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().MAX_D1D, "");
340 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().MAX_Q1D, "");
341
342 const auto B = Reshape(b.Read(), Q1D, D1D);
343 const auto G = Reshape(g.Read(), Q1D, D1D);
344 // note the different shape for D, this is a (symmetric) matrix so we only
345 // store necessary entries
346 MFEM_VERIFY(d.Size() == Q1D*Q1D*4*2*NE, "");
347 const auto D = Reshape(d.Read(), Q1D*Q1D, /*3*/4, 2, NE);
348 auto Y = Reshape(y.ReadWrite(), D1D, D1D, 2, NE);
349
350 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
351 {
352 const int D1D = T_D1D ? T_D1D : d1d;
353 const int Q1D = T_Q1D ? T_Q1D : q1d;
354 constexpr int MD1 = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
355 constexpr int MQ1 = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
356 // gradphi \cdot Q \gradphi has four terms
357 real_t QD0[MQ1][MD1];
358 real_t QD1[MQ1][MD1];
359 real_t QD2[MQ1][MD1];
360 for (int qx = 0; qx < Q1D; ++qx)
361 {
362 for (int dy = 0; dy < D1D; ++dy)
363 {
364 QD0[qx][dy] = 0.0;
365 QD1[qx][dy] = 0.0;
366 QD2[qx][dy] = 0.0;
367 for (int qy = 0; qy < Q1D; ++qy)
368 {
369 const int q = qx + qy * Q1D;
370 const real_t D0 = D(q,0,0,e);
371 const real_t D1 = D(q,1,0,e);
372 const real_t D2 = D(q,3/*2*/,0,e); // size from 3 (symmetric) to 4 (dims x dims)
373 QD0[qx][dy] += B(qy, dy) * B(qy, dy) * D0;
374 QD1[qx][dy] += B(qy, dy) * G(qy, dy) * D1;
375 QD2[qx][dy] += G(qy, dy) * G(qy, dy) * D2;
376 }
377 }
378 }
379 for (int dy = 0; dy < D1D; ++dy)
380 {
381 for (int dx = 0; dx < D1D; ++dx)
382 {
383 real_t temp = 0.0;
384 for (int qx = 0; qx < Q1D; ++qx)
385 {
386 temp += G(qx, dx) * G(qx, dx) * QD0[qx][dy];
387 temp += G(qx, dx) * B(qx, dx) * QD1[qx][dy];
388 temp += B(qx, dx) * G(qx, dx) * QD1[qx][dy];
389 temp += B(qx, dx) * B(qx, dx) * QD2[qx][dy];
390 }
391 Y(dx,dy,0,e) += temp;
392 Y(dx,dy,1,e) += temp;
393 }
394 }
395 });
396}
397
398template<int T_D1D = 0, int T_Q1D = 0>
399static void PAVectorDiffusionDiagonal3D(const int NE,
400 const Array<real_t> &b,
401 const Array<real_t> &g,
402 const Vector &d,
403 Vector &y,
404 const int d1d = 0,
405 const int q1d = 0)
406{
407 constexpr int DIM = 3;
408 const int D1D = T_D1D ? T_D1D : d1d;
409 const int Q1D = T_Q1D ? T_Q1D : q1d;
410 const int max_q1d = T_Q1D ? T_Q1D : DeviceDofQuadLimits::Get().MAX_Q1D;
411 const int max_d1d = T_D1D ? T_D1D : DeviceDofQuadLimits::Get().MAX_D1D;
412 MFEM_VERIFY(D1D <= max_d1d, "");
413 MFEM_VERIFY(Q1D <= max_q1d, "");
414 auto B = Reshape(b.Read(), Q1D, D1D);
415 auto G = Reshape(g.Read(), Q1D, D1D);
416 MFEM_VERIFY(d.Size() == Q1D*Q1D*Q1D*9*3*NE, "");
417 auto Q = Reshape(d.Read(), Q1D*Q1D*Q1D, 9/*PA_SIZE:dims*dims*/, 3/*VDIM*/, NE);
418 auto Y = Reshape(y.ReadWrite(), D1D, D1D, D1D, 3, NE);
419 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
420 {
421 const int D1D = T_D1D ? T_D1D : d1d;
422 const int Q1D = T_Q1D ? T_Q1D : q1d;
423 constexpr int MD1 = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
424 constexpr int MQ1 = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
425 real_t QQD[MQ1][MQ1][MD1];
426 real_t QDD[MQ1][MD1][MD1];
427 for (int i = 0; i < DIM; ++i)
428 {
429 for (int j = 0; j < DIM; ++j)
430 {
431 // first tensor contraction, along z direction
432 for (int qx = 0; qx < Q1D; ++qx)
433 {
434 for (int qy = 0; qy < Q1D; ++qy)
435 {
436 for (int dz = 0; dz < D1D; ++dz)
437 {
438 QQD[qx][qy][dz] = 0.0;
439 for (int qz = 0; qz < Q1D; ++qz)
440 {
441 const int q = qx + (qy + qz * Q1D) * Q1D;
442 const int k = j >= i ?
443 3 - (3-i)*(2-i)/2 + j:
444 3 - (3-j)*(2-j)/2 + i;
445 // using 6 symmetric values
446 const real_t O = Q(q,k,0,e);
447 const real_t Bz = B(qz,dz);
448 const real_t Gz = G(qz,dz);
449 const real_t L = i==2 ? Gz : Bz;
450 const real_t R = j==2 ? Gz : Bz;
451 QQD[qx][qy][dz] += L * O * R;
452 }
453 }
454 }
455 }
456 // second tensor contraction, along y direction
457 for (int qx = 0; qx < Q1D; ++qx)
458 {
459 for (int dz = 0; dz < D1D; ++dz)
460 {
461 for (int dy = 0; dy < D1D; ++dy)
462 {
463 QDD[qx][dy][dz] = 0.0;
464 for (int qy = 0; qy < Q1D; ++qy)
465 {
466 const real_t By = B(qy,dy);
467 const real_t Gy = G(qy,dy);
468 const real_t L = i==1 ? Gy : By;
469 const real_t R = j==1 ? Gy : By;
470 QDD[qx][dy][dz] += L * QQD[qx][qy][dz] * R;
471 }
472 }
473 }
474 }
475 // third tensor contraction, along x direction
476 for (int dz = 0; dz < D1D; ++dz)
477 {
478 for (int dy = 0; dy < D1D; ++dy)
479 {
480 for (int dx = 0; dx < D1D; ++dx)
481 {
482 real_t temp = 0.0;
483 for (int qx = 0; qx < Q1D; ++qx)
484 {
485 const real_t Bx = B(qx,dx);
486 const real_t Gx = G(qx,dx);
487 const real_t L = i==0 ? Gx : Bx;
488 const real_t R = j==0 ? Gx : Bx;
489 temp += L * QDD[qx][dy][dz] * R;
490 }
491 Y(dx, dy, dz, 0, e) += temp;
492 Y(dx, dy, dz, 1, e) += temp;
493 Y(dx, dy, dz, 2, e) += temp;
494 }
495 }
496 }
497 }
498 }
499 });
500}
501
502static void PAVectorDiffusionAssembleDiagonal(const int dim,
503 const int D1D,
504 const int Q1D,
505 const int NE,
506 const Array<real_t> &B,
507 const Array<real_t> &G,
508 const Vector &op,
509 Vector &y)
510{
511 if (dim == 2)
512 {
513 return PAVectorDiffusionDiagonal2D(NE, B, G, op, y, D1D, Q1D);
514 }
515 else if (dim == 3)
516 {
517 return PAVectorDiffusionDiagonal3D(NE, B, G, op, y, D1D, Q1D);
518 }
519 MFEM_ABORT("Dimension not implemented.");
520}
521
523{
524 if (DeviceCanUseCeed())
525 {
526 ceedOp->GetDiagonal(diag);
527 }
528 else
529 {
530 MFEM_VERIFY(!VQ && !MQ, "VQ and MQ not supported.");
531 PAVectorDiffusionAssembleDiagonal(dim, dofs1D, quad1D, ne,
532 maps->B, maps->G,
533 pa_data, diag);
534 }
535}
536
537/*
538// PA Diffusion Apply kernel
539void VectorDiffusionIntegrator::AddMultPA(const Vector &x, Vector &y) const
540{
541 if (DeviceCanUseCeed())
542 {
543 ceedOp->AddMult(x, y);
544 }
545 else
546 {
547 const int D1D = dofs1D;
548 const int Q1D = quad1D;
549 const Array<real_t> &B = maps->B;
550 const Array<real_t> &G = maps->G;
551 const Array<real_t> &Bt = maps->Bt;
552 const Array<real_t> &Gt = maps->Gt;
553 const Vector &D = pa_data;
554 ApplyPAKernels::Run(dim, sdim, D1D, Q1D, ne, B, G, Bt, Gt, D, x, y, D1D,
555 Q1D, sdim);
556 }
557}
558
559/// \cond DO_NOT_DOCUMENT
560
561VectorDiffusionIntegrator::ApplyKernelType
562VectorDiffusionIntegrator::ApplyPAKernels::Fallback(int DIM, int, int, int)
563{
564 if (DIM == 2) { return internal::PAVectorDiffusionApply2D; }
565 else if (DIM == 3) { return internal::PAVectorDiffusionApply3D; }
566 else { MFEM_ABORT(""); }
567}
568
569VectorDiffusionIntegrator::Kernels::Kernels()
570{
571 VectorDiffusionIntegrator::AddSpecialization<2, 3, 2, 2>();
572 VectorDiffusionIntegrator::AddSpecialization<2, 3, 3, 3>();
573 VectorDiffusionIntegrator::AddSpecialization<2, 3, 4, 4>();
574 VectorDiffusionIntegrator::AddSpecialization<2, 3, 5, 5>();
575}
576
577/// \endcond DO_NOT_DOCUMENT
578*/
579
580} // namespace mfem
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:381
Abstract base class BilinearFormIntegrator.
Class to represent a coefficient evaluated at quadrature points.
void SetConstant(real_t constant)
Set this vector to the given constant.
void Project(Coefficient &coeff)
Evaluate the given Coefficient at the quadrature points defined by qs.
int GetVDim() const
Return the number of values per quadrature point.
void ProjectTranspose(MatrixCoefficient &coeff)
Project the transpose of coeff.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:277
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe)
Array< real_t > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition fe_base.hpp:214
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:165
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition fe_base.hpp:193
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:178
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:182
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:678
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:869
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:644
int GetVDim() const
Returns the vector dimension of the finite element space.
Definition fespace.hpp:819
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3860
Abstract class for all finite elements.
Definition fe_base.hpp:247
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.cpp:375
Vector J
Jacobians of the element transformations at all quadrature points.
Definition mesh.hpp:3115
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
const IntegrationRule * IntRule
Base class for Matrix Coefficients that optionally depend on time and space.
int GetVDim() const
For backward compatibility get the width of the matrix.
Mesh data type.
Definition mesh.hpp:65
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1309
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
Definition mesh.cpp:883
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:7558
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:164
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
void AssemblePA(const FiniteElementSpace &fes) override
Method defining partial assembly.
const DofToQuad * maps
Not owned.
VectorDiffusionIntegrator(const IntegrationRule *ir=nullptr)
void AddMultPA(const Vector &x, Vector &y) const override
Method for partially assembled action.
const GeometricFactors * geom
Not owned.
void AssembleDiagonalPA(Vector &diag) override
Assemble diagonal and add it to Vector diag.
Vector data type.
Definition vector.hpp:82
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:520
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:536
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:528
void GetDiagonal(mfem::Vector &diag) const
Definition operator.cpp:104
void AddMult(const mfem::Vector &x, mfem::Vector &y, const real_t a=1.0) const override
Operator application: y+=A(x) (default) or y+=a*A(x).
Definition operator.cpp:72
Represent a DiffusionIntegrator with AssemblyLevel::Partial using libCEED.
Definition diffusion.hpp:27
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
constexpr int DIM
mfem::real_t real_t
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:138
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
Definition util.cpp:33
@ FULL
Store the coefficient as a full QuadratureFunction.
void forall_2D(int N, int X, int Y, lambda &&body)
Definition forall.hpp:925
void forall_3D(int N, int X, int Y, int Z, lambda &&body)
Definition forall.hpp:937
float real_t
Definition config.hpp:46
MemoryType
Memory types supported by MFEM.
void forall(int N, lambda &&body)
Definition forall.hpp:839
static const DeviceDofQuadLimits & Get()
Return a const reference to the DeviceDofQuadLimits singleton.
Definition forall.hpp:128
int MAX_D1D
Maximum number of 1D nodal points.
Definition forall.hpp:118
int MAX_Q1D
Maximum number of 1D quadrature points.
Definition forall.hpp:119