MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
parameterspace.hpp
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#pragma once
12
13#include "../fe/fe_base.hpp"
14#include "../../fem/fespace.hpp"
15
16namespace mfem::future
17{
18
19/// Base class for parametric spaces
21{
22public:
24
25 /// @brief Get vector dimension at each point
26 ///
27 /// This is the number of components at each point in the parametric space.
28 int GetVDim() const { return vdim; }
29
30 /// Get DofToQuad information
31 const DofToQuad& GetDofToQuad() const { return dtq; }
32
33 /// Get total size of the space (T-vector size)
34 ///
35 /// returns the true size vsize of the space
36 virtual int GetTrueVSize() const = 0;
37
38 /// Get local vector size (L-vector size)
39 ///
40 /// returns the local size of the space
41 virtual int GetVSize() const = 0;
42
43 /// Get spatial dimension
44 ///
45 /// returns always 1.
46 int Dimension() const
47 {
48 return 1;
49 }
50
51 /// @brief Get T-vector to L-vector transformation
52 ///
53 /// returns identity by default that is lazy evaluated.
54 virtual const Operator* GetProlongationMatrix() const
55 {
56 if (!prolongation)
57 {
59 }
60 return prolongation.get();
61 }
62
63 /// @brief Get L-vector to E-vector transformation
64 /// @note This is a mock call to replicate interface of FiniteElementSpace.
65 /// It should not be used by a user.
66 ///
67 /// returns identity by default that is lazy evaluated.
69 {
70 if (!elem_restr)
71 {
73 }
74 return elem_restr.get();
75 }
76
77protected:
78 int vdim;
80 mutable std::unique_ptr<Operator> prolongation;
81 mutable std::unique_ptr<Operator> elem_restr;
82};
83
84/// @brief Uniform parameter space
86{
87public:
88 /// @brief Constructor for a uniform parameter space
89 ///
90 /// @param mesh The mesh to determine dimension and number of elements.
91 /// @param ir The integration rule to determine the number of quadrature points.
92 /// @param vdim The vector dimension at each point.
93 /// @param used_in_tensor_product If true, the number of quadrature points is
94 /// calculated as the nth root of the number of points in the integration rule,
95 /// where n is the mesh dimension. If false, the number of quadrature points is
96 /// taken directly from the integration rule.
98 bool used_in_tensor_product = true) :
100 {
101 // Setup DofToQuad information
102 dtq.nqpt = (int)floor(std::pow(ir.GetNPoints(), 1.0 / mesh.Dimension()) + 0.5);
103 dtq.ndof = dtq.nqpt;
104 dtq.mode = used_in_tensor_product ? DofToQuad::TENSOR : DofToQuad::FULL;
105
106 // Calculate sizes
107 const int num_qp = used_in_tensor_product ?
108 static_cast<int>(std::pow(dtq.nqpt, mesh.Dimension())) :
109 ir.GetNPoints();
110
111 tsize = vdim * num_qp * mesh.GetNE();
112 lsize = tsize;
113 }
114
115 int GetTrueVSize() const override
116 {
117 return tsize;
118 }
119
120 int GetVSize() const override
121 {
122 return lsize;
123 }
124
125private:
126 /// T-vector size
127 int tsize;
128
129 /// L-vector size
130 int lsize;
131};
132
134{
135public:
140
141 /// @brief Get the ParameterSpace
143 {
144 return space;
145 }
146
147 using Vector::operator=;
148
149private:
150 /// the parametric space
151 ParameterSpace &space;
152};
153
154} // namespace mfem::future
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition fe_base.hpp:141
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
Definition fe_base.hpp:174
@ FULL
Full multidimensional representation which does not use tensor product structure. The ordering of the...
Definition fe_base.hpp:158
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:165
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
Identity Operator I: x -> x.
Definition operator.hpp:815
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
Mesh data type.
Definition mesh.hpp:65
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
Abstract operator.
Definition operator.hpp:25
Vector data type.
Definition vector.hpp:82
ParameterFunction(ParameterSpace &space)
const ParameterSpace & GetParameterSpace() const
Get the ParameterSpace.
Base class for parametric spaces.
virtual const Operator * GetElementRestriction(ElementDofOrdering o) const
Get L-vector to E-vector transformation.
virtual int GetVSize() const =0
std::unique_ptr< Operator > prolongation
const DofToQuad & GetDofToQuad() const
Get DofToQuad information.
virtual const Operator * GetProlongationMatrix() const
Get T-vector to L-vector transformation.
virtual int GetTrueVSize() const =0
int GetVDim() const
Get vector dimension at each point.
std::unique_ptr< Operator > elem_restr
UniformParameterSpace(Mesh &mesh, const IntegrationRule &ir, int vdim, bool used_in_tensor_product=true)
Constructor for a uniform parameter space.
string space
int GetTrueVSize(const FieldDescriptor &f)
Get the true dof size of a field descriptor.
Definition util.hpp:829
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:47