MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
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. More...
 
 GridFunction (FiniteElementSpace *f)
 Construct a GridFunction associated with the FiniteElementSpace *f. More...
 
 GridFunction (FiniteElementSpace *f, double *data)
 Construct a GridFunction using previously allocated array data. More...
 
 GridFunction (Mesh *m, std::istream &input)
 Construct a GridFunction on the given Mesh, using the data from input. More...
 
 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. More...
 
void MakeOwner (FiniteElementCollection *fec_)
 Make the GridFunction the owner of fec and fes. More...
 
FiniteElementCollectionOwnFEC ()
 
int VectorDim () const
 
const VectorGetTrueVector () const
 Read only access to the (optional) internal true-dof Vector. More...
 
VectorGetTrueVector ()
 Read and write access to the (optional) internal true-dof Vector. More...
 
void GetTrueDofs (Vector &tv) const
 Extract the true-dofs from the GridFunction. More...
 
void SetTrueVector ()
 Shortcut for calling GetTrueDofs() with GetTrueVector() as argument. More...
 
virtual void SetFromTrueDofs (const Vector &tv)
 Set the GridFunction from the given true-dof vector. More...
 
void SetFromTrueVector ()
 Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument. More...
 
void GetNodalValues (int i, Array< double > &nval, int vdim=1) const
 Returns the values in the vertices of i'th element for dimension vdim. More...
 
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. More...
 
void GetNodalValues (Vector &nval, int vdim=1) const
 Return the values as a vector on mesh vertices for dimension vdim. More...
 
void GetVectorFieldNodalValues (Vector &val, int comp) const
 
void ProjectVectorFieldOn (GridFunction &vec_field, int comp=0)
 
void GetDerivative (int comp, int der_comp, GridFunction &der)
 
double GetDivergence (ElementTransformation &tr) const
 
void GetCurl (ElementTransformation &tr, Vector &curl) const
 
void GetGradient (ElementTransformation &tr, Vector &grad) const
 
void GetGradients (ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
 
void GetGradients (const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
 
void GetVectorGradient (ElementTransformation &tr, DenseMatrix &grad) const
 
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, double min_=0.0, double 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. More...
 
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). More...
 
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. More...
 
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). More...
 
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. More...
 
void ProjectCoefficient (VectorCoefficient &vcoeff, int attribute)
 Project vcoeff VectorCoefficient to this GridFunction, only projecting onto elements with the given attribute. More...
 
void ProjectCoefficient (Coefficient *coeff[])
 Analogous to the version with argument vcoeff VectorCoefficient but using an array of scalar coefficients for each component. More...
 
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. More...
 
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. More...
 
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. More...
 
