MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::GridFunction Class Reference

Class for grid function - Vector with associated FE space. More...

#include <gridfunc.hpp>

Inheritance diagram for mfem::GridFunction:
[legend]
Collaboration diagram for mfem::GridFunction:
[legend]

Public Types

enum  AvgType { ARITHMETIC , HARMONIC }
 

Public Member Functions

 GridFunction ()
 
 GridFunction (const GridFunction &orig)
 Copy constructor. The internal true-dof vector t_vec is not copied.
 
 GridFunction (FiniteElementSpace *f)
 Construct a GridFunction associated with the FiniteElementSpace *f.
 
 GridFunction (FiniteElementSpace *f, real_t *data)
 Construct a GridFunction using previously allocated array data.
 
 GridFunction (FiniteElementSpace *f, Vector &base, int base_offset=0)
 Construct a GridFunction using previously allocated Vector base starting at the given offset, base_offset.
 
 GridFunction (Mesh *m, std::istream &input)
 Construct a GridFunction on the given Mesh, using the data from input.
 
 GridFunction (Mesh *m, GridFunction *gf_array[], int num_pieces)
 
GridFunctionoperator= (const GridFunction &rhs)
 Copy assignment. Only the data of the base class Vector is copied.
 
void MakeOwner (FiniteElementCollection *fec_)
 Make the GridFunction the owner of fec and fes.
 
FiniteElementCollectionOwnFEC ()
 
int VectorDim () const
 
int CurlDim () const
 
const VectorGetTrueVector () const
 Read only access to the (optional) internal true-dof Vector.
 
VectorGetTrueVector ()
 Read and write access to the (optional) internal true-dof Vector.
 
void GetTrueDofs (Vector &tv) const
 Extract the true-dofs from the GridFunction.
 
void SetTrueVector ()
 Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
 
virtual void SetFromTrueDofs (const Vector &tv)
 Set the GridFunction from the given true-dof vector.
 
void SetFromTrueVector ()
 Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
 
void GetNodalValues (int i, Array< real_t > &nval, int vdim=1) const
 Returns the values in the vertices of i'th element for dimension vdim.
 
void GetLaplacians (int i, const IntegrationRule &ir, Vector &laps, int vdim=1) const
 
void GetLaplacians (int i, const IntegrationRule &ir, Vector &laps, DenseMatrix &tr, int vdim=1) const
 
void GetHessians (int i, const IntegrationRule &ir, DenseMatrix &hess, int vdim=1) const
 
void GetHessians (int i, const IntegrationRule &ir, DenseMatrix &hess, DenseMatrix &tr, int vdim=1) const
 
void GetValuesFrom (const GridFunction &orig_func)
 
void GetBdrValuesFrom (const GridFunction &orig_func)
 
void GetVectorFieldValues (int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
 
void ReorderByNodes ()
 For a vector grid function, makes sure that the ordering is byNODES.
 
void GetNodalValues (Vector &nval, int vdim=1) const
 Return the values as a vector on mesh vertices for dimension vdim.
 
void GetVectorFieldNodalValues (Vector &val, int comp) const
 
void ProjectVectorFieldOn (GridFunction &vec_field, int comp=0)
 
void GetDerivative (int comp, int der_comp, GridFunction &der)
 Compute a certain derivative of a function's component. Derivatives of the function are computed at the DOF locations of der, and averaged over overlapping DOFs. Thus this function projects the derivative to the FiniteElementSpace of der.
 
real_t GetDivergence (ElementTransformation &tr) const
 
void GetCurl (ElementTransformation &tr, Vector &curl) const
 
void GetGradient (ElementTransformation &tr, Vector &grad) const
 Gradient of a scalar function at a quadrature point.
 
void GetGradients (ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
 Extension of GetGradient(...) for a collection of IntegrationPoints.
 
void GetGradients (const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
 Extension of GetGradient(...) for a collection of IntegrationPoints.
 
void GetVectorGradient (ElementTransformation &tr, DenseMatrix &grad) const
 Compute the vector gradient with respect to the physical element variable.
 
void GetVectorGradientHat (ElementTransformation &T, DenseMatrix &gh) const
 Compute the vector gradient with respect to the reference element variable.
 
void GetElementAverages (GridFunction &avgs) const
 
virtual void GetElementDofValues (int el, Vector &dof_vals) const
 
void ImposeBounds (int i, const Vector &weights, const Vector &lo_, const Vector &hi_)
 
void ImposeBounds (int i, const Vector &weights, real_t min_=0.0, real_t max_=infinity())
 
void RestrictConforming ()
 
void ProjectGridFunction (const GridFunction &src)
 Project the src GridFunction to this GridFunction, both of which must be on the same mesh.
 
virtual void ProjectCoefficient (Coefficient &coeff)
 Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of the FiniteElementSpace fes. Note that this is usually interpolation at the degrees of freedom in each element (not L2 projection).
 
void ProjectCoefficient (Coefficient &coeff, Array< int > &dofs, int vd=0)
 Project coeff Coefficient to this GridFunction, using one element for each degree of freedom in dofs and nodal interpolation on that element.
 
void ProjectCoefficient (VectorCoefficient &vcoeff)
 Project vcoeff VectorCoefficient to this GridFunction. The projection computation depends on the choice of the FiniteElementSpace fes. Note that this is usually interpolation at the degrees of freedom in each element (not L2 projection).
 
void ProjectCoefficient (VectorCoefficient &vcoeff, Array< int > &dofs)
 Project vcoeff VectorCoefficient to this GridFunction, using one element for each degree of freedom in dofs and nodal interpolation on that element.
 
void ProjectCoefficient (VectorCoefficient &vcoeff, int attribute)
 Project vcoeff VectorCoefficient to this GridFunction, only projecting onto elements with the given attribute.
 
void ProjectCoefficient (Coefficient *coeff[])
 Analogous to the version with argument vcoeff VectorCoefficient but using an array of scalar coefficients for each component.
 
virtual void ProjectDiscCoefficient (VectorCoefficient &coeff)
 Project a discontinuous vector coefficient as a grid function on a continuous finite element space. The values in shared dofs are determined from the element with maximal attribute.
 
virtual void ProjectDiscCoefficient (Coefficient &coeff, AvgType type)
 Projects a discontinuous coefficient so that the values in shared vdofs are computed by taking an average of the possible values.
 
virtual void ProjectDiscCoefficient (VectorCoefficient &coeff, AvgType type)
 Projects a discontinuous vector coefficient so that the values in shared vdofs are computed by taking an average of the possible values.
 
virtual void CountElementsPerVDof (Array< int > &elem_per_vdof) const
 For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteElementSpace::GetElementVDofs().
 
void ProjectBdrCoefficient (Coefficient &coeff, const Array< int > &attr)
 Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the boundary attributes marked in the attr array.
 
virtual void ProjectBdrCoefficient (VectorCoefficient &vcoeff, const Array< int > &attr)
 Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with the boundary attributes marked in the attr array.
 
virtual void ProjectBdrCoefficient (Coefficient *coeff[], const Array< int > &attr)
 Project a set of Coefficients on the components of the GridFunction, modifying only DOFs on the boundary associated with the boundary attributed marked in the attr array.
 
void ProjectBdrCoefficientNormal (VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
 
virtual void ProjectBdrCoefficientTangent (VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
 Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attributes that are marked in bdr_attr are projected. Assumes ND-type VectorFE GridFunction.
 
virtual real_t ComputeL2Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
 
virtual real_t ComputeElementGradError (int ielem, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
 Returns ||grad u_ex - grad u_h||_L2 in element ielem for H1 or L2 elements.
 
virtual real_t ComputeL2Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
 Returns ||u_ex - u_h||_L2 for H1 or L2 elements.
 
virtual real_t ComputeL2Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
 
virtual real_t ComputeGradError (VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
 Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
 
virtual real_t ComputeCurlError (VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
 Returns ||curl u_ex - curl u_h||_L2 for ND elements.
 
virtual real_t ComputeDivError (Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
 Returns ||div u_ex - div u_h||_L2 for RT elements.
 
virtual real_t ComputeDGFaceJumpError (Coefficient *exsol, Coefficient *ell_coeff, class JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
 
MFEM_DEPRECATED real_t ComputeDGFaceJumpError (Coefficient *exsol, Coefficient *ell_coeff, real_t Nu, const IntegrationRule *irs[]=NULL) const
 Returns the Face Jumps error for L2 elements, with 1/h scaling.
 
virtual real_t ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t Nu, int norm_type) const
 
virtual real_t ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeHDivError (VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
 Returns the error measured in H(div)-norm for RT elements.
 
virtual real_t ComputeHCurlError (VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
 Returns the error measured in H(curl)-norm for ND elements.
 
virtual real_t ComputeMaxError (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeMaxError (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeMaxError (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeL1Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeL1Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeW11Error (Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, const Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeL1Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeLpError (const real_t p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
 
virtual void ComputeElementLpErrors (const real_t p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementL1Errors (Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementL2Errors (Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementMaxErrors (Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 
virtual real_t ComputeLpError (const real_t p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementLpErrors (const real_t p, VectorCoefficient &exsol, Vector &error, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementL1Errors (VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementL2Errors (VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementMaxErrors (VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeFlux (BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
 
GridFunctionoperator= (real_t value)
 Redefine '=' for GridFunction = constant.
 
GridFunctionoperator= (const Vector &v)
 Copy the data from v.
 
virtual void Update ()
 Transform by the Space UpdateMatrix (e.g., on Mesh change).
 
long GetSequence () const
 
FiniteElementSpaceFESpace ()
 
const FiniteElementSpaceFESpace () const
 
virtual void SetSpace (FiniteElementSpace *f)
 Associate a new FiniteElementSpace with the GridFunction.
 
virtual void MakeRef (FiniteElementSpace *f, real_t *v)
 Make the GridFunction reference external data on a new FiniteElementSpace.
 
virtual void MakeRef (FiniteElementSpace *f, Vector &v, int v_offset)
 Make the GridFunction reference external data on a new FiniteElementSpace.
 
void MakeTRef (FiniteElementSpace *f, real_t *tv)
 Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
 
void MakeTRef (FiniteElementSpace *f, Vector &tv, int tv_offset)
 Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
 
virtual void Save (std::ostream &out) const
 Save the GridFunction to an output stream.
 
virtual void Save (const char *fname, int precision=16) const
 
virtual void Save (adios2stream &out, const std::string &variable_name, const adios2stream::data_type type=adios2stream::data_type::point_data) const
 Save the GridFunction to a binary output stream using adios2 bp format.
 
void SaveVTK (std::ostream &out, const std::string &field_name, int ref)
 Write the GridFunction in VTK format. Note that Mesh::PrintVTK must be called first. The parameter ref > 0 must match the one used in Mesh::PrintVTK.
 
void SaveSTL (std::ostream &out, int TimesToRefine=1)
 Write the GridFunction in STL format. Note that the mesh dimension must be 2 and that quad elements will be broken into two triangles.
 
virtual ~GridFunction ()
 Destroys grid function.
 
void MakeRef (Vector &base, int offset, int size)
 Reset the Vector to be a reference to a sub-vector of base.
 
void MakeRef (Vector &base, int offset)
 Reset the Vector to be a reference to a sub-vector of base without changing its current size.
 
Element index Get Value Methods

These methods take an element index and return the interpolated value of the field at a given reference point within the element.

Warning
These methods retrieve and use the ElementTransformation object from the mfem::Mesh. This can alter the state of the element transformation object and can also lead to unexpected results when the ElementTransformation object is already in use such as when these methods are called from within an integration loop. Consider using GetValue(ElementTransformation &T, ...) instead.
virtual real_t GetValue (int i, const IntegrationPoint &ip, int vdim=1) const
 
virtual void GetVectorValue (int i, const IntegrationPoint &ip, Vector &val) const
 
Element Index Get Values Methods

These are convenience methods for repeatedly calling GetValue for multiple points within a given element. The GetValues methods are optimized and should perform better than repeatedly calling GetValue. The GetVectorValues method simply calls GetVectorValue repeatedly.

Warning
These methods retrieve and use the ElementTransformation object from the mfem::Mesh. This can alter the state of the element transformation object and can also lead to unexpected results when the ElementTransformation object is already in use such as when these methods are called from within an integration loop. Consider using GetValues(ElementTransformation &T, ...) instead.
void GetValues (int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
 
void GetValues (int i, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
 
void GetVectorValues (int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
 
ElementTransformation Get Value Methods

These member functions are designed for use within GridFunctionCoefficient objects. These can be used with ElementTransformation objects coming from either Mesh::GetElementTransformation() or Mesh::GetBdrElementTransformation().

Note
These methods do not reset the ElementTransformation object so they should be safe to use within integration loops or other contexts where the ElementTransformation is already in use.
virtual real_t GetValue (ElementTransformation &T, const IntegrationPoint &ip, int comp=0, Vector *tr=NULL) const
 
virtual void GetVectorValue (ElementTransformation &T, const IntegrationPoint &ip, Vector &val, Vector *tr=NULL) const
 
ElementTransformation Get Values Methods

These are convenience methods for repeatedly calling GetValue for multiple points within a given element. They work by calling either the ElementTransformation or FaceElementTransformations versions described above. Consequently, these methods should not be expected to run faster than calling the above methods in an external loop.

Note
These methods do not reset the ElementTransformation object so they should be safe to use within integration loops or other contexts where the ElementTransformation is already in use.
These methods can also be used with FaceElementTransformations objects.
void GetValues (ElementTransformation &T, const IntegrationRule &ir, Vector &vals, int comp=0, DenseMatrix *tr=NULL) const
 
void GetVectorValues (ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix *tr=NULL) const
 
Face Index Get Values Methods

These methods are designed to work with Discontinuous Galerkin basis functions. They compute field values on the interface between elements, or on boundary elements, by interpolating the field in a neighboring element. The side argument indices which neighboring element should be used: 0, 1, or 2 (automatically chosen).

Warning
These methods retrieve and use the FaceElementTransformations object from the mfem::Mesh. This can alter the state of the face element transformations object and can also lead to unexpected results when the FaceElementTransformations object is already in use such as when these methods are called from within an integration loop. Consider using GetValues(ElementTransformation &T, ...) instead.
int GetFaceValues (int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
 
int GetFaceVectorValues (int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
 
- Public Member Functions inherited from mfem::Vector
 Vector ()
 
 Vector (const Vector &)
 Copy constructor. Allocates a new data array and copies the data.
 
 Vector (Vector &&v)
 Move constructor. "Steals" data from its argument.
 
 Vector (int s)
 Creates vector of size s.
 
 Vector (real_t *data_, int size_)
 Creates a vector referencing an array of doubles, owned by someone else.
 
 Vector (Vector &base, int base_offset, int size_)
 Create a Vector referencing a sub-vector of the Vector base starting at the given offset, base_offset, and size size_.
 
 Vector (int size_, MemoryType mt)
 Create a Vector of size size_ using MemoryType mt.
 
 Vector (int size_, MemoryType h_mt, MemoryType d_mt)
 Create a Vector of size size_ using host MemoryType h_mt and device MemoryType d_mt.
 
template<int N, typename T = real_t>
 Vector (const T(&values)[N])
 Create a vector using a braced initializer list.
 
virtual void UseDevice (bool use_dev) const
 Enable execution of Vector operations using the mfem::Device.
 
virtual bool UseDevice () const
 Return the device flag of the Memory object used by the Vector.
 
void Load (std::istream **in, int np, int *dim)
 Reads a vector from multiple files.
 
void Load (std::istream &in, int Size)
 Load a vector from an input stream.
 
void Load (std::istream &in)
 Load a vector from an input stream, reading the size from the stream.
 
void SetSize (int s)
 Resize the vector to size s.
 
void SetSize (int s, MemoryType mt)
 Resize the vector to size s using MemoryType mt.
 
void SetSize (int s, const Vector &v)
 Resize the vector to size s using the MemoryType of v.
 
void SetData (real_t *d)
 
void SetDataAndSize (real_t *d, int s)
 Set the Vector data and size.
 
void NewDataAndSize (real_t *d, int s)
 Set the Vector data and size, deleting the old data, if owned.
 
void NewMemoryAndSize (const Memory< real_t > &mem, int s, bool own_mem)
 Reset the Vector to use the given external Memory mem and size s.
 
void MakeRef (Vector &base, int offset, int size)
 Reset the Vector to be a reference to a sub-vector of base.
 
void MakeRef (Vector &base, int offset)
 Reset the Vector to be a reference to a sub-vector of base without changing its current size.
 
void MakeDataOwner () const
 Set the Vector data (host pointer) ownership flag.
 
void Destroy ()
 Destroy a vector.
 
void DeleteDevice (bool copy_to_host=true)
 Delete the device pointer, if owned. If copy_to_host is true and the data is valid only on device, move it to host before deleting. Invalidates the device memory.
 
int Size () const
 Returns the size of the vector.
 
int Capacity () const
 Return the size of the currently allocated data array.
 
real_tGetData () const
 Return a pointer to the beginning of the Vector data.
 
MFEM_DEPRECATED operator real_t * ()
 Conversion to double *. Deprecated.
 
MFEM_DEPRECATED operator const real_t * () const
 Conversion to const double *. Deprecated.
 
real_tbegin ()
 STL-like begin.
 
real_tend ()
 STL-like end.
 
const real_tbegin () const
 STL-like begin (const version).
 
const real_tend () const
 STL-like end (const version).
 
Memory< real_t > & GetMemory ()
 Return a reference to the Memory object used by the Vector.
 
const Memory< real_t > & GetMemory () const
 Return a reference to the Memory object used by the Vector, const version.
 
void SyncMemory (const Vector &v) const
 Update the memory location of the vector to match v.
 
void SyncAliasMemory (const Vector &v) const
 Update the alias memory location of the vector to match v.
 
bool OwnsData () const
 Read the Vector data (host pointer) ownership flag.
 
void StealData (real_t **p)
 Changes the ownership of the data; after the call the Vector is empty.
 
real_tStealData ()
 Changes the ownership of the data; after the call the Vector is empty.
 
real_tElem (int i)
 Access Vector entries. Index i = 0 .. size-1.
 
const real_tElem (int i) const
 Read only access to Vector entries. Index i = 0 .. size-1.
 
real_toperator() (int i)
 Access Vector entries using () for 0-based indexing.
 
const real_toperator() (int i) const
 Read only access to Vector entries using () for 0-based indexing.
 
real_toperator[] (int i)
 Access Vector entries using [] for 0-based indexing.
 
const real_toperator[] (int i) const
 Read only access to Vector entries using [] for 0-based indexing.
 
real_t operator* (const real_t *) const
 Dot product with a double * array.
 
real_t operator* (const Vector &v) const
 Return the inner-product.
 
Vectoroperator= (const real_t *v)
 Copy Size() entries from v.
 
Vectoroperator= (const Vector &v)
 Copy assignment.
 
Vectoroperator= (Vector &&v)
 Move assignment.
 
Vectoroperator= (real_t value)
 Redefine '=' for vector = constant.
 
Vectoroperator*= (real_t c)
 
Vectoroperator*= (const Vector &v)
 Component-wise scaling: (*this)(i) *= v(i)
 
Vectoroperator/= (real_t c)
 
Vectoroperator/= (const Vector &v)
 Component-wise division: (*this)(i) /= v(i)
 
Vectoroperator-= (real_t c)
 
Vectoroperator-= (const Vector &v)
 
Vectoroperator+= (real_t c)
 
Vectoroperator+= (const Vector &v)
 
VectorAdd (const real_t a, const Vector &Va)
 (*this) += a * Va
 
VectorSet (const real_t a, const Vector &x)
 (*this) = a * x
 
void SetVector (const Vector &v, int offset)
 
void AddSubVector (const Vector &v, int offset)
 
void Neg ()
 (*this) = -(*this)
 
void Reciprocal ()
 (*this)(i) = 1.0 / (*this)(i)
 
void Swap (Vector &other)
 Swap the contents of two Vectors.
 
void cross3D (const Vector &vin, Vector &vout) const
 
void median (const Vector &lo, const Vector &hi)
 v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi.
 
void GetSubVector (const Array< int > &dofs, Vector &elemvect) const
 Extract entries listed in dofs to the output Vector elemvect.
 
void GetSubVector (const Array< int > &dofs, real_t *elem_data) const
 Extract entries listed in dofs to the output array elem_data.
 
void SetSubVector (const Array< int > &dofs, const real_t value)
 Set the entries listed in dofs to the given value.
 
void SetSubVector (const Array< int > &dofs, const Vector &elemvect)
 Set the entries listed in dofs to the values given in the elemvect Vector. Negative dof values cause the -dof-1 position in this Vector to receive the -val from elemvect.
 
void SetSubVector (const Array< int > &dofs, real_t *elem_data)
 Set the entries listed in dofs to the values given the , elem_data array. Negative dof values cause the -dof-1 position in this Vector to receive the -val from elem_data.
 
void AddElementVector (const Array< int > &dofs, const Vector &elemvect)
 Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof-1 position in this Vector to add the -val from elemvect.
 
void AddElementVector (const Array< int > &dofs, real_t *elem_data)
 Add elements of the elem_data array to the entries listed in dofs. Negative dof values cause the -dof-1 position in this Vector to add the -val from elem_data.
 
void AddElementVector (const Array< int > &dofs, const real_t a, const Vector &elemvect)
 Add times the elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof-1 position in this Vector to add the -a*val from elemvect.
 
void SetSubVectorComplement (const Array< int > &dofs, const real_t val)
 Set all vector entries NOT in the dofs Array to the given val.
 
void Print (std::ostream &out=mfem::out, int width=8) const
 Prints vector to stream out.
 
void Print (adios2stream &out, const std::string &variable_name) const
 
void Print_HYPRE (std::ostream &out) const
 Prints vector to stream out in HYPRE_Vector format.
 
void PrintHash (std::ostream &out) const
 Print the Vector size and hash of its data.
 
void Randomize (int seed=0)
 Set random values in the vector.
 
real_t Norml2 () const
 Returns the l2 norm of the vector.
 
real_t Normlinf () const
 Returns the l_infinity norm of the vector.
 
real_t Norml1 () const
 Returns the l_1 norm of the vector.
 
real_t Normlp (real_t p) const
 Returns the l_p norm of the vector.
 
real_t Max () const
 Returns the maximal element of the vector.
 
real_t Min () const
 Returns the minimal element of the vector.
 
real_t Sum () const
 Return the sum of the vector entries.
 
real_t DistanceSquaredTo (const real_t *p) const
 Compute the square of the Euclidean distance to another vector.
 
real_t DistanceSquaredTo (const Vector &p) const
 Compute the square of the Euclidean distance to another vector.
 
real_t DistanceTo (const real_t *p) const
 Compute the Euclidean distance to another vector.
 
real_t DistanceTo (const Vector &p) const
 Compute the Euclidean distance to another vector.
 
int CheckFinite () const
 Count the number of entries in the Vector for which isfinite is false, i.e. the entry is a NaN or +/-Inf.
 
virtual ~Vector ()
 Destroys vector.
 
virtual const real_tRead (bool on_dev=true) const
 Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
 
virtual const real_tHostRead () const
 Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
 
virtual real_tWrite (bool on_dev=true)
 Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
 
virtual real_tHostWrite ()
 Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
 
virtual real_tReadWrite (bool on_dev=true)
 Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
 
virtual real_tHostReadWrite ()
 Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
 

Protected Member Functions

void SaveSTLTri (std::ostream &out, real_t p1[], real_t p2[], real_t p3[])
 
void ProjectDeltaCoefficient (DeltaCoefficient &delta_coeff, real_t &integral)
 
void SumFluxAndCount (BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
 
void ProjectDiscCoefficient (VectorCoefficient &coeff, Array< int > &dof_attr)
 
void LegacyNCReorder ()
 Loading helper.
 
void Destroy ()
 
void AccumulateAndCountZones (Coefficient &coeff, AvgType type, Array< int > &zones_per_vdof)
 Accumulates (depending on type) the values of coeff at all shared vdofs and counts in how many zones each vdof appears.
 
void AccumulateAndCountZones (VectorCoefficient &vcoeff, AvgType type, Array< int > &zones_per_vdof)
 Accumulates (depending on type) the values of vcoeff at all shared vdofs and counts in how many zones each vdof appears.
 
void AccumulateAndCountDerivativeValues (int comp, int der_comp, GridFunction &der, Array< int > &zones_per_dof)
 Used for the serial and parallel implementations of the GetDerivative() method; see its documentation.
 
void AccumulateAndCountBdrValues (Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr, Array< int > &values_counter)
 
void AccumulateAndCountBdrTangentValues (VectorCoefficient &vcoeff, const Array< int > &bdr_attr, Array< int > &values_counter)
 
void ComputeMeans (AvgType type, Array< int > &zones_per_vdof)
 

Protected Attributes

FiniteElementSpacefes
 FE space on which the grid function lives. Owned if fec is not NULL.
 
FiniteElementCollectionfec
 Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner().
 
long fes_sequence
 
Vector t_vec
 
- Protected Attributes inherited from mfem::Vector
Memory< real_tdata
 
int size
 

Detailed Description

Class for grid function - Vector with associated FE space.

Definition at line 30 of file gridfunc.hpp.

Member Enumeration Documentation

◆ AvgType

Enumerator
ARITHMETIC 
HARMONIC 

Definition at line 422 of file gridfunc.hpp.

Constructor & Destructor Documentation

◆ GridFunction() [1/7]

mfem::GridFunction::GridFunction ( )
inline

Definition at line 75 of file gridfunc.hpp.

◆ GridFunction() [2/7]

mfem::GridFunction::GridFunction ( const GridFunction & orig)
inline

Copy constructor. The internal true-dof vector t_vec is not copied.

Definition at line 78 of file gridfunc.hpp.

◆ GridFunction() [3/7]

mfem::GridFunction::GridFunction ( FiniteElementSpace * f)
inline

Construct a GridFunction associated with the FiniteElementSpace *f.

Definition at line 83 of file gridfunc.hpp.

◆ GridFunction() [4/7]

mfem::GridFunction::GridFunction ( FiniteElementSpace * f,
real_t * data )
inline

Construct a GridFunction using previously allocated array data.

The GridFunction does not assume ownership of data which is assumed to be of size at least f->GetVSize(). Similar to the Vector constructor for externally allocated array, the pointer data can be NULL. The data array can be replaced later using the method SetData().

Definition at line 92 of file gridfunc.hpp.

◆ GridFunction() [5/7]

mfem::GridFunction::GridFunction ( FiniteElementSpace * f,
Vector & base,
int base_offset = 0 )
inline

Construct a GridFunction using previously allocated Vector base starting at the given offset, base_offset.

Definition at line 98 of file gridfunc.hpp.

◆ GridFunction() [6/7]

mfem::GridFunction::GridFunction ( Mesh * m,
std::istream & input )

Construct a GridFunction on the given Mesh, using the data from input.

The content of input should be in the format created by the method Save(). The reconstructed FiniteElementSpace and FiniteElementCollection are owned by the GridFunction.

Definition at line 35 of file gridfunc.cpp.

◆ GridFunction() [7/7]

mfem::GridFunction::GridFunction ( Mesh * m,
GridFunction * gf_array[],
int num_pieces )

Definition at line 76 of file gridfunc.cpp.

◆ ~GridFunction()

virtual mfem::GridFunction::~GridFunction ( )
inlinevirtual

Destroys grid function.

Definition at line 764 of file gridfunc.hpp.

Member Function Documentation

◆ AccumulateAndCountBdrTangentValues()

void mfem::GridFunction::AccumulateAndCountBdrTangentValues ( VectorCoefficient & vcoeff,
const Array< int > & bdr_attr,
Array< int > & values_counter )
protected

Definition at line 2202 of file gridfunc.cpp.

◆ AccumulateAndCountBdrValues()

void mfem::GridFunction::AccumulateAndCountBdrValues ( Coefficient * coeff[],
VectorCoefficient * vcoeff,
const Array< int > & attr,
Array< int > & values_counter )
protected

Definition at line 2060 of file gridfunc.cpp.

◆ AccumulateAndCountDerivativeValues()

void mfem::GridFunction::AccumulateAndCountDerivativeValues ( int comp,
int der_comp,
GridFunction & der,
Array< int > & zones_per_dof )
protected

Used for the serial and parallel implementations of the GetDerivative() method; see its documentation.

Definition at line 1324 of file gridfunc.cpp.

◆ AccumulateAndCountZones() [1/2]

void mfem::GridFunction::AccumulateAndCountZones ( Coefficient & coeff,
AvgType type,
Array< int > & zones_per_vdof )
protected

Accumulates (depending on type) the values of coeff at all shared vdofs and counts in how many zones each vdof appears.

Definition at line 1964 of file gridfunc.cpp.

◆ AccumulateAndCountZones() [2/2]

void mfem::GridFunction::AccumulateAndCountZones ( VectorCoefficient & vcoeff,
AvgType type,
Array< int > & zones_per_vdof )
protected

Accumulates (depending on type) the values of vcoeff at all shared vdofs and counts in how many zones each vdof appears.

Definition at line 2005 of file gridfunc.cpp.

◆ ComputeCurlError()

real_t mfem::GridFunction::ComputeCurlError ( VectorCoefficient * excurl,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns ||curl u_ex - curl u_h||_L2 for ND elements.

Reimplemented in mfem::ParGridFunction.

Definition at line 2892 of file gridfunc.cpp.

◆ ComputeDGFaceJumpError() [1/2]

real_t mfem::GridFunction::ComputeDGFaceJumpError ( Coefficient * exsol,
Coefficient * ell_coeff,
class JumpScaling jump_scaling,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns the Face Jumps error for L2 elements. The error can be weighted by a constant nu, by nu/h, or nu*p^2/h, depending on the value of jump_scaling.

Reimplemented in mfem::ParGridFunction.

Definition at line 2969 of file gridfunc.cpp.

◆ ComputeDGFaceJumpError() [2/2]

real_t mfem::GridFunction::ComputeDGFaceJumpError ( Coefficient * exsol,
Coefficient * ell_coeff,
real_t Nu,
const IntegrationRule * irs[] = NULL ) const

Returns the Face Jumps error for L2 elements, with 1/h scaling.

Definition at line 3083 of file gridfunc.cpp.

◆ ComputeDivError()

real_t mfem::GridFunction::ComputeDivError ( Coefficient * exdiv,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns ||div u_ex - div u_h||_L2 for RT elements.

Reimplemented in mfem::ParGridFunction.

Definition at line 2933 of file gridfunc.cpp.

◆ ComputeElementGradError()

real_t mfem::GridFunction::ComputeElementGradError ( int ielem,
VectorCoefficient * exgrad,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns ||grad u_ex - grad u_h||_L2 in element ielem for H1 or L2 elements.

Definition at line 2814 of file gridfunc.cpp.

◆ ComputeElementL1Errors() [1/2]

virtual void mfem::GridFunction::ComputeElementL1Errors ( Coefficient & exsol,
Vector & error,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Definition at line 623 of file gridfunc.hpp.

◆ ComputeElementL1Errors() [2/2]

virtual void mfem::GridFunction::ComputeElementL1Errors ( VectorCoefficient & exsol,
Vector & error,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Definition at line 659 of file gridfunc.hpp.

◆ ComputeElementL2Errors() [1/2]

virtual void mfem::GridFunction::ComputeElementL2Errors ( Coefficient & exsol,
Vector & error,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Definition at line 629 of file gridfunc.hpp.

◆ ComputeElementL2Errors() [2/2]

virtual void mfem::GridFunction::ComputeElementL2Errors ( VectorCoefficient & exsol,
Vector & error,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Definition at line 665 of file gridfunc.hpp.

◆ ComputeElementLpErrors() [1/2]

void mfem::GridFunction::ComputeElementLpErrors ( const real_t p,
Coefficient & exsol,
Vector & error,
Coefficient * weight = NULL,
const IntegrationRule * irs[] = NULL ) const
virtual

Compute the Lp error in each element of the mesh and store the results in the Vector error. The result should be of length number of elements, for example an L2 GridFunction of order zero using map type VALUE.

Definition at line 3365 of file gridfunc.cpp.

◆ ComputeElementLpErrors() [2/2]

void mfem::GridFunction::ComputeElementLpErrors ( const real_t p,
VectorCoefficient & exsol,
Vector & error,
Coefficient * weight = NULL,
VectorCoefficient * v_weight = NULL,
const IntegrationRule * irs[] = NULL ) const
virtual

Compute the Lp error in each element of the mesh and store the results in the Vector @ error. The result should be of length number of elements, for example an L2 GridFunction of order zero using map type VALUE.

Definition at line 3522 of file gridfunc.cpp.

◆ ComputeElementMaxErrors() [1/2]

virtual void mfem::GridFunction::ComputeElementMaxErrors ( Coefficient & exsol,
Vector & error,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Definition at line 635 of file gridfunc.hpp.

◆ ComputeElementMaxErrors() [2/2]

virtual void mfem::GridFunction::ComputeElementMaxErrors ( VectorCoefficient & exsol,
Vector & error,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Definition at line 671 of file gridfunc.hpp.

◆ ComputeFlux()

void mfem::GridFunction::ComputeFlux ( BilinearFormIntegrator & blfi,
GridFunction & flux,
bool wcoef = true,
int subdomain = -1 )
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 308 of file gridfunc.cpp.

◆ ComputeGradError()

real_t mfem::GridFunction::ComputeGradError ( VectorCoefficient * exgrad,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.

Reimplemented in mfem::ParGridFunction.

Definition at line 2852 of file gridfunc.cpp.

◆ ComputeH1Error() [1/2]

real_t mfem::GridFunction::ComputeH1Error ( Coefficient * exsol,
VectorCoefficient * exgrad,
Coefficient * ell_coef,
real_t Nu,
int norm_type ) const
virtual

This method is kept for backward compatibility.

Returns either the H1-seminorm, or the DG face jumps error, or both depending on norm_type = 1, 2, 3. Additional arguments for the DG face jumps norm: ell_coeff: mesh-depended coefficient (weight) Nu: scalar constant weight

Reimplemented in mfem::ParGridFunction.

Definition at line 3092 of file gridfunc.cpp.

◆ ComputeH1Error() [2/2]

real_t mfem::GridFunction::ComputeH1Error ( Coefficient * exsol,
VectorCoefficient * exgrad,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns the error measured in H1-norm for H1 elements or in "broken" H1-norm for L2 elements

Reimplemented in mfem::ParGridFunction.

Definition at line 3109 of file gridfunc.cpp.

◆ ComputeHCurlError()

real_t mfem::GridFunction::ComputeHCurlError ( VectorCoefficient * exsol,
VectorCoefficient * excurl,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns the error measured in H(curl)-norm for ND elements.

Reimplemented in mfem::ParGridFunction.

Definition at line 3127 of file gridfunc.cpp.

◆ ComputeHDivError()

real_t mfem::GridFunction::ComputeHDivError ( VectorCoefficient * exsol,
Coefficient * exdiv,
const IntegrationRule * irs[] = NULL ) const
virtual

Returns the error measured in H(div)-norm for RT elements.

Reimplemented in mfem::ParGridFunction.

Definition at line 3118 of file gridfunc.cpp.

◆ ComputeL1Error() [1/3]

virtual real_t mfem::GridFunction::ComputeL1Error ( Coefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Reimplemented in mfem::ParGridFunction.

Definition at line 594 of file gridfunc.hpp.

◆ ComputeL1Error() [2/3]

virtual real_t mfem::GridFunction::ComputeL1Error ( Coefficient * exsol[],
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Reimplemented in mfem::ParGridFunction.

Definition at line 590 of file gridfunc.hpp.

◆ ComputeL1Error() [3/3]

virtual real_t mfem::GridFunction::ComputeL1Error ( VectorCoefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Reimplemented in mfem::ParGridFunction.

Definition at line 602 of file gridfunc.hpp.

◆ ComputeL2Error() [1/3]

virtual real_t mfem::GridFunction::ComputeL2Error ( Coefficient & exsol,
const IntegrationRule * irs[] = NULL,
const Array< int > * elems = NULL ) const
inlinevirtual

Returns ||u_ex - u_h||_L2 for H1 or L2 elements.

Reimplemented in mfem::ParGridFunction.

Definition at line 513 of file gridfunc.hpp.

◆ ComputeL2Error() [2/3]

real_t mfem::GridFunction::ComputeL2Error ( Coefficient * exsol[],
const IntegrationRule * irs[] = NULL,
const Array< int > * elems = NULL ) const
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 2718 of file gridfunc.cpp.

◆ ComputeL2Error() [3/3]

real_t mfem::GridFunction::ComputeL2Error ( VectorCoefficient & exsol,
const IntegrationRule * irs[] = NULL,
const Array< int > * elems = NULL ) const
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 2773 of file gridfunc.cpp.

◆ ComputeLpError() [1/2]

real_t mfem::GridFunction::ComputeLpError ( const real_t p,
Coefficient & exsol,
Coefficient * weight = NULL,
const IntegrationRule * irs[] = NULL,
const Array< int > * elems = NULL ) const
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 3298 of file gridfunc.cpp.

◆ ComputeLpError() [2/2]

real_t mfem::GridFunction::ComputeLpError ( const real_t p,
VectorCoefficient & exsol,
Coefficient * weight = NULL,
VectorCoefficient * v_weight = NULL,
const IntegrationRule * irs[] = NULL ) const
virtual

When given a vector weight, compute the pointwise (scalar) error as the dot product of the vector error with the vector weight. Otherwise, the scalar error is the l_2 norm of the vector error.

Reimplemented in mfem::ParGridFunction.

Definition at line 3431 of file gridfunc.cpp.

◆ ComputeMaxError() [1/3]

virtual real_t mfem::GridFunction::ComputeMaxError ( Coefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Reimplemented in mfem::ParGridFunction.

Definition at line 575 of file gridfunc.hpp.

◆ ComputeMaxError() [2/3]

real_t mfem::GridFunction::ComputeMaxError ( Coefficient * exsol[],
const IntegrationRule * irs[] = NULL ) const
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 3136 of file gridfunc.cpp.

◆ ComputeMaxError() [3/3]

virtual real_t mfem::GridFunction::ComputeMaxError ( VectorCoefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlinevirtual

Reimplemented in mfem::ParGridFunction.

Definition at line 584 of file gridfunc.hpp.

◆ ComputeMeans()

void mfem::GridFunction::ComputeMeans ( AvgType type,
Array< int > & zones_per_vdof )
protected

Definition at line 2254 of file gridfunc.cpp.

◆ ComputeW11Error()

real_t mfem::GridFunction::ComputeW11Error ( Coefficient * exsol,
VectorCoefficient * exgrad,
int norm_type,
const Array< int > * elems = NULL,
const IntegrationRule * irs[] = NULL ) const
virtual

Definition at line 3192 of file gridfunc.cpp.

◆ CountElementsPerVDof()

void mfem::GridFunction::CountElementsPerVDof ( Array< int > & elem_per_vdof) const
virtual

For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteElementSpace::GetElementVDofs().

Reimplemented in mfem::ParGridFunction.

Definition at line 1947 of file gridfunc.cpp.

◆ CurlDim()

int mfem::GridFunction::CurlDim ( ) const

Definition at line 346 of file gridfunc.cpp.

◆ Destroy()

void mfem::GridFunction::Destroy ( )
protected

Definition at line 154 of file gridfunc.cpp.

◆ FESpace() [1/2]

FiniteElementSpace * mfem::GridFunction::FESpace ( )
inline

Definition at line 696 of file gridfunc.hpp.

◆ FESpace() [2/2]

const FiniteElementSpace * mfem::GridFunction::FESpace ( ) const
inline

Definition at line 697 of file gridfunc.hpp.

◆ GetBdrValuesFrom()

void mfem::GridFunction::GetBdrValuesFrom ( const GridFunction & orig_func)

Definition at line 1150 of file gridfunc.cpp.

◆ GetCurl()

void mfem::GridFunction::GetCurl ( ElementTransformation & tr,
Vector & curl ) const

Definition at line 1490 of file gridfunc.cpp.

◆ GetDerivative()

void mfem::GridFunction::GetDerivative ( int comp,
int der_comp,
GridFunction & der )

Compute a certain derivative of a function's component. Derivatives of the function are computed at the DOF locations of der, and averaged over overlapping DOFs. Thus this function projects the derivative to the FiniteElementSpace of der.

Parameters
[in]compIndex of the function's component to be differentiated. The index is 1-based, i.e., use 1 for scalar functions.
[in]der_compUse 0/1/2 for derivatives in x/y/z directions.
[out]derThe resulting derivative (scalar function). The FiniteElementSpace of this function must be set before the call.

Definition at line 1377 of file gridfunc.cpp.

◆ GetDivergence()

real_t mfem::GridFunction::GetDivergence ( ElementTransformation & tr) const

Definition at line 1404 of file gridfunc.cpp.

◆ GetElementAverages()

void mfem::GridFunction::GetElementAverages ( GridFunction & avgs) const

Compute \( (\int_{\Omega} (*this) \psi_i)/(\int_{\Omega} \psi_i) \), where \( \psi_i \) are the basis functions for the FE space of avgs. Both FE spaces should be scalar and on the same mesh.

Definition at line 1729 of file gridfunc.cpp.

◆ GetElementDofValues()

void mfem::GridFunction::GetElementDofValues ( int el,
Vector & dof_vals ) const
virtual

Sets the output vector dof_vals to the values of the degrees of freedom of element el.

Reimplemented in mfem::ParGridFunction.

Definition at line 1769 of file gridfunc.cpp.

◆ GetFaceValues()

int mfem::GridFunction::GetFaceValues ( int i,
int side,
const IntegrationRule & ir,
Vector & vals,
DenseMatrix & tr,
int vdim = 1 ) const

Compute a collection of scalar values from within the face indicated by the index i.

Definition at line 664 of file gridfunc.cpp.

◆ GetFaceVectorValues()

int mfem::GridFunction::GetFaceVectorValues ( int i,
int side,
const IntegrationRule & ir,
DenseMatrix & vals,
DenseMatrix & tr ) const

Compute a collection of vector values from within the face indicated by the index i.

Definition at line 1060 of file gridfunc.cpp.

◆ GetGradient()

void mfem::GridFunction::GetGradient ( ElementTransformation & tr,
Vector & grad ) const

Gradient of a scalar function at a quadrature point.

Note
It is assumed that the IntegrationPoint of interest has been specified by ElementTransformation::SetIntPoint() before calling GetGradient().
Can be used from a ParGridFunction when tr is an ElementTransformation of a face-neighbor element and face-neighbor data has been exchanged.

Definition at line 1582 of file gridfunc.cpp.

◆ GetGradients() [1/2]

void mfem::GridFunction::GetGradients ( const int elem,
const IntegrationRule & ir,
DenseMatrix & grad ) const
inline

Extension of GetGradient(...) for a collection of IntegrationPoints.

Definition at line 346 of file gridfunc.hpp.

◆ GetGradients() [2/2]

void mfem::GridFunction::GetGradients ( ElementTransformation & tr,
const IntegrationRule & ir,
DenseMatrix & grad ) const

Extension of GetGradient(...) for a collection of IntegrationPoints.

Definition at line 1647 of file gridfunc.cpp.

◆ GetHessians() [1/2]

void mfem::GridFunction::GetHessians ( int i,
const IntegrationRule & ir,
DenseMatrix & hess,
DenseMatrix & tr,
int vdim = 1 ) const

Definition at line 651 of file gridfunc.cpp.

◆ GetHessians() [2/2]

void mfem::GridFunction::GetHessians ( int i,
const IntegrationRule & ir,
DenseMatrix & hess,
int vdim = 1 ) const

Definition at line 608 of file gridfunc.cpp.

◆ GetLaplacians() [1/2]

void mfem::GridFunction::GetLaplacians ( int i,
const IntegrationRule & ir,
Vector & laps,
DenseMatrix & tr,
int vdim = 1 ) const

Definition at line 596 of file gridfunc.cpp.

◆ GetLaplacians() [2/2]

void mfem::GridFunction::GetLaplacians ( int i,
const IntegrationRule & ir,
Vector & laps,
int vdim = 1 ) const

Definition at line 569 of file gridfunc.cpp.

◆ GetNodalValues() [1/2]

void mfem::GridFunction::GetNodalValues ( int i,
Array< real_t > & nval,
int vdim = 1 ) const

Returns the values in the vertices of i'th element for dimension vdim.

Definition at line 395 of file gridfunc.cpp.

◆ GetNodalValues() [2/2]

void mfem::GridFunction::GetNodalValues ( Vector & nval,
int vdim = 1 ) const

Return the values as a vector on mesh vertices for dimension vdim.

Definition at line 1920 of file gridfunc.cpp.

◆ GetSequence()

long mfem::GridFunction::GetSequence ( ) const
inline

Return update counter, similar to Mesh::GetSequence(). Used to check if it is up to date with the space.

Definition at line 694 of file gridfunc.hpp.

◆ GetTrueDofs()

void mfem::GridFunction::GetTrueDofs ( Vector & tv) const

Extract the true-dofs from the GridFunction.

Definition at line 366 of file gridfunc.cpp.

◆ GetTrueVector() [1/2]

Vector & mfem::GridFunction::GetTrueVector ( )
inline

Read and write access to the (optional) internal true-dof Vector.

Note that t_vec is set if it is not allocated or set already.

Definition at line 137 of file gridfunc.hpp.

◆ GetTrueVector() [2/2]

const Vector & mfem::GridFunction::GetTrueVector ( ) const
inline

Read only access to the (optional) internal true-dof Vector.

Definition at line 130 of file gridfunc.hpp.

◆ GetValue() [1/2]

real_t mfem::GridFunction::GetValue ( ElementTransformation & T,
const IntegrationPoint & ip,
int comp = 0,
Vector * tr = NULL ) const
virtual

Return a scalar value from within the element indicated by the ElementTransformation Object.

Reimplemented in mfem::ParGridFunction.

Definition at line 723 of file gridfunc.cpp.

◆ GetValue() [2/2]

real_t mfem::GridFunction::GetValue ( int i,
const IntegrationPoint & ip,
int vdim = 1 ) const
virtual

Return a scalar value from within the given element.

Reimplemented in mfem::ParGridFunction.

Definition at line 449 of file gridfunc.cpp.

◆ GetValues() [1/3]

void mfem::GridFunction::GetValues ( ElementTransformation & T,
const IntegrationRule & ir,
Vector & vals,
int comp = 0,
DenseMatrix * tr = NULL ) const

Compute a collection of scalar values from within the element indicated by the ElementTransformation object.

Definition at line 838 of file gridfunc.cpp.

◆ GetValues() [2/3]

void mfem::GridFunction::GetValues ( int i,
const IntegrationRule & ir,
Vector & vals,
DenseMatrix & tr,
int vdim = 1 ) const

Compute a collection of vector values from within the element indicated by the index i.

Definition at line 558 of file gridfunc.cpp.

◆ GetValues() [3/3]

void mfem::GridFunction::GetValues ( int i,
const IntegrationRule & ir,
Vector & vals,
int vdim = 1 ) const

Compute a collection of scalar values from within the element indicated by the index i.

Definition at line 521 of file gridfunc.cpp.

◆ GetValuesFrom()

void mfem::GridFunction::GetValuesFrom ( const GridFunction & orig_func)

Definition at line 1104 of file gridfunc.cpp.

◆ GetVectorFieldNodalValues()

void mfem::GridFunction::GetVectorFieldNodalValues ( Vector & val,
int comp ) const

Definition at line 1253 of file gridfunc.cpp.

◆ GetVectorFieldValues()

void mfem::GridFunction::GetVectorFieldValues ( int i,
const IntegrationRule & ir,
DenseMatrix & vals,
DenseMatrix & tr,
int comp = 0 ) const

Definition at line 1188 of file gridfunc.cpp.

◆ GetVectorGradient()

void mfem::GridFunction::GetVectorGradient ( ElementTransformation & tr,
DenseMatrix & grad ) const

Compute the vector gradient with respect to the physical element variable.

Definition at line 1671 of file gridfunc.cpp.

◆ GetVectorGradientHat()

void mfem::GridFunction::GetVectorGradientHat ( ElementTransformation & T,
DenseMatrix & gh ) const

Compute the vector gradient with respect to the reference element variable.

Definition at line 1388 of file gridfunc.cpp.

◆ GetVectorValue() [1/2]

void mfem::GridFunction::GetVectorValue ( ElementTransformation & T,
const IntegrationPoint & ip,
Vector & val,
Vector * tr = NULL ) const
virtual

Return a vector value from within the element indicated by the ElementTransformation Object.

Reimplemented in mfem::ParGridFunction.

Definition at line 858 of file gridfunc.cpp.

◆ GetVectorValue() [2/2]

void mfem::GridFunction::GetVectorValue ( int i,
const IntegrationPoint & ip,
Vector & val ) const
virtual

Return a vector value from within the given element.

Reimplemented in mfem::ParGridFunction.

Definition at line 476 of file gridfunc.cpp.

◆ GetVectorValues() [1/2]

void mfem::GridFunction::GetVectorValues ( ElementTransformation & T,
const IntegrationRule & ir,
DenseMatrix & vals,
DenseMatrix * tr = NULL ) const

Compute a collection of vector values from within the element indicated by the ElementTransformation object.

Definition at line 998 of file gridfunc.cpp.

◆ GetVectorValues() [2/2]

void mfem::GridFunction::GetVectorValues ( int i,
const IntegrationRule & ir,
DenseMatrix & vals,
DenseMatrix & tr ) const

Definition at line 714 of file gridfunc.cpp.

◆ ImposeBounds() [1/2]

void mfem::GridFunction::ImposeBounds ( int i,
const Vector & weights,
const Vector & lo_,
const Vector & hi_ )

Impose the given bounds on the function's DOFs while preserving its local integral (described in terms of the given weights) on the i'th element through SLBPQ optimization. Intended to be used for discontinuous FE functions.

Definition at line 1832 of file gridfunc.cpp.

◆ ImposeBounds() [2/2]

void mfem::GridFunction::ImposeBounds ( int i,
const Vector & weights,
real_t min_ = 0.0,
real_t max_ = infinity() )

Definition at line 1868 of file gridfunc.cpp.

◆ LegacyNCReorder()

void mfem::GridFunction::LegacyNCReorder ( )
protected

Loading helper.

Definition at line 3868 of file gridfunc.cpp.

◆ MakeOwner()

void mfem::GridFunction::MakeOwner ( FiniteElementCollection * fec_)
inline

Make the GridFunction the owner of fec and fes.

If the new FiniteElementCollection, fec_, is NULL, ownership of fec and fes is taken away.

Definition at line 122 of file gridfunc.hpp.

◆ MakeRef() [1/4]

void mfem::GridFunction::MakeRef ( FiniteElementSpace * f,
real_t * v )
virtual

Make the GridFunction reference external data on a new FiniteElementSpace.

This method changes the FiniteElementSpace associated with the GridFunction and sets the pointer v as external data in the GridFunction.

Reimplemented in mfem::ParGridFunction.

Definition at line 203 of file gridfunc.cpp.

◆ MakeRef() [2/4]

void mfem::GridFunction::MakeRef ( FiniteElementSpace * f,
Vector & v,
int v_offset )
virtual

Make the GridFunction reference external data on a new FiniteElementSpace.

This method changes the FiniteElementSpace associated with the GridFunction and sets the data of the Vector v (plus the v_offset) as external data in the GridFunction.

Note
This version of the method will also perform bounds checks when the build option MFEM_DEBUG is enabled.

Reimplemented in mfem::ParGridFunction.

Definition at line 211 of file gridfunc.cpp.

◆ MakeRef() [3/4]

void mfem::Vector::MakeRef ( Vector & base,
int offset )
inline

Reset the Vector to be a reference to a sub-vector of base without changing its current size.

Definition at line 203 of file vector.hpp.

◆ MakeRef() [4/4]

void mfem::Vector::MakeRef ( Vector & base,
int offset,
int size )
inline

Reset the Vector to be a reference to a sub-vector of base.

Definition at line 199 of file vector.hpp.

◆ MakeTRef() [1/2]

void mfem::GridFunction::MakeTRef ( FiniteElementSpace * f,
real_t * tv )

Associate a new FiniteElementSpace and new true-dof data with the GridFunction.

Definition at line 221 of file gridfunc.cpp.

◆ MakeTRef() [2/2]

void mfem::GridFunction::MakeTRef ( FiniteElementSpace * f,
Vector & tv,
int tv_offset )

Associate a new FiniteElementSpace and new true-dof data with the GridFunction.

  • If the prolongation matrix of f is trivial (i.e. its method FiniteElementSpace::GetProlongationMatrix() returns NULL), this method calls MakeRef() with the same arguments.
  • Otherwise, this method calls SetSpace() with argument f.
  • The internal true-dof vector is set to reference the sub-vector of tv starting at the offset tv_offset.

Definition at line 235 of file gridfunc.cpp.

◆ operator=() [1/3]

GridFunction & mfem::GridFunction::operator= ( const GridFunction & rhs)
inline

Copy assignment. Only the data of the base class Vector is copied.

It is assumed that this object and rhs use FiniteElementSpaces that have the same size.

Note
Defining this method overwrites the implicitly defined copy assignment operator.

Definition at line 116 of file gridfunc.hpp.

◆ operator=() [2/3]

GridFunction & mfem::GridFunction::operator= ( const Vector & v)

Copy the data from v.

The size of v must be equal to the size of the associated FiniteElementSpace fes.

Definition at line 3621 of file gridfunc.cpp.

◆ operator=() [3/3]

GridFunction & mfem::GridFunction::operator= ( real_t value)

Redefine '=' for GridFunction = constant.

Definition at line 3615 of file gridfunc.cpp.

◆ OwnFEC()

FiniteElementCollection * mfem::GridFunction::OwnFEC ( )
inline

Definition at line 124 of file gridfunc.hpp.

◆ ProjectBdrCoefficient() [1/3]

void mfem::GridFunction::ProjectBdrCoefficient ( Coefficient & coeff,
const Array< int > & attr )
inline

Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the boundary attributes marked in the attr array.

Definition at line 469 of file gridfunc.hpp.

◆ ProjectBdrCoefficient() [2/3]

void mfem::GridFunction::ProjectBdrCoefficient ( Coefficient * coeff[],
const Array< int > & attr )
virtual

Project a set of Coefficients on the components of the GridFunction, modifying only DOFs on the boundary associated with the boundary attributed marked in the attr array.

If a Coefficient pointer in the array coeff is NULL, that component will not be touched.

Reimplemented in mfem::ParGridFunction.

Definition at line 2596 of file gridfunc.cpp.

◆ ProjectBdrCoefficient() [3/3]

void mfem::GridFunction::ProjectBdrCoefficient ( VectorCoefficient & vcoeff,
const Array< int > & attr )
virtual

Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with the boundary attributes marked in the attr array.

Reimplemented in mfem::ParGridFunction.

Definition at line 2578 of file gridfunc.cpp.

◆ ProjectBdrCoefficientNormal()

void mfem::GridFunction::ProjectBdrCoefficientNormal ( VectorCoefficient & vcoeff,
const Array< int > & bdr_attr )

Project the normal component of the given VectorCoefficient on the boundary. Only boundary attributes that are marked in 'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction.

Definition at line 2626 of file gridfunc.cpp.

◆ ProjectBdrCoefficientTangent()

void mfem::GridFunction::ProjectBdrCoefficientTangent ( VectorCoefficient & vcoeff,
const Array< int > & bdr_attr )
virtual

Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attributes that are marked in bdr_attr are projected. Assumes ND-type VectorFE GridFunction.

Reimplemented in mfem::ParGridFunction.

Definition at line 2701 of file gridfunc.cpp.

◆ ProjectCoefficient() [1/6]

void mfem::GridFunction::ProjectCoefficient ( Coefficient & coeff)
virtual

Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of the FiniteElementSpace fes. Note that this is usually interpolation at the degrees of freedom in each element (not L2 projection).

Reimplemented in mfem::ParGridFunction.

Definition at line 2346 of file gridfunc.cpp.

◆ ProjectCoefficient() [2/6]

void mfem::GridFunction::ProjectCoefficient ( Coefficient & coeff,
Array< int > & dofs,
int vd = 0 )

Project coeff Coefficient to this GridFunction, using one element for each degree of freedom in dofs and nodal interpolation on that element.

Definition at line 2378 of file gridfunc.cpp.

◆ ProjectCoefficient() [3/6]

void mfem::GridFunction::ProjectCoefficient ( Coefficient * coeff[])

Analogous to the version with argument vcoeff VectorCoefficient but using an array of scalar coefficients for each component.

Definition at line 2483 of file gridfunc.cpp.

◆ ProjectCoefficient() [4/6]

void mfem::GridFunction::ProjectCoefficient ( VectorCoefficient & vcoeff)

Project vcoeff VectorCoefficient to this GridFunction. The projection computation depends on the choice of the FiniteElementSpace fes. Note that this is usually interpolation at the degrees of freedom in each element (not L2 projection).

Definition at line 2404 of file gridfunc.cpp.

◆ ProjectCoefficient() [5/6]

void mfem::GridFunction::ProjectCoefficient ( VectorCoefficient & vcoeff,
Array< int > & dofs )

Project vcoeff VectorCoefficient to this GridFunction, using one element for each degree of freedom in dofs and nodal interpolation on that element.

Definition at line 2425 of file gridfunc.cpp.

◆ ProjectCoefficient() [6/6]

void mfem::GridFunction::ProjectCoefficient ( VectorCoefficient & vcoeff,
int attribute )

Project vcoeff VectorCoefficient to this GridFunction, only projecting onto elements with the given attribute.

Definition at line 2457 of file gridfunc.cpp.

◆ ProjectDeltaCoefficient()

void mfem::GridFunction::ProjectDeltaCoefficient ( DeltaCoefficient & delta_coeff,
real_t & integral )
protected

Definition at line 2279 of file gridfunc.cpp.

◆ ProjectDiscCoefficient() [1/4]

void mfem::GridFunction::ProjectDiscCoefficient ( Coefficient & coeff,
AvgType type )
virtual

Projects a discontinuous coefficient so that the values in shared vdofs are computed by taking an average of the possible values.

Reimplemented in mfem::ParGridFunction.

Definition at line 2558 of file gridfunc.cpp.

◆ ProjectDiscCoefficient() [2/4]

void mfem::GridFunction::ProjectDiscCoefficient ( VectorCoefficient & coeff)
virtual

Project a discontinuous vector coefficient as a grid function on a continuous finite element space. The values in shared dofs are determined from the element with maximal attribute.

Reimplemented in mfem::ParGridFunction.

Definition at line 2552 of file gridfunc.cpp.

◆ ProjectDiscCoefficient() [3/4]

void mfem::GridFunction::ProjectDiscCoefficient ( VectorCoefficient & coeff,
Array< int > & dof_attr )
protected

Project a discontinuous vector coefficient in a continuous space and return in dof_attr the maximal attribute of the elements containing each degree of freedom.

Definition at line 2520 of file gridfunc.cpp.

◆ ProjectDiscCoefficient() [4/4]

void mfem::GridFunction::ProjectDiscCoefficient ( VectorCoefficient & coeff,
AvgType type )
virtual

Projects a discontinuous vector coefficient so that the values in shared vdofs are computed by taking an average of the possible values.

Reimplemented in mfem::ParGridFunction.

Definition at line 2569 of file gridfunc.cpp.

◆ ProjectGridFunction()

void mfem::GridFunction::ProjectGridFunction ( const GridFunction & src)

Project the src GridFunction to this GridFunction, both of which must be on the same mesh.

The current implementation assumes that all elements use the same projection matrix.

Definition at line 1780 of file gridfunc.cpp.

◆ ProjectVectorFieldOn()

void mfem::GridFunction::ProjectVectorFieldOn ( GridFunction & vec_field,
int comp = 0 )

Definition at line 1284 of file gridfunc.cpp.

◆ ReorderByNodes()

void mfem::GridFunction::ReorderByNodes ( )

For a vector grid function, makes sure that the ordering is byNODES.

Definition at line 1226 of file gridfunc.cpp.

◆ RestrictConforming()

void mfem::GridFunction::RestrictConforming ( )

On a non-conforming mesh, make sure the function lies in the conforming space by multiplying with R and then with P, the conforming restriction and prolongation matrices of the space, respectively.

Definition at line 1907 of file gridfunc.cpp.

◆ Save() [1/3]

void mfem::GridFunction::Save ( adios2stream & out,
const std::string & variable_name,
const adios2stream::data_type type = adios2stream::data_type::point_data ) const
virtual

Save the GridFunction to a binary output stream using adios2 bp format.

Reimplemented in mfem::ParGridFunction.

Definition at line 3661 of file gridfunc.cpp.

◆ Save() [2/3]

void mfem::GridFunction::Save ( const char * fname,
int precision = 16 ) const
virtual

Save the GridFunction to a file. The given precision will be used for ASCII output.

Reimplemented in mfem::ParGridFunction.

Definition at line 3653 of file gridfunc.cpp.

◆ Save() [3/3]

void mfem::GridFunction::Save ( std::ostream & out) const
virtual

Save the GridFunction to an output stream.

Reimplemented in mfem::ParGridFunction.

Definition at line 3628 of file gridfunc.cpp.

◆ SaveSTL()

void mfem::GridFunction::SaveSTL ( std::ostream & out,
int TimesToRefine = 1 )

Write the GridFunction in STL format. Note that the mesh dimension must be 2 and that quad elements will be broken into two triangles.

Definition at line 3768 of file gridfunc.cpp.

◆ SaveSTLTri()

void mfem::GridFunction::SaveSTLTri ( std::ostream & out,
real_t p1[],
real_t p2[],
real_t p3[] )
protected

Definition at line 3748 of file gridfunc.cpp.

◆ SaveVTK()

void mfem::GridFunction::SaveVTK ( std::ostream & out,
const std::string & field_name,
int ref )

Write the GridFunction in VTK format. Note that Mesh::PrintVTK must be called first. The parameter ref > 0 must match the one used in Mesh::PrintVTK.

Definition at line 3669 of file gridfunc.cpp.

◆ SetFromTrueDofs()

void mfem::GridFunction::SetFromTrueDofs ( const Vector & tv)
virtual

Set the GridFunction from the given true-dof vector.

Reimplemented in mfem::ParGridFunction.

Definition at line 381 of file gridfunc.cpp.

◆ SetFromTrueVector()

void mfem::GridFunction::SetFromTrueVector ( )
inline

Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.

Definition at line 150 of file gridfunc.hpp.

◆ SetSpace()

void mfem::GridFunction::SetSpace ( FiniteElementSpace * f)
virtual

Associate a new FiniteElementSpace with the GridFunction.

The GridFunction is resized using the SetSize() method.

Reimplemented in mfem::ParGridFunction.

Definition at line 195 of file gridfunc.cpp.

◆ SetTrueVector()

void mfem::GridFunction::SetTrueVector ( )
inline

Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.

Definition at line 144 of file gridfunc.hpp.

◆ SumFluxAndCount()

void mfem::GridFunction::SumFluxAndCount ( BilinearFormIntegrator & blfi,
GridFunction & flux,
Array< int > & counts,
bool wcoef,
int subdomain )
protected

Definition at line 251 of file gridfunc.cpp.

◆ Update()

void mfem::GridFunction::Update ( )
virtual

Transform by the Space UpdateMatrix (e.g., on Mesh change).

Reimplemented in mfem::ParGridFunction.

Definition at line 164 of file gridfunc.cpp.

◆ VectorDim()

int mfem::GridFunction::VectorDim ( ) const

Definition at line 323 of file gridfunc.cpp.

Member Data Documentation

◆ fec

FiniteElementCollection* mfem::GridFunction::fec
protected

Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner().

If not NULL, this pointer is owned by the GridFunction.

Definition at line 40 of file gridfunc.hpp.

◆ fes

FiniteElementSpace* mfem::GridFunction::fes
protected

FE space on which the grid function lives. Owned if fec is not NULL.

Definition at line 34 of file gridfunc.hpp.

◆ fes_sequence

long mfem::GridFunction::fes_sequence
protected

Definition at line 42 of file gridfunc.hpp.

◆ t_vec

Vector mfem::GridFunction::t_vec
protected

Optional, internal true-dof vector: if the FiniteElementSpace fes has a non-trivial (i.e. not NULL) prolongation operator, this Vector may hold associated true-dof values - either owned or external.

Definition at line 47 of file gridfunc.hpp.


The documentation for this class was generated from the following files: