MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
linearform.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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#ifndef MFEM_LINEARFORM
13#define MFEM_LINEARFORM
14
15#include "../config/config.hpp"
16#include "lininteg.hpp"
17#include "linearform_ext.hpp"
18#include "gridfunc.hpp"
19
20namespace mfem
21{
22
23/// Vector with associated FE space and LinearFormIntegrators.
24class LinearForm : public Vector
25{
27
28protected:
29 /// FE space on which the LinearForm lives. Not owned.
31
32 /** @brief Extension for supporting different assembly levels. */
34
35 /// @brief Should we use the device-compatible fast assembly algorithm (false
36 /// by default)
37 bool fast_assembly = false;
38
39 /** @brief Indicates the LinearFormIntegrator%s stored in #domain_integs,
40 #domain_delta_integs, #boundary_integs, and #boundary_face_integs are
41 owned by another LinearForm. */
42 int extern_lfs = 0;
43
44 /// Set of Domain Integrators to be applied.
46 /// Element attribute marker (should be of length mesh->attributes.Max() or
47 /// 0 if mesh->attributes is empty)
48 /// Includes all by default.
49 /// 0 - ignore attribute
50 /// 1 - include attribute
52
53 /// Separate array for integrators with delta function coefficients.
55
56 /// Set of Boundary Integrators to be applied.
58 /// Entries are not owned.
60
61 /// Set of Boundary Face Integrators to be applied.
64
65 /// Set of Internal Face Integrators to be applied.
67
68 /// The element ids where the centers of the delta functions lie
70
71 /// The reference coordinates where the centers of the delta functions lie
73
74 /// If true, the delta locations are not (re)computed during assembly.
76 { return (domain_delta_integs_elem_id.Size() != 0); }
77
78 /// Force (re)computation of delta locations.
80
81private:
82 /// Copy construction is not supported; body is undefined.
83 LinearForm(const LinearForm &);
84
85public:
86 /// Creates linear form associated with FE space @a *f.
87 /** The pointer @a f is not owned by the newly constructed object. */
89 { fes = f; UseDevice(true); }
90
91 /** @brief Create a LinearForm on the FiniteElementSpace @a f, using the
92 same integrators as the LinearForm @a lf.
93
94 The pointer @a f is not owned by the newly constructed object.
95
96 The integrators in @a lf are copied as pointers and they are not owned by
97 the newly constructed LinearForm. */
99
100 /// Create an empty LinearForm without an associated FiniteElementSpace.
101 /** The associated FiniteElementSpace can be set later using one of the
102 methods: Update(FiniteElementSpace *) or
103 Update(FiniteElementSpace *, Vector &, int). */
105 { fes = NULL; UseDevice(true); }
106
107 /// Construct a LinearForm using previously allocated array @a data.
108 /** The LinearForm does not assume ownership of @a data which is assumed to
109 be of size at least `f->GetVSize()`. Similar to the Vector constructor
110 for externally allocated array, the pointer @a data can be NULL. The data
111 array can be replaced later using the method SetData(). */
113 { fes = f; }
114
115 /// Copy assignment. Only the data of the base class Vector is copied.
116 /** It is assumed that this object and @a rhs use FiniteElementSpace%s that
117 have the same size.
118
119 @note Defining this method overwrites the implicitly defined copy
120 assignment operator. */
122 { return operator=((const Vector &)rhs); }
123
124 /// (DEPRECATED) Return the FE space associated with the LinearForm.
125 /** @deprecated Use FESpace() instead. */
126 MFEM_DEPRECATED FiniteElementSpace *GetFES() { return fes; }
127
128 /// Read+write access to the associated FiniteElementSpace.
130 /// Read-only access to the associated FiniteElementSpace.
131 const FiniteElementSpace *FESpace() const { return fes; }
132
133 /// Adds new Domain Integrator. Assumes ownership of @a lfi.
135 /// Adds new Domain Integrator restricted to certain elements specified by
136 /// the @a elem_attr_marker.
138 Array<int> &elem_marker);
139
140 /// Adds new Boundary Integrator. Assumes ownership of @a lfi.
142
143 /** @brief Add new Boundary Integrator, restricted to the given boundary
144 attributes.
145
146 Assumes ownership of @a lfi. The array @a bdr_attr_marker is stored
147 internally as a pointer to the given Array<int> object. */
149 Array<int> &bdr_attr_marker);
150
151 /// Adds new Boundary Face Integrator. Assumes ownership of @a lfi.
153
154 /** @brief Add new Boundary Face Integrator, restricted to the given boundary
155 attributes.
156
157 Assumes ownership of @a lfi. The array @a bdr_attr_marker is stored
158 internally as a pointer to the given Array<int> object. */
160 Array<int> &bdr_attr_marker);
161
162 /// Adds new Interior Face Integrator. Assumes ownership of @a lfi.
164
165 /** @brief Access all integrators added with AddDomainIntegrator() which are
166 not DeltaLFIntegrator%s or they are DeltaLFIntegrator%s with non-delta
167 coefficients. */
169
170 /** @brief Access all domain markers added with AddDomainIntegrator().
171 If no marker was specified when the integrator was added, the
172 corresponding pointer (to Array<int>) will be NULL. */
174
175 /** @brief Access all integrators added with AddDomainIntegrator() which are
176 DeltaLFIntegrator%s with delta coefficients. */
178
179 /// Access all integrators added with AddBoundaryIntegrator().
181
182 /// Access all integrators added with AddBdrFaceIntegrator().
184
185 /// Access all integrators added with AddInteriorFaceIntegrator().
187
188 /** @brief Access all boundary markers added with AddBdrFaceIntegrator().
189 If no marker was specified when the integrator was added, the
190 corresponding pointer (to Array<int>) will be NULL. */
192
193 /// @brief Which assembly algorithm to use: the new device-compatible fast
194 /// assembly (true), or the legacy CPU-only algorithm (false).
195 /** If not set, the default value is false. If used, this method must be
196 called before assembly. */
197 void UseFastAssembly(bool use_fa);
198
199 /// Assembles the linear form i.e. sums over all domain/bdr integrators.
200 /** When @ref UseFastAssembly "UseFastAssembly(true)" has been called and the
201 linearform assembly is compatible with device execution, it will be
202 executed on the device. */
203 void Assemble();
204
205 /// Return true if assembly on device is supported, false otherwise.
206 virtual bool SupportsDevice() const;
207
208 /// Assembles delta functions of the linear form
209 void AssembleDelta();
210
211 /// Update the object according to the associated FE space #fes.
212 /** This method should be called when the associated FE space #fes has been
213 updated, e.g. after its associated Mesh object has been refined.
214
215 @note This method does not perform assembly. */
216 void Update();
217
218 /// Associate a new FE space, @a *f, with this object and Update() it. */
220
221 /** @brief Associate a new FE space, @a *f, with this object and use the data
222 of @a v, offset by @a v_offset, to initialize this object's Vector::data.
223
224 @note This method does not perform assembly. */
225 void Update(FiniteElementSpace *f, Vector &v, int v_offset);
226
227 /** @brief Make the LinearForm reference external data on a new
228 FiniteElementSpace. */
229 /** This method changes the FiniteElementSpace associated with the LinearForm
230 @a *f and sets the data of the Vector @a v (plus the @a v_offset) as
231 external data in the LinearForm.
232
233 @note This version of the method will also perform bounds checks when the
234 build option MFEM_DEBUG is enabled. */
235 virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
236
237 /// Return the action of the LinearForm as a linear mapping.
238 /** Linear forms are linear functionals which map GridFunctions to
239 the real numbers. This method performs this mapping which in
240 this case is equivalent as an inner product of the LinearForm
241 and GridFunction. */
242 real_t operator()(const GridFunction &gf) const { return (*this)*gf; }
243
244 /// Redefine '=' for LinearForm = constant.
246
247 /// Copy the data from @a v.
248 /** The size of @a v must be equal to the size of the associated
249 FiniteElementSpace #fes. */
250 LinearForm &operator=(const Vector &v);
251
252 /// Destroys linear form.
253 ~LinearForm();
254};
255
256}
257
258#endif
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
Class extending the LinearForm class to support assembly on devices.
Abstract base class LinearFormIntegrator.
Definition lininteg.hpp:25
Vector with associated FE space and LinearFormIntegrators.
void AddInteriorFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Interior Face Integrator. Assumes ownership of lfi.
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Array< Array< int > * > boundary_face_integs_marker
Entries not owned.
FiniteElementSpace * FESpace()
Read+write access to the associated FiniteElementSpace.
void AssembleDelta()
Assembles delta functions of the linear form.
virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset)
Make the LinearForm reference external data on a new FiniteElementSpace.
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Array< Array< int > * > * GetFLFI_Marker()
Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the in...
MFEM_DEPRECATED FiniteElementSpace * GetFES()
(DEPRECATED) Return the FE space associated with the LinearForm.
LinearForm()
Create an empty LinearForm without an associated FiniteElementSpace.
Array< LinearFormIntegrator * > boundary_face_integs
Set of Boundary Face Integrators to be applied.
LinearForm(FiniteElementSpace *f)
Creates linear form associated with FE space *f.
void UseFastAssembly(bool use_fa)
Which assembly algorithm to use: the new device-compatible fast assembly (true), or the legacy CPU-on...
Array< int > domain_delta_integs_elem_id
The element ids where the centers of the delta functions lie.
virtual bool SupportsDevice() const
Return true if assembly on device is supported, false otherwise.
LinearForm & operator=(const LinearForm &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Array< Array< int > * > boundary_integs_marker
Entries are not owned.
bool HaveDeltaLocations()
If true, the delta locations are not (re)computed during assembly.
FiniteElementSpace * fes
FE space on which the LinearForm lives. Not owned.
LinearForm(FiniteElementSpace *f, real_t *data)
Construct a LinearForm using previously allocated array data.
Array< Array< int > * > domain_integs_marker
Array< LinearFormIntegrator * > interior_face_integs
Set of Internal Face Integrators to be applied.
real_t operator()(const GridFunction &gf) const
Return the action of the LinearForm as a linear mapping.
LinearFormExtension * ext
Extension for supporting different assembly levels.
void ResetDeltaLocations()
Force (re)computation of delta locations.
void Update()
Update the object according to the associated FE space fes.
Array< LinearFormIntegrator * > * GetFLFI()
Access all integrators added with AddBdrFaceIntegrator().
Array< DeltaLFIntegrator * > domain_delta_integs
Separate array for integrators with delta function coefficients.
Array< LinearFormIntegrator * > domain_integs
Set of Domain Integrators to be applied.
Array< LinearFormIntegrator * > * GetDLFI()
Access all integrators added with AddDomainIntegrator() which are not DeltaLFIntegrators or they are ...
Array< LinearFormIntegrator * > * GetBLFI()
Access all integrators added with AddBoundaryIntegrator().
void Update(FiniteElementSpace *f)
Associate a new FE space, *f, with this object and Update() it. *‍/.
~LinearForm()
Destroys linear form.
Array< LinearFormIntegrator * > * GetIFLFI()
Access all integrators added with AddInteriorFaceIntegrator().
Array< IntegrationPoint > domain_delta_integs_ip
The reference coordinates where the centers of the delta functions lie.
Array< Array< int > * > * GetDLFI_Marker()
Access all domain markers added with AddDomainIntegrator(). If no marker was specified when the integ...
const FiniteElementSpace * FESpace() const
Read-only access to the associated FiniteElementSpace.
Array< LinearFormIntegrator * > boundary_integs
Set of Boundary Integrators to be applied.
Array< DeltaLFIntegrator * > * GetDLFI_Delta()
Access all integrators added with AddDomainIntegrator() which are DeltaLFIntegrators with delta coeff...
int extern_lfs
Indicates the LinearFormIntegrators stored in domain_integs, domain_delta_integs, boundary_integs,...
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
bool fast_assembly
Should we use the device-compatible fast assembly algorithm (false by default)
Vector data type.
Definition vector.hpp:80
Memory< real_t > data
Definition vector.hpp:83
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition vector.hpp:139
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30