void ProjectBdrCoefficient (Coefficient &coeff, 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. More...
 
virtual void ProjectBdrCoefficient (VectorCoefficient &vcoeff, 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. More...
 
virtual void ProjectBdrCoefficient (Coefficient *coeff[], 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. More...
 
void ProjectBdrCoefficientNormal (VectorCoefficient &vcoeff, Array< int > &bdr_attr)
 
virtual void ProjectBdrCoefficientTangent (VectorCoefficient &vcoeff, 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. More...
 
virtual double ComputeL2Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeL2Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeL2Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
 
virtual double ComputeGradError (VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
 Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements. More...
 
virtual double ComputeCurlError (VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
 Returns ||curl u_ex - curl u_h||_L2 for ND elements. More...
 
virtual double ComputeDivError (Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
 Returns ||div u_ex - div u_h||_L2 for RT elements. More...
 
virtual double ComputeDGFaceJumpError (Coefficient *exsol, Coefficient *ell_coeff, class JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
 
MFEM_DEPRECATED double ComputeDGFaceJumpError (Coefficient *exsol, Coefficient *ell_coeff, double Nu, const IntegrationRule *irs[]=NULL) const
 Returns the Face Jumps error for L2 elements, with 1/h scaling. More...
 
virtual double ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
 
virtual double ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeHDivError (VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
 Returns the error measured in H(div)-norm for RT elements. More...
 
virtual double ComputeHCurlError (VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
 Returns the error measured in H(curl)-norm for ND elements. More...
 
virtual double ComputeMaxError (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeMaxError (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeMaxError (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeL1Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeW11Error (Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeL1Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeLpError (const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementLpErrors (const double 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 double ComputeLpError (const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
 
virtual void ComputeElementLpErrors (const double 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= (double value)
 Redefine '=' for GridFunction = constant. More...
 
GridFunctionoperator= (const Vector &v)
 Copy the data from v. More...
 
virtual void Update ()
 Transform by the Space UpdateMatrix (e.g., on Mesh change). More...
 
FiniteElementSpaceFESpace ()
 
const FiniteElementSpaceFESpace () const
 
virtual void SetSpace (FiniteElementSpace *f)
 Associate a new FiniteElementSpace with the GridFunction. More...
 
virtual void MakeRef (FiniteElementSpace *f, double *v)
 Make the GridFunction reference external data on a new FiniteElementSpace. More...
 
virtual void MakeRef (FiniteElementSpace *f, Vector &v, int v_offset)
 Make the GridFunction reference external data on a new FiniteElementSpace. More...
 
void MakeTRef (FiniteElementSpace *f, double *tv)
 Associate a new FiniteElementSpace and new true-dof data with the GridFunction. More...
 
void MakeTRef (FiniteElementSpace *f, Vector &tv, int tv_offset)
 Associate a new FiniteElementSpace and new true-dof data with the GridFunction. More...
 
virtual void Save (std::ostream &out) const
 Save the GridFunction to an output stream. More...
 
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. More...
 
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. More...
 
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. More...
 
virtual ~GridFunction ()
 Destroys grid function. More...
 
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 double 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 double 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 ()
 Default constructor for Vector. Sets size = 0 and data = NULL. More...
 
 Vector (const Vector &)
 Copy constructor. Allocates a new data array and copies the data. More...
 
 Vector (int s)
 Creates vector of size s. More...
 
 Vector (double *data_, int size_)
 Creates a vector referencing an array of doubles, owned by someone else. More...
 
 Vector (int size_, MemoryType mt)
 Create a Vector of size size_ using MemoryType mt. More...
 
 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. More...
 
template<int N>
 Vector (const double(&values)[N])
 Create a vector using a braced initializer list. More...
 
virtual void UseDevice (bool use_dev) const
 Enable execution of Vector operations using the mfem::Device. More...
 
virtual bool UseDevice () const
 Return the device flag of the Memory object used by the Vector. More...
 
void Load (std::istream **in, int np, int *dim)
 Reads a vector from multiple files. More...
 
void Load (std::istream &in, int Size)
 Load a vector from an input stream. More...
 
void Load (std::istream &in)
 Load a vector from an input stream, reading the size from the stream. More...
 
void SetSize (int s)
 Resize the vector to size s. More...
 
void SetSize (int s, MemoryType mt)
 Resize the vector to size s using MemoryType mt. More...
 
void SetSize (int s, Vector &v)
 Resize the vector to size s using the MemoryType of v. More...
 
void SetData (double *d)
 
void SetDataAndSize (double *d, int s)
 Set the Vector data and size. More...
 
void NewDataAndSize (double *d, int s)
 Set the Vector data and size, deleting the old data, if owned. More...
 
void NewMemoryAndSize (const Memory< double > &mem, int s, bool own_mem)
 Reset the Vector to use the given external Memory mem and size s. More...
 
void MakeRef (Vector &base, int offset, int size)
 Reset the Vector to be a reference to a sub-vector of base. More...
 
void MakeRef (Vector &base, int offset)
 Reset the Vector to be a reference to a sub-vector of base without changing its current size. More...
 
void MakeDataOwner () const
 Set the Vector data (host pointer) ownership flag. More...
 
void Destroy ()
 Destroy a vector. More...
 
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. More...
 
int Size () const
 Returns the size of the vector. More...
 
int Capacity () const
 Return the size of the currently allocated data array. More...
 
double * GetData () const
 Return a pointer to the beginning of the Vector data. More...
 
 operator double * ()
 Conversion to double *. More...
 
 operator const double * () const
 Conversion to const double *. More...
 
double * begin ()
 STL-like begin. More...
 
double * end ()
 STL-like end. More...
 
const double * begin () const
 STL-like begin (const version). More...
 
const double * end () const
 STL-like end (const version). More...
 
Memory< double > & GetMemory ()
 Return a reference to the Memory object used by the Vector. More...
 
const Memory< double > & GetMemory () const
 Return a reference to the Memory object used by the Vector, const version. More...
 
void SyncMemory (const Vector &v) const
 Update the memory location of the vector to match v. More...
 
void SyncAliasMemory (const Vector &v) const
 Update the alias memory location of the vector to match v. More...
 
bool OwnsData () const
 Read the Vector data (host pointer) ownership flag. More...
 
void StealData (double **p)
 Changes the ownership of the data; after the call the Vector is empty. More...
 
double * StealData ()
 Changes the ownership of the data; after the call the Vector is empty. More...
 
double & Elem (int i)
 Access Vector entries. Index i = 0 .. size-1. More...
 
const double & Elem (int i) const
 Read only access to Vector entries. Index i = 0 .. size-1. More...
 
double & operator() (int i)
 Access Vector entries using () for 0-based indexing. More...
 
const double & operator() (int i) const
 Read only access to Vector entries using () for 0-based indexing. More...
 
double operator* (const double *) const
 Dot product with a double * array. More...
 
double operator* (const Vector &v) const
 Return the inner-product. More...
 
Vectoroperator= (const double *v)
 Copy Size() entries from v. More...
 
Vectoroperator= (const Vector &v)
 Copy assignment. More...
 
Vectoroperator= (double value)
 Redefine '=' for vector = constant. More...
 
Vectoroperator*= (double c)
 
Vectoroperator*= (const Vector &v)
 Component-wise scaling: (*this)(i) *= v(i) More...
 
Vectoroperator/= (double c)
 
Vectoroperator/= (const Vector &v)
 Component-wise division: (*this)(i) /= v(i) More...
 
Vectoroperator-= (double c)
 
Vectoroperator-= (const Vector &v)
 
Vectoroperator+= (double c)
 
Vectoroperator+= (const Vector &v)
 
VectorAdd (const double a, const Vector &Va)
 (*this) += a * Va More...
 
VectorSet (const double a, const Vector &x)
 (*this) = a * x More...
 
void SetVector (const Vector &v, int offset)
 
void Neg ()
 (*this) = -(*this) More...
 
void Swap (Vector &other)
 Swap the contents of two Vectors. More...
 
void median (const Vector &lo, const Vector &hi)
 v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi. More...
 
void GetSubVector (const Array< int > &dofs, Vector &elemvect) const
 Extract entries listed in dofs to the output Vector elemvect. More...
 
void GetSubVector (const Array< int > &dofs, double *elem_data) const
 Extract entries listed in dofs to the output array elem_data. More...
 
void SetSubVector (const Array< int > &dofs, const double value)
 Set the entries listed in dofs to the given value. More...
 
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. More...
 
void SetSubVector (const Array< int > &dofs, double *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. More...
 
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. More...
 
void AddElementVector (const Array< int > &dofs, double *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. More...
 
void AddElementVector (const Array< int > &dofs, const double 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. More...
 
void SetSubVectorComplement (const Array< int > &dofs, const double val)
 Set all vector entries NOT in the dofs Array to the given val. More...
 
void Print (std::ostream &out=mfem::out, int width=8) const
 Prints vector to stream out. More...
 
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. More...
 
void PrintHash (std::ostream &out) const
 Print the Vector size and hash of its data. More...
 
void Randomize (int seed=0)
 Set random values in the vector. More...
 
double Norml2 () const
 Returns the l2 norm of the vector. More...
 
double Normlinf () const
 Returns the l_infinity norm of the vector. More...
 
double Norml1 () const
 Returns the l_1 norm of the vector. More...
 
double Normlp (double p) const
 Returns the l_p norm of the vector. More...
 
double Max () const
 Returns the maximal element of the vector. More...
 
double Min () const
 Returns the minimal element of the vector. More...
 
double Sum () const
 Return the sum of the vector entries. More...
 
double DistanceSquaredTo (const double *p) const
 Compute the square of the Euclidean distance to another vector. More...
 
double DistanceTo (const double *p) const
 Compute the Euclidean distance to another vector. More...
 
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. More...
 
virtual ~Vector ()
 Destroys vector. More...
 
virtual const double * Read (bool on_dev=true) const
 Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev). More...
 
virtual const double * HostRead () const
 Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false). More...
 
virtual double * Write (bool on_dev=true)
 Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev). More...
 
virtual double * HostWrite ()
 Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false). More...
 
virtual double * ReadWrite (bool on_dev=true)
 Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev). More...
 
virtual double * HostReadWrite ()
 Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false). More...
 
MFEM_DEPRECATED Vector (N_Vector nv)
 (DEPRECATED) Construct a wrapper Vector from SUNDIALS N_Vector. More...
 
virtual MFEM_DEPRECATED N_Vector ToNVector ()
 (DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. More...
 
virtual MFEM_DEPRECATED void ToNVector (N_Vector &nv, long global_length=0)
 Update an existing wrapper SUNDIALS N_Vector to point to this Vector. More...
 

Protected Member Functions

void SaveSTLTri (std::ostream &out, double p1[], double p2[], double p3[])
 
void GetVectorGradientHat (ElementTransformation &T, DenseMatrix &gh) const
 
void ProjectDeltaCoefficient (DeltaCoefficient &delta_coeff, double &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. More...
 
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. More...
 
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. More...
 
void AccumulateAndCountBdrValues (Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr, Array< int > &values_counter)
 
void AccumulateAndCountBdrTangentValues (VectorCoefficient &vcoeff, 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. More...
 
FiniteElementCollectionfec
 Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner(). More...
 
long fes_sequence
 
Vector t_vec
 
- Protected Attributes inherited from mfem::Vector
Memory< double > data
 
int size
 

Detailed Description

Class for grid function - Vector with associated FE space.

Definition at line 30 of file gridfunc.hpp.

Member Enumeration Documentation

Enumerator
ARITHMETIC 
HARMONIC 

Definition at line 388 of file gridfunc.hpp.

Constructor & Destructor Documentation

mfem::GridFunction::GridFunction ( )
inline

Definition at line 77 of file gridfunc.hpp.

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

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

Definition at line 80 of file gridfunc.hpp.

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

Construct a GridFunction associated with the FiniteElementSpace *f.

Definition at line 85 of file gridfunc.hpp.

mfem::GridFunction::GridFunction ( FiniteElementSpace f,
double *  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 94 of file gridfunc.hpp.

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.

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

Definition at line 76 of file gridfunc.cpp.

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

Destroys grid function.

Definition at line 697 of file gridfunc.hpp.

Member Function Documentation

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

Definition at line 2139 of file gridfunc.cpp.

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

Definition at line 1997 of file gridfunc.cpp.

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 1901 of file gridfunc.cpp.

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 1942 of file gridfunc.cpp.

double 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 2749 of file gridfunc.cpp.

double 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 2827 of file gridfunc.cpp.

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

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

Definition at line 2941 of file gridfunc.cpp.

double 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 2791 of file gridfunc.cpp.

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

Definition at line 560 of file gridfunc.hpp.

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

Definition at line 596 of file gridfunc.hpp.

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

Definition at line 566 of file gridfunc.hpp.

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

Definition at line 602 of file gridfunc.hpp.

void mfem::GridFunction::ComputeElementLpErrors ( const double  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 3221 of file gridfunc.cpp.

void mfem::GridFunction::ComputeElementLpErrors ( const double  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 3378 of file gridfunc.cpp.

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

Definition at line 572 of file gridfunc.hpp.

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

Definition at line 608 of file gridfunc.hpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 296 of file gridfunc.cpp.

double 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 2709 of file gridfunc.cpp.

double mfem::GridFunction::ComputeH1Error ( Coefficient exsol,
VectorCoefficient exgrad,
Coefficient ell_coef,
double  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 2950 of file gridfunc.cpp.

double 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 2967 of file gridfunc.cpp.

double 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 2985 of file gridfunc.cpp.

double 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 2976 of file gridfunc.cpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 535 of file gridfunc.hpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 543 of file gridfunc.hpp.

virtual double mfem::GridFunction::ComputeL2Error ( Coefficient exsol,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

Reimplemented in mfem::ParGridFunction.

Definition at line 456 of file gridfunc.hpp.

double mfem::GridFunction::ComputeL2Error ( Coefficient exsol[],
const IntegrationRule irs[] = NULL 
) const
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 2615 of file gridfunc.cpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 2668 of file gridfunc.cpp.

double mfem::GridFunction::ComputeLpError ( const double  p,
Coefficient exsol,
Coefficient weight = NULL,
const IntegrationRule irs[] = NULL 
) const
virtual

Reimplemented in mfem::ParGridFunction.

Definition at line 3156 of file gridfunc.cpp.

double mfem::GridFunction::ComputeLpError ( const double  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 3287 of file gridfunc.cpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 520 of file gridfunc.hpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 2994 of file gridfunc.cpp.

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

Reimplemented in mfem::ParGridFunction.

Definition at line 529 of file gridfunc.hpp.

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

Definition at line 2190 of file gridfunc.cpp.

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

Definition at line 3050 of file gridfunc.cpp.

void mfem::GridFunction::Destroy ( )
protected

Definition at line 154 of file gridfunc.cpp.

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

Definition at line 629 of file gridfunc.hpp.

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

Definition at line 630 of file gridfunc.hpp.

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

Definition at line 1124 of file gridfunc.cpp.

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

Definition at line 1466 of file gridfunc.cpp.

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

Definition at line 1301 of file gridfunc.cpp.

double mfem::GridFunction::GetDivergence ( ElementTransformation tr) const

Definition at line 1380 of file gridfunc.cpp.

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 1723 of file gridfunc.cpp.

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 1753 of file gridfunc.cpp.

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 605 of file gridfunc.cpp.

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 1044 of file gridfunc.cpp.

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

Definition at line 1568 of file gridfunc.cpp.

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

Definition at line 1636 of file gridfunc.cpp.

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

Definition at line 318 of file gridfunc.hpp.

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

Definition at line 549 of file gridfunc.cpp.

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

Definition at line 592 of file gridfunc.cpp.

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

Definition at line 510 of file gridfunc.cpp.

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

Definition at line 537 of file gridfunc.cpp.

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

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

Definition at line 361 of file gridfunc.cpp.

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 1875 of file gridfunc.cpp.

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

Extract the true-dofs from the GridFunction.

Definition at line 332 of file gridfunc.cpp.

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

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

Note that the returned Vector may be empty, if not previously allocated or set.

Definition at line 127 of file gridfunc.hpp.

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

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

Note that the returned Vector may be empty, if not previously allocated or set.

Definition at line 131 of file gridfunc.hpp.

double 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 402 of file gridfunc.cpp.

double 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 714 of file gridfunc.cpp.

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 466 of file gridfunc.cpp.

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 499 of file gridfunc.cpp.

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 832 of file gridfunc.cpp.

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

Definition at line 1087 of file gridfunc.cpp.

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

Definition at line 1230 of file gridfunc.cpp.

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

Definition at line 1161 of file gridfunc.cpp.

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

Definition at line 1661 of file gridfunc.cpp.

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

Definition at line 1361 of file gridfunc.cpp.

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 425 of file gridfunc.cpp.

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 852 of file gridfunc.cpp.

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

Definition at line 655 of file gridfunc.cpp.

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 987 of file gridfunc.cpp.

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 1804 of file gridfunc.cpp.

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

Definition at line 1831 of file gridfunc.cpp.

void mfem::GridFunction::LegacyNCReorder ( )
protected

Loading helper.

Definition at line 3724 of file gridfunc.cpp.

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 118 of file gridfunc.hpp.

void mfem::GridFunction::MakeRef ( FiniteElementSpace f,
double *  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 201 of file gridfunc.cpp.

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 209 of file gridfunc.cpp.

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

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

Definition at line 219 of file gridfunc.cpp.

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 233 of file gridfunc.cpp.

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 112 of file gridfunc.hpp.

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

Redefine '=' for GridFunction = constant.

Definition at line 3471 of file gridfunc.cpp.

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 3477 of file gridfunc.cpp.

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

Definition at line 120 of file gridfunc.hpp.

void mfem::GridFunction::ProjectBdrCoefficient ( Coefficient coeff,
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 424 of file gridfunc.hpp.

void mfem::GridFunction::ProjectBdrCoefficient ( VectorCoefficient vcoeff,
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 2491 of file gridfunc.cpp.

void mfem::GridFunction::ProjectBdrCoefficient ( Coefficient coeff[],
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 2509 of file gridfunc.cpp.

void mfem::GridFunction::ProjectBdrCoefficientNormal ( VectorCoefficient vcoeff,
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 2527 of file gridfunc.cpp.

void mfem::GridFunction::ProjectBdrCoefficientTangent ( VectorCoefficient vcoeff,
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 2598 of file gridfunc.cpp.

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 2278 of file gridfunc.cpp.

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 2305 of file gridfunc.cpp.

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 2331 of file gridfunc.cpp.

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 2346 of file gridfunc.cpp.

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 2378 of file gridfunc.cpp.

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 2398 of file gridfunc.cpp.

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

Definition at line 2215 of file gridfunc.cpp.

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 2433 of file gridfunc.cpp.

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 2465 of file gridfunc.cpp.

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 2471 of file gridfunc.cpp.

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 2482 of file gridfunc.cpp.

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 1760 of file gridfunc.cpp.

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

Definition at line 1261 of file gridfunc.cpp.

void mfem::GridFunction::ReorderByNodes ( )

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

Definition at line 1203 of file gridfunc.cpp.

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 1862 of file gridfunc.cpp.

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

Save the GridFunction to an output stream.

Reimplemented in mfem::ParGridFunction.

Definition at line 3484 of file gridfunc.cpp.

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 3509 of file gridfunc.cpp.

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 3517 of file gridfunc.cpp.

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 3624 of file gridfunc.cpp.

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

Definition at line 3604 of file gridfunc.cpp.

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 3525 of file gridfunc.cpp.

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

Set the GridFunction from the given true-dof vector.

Reimplemented in mfem::ParGridFunction.

Definition at line 347 of file gridfunc.cpp.

void mfem::GridFunction::SetFromTrueVector ( )
inline

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

Definition at line 143 of file gridfunc.hpp.

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 193 of file gridfunc.cpp.

void mfem::GridFunction::SetTrueVector ( )
inline

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

Definition at line 137 of file gridfunc.hpp.

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

Definition at line 249 of file gridfunc.cpp.

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.

int mfem::GridFunction::VectorDim ( ) const

Definition at line 311 of file gridfunc.cpp.

Member Data Documentation

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.

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.

long mfem::GridFunction::fes_sequence
protected

Definition at line 42 of file gridfunc.hpp.

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: