MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_mixedcurl_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"
13#include "../gridfunc.hpp"
14#include "../qfunction.hpp"
17
18namespace mfem
19{
20
22 const FiniteElementSpace &test_fes)
23{
24 // Assumes tensor-product elements
25 Mesh *mesh = trial_fes.GetMesh();
26 const FiniteElement *fel = trial_fes.GetTypicalFE(); // In H(curl)
27 const FiniteElement *eltest = test_fes.GetTypicalFE(); // In scalar space
28
30 dynamic_cast<const VectorTensorFiniteElement*>(fel);
31 MFEM_VERIFY(el != NULL, "Only VectorTensorFiniteElement is supported!");
32
34 {
35 MFEM_ABORT("Unknown kernel.");
36 }
37
38 const IntegrationRule *ir
39 = IntRule ? IntRule : &MassIntegrator::GetRule(*eltest, *eltest,
41
42 const int dims = el->GetDim();
43 MFEM_VERIFY(dims == 2, "");
44
45 const int nq = ir->GetNPoints();
46 dim = mesh->Dimension();
47 MFEM_VERIFY(dim == 2, "");
48
49 ne = test_fes.GetNE();
52 dofs1D = mapsC->ndof;
53 quad1D = mapsC->nqpt;
54
55 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
56
57 if (el->GetOrder() == eltest->GetOrder())
58 {
60 }
61 else
62 {
63 dofs1Dtest = dofs1D - 1;
64 }
65
67
68 QuadratureSpace qs(*mesh, *ir);
70
71 if (dim == 2)
72 {
73 internal::PAHcurlL2Setup2D(quad1D, ne, ir->GetWeights(), coeff, pa_data);
74 }
75 else
76 {
77 MFEM_ABORT("Unsupported dimension!");
78 }
79}
80
82{
83 if (dim == 2)
84 {
85 internal::PAHcurlL2Apply2D(dofs1D, dofs1Dtest, quad1D, ne, mapsO->B,
87 x, y);
88 }
89 else
90 {
91 MFEM_ABORT("Unsupported dimension!");
92 }
93}
94
96 Vector &y) const
97{
98 if (dim == 2)
99 {
100 internal::PAHcurlL2ApplyTranspose2D(dofs1D, dofs1Dtest, quad1D, ne, mapsO->B,
101 mapsO->Bt, mapsC->B, mapsC->Gt, pa_data,
102 x, y);
103 }
104 else
105 {
106 MFEM_ABORT("Unsupported dimension!");
107 }
108}
109
111 const FiniteElementSpace &test_fes)
112{
113 // Assumes tensor-product elements, with vector test and trial spaces.
114 Mesh *mesh = trial_fes.GetMesh();
115 const FiniteElement *trial_fel = trial_fes.GetTypicalFE();
116 const FiniteElement *test_fel = test_fes.GetTypicalFE();
117
118 const VectorTensorFiniteElement *trial_el =
119 dynamic_cast<const VectorTensorFiniteElement*>(trial_fel);
120 MFEM_VERIFY(trial_el != NULL, "Only VectorTensorFiniteElement is supported!");
121
122 const VectorTensorFiniteElement *test_el =
123 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
124 MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
125
126 const IntegrationRule *ir
127 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
129 const int dims = trial_el->GetDim();
130 MFEM_VERIFY(dims == 3, "");
131
132 const int nq = ir->GetNPoints();
133 dim = mesh->Dimension();
134 MFEM_VERIFY(dim == 3, "");
135
136 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
137
138 ne = trial_fes.GetNE();
140 mapsC = &trial_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
141 mapsO = &trial_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
142 mapsCtest = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
143 mapsOtest = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
144 dofs1D = mapsC->ndof;
145 quad1D = mapsC->nqpt;
146 dofs1Dtest = mapsCtest->ndof;
147
148 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
149
150 testType = test_el->GetDerivType();
151 trialType = trial_el->GetDerivType();
152
153 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
154 coeffDim = (DQ ? 3 : 1);
155
156 const bool curlSpaces = (testType == mfem::FiniteElement::CURL &&
157 trialType == mfem::FiniteElement::CURL);
158
159 const int ndata = curlSpaces ? (coeffDim == 1 ? 1 : 9) : symmDims;
160 pa_data.SetSize(ndata * nq * ne, Device::GetMemoryType());
161
162 QuadratureSpace qs(*mesh, *ir);
164 if (Q) { coeff.Project(*Q); }
165 else if (DQ) { coeff.Project(*DQ); }
166 else { coeff.SetConstant(1.0); }
167
168 if (testType == mfem::FiniteElement::CURL &&
169 trialType == mfem::FiniteElement::CURL && dim == 3)
170 {
171 if (coeffDim == 1)
172 {
173 internal::PAHcurlL2Setup3D(nq, coeffDim, ne, ir->GetWeights(), coeff, pa_data);
174 }
175 else
176 {
177 internal::PAHcurlHdivMassSetup3D(quad1D, coeffDim, ne, false, ir->GetWeights(),
178 geom->J, coeff, pa_data);
179 }
180 }
181 else if (testType == mfem::FiniteElement::DIV &&
182 trialType == mfem::FiniteElement::CURL && dim == 3 &&
183 test_fel->GetOrder() == trial_fel->GetOrder())
184 {
185 internal::PACurlCurlSetup3D(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
186 coeff, pa_data);
187 }
188 else
189 {
190 MFEM_ABORT("Unknown kernel.");
191 }
192}
193
195{
196 if (testType == mfem::FiniteElement::CURL &&
197 trialType == mfem::FiniteElement::CURL && dim == 3)
198 {
199 const int ndata = coeffDim == 1 ? 1 : 9;
200
202 {
203 const int ID = (dofs1D << 4) | quad1D;
204 switch (ID)
205 {
206 case 0x23:
207 return internal::SmemPAHcurlL2Apply3D<2,3>(
208 dofs1D, quad1D, ndata, ne,
209 mapsO->B, mapsC->B, mapsC->G,
210 pa_data, x, y);
211 case 0x34:
212 return internal::SmemPAHcurlL2Apply3D<3,4>(
213 dofs1D, quad1D, ndata, ne,
214 mapsO->B, mapsC->B, mapsC->G,
215 pa_data, x, y);
216 case 0x45:
217 return internal::SmemPAHcurlL2Apply3D<4,5>(
218 dofs1D, quad1D, ndata, ne,
219 mapsO->B, mapsC->B, mapsC->G,
220 pa_data, x, y);
221 case 0x56:
222 return internal::SmemPAHcurlL2Apply3D<5,6>(
223 dofs1D, quad1D, ndata, ne,
224 mapsO->B, mapsC->B, mapsC->G,
225 pa_data, x, y);
226 default:
227 return internal::SmemPAHcurlL2Apply3D(
228 dofs1D, quad1D, ndata, ne,
229 mapsO->B, mapsC->B, mapsC->G,
230 pa_data, x, y);
231 }
232 }
233 else
234 {
235 internal::PAHcurlL2Apply3D(dofs1D, quad1D, ndata, ne, mapsO->B, mapsC->B,
236 mapsO->Bt, mapsC->Bt, mapsC->G, pa_data, x, y);
237 }
238 }
239 else if (testType == mfem::FiniteElement::DIV &&
240 trialType == mfem::FiniteElement::CURL && dim == 3)
241 {
242 internal::PAHcurlHdivApply3D(dofs1D, dofs1Dtest, quad1D, ne, mapsO->B,
243 mapsC->B, mapsOtest->Bt, mapsCtest->Bt, mapsC->G,
244 pa_data, x, y);
245 }
246 else
247 {
248 MFEM_ABORT("Unsupported dimension or space!");
249 }
250}
251
253 Vector &y) const
254{
255 if (testType == mfem::FiniteElement::DIV &&
256 trialType == mfem::FiniteElement::CURL && dim == 3)
257 {
258 internal::PAHcurlHdivApplyTranspose3D(dofs1D, dofs1Dtest, quad1D, ne, mapsO->B,
259 mapsC->B, mapsOtest->Bt, mapsCtest->Bt,
260 mapsC->Gt, pa_data, x, y);
261 }
262 else
263 {
264 MFEM_ABORT("Unsupported dimension or space!");
265 }
266}
267
269 &trial_fes,
270 const FiniteElementSpace &test_fes)
271{
272 // Assumes tensor-product elements, with vector test and trial spaces.
273 Mesh *mesh = trial_fes.GetMesh();
274 const FiniteElement *trial_fel = trial_fes.GetTypicalFE();
275 const FiniteElement *test_fel = test_fes.GetTypicalFE();
276
277 const VectorTensorFiniteElement *trial_el =
278 dynamic_cast<const VectorTensorFiniteElement*>(trial_fel);
279 MFEM_VERIFY(trial_el != NULL, "Only VectorTensorFiniteElement is supported!");
280
281 const VectorTensorFiniteElement *test_el =
282 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
283 MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
284
285 const IntegrationRule *ir
286 = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
288 const int dims = trial_el->GetDim();
289 MFEM_VERIFY(dims == 3, "");
290
291 const int nq = ir->GetNPoints();
292 dim = mesh->Dimension();
293 MFEM_VERIFY(dim == 3, "");
294
295 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
296
297 ne = trial_fes.GetNE();
299 mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
300 mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
301 dofs1D = mapsC->ndof;
302 quad1D = mapsC->nqpt;
303
304 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
305
306 testType = test_el->GetDerivType();
307 trialType = trial_el->GetDerivType();
308
309 const bool curlSpaces = (testType == mfem::FiniteElement::CURL &&
310 trialType == mfem::FiniteElement::CURL);
311
312 const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
313
314 coeffDim = DQ ? 3 : 1;
315 const int ndata = curlSpaces ? (DQ ? 9 : 1) : symmDims;
316
317 pa_data.SetSize(ndata * nq * ne, Device::GetMemoryType());
318
319 QuadratureSpace qs(*mesh, *ir);
321 if (Q) { coeff.Project(*Q); }
322 else if (DQ) { coeff.Project(*DQ); }
323 else { coeff.SetConstant(1.0); }
324
325 if (trialType == mfem::FiniteElement::CURL && dim == 3)
326 {
327 if (coeffDim == 1)
328 {
329 internal::PAHcurlL2Setup3D(nq, coeffDim, ne, ir->GetWeights(), coeff, pa_data);
330 }
331 else
332 {
333 internal::PAHcurlHdivMassSetup3D(quad1D, coeffDim, ne, false, ir->GetWeights(),
334 geom->J, coeff, pa_data);
335 }
336 }
337 else if (trialType == mfem::FiniteElement::DIV && dim == 3 &&
338 test_el->GetOrder() == trial_el->GetOrder())
339 {
340 internal::PACurlCurlSetup3D(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
341 coeff, pa_data);
342 }
343 else
344 {
345 MFEM_ABORT("Unknown kernel.");
346 }
347}
348
350{
351 if (testType == mfem::FiniteElement::CURL &&
352 trialType == mfem::FiniteElement::CURL && dim == 3)
353 {
354 const int ndata = coeffDim == 1 ? 1 : 9;
356 {
357 const int ID = (dofs1D << 4) | quad1D;
358 switch (ID)
359 {
360 case 0x23:
361 return internal::SmemPAHcurlL2ApplyTranspose3D<2,3>(
362 dofs1D, quad1D, ndata,
363 ne, mapsO->B, mapsC->B,
364 mapsC->G, pa_data, x, y);
365 case 0x34:
366 return internal::SmemPAHcurlL2ApplyTranspose3D<3,4>(
367 dofs1D, quad1D, ndata,
368 ne, mapsO->B, mapsC->B,
369 mapsC->G, pa_data, x, y);
370 case 0x45:
371 return internal::SmemPAHcurlL2ApplyTranspose3D<4,5>(
372 dofs1D, quad1D, ndata,
373 ne, mapsO->B, mapsC->B,
374 mapsC->G, pa_data, x, y);
375 case 0x56:
376 return internal::SmemPAHcurlL2ApplyTranspose3D<5,6>(
377 dofs1D, quad1D, ndata,
378 ne, mapsO->B, mapsC->B,
379 mapsC->G, pa_data, x, y);
380 default:
381 return internal::SmemPAHcurlL2ApplyTranspose3D(
382 dofs1D, quad1D, ndata, ne,
383 mapsO->B, mapsC->B,
384 mapsC->G, pa_data, x, y);
385 }
386 }
387 else
388 {
389 internal::PAHcurlL2ApplyTranspose3D(dofs1D, quad1D, ndata, ne, mapsO->B,
390 mapsC->B, mapsO->Bt, mapsC->Bt, mapsC->Gt,
391 pa_data, x, y);
392 }
393 }
394 else if (testType == mfem::FiniteElement::CURL &&
395 trialType == mfem::FiniteElement::DIV && dim == 3)
396 {
397 internal::PAHcurlHdivApplyTranspose3D(dofs1D, dofs1D, quad1D, ne, mapsO->B,
398 mapsC->B, mapsO->Bt, mapsC->Bt,
399 mapsC->Gt, pa_data, x, y);
400 }
401 else
402 {
403 MFEM_ABORT("Unsupported dimension or space!");
404 }
405}
406
408 Vector &y) const
409{
410 if (testType == mfem::FiniteElement::CURL &&
411 trialType == mfem::FiniteElement::DIV && dim == 3)
412 {
413 internal::PAHcurlHdivApply3D(dofs1D, dofs1D, quad1D, ne, mapsO->B,
414 mapsC->B, mapsO->Bt, mapsC->Bt, mapsC->G,
415 pa_data, x, y);
416 }
417 else
418 {
419 MFEM_ABORT("Unsupported dimension or space!");
420 }
421}
422
423} // namespace mfem
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.
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition device.hpp:278
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition device.hpp:259
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
Array< real_t > Gt
Transpose of G.
Definition fe_base.hpp:221
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
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:891
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
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:338
int GetDerivType() const
Returns the FiniteElement::DerivType of the element describing the spatial derivative method implemen...
Definition fe_base.hpp:365
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:321
@ DIV
Implements CalcDivShape methods.
Definition fe_base.hpp:306
@ CURL
Implements CalcCurlShape methods.
Definition fe_base.hpp:307
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 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
const DofToQuad * mapsO
Not owned. DOF-to-quad map, open.
const DofToQuad * mapsC
Not owned. DOF-to-quad map, closed.
void AddMultTransposePA(const Vector &x, Vector &y) const override
Method for partially assembled transposed action.
void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
void AddMultPA(const Vector &, Vector &) const override
Method for partially assembled action.
void AddMultTransposePA(const Vector &, Vector &) const override
Method for partially assembled transposed action.
void AddMultPA(const Vector &, Vector &) const override
Method for partially assembled action.
void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
DiagonalMatrixCoefficient * DQ
void AddMultPA(const Vector &, Vector &) const override
Method for partially assembled action.
void AddMultTransposePA(const Vector &, Vector &) const override
Method for partially assembled transposed action.
void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:120
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
Definition fe_base.hpp:1362
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const override
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.hpp:1354
Vector data type.
Definition vector.hpp:82
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
@ FULL
Store the coefficient as a full QuadratureFunction.
@ DEVICE_MASK
Biwise-OR of all device backends.
Definition device.hpp:97