MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_vectorfediv_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"
16
17namespace mfem
18{
19
20void
22 const FiniteElementSpace &test_fes)
23{
24 // Assumes tensor-product elements, with a vector test space and
25 // scalar trial space.
26 Mesh *mesh = trial_fes.GetMesh();
27 const FiniteElement *trial_fel = trial_fes.GetTypicalFE();
28 const FiniteElement *test_fel = test_fes.GetTypicalFE();
29
30 const VectorTensorFiniteElement *trial_el =
31 dynamic_cast<const VectorTensorFiniteElement*>(trial_fel);
32 MFEM_VERIFY(trial_el != NULL, "Only VectorTensorFiniteElement is supported!");
33
34 const NodalTensorFiniteElement *test_el =
35 dynamic_cast<const NodalTensorFiniteElement*>(test_fel);
36 MFEM_VERIFY(test_el != NULL, "Only NodalTensorFiniteElement is supported!");
37
39 *trial_el, *trial_el,
41
42 const int dims = trial_el->GetDim();
43 MFEM_VERIFY(dims == 2 || dims == 3, "");
44
45 const int nq = ir->GetNPoints();
46 dim = mesh->Dimension();
47 MFEM_VERIFY(dim == 2 || dim == 3, "");
48
49 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder() + 1, "");
50
51 ne = trial_fes.GetNE();
52 mapsC = &trial_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
53 mapsO = &trial_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
54 dofs1D = mapsC->ndof;
55 quad1D = mapsC->nqpt;
56
57 L2mapsO = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
58 L2dofs1D = L2mapsO->ndof;
59
60 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
61 if (dim == 2)
62 {
63 MFEM_VERIFY(nq == quad1D * quad1D, "");
64 }
65 else
66 {
67 MFEM_VERIFY(nq == quad1D * quad1D * quad1D, "");
68 }
69
70 pa_data.SetSize(nq * ne, Device::GetMemoryType());
71
72 QuadratureSpace qs(*mesh, *ir);
74
75 if (test_el->GetMapType() == FiniteElement::INTEGRAL)
76 {
77 const GeometricFactors *geom =
79 coeff /= geom->detJ;
80 }
81
82 if (trial_el->GetDerivType() == mfem::FiniteElement::DIV && dim == 3)
83 {
84 internal::PAHdivL2Setup3D(quad1D, ne, ir->GetWeights(), coeff, pa_data);
85 }
86 else if (trial_el->GetDerivType() == mfem::FiniteElement::DIV && dim == 2)
87 {
88 internal::PAHdivL2Setup2D(quad1D, ne, ir->GetWeights(), coeff, pa_data);
89 }
90 else
91 {
92 MFEM_ABORT("Unknown kernel.");
93 }
94}
95
97 Vector &diag)
98{
99 if (dim == 3)
100 {
101 internal::PAHdivL2AssembleDiagonal_ADAt_3D(dofs1D, quad1D, L2dofs1D, ne,
102 L2mapsO->B,
103 mapsC->Gt, mapsO->Bt, pa_data, D, diag);
104 }
105 else if (dim == 2)
106 {
107 internal::PAHdivL2AssembleDiagonal_ADAt_2D(dofs1D, quad1D, L2dofs1D, ne,
108 L2mapsO->B,
109 mapsC->Gt, mapsO->Bt, pa_data, D, diag);
110 }
111 else
112 {
113 MFEM_ABORT("Unsupported dimension!");
114 }
115}
116
118{
119 if (dim == 3)
120 {
121 internal::PAHdivL2Apply3D(dofs1D, quad1D, L2dofs1D, ne, mapsO->B, mapsC->G,
122 L2mapsO->Bt, pa_data, x, y);
123 }
124 else if (dim == 2)
125 {
126 internal::PAHdivL2Apply2D(dofs1D, quad1D, L2dofs1D, ne, mapsO->B, mapsC->G,
127 L2mapsO->Bt, pa_data, x, y);
128 }
129 else
130 {
131 MFEM_ABORT("Unsupported dimension!");
132 }
133}
134
136 Vector &y) const
137{
138 if (dim == 3)
139 {
140 internal::PAHdivL2ApplyTranspose3D(dofs1D, quad1D, L2dofs1D, ne, L2mapsO->B,
141 mapsC->Gt, mapsO->Bt, pa_data, x, y);
142 }
143 else if (dim == 2)
144 {
145 internal::PAHdivL2ApplyTranspose2D(dofs1D, quad1D, L2dofs1D, ne, L2mapsO->B,
146 mapsC->Gt, mapsO->Bt, pa_data, x, y);
147 }
148 else
149 {
150 MFEM_ABORT("Unsupported dimension!");
151 }
152}
153
154} // namespace mfem
Class to represent a coefficient evaluated at quadrature points.
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition device.hpp:278
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
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition fe_base.hpp:360
@ DIV
Implements CalcDivShape methods.
Definition fe_base.hpp:306
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition mesh.hpp:2937
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition mesh.hpp:2982
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 & 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.cpp:2672
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:120
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
void AssembleDiagonalPA_ADAt(const Vector &D, Vector &diag) override
Assemble diagonal of ( is this integrator) and add it to diag.
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.