MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_curlcurl_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 "../qfunction.hpp"
14
15namespace mfem
16{
17
19{
20 // Assumes tensor-product elements
21 Mesh *mesh = fes.GetMesh();
22 const FiniteElement *fel = fes.GetTypicalFE();
23
25 dynamic_cast<const VectorTensorFiniteElement*>(fel);
26 MFEM_VERIFY(el != NULL, "Only VectorTensorFiniteElement is supported!");
27
28 const IntegrationRule *ir
31
32 const int dims = el->GetDim();
33 MFEM_VERIFY(dims == 2 || dims == 3, "");
34
35 nq = ir->GetNPoints();
36 dim = mesh->Dimension();
37 MFEM_VERIFY(dim == 2 || dim == 3, "");
38
39 ne = fes.GetNE();
43 dofs1D = mapsC->ndof;
44 quad1D = mapsC->nqpt;
45
46 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
47
48 QuadratureSpace qs(*mesh, *ir);
50 if (Q) { coeff.Project(*Q); }
51 else if (MQ) { coeff.ProjectTranspose(*MQ); }
52 else if (DQ) { coeff.Project(*DQ); }
53 else { coeff.SetConstant(1.0); }
54
55 const int coeff_dim = coeff.GetVDim();
56 symmetric = (coeff_dim != dim*dim);
57 const int sym_dims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
58 const int ndata = (dim == 2) ? 1 : (symmetric ? sym_dims : dim*dim);
60
62 {
63 MFEM_ABORT("Unknown kernel.");
64 }
65
66 if (dim == 3)
67 {
68 internal::PACurlCurlSetup3D(quad1D, coeff_dim, ne, ir->GetWeights(), geom->J,
69 coeff, pa_data);
70 }
71 else
72 {
73 internal::PACurlCurlSetup2D(quad1D, ne, ir->GetWeights(), geom->J, coeff,
74 pa_data);
75 }
76}
77
79{
80 if (dim == 3)
81 {
83 {
84 const int ID = (dofs1D << 4) | quad1D;
85 switch (ID)
86 {
87 case 0x23:
88 return internal::SmemPACurlCurlAssembleDiagonal3D<2,3>(
89 dofs1D,
90 quad1D,
92 mapsO->B, mapsC->B,
93 mapsO->G, mapsC->G,
94 pa_data, diag);
95 case 0x34:
96 return internal::SmemPACurlCurlAssembleDiagonal3D<3,4>(
97 dofs1D,
98 quad1D,
100 mapsO->B, mapsC->B,
101 mapsO->G, mapsC->G,
102 pa_data, diag);
103 case 0x45:
104 return internal::SmemPACurlCurlAssembleDiagonal3D<4,5>(
105 dofs1D,
106 quad1D,
107 symmetric, ne,
108 mapsO->B, mapsC->B,
109 mapsO->G, mapsC->G,
110 pa_data, diag);
111 case 0x56:
112 return internal::SmemPACurlCurlAssembleDiagonal3D<5,6>(
113 dofs1D,
114 quad1D,
115 symmetric, ne,
116 mapsO->B, mapsC->B,
117 mapsO->G, mapsC->G,
118 pa_data, diag);
119 default:
120 return internal::SmemPACurlCurlAssembleDiagonal3D(
121 dofs1D, quad1D,
122 symmetric, ne,
123 mapsO->B, mapsC->B,
124 mapsO->G, mapsC->G,
125 pa_data, diag);
126 }
127 }
128 else
129 {
130 internal::PACurlCurlAssembleDiagonal3D(dofs1D, quad1D, symmetric, ne,
131 mapsO->B, mapsC->B,
132 mapsO->G, mapsC->G,
133 pa_data, diag);
134 }
135 }
136 else if (dim == 2)
137 {
138 internal::PACurlCurlAssembleDiagonal2D(dofs1D, quad1D, ne,
139 mapsO->B, mapsC->G, pa_data, diag);
140 }
141 else
142 {
143 MFEM_ABORT("Unsupported dimension!");
144 }
145}
146
148{
149 if (dim == 3)
150 {
152 {
153 const int ID = (dofs1D << 4) | quad1D;
154 switch (ID)
155 {
156 case 0x23:
157 return internal::SmemPACurlCurlApply3D<2,3>(
158 dofs1D, quad1D,
159 symmetric, ne,
160 mapsO->B, mapsC->B, mapsO->Bt, mapsC->Bt,
161 mapsC->G, mapsC->Gt, pa_data, x, y);
162 case 0x34:
163 return internal::SmemPACurlCurlApply3D<3,4>(
164 dofs1D, quad1D,
165 symmetric, ne,
166 mapsO->B, mapsC->B, mapsO->Bt, mapsC->Bt,
167 mapsC->G, mapsC->Gt, pa_data, x, y);
168 case 0x45:
169 return internal::SmemPACurlCurlApply3D<4,5>(
170 dofs1D, quad1D,
171 symmetric, ne,
172 mapsO->B, mapsC->B, mapsO->Bt, mapsC->Bt,
173 mapsC->G, mapsC->Gt, pa_data, x, y);
174 case 0x56:
175 return internal::SmemPACurlCurlApply3D<5,6>(
176 dofs1D, quad1D,
177 symmetric, ne,
178 mapsO->B, mapsC->B, mapsO->Bt, mapsC->Bt,
179 mapsC->G, mapsC->Gt, pa_data, x, y);
180 default:
181 return internal::SmemPACurlCurlApply3D(
183 mapsO->B, mapsC->B, mapsO->Bt, mapsC->Bt,
184 mapsC->G, mapsC->Gt, pa_data, x, y);
185 }
186 }
187 else
188 {
189 internal::PACurlCurlApply3D(dofs1D, quad1D, symmetric, ne, mapsO->B, mapsC->B,
190 mapsO->Bt, mapsC->Bt, mapsC->G, mapsC->Gt,
191 pa_data, x, y);
192 }
193 }
194 else if (dim == 2)
195 {
196 internal::PACurlCurlApply2D(dofs1D, quad1D, ne, mapsO->B, mapsO->Bt,
197 mapsC->G, mapsC->Gt, pa_data, x, y);
198 }
199 else
200 {
201 MFEM_ABORT("Unsupported dimension!");
202 }
203}
204
205} // 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.
int GetVDim() const
Return the number of values per quadrature point.
void ProjectTranspose(MatrixCoefficient &coeff)
Project the transpose of coeff.
const GeometricFactors * geom
Not owned.
bool symmetric
False if using a nonsymmetric matrix coefficient.
MatrixCoefficient * MQ
const DofToQuad * mapsC
Not owned. DOF-to-quad map, closed.
void AddMultPA(const Vector &x, Vector &y) const override
Method for partially assembled action.
void AssembleDiagonalPA(Vector &diag) override
Assemble diagonal and add it to Vector diag.
void AssemblePA(const FiniteElementSpace &fes) override
Method defining partial assembly.
const DofToQuad * mapsO
Not owned. DOF-to-quad map, open.
DiagonalMatrixCoefficient * DQ
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 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
@ 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
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
@ SYMMETRIC
Store the triangular part of symmetric matrices.
@ DEVICE_MASK
Biwise-OR of all device backends.
Definition device.hpp:97