MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_divdiv_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
21{
22 // Assumes tensor-product elements
23 Mesh *mesh = fes.GetMesh();
24 const FiniteElement *fel = fes.GetTypicalFE();
25
27 dynamic_cast<const VectorTensorFiniteElement*>(fel);
28 MFEM_VERIFY(el != NULL, "Only VectorTensorFiniteElement is supported!");
29
31 (*el, *el, *mesh->GetTypicalElementTransformation());
32
33 const int dims = el->GetDim();
34 MFEM_VERIFY(dims == 2 || dims == 3, "");
35
36 const int nq = ir->GetNPoints();
37 dim = mesh->Dimension();
38 MFEM_VERIFY(dim == 2 || dim == 3, "");
39
40 ne = fes.GetNE();
42 mapsC = &el->GetDofToQuad(*ir, DofToQuad::TENSOR);
43 mapsO = &el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
44 dofs1D = mapsC->ndof;
45 quad1D = mapsC->nqpt;
46
47 MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
48
49 pa_data.SetSize(nq * ne, Device::GetMemoryType());
50
51 QuadratureSpace qs(*mesh, *ir);
53
54 if (el->GetDerivType() == mfem::FiniteElement::DIV && dim == 3)
55 {
56 internal::PADivDivSetup3D(quad1D, ne, ir->GetWeights(), geom->J, coeff,
57 pa_data);
58 }
59 else if (el->GetDerivType() == mfem::FiniteElement::DIV && dim == 2)
60 {
61 internal::PADivDivSetup2D(quad1D, ne, ir->GetWeights(), geom->J, coeff,
62 pa_data);
63 }
64 else
65 {
66 MFEM_ABORT("Unknown kernel.");
67 }
68}
69
71{
72 if (dim == 3)
73 {
74 internal::PADivDivAssembleDiagonal3D(dofs1D, quad1D, ne,
75 mapsO->B, mapsC->G, pa_data, diag);
76 }
77 else
78 {
79 internal::PADivDivAssembleDiagonal2D(dofs1D, quad1D, ne,
80 mapsO->B, mapsC->G, pa_data, diag);
81 }
82}
83
85{
86 if (dim == 3)
87 internal::PADivDivApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->G,
88 mapsO->Bt, mapsC->Gt, pa_data, x, y);
89 else if (dim == 2)
90 internal::PADivDivApply2D(dofs1D, quad1D, ne, mapsO->B, mapsC->G,
91 mapsO->Bt, mapsC->Gt, pa_data, x, y);
92 else
93 {
94 MFEM_ABORT("Unsupported dimension!");
95 }
96}
97
98} // namespace mfem
Class to represent a coefficient evaluated at quadrature points.
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition device.hpp:278
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.
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
@ DIV
Implements CalcDivShape methods.
Definition fe_base.hpp:306
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
@ FULL
Store the coefficient as a full QuadratureFunction.