MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_vecmass_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
13#include "../bilininteg.hpp"
14#include "../gridfunc.hpp"
16
17namespace mfem
18{
19
21{
22 // Assuming the same element type
23 Mesh *mesh = fes.GetMesh();
24 const FiniteElement &el = *fes.GetTypicalFE();
26 const IntegrationRule *ir
27 = IntRule ? IntRule : &MassIntegrator::GetRule(el, el, *T);
28 if (DeviceCanUseCeed())
29 {
30 delete ceedOp;
31 const bool mixed = mesh->GetNumGeometries(mesh->Dimension()) > 1 ||
32 fes.IsVariableOrder();
33 if (mixed)
34 {
35 ceedOp = new ceed::MixedPAMassIntegrator(*this, fes, Q);
36 }
37 else
38 {
39 ceedOp = new ceed::PAMassIntegrator(fes, *ir, Q);
40 }
41 return;
42 }
43 dim = mesh->Dimension();
44 ne = fes.GetMesh()->GetNE();
45 nq = ir->GetNPoints();
49 dofs1D = maps->ndof;
50 quad1D = maps->nqpt;
52 real_t coeff = 1.0;
53 if (Q)
54 {
55 ConstantCoefficient *cQ = dynamic_cast<ConstantCoefficient*>(Q);
56 MFEM_VERIFY(cQ != NULL, "Only ConstantCoefficient is supported.");
57 coeff = cQ->constant;
58 }
59 if (!(dim == 2 || dim == 3))
60 {
61 MFEM_ABORT("Dimension not supported.");
62 }
63 if (dim == 2)
64 {
65 const real_t constant = coeff;
66 const int NE = ne;
67 const int NQ = nq;
68 auto w = ir->GetWeights().Read();
69 auto J = Reshape(geom->J.Read(), NQ,2,2,NE);
70 auto v = Reshape(pa_data.Write(), NQ, NE);
71 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
72 {
73 for (int q = 0; q < NQ; ++q)
74 {
75 const real_t J11 = J(q,0,0,e);
76 const real_t J12 = J(q,1,0,e);
77 const real_t J21 = J(q,0,1,e);
78 const real_t J22 = J(q,1,1,e);
79 const real_t detJ = (J11*J22)-(J21*J12);
80 v(q,e) = w[q] * constant * detJ;
81 }
82 });
83 }
84 if (dim == 3)
85 {
86 const real_t constant = coeff;
87 const int NE = ne;
88 const int NQ = nq;
89 auto W = ir->GetWeights().Read();
90 auto J = Reshape(geom->J.Read(), NQ,3,3,NE);
91 auto v = Reshape(pa_data.Write(), NQ,NE);
92 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
93 {
94 for (int q = 0; q < NQ; ++q)
95 {
96 const real_t J11 = J(q,0,0,e), J12 = J(q,0,1,e), J13 = J(q,0,2,e);
97 const real_t J21 = J(q,1,0,e), J22 = J(q,1,1,e), J23 = J(q,1,2,e);
98 const real_t J31 = J(q,2,0,e), J32 = J(q,2,1,e), J33 = J(q,2,2,e);
99 const real_t detJ = J11 * (J22 * J33 - J32 * J23) -
100 J21 * (J12 * J33 - J32 * J13) +
101 J31 * (J12 * J23 - J22 * J13);
102 v(q,e) = W[q] * constant * detJ;
103 }
104 });
105 }
106}
107
108template<const int T_D1D = 0, const int T_Q1D = 0>
109static void PAVectorMassAssembleDiagonal2D(const int NE,
110 const Array<real_t> &B_,
111 const Array<real_t> &Bt_,
112 const Vector &op_,
113 Vector &diag_,
114 const int d1d = 0,
115 const int q1d = 0)
116{
117 const int D1D = T_D1D ? T_D1D : d1d;
118 const int Q1D = T_Q1D ? T_Q1D : q1d;
119 constexpr int VDIM = 2;
120 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().MAX_D1D, "");
121 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().MAX_Q1D, "");
122 auto B = Reshape(B_.Read(), Q1D, D1D);
123 auto op = Reshape(op_.Read(), Q1D, Q1D, NE);
124 auto y = Reshape(diag_.ReadWrite(), D1D, D1D, VDIM, NE);
125 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
126 {
127 const int D1D = T_D1D ? T_D1D : d1d;
128 const int Q1D = T_Q1D ? T_Q1D : q1d;
129 constexpr int max_D1D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
130 constexpr int max_Q1D = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
131
132 real_t temp[max_Q1D][max_D1D];
133 for (int qx = 0; qx < Q1D; ++qx)
134 {
135 for (int dy = 0; dy < D1D; ++dy)
136 {
137 temp[qx][dy] = 0.0;
138 for (int qy = 0; qy < Q1D; ++qy)
139 {
140 temp[qx][dy] += B(qy, dy) * B(qy, dy) * op(qx, qy, e);
141 }
142 }
143 }
144 for (int dy = 0; dy < D1D; ++dy)
145 {
146 for (int dx = 0; dx < D1D; ++dx)
147 {
148 real_t temp1 = 0.0;
149 for (int qx = 0; qx < Q1D; ++qx)
150 {
151 temp1 += B(qx, dx) * B(qx, dx) * temp[qx][dy];
152 }
153 y(dx, dy, 0, e) = temp1;
154 y(dx, dy, 1, e) = temp1;
155 }
156 }
157 });
158}
159
160template<const int T_D1D = 0, const int T_Q1D = 0>
161static void PAVectorMassAssembleDiagonal3D(const int NE,
162 const Array<real_t> &B_,
163 const Array<real_t> &Bt_,
164 const Vector &op_,
165 Vector &diag_,
166 const int d1d = 0,
167 const int q1d = 0)
168{
169 const int D1D = T_D1D ? T_D1D : d1d;
170 const int Q1D = T_Q1D ? T_Q1D : q1d;
171 constexpr int VDIM = 3;
172 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().MAX_D1D, "");
173 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().MAX_Q1D, "");
174 auto B = Reshape(B_.Read(), Q1D, D1D);
175 auto op = Reshape(op_.Read(), Q1D, Q1D, Q1D, NE);
176 auto y = Reshape(diag_.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
177 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
178 {
179 const int D1D = T_D1D ? T_D1D : d1d; // nvcc workaround
180 const int Q1D = T_Q1D ? T_Q1D : q1d;
181 // the following variables are evaluated at compile time
182 constexpr int max_D1D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
183 constexpr int max_Q1D = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
184
185 real_t temp[max_Q1D][max_Q1D][max_D1D];
186 for (int qx = 0; qx < Q1D; ++qx)
187 {
188 for (int qy = 0; qy < Q1D; ++qy)
189 {
190 for (int dz = 0; dz < D1D; ++dz)
191 {
192 temp[qx][qy][dz] = 0.0;
193 for (int qz = 0; qz < Q1D; ++qz)
194 {
195 temp[qx][qy][dz] += B(qz, dz) * B(qz, dz) * op(qx, qy, qz, e);
196 }
197 }
198 }
199 }
200 real_t temp2[max_Q1D][max_D1D][max_D1D];
201 for (int qx = 0; qx < Q1D; ++qx)
202 {
203 for (int dz = 0; dz < D1D; ++dz)
204 {
205 for (int dy = 0; dy < D1D; ++dy)
206 {
207 temp2[qx][dy][dz] = 0.0;
208 for (int qy = 0; qy < Q1D; ++qy)
209 {
210 temp2[qx][dy][dz] += B(qy, dy) * B(qy, dy) * temp[qx][qy][dz];
211 }
212 }
213 }
214 }
215 for (int dz = 0; dz < D1D; ++dz)
216 {
217 for (int dy = 0; dy < D1D; ++dy)
218 {
219 for (int dx = 0; dx < D1D; ++dx)
220 {
221 real_t temp3 = 0.0;
222 for (int qx = 0; qx < Q1D; ++qx)
223 {
224 temp3 += B(qx, dx) * B(qx, dx)
225 * temp2[qx][dy][dz];
226 }
227 y(dx, dy, dz, 0, e) = temp3;
228 y(dx, dy, dz, 1, e) = temp3;
229 y(dx, dy, dz, 2, e) = temp3;
230 }
231 }
232 }
233 });
234}
235
236static void PAVectorMassAssembleDiagonal(const int dim,
237 const int D1D,
238 const int Q1D,
239 const int NE,
240 const Array<real_t> &B,
241 const Array<real_t> &Bt,
242 const Vector &op,
243 Vector &y)
244{
245 if (dim == 2)
246 {
247 return PAVectorMassAssembleDiagonal2D(NE, B, Bt, op, y, D1D, Q1D);
248 }
249 else if (dim == 3)
250 {
251 return PAVectorMassAssembleDiagonal3D(NE, B, Bt, op, y, D1D, Q1D);
252 }
253 MFEM_ABORT("Dimension not implemented.");
254}
255
257{
258 if (DeviceCanUseCeed())
259 {
260 ceedOp->GetDiagonal(diag);
261 }
262 else
263 {
264 PAVectorMassAssembleDiagonal(dim, dofs1D, quad1D, ne,
265 maps->B, maps->Bt,
266 pa_data, diag);
267 }
268}
269
270template<const int T_D1D = 0, const int T_Q1D = 0>
271static void PAVectorMassApply2D(const int NE,
272 const Array<real_t> &B_,
273 const Array<real_t> &Bt_,
274 const Vector &op_,
275 const Vector &x_,
276 Vector &y_,
277 const int d1d = 0,
278 const int q1d = 0)
279{
280 const int D1D = T_D1D ? T_D1D : d1d;
281 const int Q1D = T_Q1D ? T_Q1D : q1d;
282 constexpr int VDIM = 2;
283 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().MAX_D1D, "");
284 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().MAX_Q1D, "");
285 auto B = Reshape(B_.Read(), Q1D, D1D);
286 auto Bt = Reshape(Bt_.Read(), D1D, Q1D);
287 auto op = Reshape(op_.Read(), Q1D, Q1D, NE);
288 auto x = Reshape(x_.Read(), D1D, D1D, VDIM, NE);
289 auto y = Reshape(y_.ReadWrite(), D1D, D1D, VDIM, NE);
290 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
291 {
292 const int D1D = T_D1D ? T_D1D : d1d; // nvcc workaround
293 const int Q1D = T_Q1D ? T_Q1D : q1d;
294 // the following variables are evaluated at compile time
295 constexpr int max_D1D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
296 constexpr int max_Q1D = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
297 real_t sol_xy[max_Q1D][max_Q1D];
298 for (int c = 0; c < VDIM; ++c)
299 {
300 for (int qy = 0; qy < Q1D; ++qy)
301 {
302 for (int qx = 0; qx < Q1D; ++qx)
303 {
304 sol_xy[qy][qx] = 0.0;
305 }
306 }
307 for (int dy = 0; dy < D1D; ++dy)
308 {
309 real_t sol_x[max_Q1D];
310 for (int qy = 0; qy < Q1D; ++qy)
311 {
312 sol_x[qy] = 0.0;
313 }
314 for (int dx = 0; dx < D1D; ++dx)
315 {
316 const real_t s = x(dx,dy,c,e);
317 for (int qx = 0; qx < Q1D; ++qx)
318 {
319 sol_x[qx] += B(qx,dx)* s;
320 }
321 }
322 for (int qy = 0; qy < Q1D; ++qy)
323 {
324 const real_t d2q = B(qy,dy);
325 for (int qx = 0; qx < Q1D; ++qx)
326 {
327 sol_xy[qy][qx] += d2q * sol_x[qx];
328 }
329 }
330 }
331 for (int qy = 0; qy < Q1D; ++qy)
332 {
333 for (int qx = 0; qx < Q1D; ++qx)
334 {
335 sol_xy[qy][qx] *= op(qx,qy,e);
336 }
337 }
338 for (int qy = 0; qy < Q1D; ++qy)
339 {
340 real_t sol_x[max_D1D];
341 for (int dx = 0; dx < D1D; ++dx)
342 {
343 sol_x[dx] = 0.0;
344 }
345 for (int qx = 0; qx < Q1D; ++qx)
346 {
347 const real_t s = sol_xy[qy][qx];
348 for (int dx = 0; dx < D1D; ++dx)
349 {
350 sol_x[dx] += Bt(dx,qx) * s;
351 }
352 }
353 for (int dy = 0; dy < D1D; ++dy)
354 {
355 const real_t q2d = Bt(dy,qy);
356 for (int dx = 0; dx < D1D; ++dx)
357 {
358 y(dx,dy,c,e) += q2d * sol_x[dx];
359 }
360 }
361 }
362 }
363 });
364}
365
366template<const int T_D1D = 0, const int T_Q1D = 0>
367static void PAVectorMassApply3D(const int NE,
368 const Array<real_t> &B_,
369 const Array<real_t> &Bt_,
370 const Vector &op_,
371 const Vector &x_,
372 Vector &y_,
373 const int d1d = 0,
374 const int q1d = 0)
375{
376 const int D1D = T_D1D ? T_D1D : d1d;
377 const int Q1D = T_Q1D ? T_Q1D : q1d;
378 constexpr int VDIM = 3;
379 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().MAX_D1D, "");
380 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().MAX_Q1D, "");
381 auto B = Reshape(B_.Read(), Q1D, D1D);
382 auto Bt = Reshape(Bt_.Read(), D1D, Q1D);
383 auto op = Reshape(op_.Read(), Q1D, Q1D, Q1D, NE);
384 auto x = Reshape(x_.Read(), D1D, D1D, D1D, VDIM, NE);
385 auto y = Reshape(y_.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
386 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
387 {
388 const int D1D = T_D1D ? T_D1D : d1d;
389 const int Q1D = T_Q1D ? T_Q1D : q1d;
390 constexpr int max_D1D = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
391 constexpr int max_Q1D = T_Q1D ? T_Q1D : DofQuadLimits::MAX_Q1D;
392 real_t sol_xyz[max_Q1D][max_Q1D][max_Q1D];
393 for (int c = 0; c < VDIM; ++ c)
394 {
395 for (int qz = 0; qz < Q1D; ++qz)
396 {
397 for (int qy = 0; qy < Q1D; ++qy)
398 {
399 for (int qx = 0; qx < Q1D; ++qx)
400 {
401 sol_xyz[qz][qy][qx] = 0.0;
402 }
403 }
404 }
405 for (int dz = 0; dz < D1D; ++dz)
406 {
407 real_t sol_xy[max_Q1D][max_Q1D];
408 for (int qy = 0; qy < Q1D; ++qy)
409 {
410 for (int qx = 0; qx < Q1D; ++qx)
411 {
412 sol_xy[qy][qx] = 0.0;
413 }
414 }
415 for (int dy = 0; dy < D1D; ++dy)
416 {
417 real_t sol_x[max_Q1D];
418 for (int qx = 0; qx < Q1D; ++qx)
419 {
420 sol_x[qx] = 0;
421 }
422 for (int dx = 0; dx < D1D; ++dx)
423 {
424 const real_t s = x(dx,dy,dz,c,e);
425 for (int qx = 0; qx < Q1D; ++qx)
426 {
427 sol_x[qx] += B(qx,dx) * s;
428 }
429 }
430 for (int qy = 0; qy < Q1D; ++qy)
431 {
432 const real_t wy = B(qy,dy);
433 for (int qx = 0; qx < Q1D; ++qx)
434 {
435 sol_xy[qy][qx] += wy * sol_x[qx];
436 }
437 }
438 }
439 for (int qz = 0; qz < Q1D; ++qz)
440 {
441 const real_t wz = B(qz,dz);
442 for (int qy = 0; qy < Q1D; ++qy)
443 {
444 for (int qx = 0; qx < Q1D; ++qx)
445 {
446 sol_xyz[qz][qy][qx] += wz * sol_xy[qy][qx];
447 }
448 }
449 }
450 }
451 for (int qz = 0; qz < Q1D; ++qz)
452 {
453 for (int qy = 0; qy < Q1D; ++qy)
454 {
455 for (int qx = 0; qx < Q1D; ++qx)
456 {
457 sol_xyz[qz][qy][qx] *= op(qx,qy,qz,e);
458 }
459 }
460 }
461 for (int qz = 0; qz < Q1D; ++qz)
462 {
463 real_t sol_xy[max_D1D][max_D1D];
464 for (int dy = 0; dy < D1D; ++dy)
465 {
466 for (int dx = 0; dx < D1D; ++dx)
467 {
468 sol_xy[dy][dx] = 0;
469 }
470 }
471 for (int qy = 0; qy < Q1D; ++qy)
472 {
473 real_t sol_x[max_D1D];
474 for (int dx = 0; dx < D1D; ++dx)
475 {
476 sol_x[dx] = 0;
477 }
478 for (int qx = 0; qx < Q1D; ++qx)
479 {
480 const real_t s = sol_xyz[qz][qy][qx];
481 for (int dx = 0; dx < D1D; ++dx)
482 {
483 sol_x[dx] += Bt(dx,qx) * s;
484 }
485 }
486 for (int dy = 0; dy < D1D; ++dy)
487 {
488 const real_t wy = Bt(dy,qy);
489 for (int dx = 0; dx < D1D; ++dx)
490 {
491 sol_xy[dy][dx] += wy * sol_x[dx];
492 }
493 }
494 }
495 for (int dz = 0; dz < D1D; ++dz)
496 {
497 const real_t wz = Bt(dz,qz);
498 for (int dy = 0; dy < D1D; ++dy)
499 {
500 for (int dx = 0; dx < D1D; ++dx)
501 {
502 y(dx,dy,dz,c,e) += wz * sol_xy[dy][dx];
503 }
504 }
505 }
506 }
507 }
508 });
509}
510
511static void PAVectorMassApply(const int dim,
512 const int D1D,
513 const int Q1D,
514 const int NE,
515 const Array<real_t> &B,
516 const Array<real_t> &Bt,
517 const Vector &op,
518 const Vector &x,
519 Vector &y)
520{
521 if (dim == 2)
522 {
523 return PAVectorMassApply2D(NE, B, Bt, op, x, y, D1D, Q1D);
524 }
525 if (dim == 3)
526 {
527 return PAVectorMassApply3D(NE, B, Bt, op, x, y, D1D, Q1D);
528 }
529 MFEM_ABORT("Unknown kernel.");
530}
531
533{
534 if (DeviceCanUseCeed())
535 {
536 ceedOp->AddMult(x, y);
537 }
538 else
539 {
540 PAVectorMassApply(dim, dofs1D, quad1D, ne, maps->B, maps->Bt, pa_data, x, y);
541 }
542}
543
544} // namespace mfem
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:337
A coefficient that is constant across space and time.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:274
@ 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
Array< real_t > Bt
Transpose of B.
Definition fe_base.hpp:199
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:709
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
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:3871
Abstract class for all finite elements.
Definition fe_base.hpp:244
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:365
Vector J
Jacobians of the element transformations at all quadrature points.
Definition mesh.hpp:2976
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
const Array< real_t > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition intrules.cpp:86
const IntegrationRule * IntRule
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &Trans)
Mesh data type.
Definition mesh.hpp:64
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
ElementTransformation * GetTypicalElementTransformation()
If the local mesh is not empty return GetElementTransformation(0); otherwise, return the identity tra...
Definition mesh.cpp:390
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:880
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:7243
void AssembleDiagonalPA(Vector &diag) override
Assemble diagonal and add it to Vector diag.
const DofToQuad * maps
Not owned.
void AddMultPA(const Vector &x, Vector &y) const override
Method for partially assembled action.
void AssemblePA(const FiniteElementSpace &fes) override
Method defining partial assembly.
const GeometricFactors * geom
Not owned.
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:494
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:510
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:502
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 MassIntegrator with AssemblyLevel::Partial using libCEED.
Definition mass.hpp:27
int dim
Definition ex24.cpp:53
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:131
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
Definition util.cpp:33
float real_t
Definition config.hpp:43
void forall(int N, lambda &&body)
Definition forall.hpp:753
static const DeviceDofQuadLimits & Get()
Return a const reference to the DeviceDofQuadLimits singleton.
Definition forall.hpp:128