MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::ParGridFunction Class Reference

Class for parallel grid function. More...

#include <pgridfunc.hpp>

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

Public Member Functions

 ParGridFunction ()
 
 ParGridFunction (const ParGridFunction &orig)
 Copy constructor. The internal vector face_nbr_data is not copied.
 
 ParGridFunction (ParFiniteElementSpace *pf)
 
 ParGridFunction (ParFiniteElementSpace *pf, real_t *data)
 Construct a ParGridFunction using previously allocated array data.
 
 ParGridFunction (ParFiniteElementSpace *pf, Vector &base, int base_offset=0)
 Construct a ParGridFunction using previously allocated Vector base starting at the given offset, base_offset.
 
 ParGridFunction (ParFiniteElementSpace *pf, GridFunction *gf)
 Construct a ParGridFunction using a GridFunction as external data.
 
 ParGridFunction (ParFiniteElementSpace *pf, HypreParVector *tv)
 Creates grid function on (all) dofs from a given vector on the true dofs, i.e. P tv.
 
 ParGridFunction (ParMesh *pmesh, const GridFunction *gf, const int *partitioning=NULL)
 Construct a local ParGridFunction from the given global GridFunction. If partitioning is NULL (default), the data from gf is NOT copied.
 
 ParGridFunction (ParMesh *pmesh, std::istream &input)
 Construct a ParGridFunction on a given ParMesh, pmesh, reading from an std::istream.
 
ParGridFunctionoperator= (const ParGridFunction &rhs)
 Copy assignment. Only the data of the base class Vector is copied.
 
ParGridFunctionoperator= (real_t value)
 Assign constant values to the ParGridFunction data.
 
ParGridFunctionoperator= (const Vector &v)
 Copy the data from a Vector to the ParGridFunction data.
 
ParFiniteElementSpaceParFESpace () const
 
void Update () override
 Transform by the Space UpdateMatrix (e.g., on Mesh change).
 
void SetSpace (FiniteElementSpace *f) override
 Associate a new FiniteElementSpace with the ParGridFunction.
 
void SetSpace (ParFiniteElementSpace *f)
 Associate a new parallel space with the ParGridFunction.
 
void MakeRef (FiniteElementSpace *f, real_t *v) override
 Make the ParGridFunction reference external data on a new FiniteElementSpace.
 
void MakeRef (ParFiniteElementSpace *f, real_t *v)
 Make the ParGridFunction reference external data on a new ParFiniteElementSpace.
 
void MakeRef (FiniteElementSpace *f, Vector &v, int v_offset) override
 Make the ParGridFunction reference external data on a new FiniteElementSpace.
 
void MakeRef (ParFiniteElementSpace *f, Vector &v, int v_offset)
 Make the ParGridFunction reference external data on a new ParFiniteElementSpace.
 
void Distribute (const Vector *tv)
 
void Distribute (const Vector &tv)
 
void AddDistribute (real_t a, const Vector *tv)
 
void AddDistribute (real_t a, const Vector &tv)
 
void SetFromTrueDofs (const Vector &tv) override
 Set the GridFunction from the given true-dof vector.
 
ParGridFunctionoperator= (const HypreParVector &tv)
 Short semantic for Distribute()
 
HypreParVectorGetTrueDofs () const
 Returns the true dofs in a new HypreParVector.
 
void ParallelAverage (Vector &tv) const
 Returns the vector averaged on the true dofs.
 
void ParallelAverage (HypreParVector &tv) const
 Returns the vector averaged on the true dofs.
 
HypreParVectorParallelAverage () const
 Returns a new vector averaged on the true dofs.
 
void ParallelProject (Vector &tv) const
 Returns the vector restricted to the true dofs.
 
void ParallelProject (HypreParVector &tv) const
 Returns the vector restricted to the true dofs.
 
HypreParVectorParallelProject () const
 Returns a new vector restricted to the true dofs.
 
void ParallelAssemble (Vector &tv) const
 Returns the vector assembled on the true dofs.
 
void ParallelAssemble (HypreParVector &tv) const
 Returns the vector assembled on the true dofs.
 
HypreParVectorParallelAssemble () const
 Returns a new vector assembled on the true dofs.
 
void ExchangeFaceNbrData ()
 
VectorFaceNbrData ()
 
const VectorFaceNbrData () const
 
real_t GetValue (int i, const IntegrationPoint &ip, int vdim=1) const override
 
real_t GetValue (ElementTransformation &T)
 
real_t GetValue (ElementTransformation &T, const IntegrationPoint &ip, int comp=0, Vector *tr=NULL) const override
 
void GetVectorValue (int i, const IntegrationPoint &ip, Vector &val) const override
 
void GetVectorValue (ElementTransformation &T, const IntegrationPoint &ip, Vector &val, Vector *tr=NULL) const override
 
void CountElementsPerVDof (Array< int > &elem_per_vdof) const override
 For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteElementSpace::GetElementVDofs().
 
void GetDerivative (int comp, int der_comp, ParGridFunction &der) const
 Parallel version of GridFunction::GetDerivative(); see its documentation.
 
void GetElementDofValues (int el, Vector &dof_vals) const override
 
void ProjectCoefficient (Coefficient &coeff) override
 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). For NURBS spaces these degrees of freedom are not available and L2 projection is resorted to as fallback.
 
void ProjectDiscCoefficient (VectorCoefficient &coeff) override
 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.
 
void ProjectDiscCoefficient (Coefficient &coeff, AvgType type) override
 Projects a discontinuous coefficient so that the values in shared vdofs are computed by taking an average of the possible values.
 
void ProjectDiscCoefficient (VectorCoefficient &vcoeff, AvgType type) override
 Projects a discontinuous vector coefficient so that the values in shared vdofs are computed by taking an average of the possible values.
 
void ProjectBdrCoefficient (VectorCoefficient &vcoeff, const Array< int > &attr) override
 Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with the boundary attributes marked in the attr array.
 
void ProjectBdrCoefficient (Coefficient *coeff[], const Array< int > &attr) override
 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 ProjectBdrCoefficientTangent (VectorCoefficient &vcoeff, const Array< int > &bdr_attr) override
 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.
 
MFEM_DEPRECATED real_t ComputeL1Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const override
 Returns ||u_ex - u_h||_L1 in parallel for H1 or L2 elements.
 
real_t ComputeL1Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const override
 Returns ||u_ex - u_h||_L1 in parallel for scalar fields.
 
real_t ComputeL1Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const override
 Returns ||u_ex - u_h||_L1 in parallel for vector fields.
 
real_t ComputeL2Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
 Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
 
real_t ComputeL2Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
 Returns ||u_ex - u_h||_L2 in parallel for scalar fields.
 
real_t ComputeL2Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
 Returns ||u_ex - u_h||_L2 in parallel for vector fields.
 
real_t ComputeGradError (VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const override
 Returns ||grad u_ex - grad u_h||_L2 in parallel for H1 or L2 elements.
 
real_t ComputeCurlError (VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const override
 Returns ||curl u_ex - curl u_h||_L2 in parallel for ND elements.
 
real_t ComputeDivError (Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const override
 Returns ||div u_ex - div u_h||_L2 in parallel for RT elements.
 
real_t ComputeDGFaceJumpError (Coefficient *exsol, Coefficient *ell_coeff, JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const override
 Returns the Face Jumps error for L2 elements.
 
real_t ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t Nu, int norm_type) const override
 
real_t ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const override
 Returns the error measured in H1-norm in parallel for H1 or L2 elements.
 
real_t ComputeHDivError (VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const override
 Returns the error measured H(div)-norm in parallel for RT elements.
 
real_t ComputeHCurlError (VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const override
 Returns the error measured H(curl)-norm in parallel for ND elements.
 
real_t ComputeMaxError (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const override
 Returns Max|u_ex - u_h| error in parallel for scalar or vector fields.
 
real_t ComputeMaxError (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const override
 Returns Max|u_ex - u_h| error in parallel for scalar fields.
 
real_t ComputeMaxError (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const override
 Returns Max|u_ex - u_h| error in parallel for vector fields.
 
real_t ComputeLpError (const real_t p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
 Returns ||u_ex - u_h||_Lp in parallel for scalar fields.
 
real_t ComputeLpError (const real_t p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const override
 Returns ||u_ex - u_h||_Lp in parallel for vector fields.
 
void ComputeFlux (BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1) override
 
void Save (std::ostream &out) const override
 
void SaveAsOne (const char *fname, int precision=16) const
 
void Save (const char *fname, int precision=16) const override
 
GridFunction GetSerialGridFunction (int save_rank, Mesh &serial_mesh) const
 Returns a GridFunction on MPI rank save_rank that does not have any duplication of vertices/nodes at processor boundaries.
 
GridFunction GetSerialGridFunction (int save_rank, FiniteElementSpace &serial_fes) const
 Returns a GridFunction on MPI rank save_rank that does not have any duplication of vertices/nodes at processor boundaries.
 
void SaveAsSerial (const char *fname, int precision=16, int save_rank=0) const
 
void Save (adios2stream &out, const std::string &variable_name, const adios2stream::data_type type=adios2stream::data_type::point_data) const override
 
void SaveAsOne (std::ostream &out=mfem::out) const
 Merge the local grid functions.
 
std::unique_ptr< ParGridFunctionProlongateToMaxOrder () const
 Return a GridFunction with the values of this, prolongated to the maximum order of all elements in the mesh.
 
virtual ~ParGridFunction ()=default
 
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 GetTrueDofs (Vector &tv) const
 Extract the true-dofs from the GridFunction.
 
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). For NURBS spaces these degrees of freedom are not available and L2 projection is resorted to as fallback.
 
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.
 
void ProjectDiscCoefficient (VectorCoefficient &coeff, Array< int > &dof_attr)
 
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.
 
- Public Member Functions inherited from mfem::GridFunction
 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_owned and fes.
 
FiniteElementCollectionOwnFEC ()
 
int VectorDim () const
 Shortcut for calling FiniteElementSpace::GetVectorDim() on the underlying fes.
 
int CurlDim () const
 Shortcut for calling FiniteElementSpace::GetCurlDim() on the underlying fes.
 
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.
 
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) const
 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
 
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.
 
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). For NURBS spaces these degrees of freedom are not available and L2 projection is resorted to as fallback.
 
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.
 
std::unique_ptr< GridFunctionProlongateToMaxOrder () const
 Return a GridFunction with the values of this, prolongated to the maximum order of all elements in the mesh.
 
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.
 
void ProjectBdrCoefficientNormal (VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
 
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.
 
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 ComputeW11Error (Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, const Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
 Returns W11 norm (or portions thereof) for H1 or L2 elements.
 
virtual void ComputeElementLpErrors (const real_t p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
 Returns ||u_ex - u_h||_Lp elementwise for H1 or L2 elements.
 
virtual void ComputeElementL1Errors (Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 Returns ||u_ex - u_h||_L1 elementwise for H1 or L2 elements.
 
virtual void ComputeElementL2Errors (Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 Returns ||u_ex - u_h||_L2 elementwise for H1 or L2 elements.
 
virtual void ComputeElementMaxErrors (Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 Returns Max|u_ex - u_h| elementwise for H1 or L2 elements.
 
virtual void ComputeElementLpErrors (const real_t p, VectorCoefficient &exsol, Vector &error, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
 Returns ||u_ex - u_h||_Lp elementwise for vector fields.
 
virtual void ComputeElementL1Errors (VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 Returns ||u_ex - u_h||_L1 elementwise for vector fields.
 
virtual void ComputeElementL2Errors (VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 Returns ||u_ex - u_h||_L2 elementwise for vector fields.
 
virtual void ComputeElementMaxErrors (VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
 Returns Max|u_ex - u_h| elementwise for vector fields.
 
GridFunctionoperator= (real_t value)
 Redefine '=' for GridFunction = constant.
 
GridFunctionoperator= (const Vector &v)
 Copy the data from v.
 
long GetSequence () const
 
FiniteElementSpaceFESpace ()
 
const FiniteElementSpaceFESpace () const
 
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.
 
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.
 
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
 
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
 
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<typename CT , int N>
 Vector (const CT(&values)[N])
 Create a vector from a statically sized C-style array of convertible type.
 
template<typename CT , typename std::enable_if< std::is_convertible< CT, real_t >::value, bool >::type = true>
 Vector (std::initializer_list< CT > values)
 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 *v) const
 
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 PrintMathematica (std::ostream &out=mfem::out) const
 Prints vector as a List for importing into Mathematica.
 
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 ProjectBdrCoefficient (Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
 
- Protected Member Functions inherited from mfem::GridFunction
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) const
 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)
 
void UpdatePRef ()
 P-refinement version of Update().
 

Protected Attributes

ParFiniteElementSpacepfes
 Points to the same object as fes.
 
Vector face_nbr_data
 Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
 
Vector send_data
 Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring processors.
 
- Protected Attributes inherited from mfem::GridFunction
FiniteElementSpacefes
 FE space on which the grid function lives. Owned if fec_owned is not NULL.
 
FiniteElementCollectionfec_owned
 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
 

Additional Inherited Members

- Public Types inherited from mfem::GridFunction
enum  AvgType { ARITHMETIC , HARMONIC }
 

Detailed Description

Class for parallel grid function.

Definition at line 49 of file pgridfunc.hpp.

Constructor & Destructor Documentation

◆ ParGridFunction() [1/9]

mfem::ParGridFunction::ParGridFunction ( )
inline

Definition at line 67 of file pgridfunc.hpp.

◆ ParGridFunction() [2/9]

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

Copy constructor. The internal vector face_nbr_data is not copied.

Definition at line 70 of file pgridfunc.hpp.

◆ ParGridFunction() [3/9]

mfem::ParGridFunction::ParGridFunction ( ParFiniteElementSpace * pf)
inline

Definition at line 73 of file pgridfunc.hpp.

◆ ParGridFunction() [4/9]

mfem::ParGridFunction::ParGridFunction ( ParFiniteElementSpace * pf,
real_t * data )
inline

Construct a ParGridFunction using previously allocated array data.

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

Definition at line 82 of file pgridfunc.hpp.

◆ ParGridFunction() [5/9]

mfem::ParGridFunction::ParGridFunction ( ParFiniteElementSpace * pf,
Vector & base,
int base_offset = 0 )
inline

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

Definition at line 87 of file pgridfunc.hpp.

◆ ParGridFunction() [6/9]

mfem::ParGridFunction::ParGridFunction ( ParFiniteElementSpace * pf,
GridFunction * gf )

Construct a ParGridFunction using a GridFunction as external data.

The parallel space *pf and the space used by *gf should match. The data from *gf is used as the local data of the ParGridFunction on each processor. The ParGridFunction does not assume ownership of the data.

Definition at line 25 of file pgridfunc.cpp.

◆ ParGridFunction() [7/9]

mfem::ParGridFunction::ParGridFunction ( ParFiniteElementSpace * pf,
HypreParVector * tv )

Creates grid function on (all) dofs from a given vector on the true dofs, i.e. P tv.

Definition at line 31 of file pgridfunc.cpp.

◆ ParGridFunction() [8/9]

mfem::ParGridFunction::ParGridFunction ( ParMesh * pmesh,
const GridFunction * gf,
const int * partitioning = NULL )

Construct a local ParGridFunction from the given global GridFunction. If partitioning is NULL (default), the data from gf is NOT copied.

Definition at line 37 of file pgridfunc.cpp.

◆ ParGridFunction() [9/9]

mfem::ParGridFunction::ParGridFunction ( ParMesh * pmesh,
std::istream & input )

Construct a ParGridFunction on a given ParMesh, pmesh, reading from an std::istream.

In the process, a ParFiniteElementSpace and a FiniteElementCollection are constructed. The new ParGridFunction assumes ownership of both.

Definition at line 81 of file pgridfunc.cpp.

◆ ~ParGridFunction()

virtual mfem::ParGridFunction::~ParGridFunction ( )
virtualdefault

Member Function Documentation

◆ AddDistribute() [1/2]

void mfem::ParGridFunction::AddDistribute ( real_t a,
const Vector & tv )
inline

Definition at line 183 of file pgridfunc.hpp.

◆ AddDistribute() [2/2]

void mfem::ParGridFunction::AddDistribute ( real_t a,
const Vector * tv )

Definition at line 148 of file pgridfunc.cpp.

◆ ComputeCurlError()

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

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

See also
GridFunction::ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 390 of file pgridfunc.hpp.

◆ ComputeDGFaceJumpError()

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

Returns the Face Jumps error for L2 elements.

Computes:

facesfjs[f]ell[f](2uexu1u2)2

Where js[f] is the jump_scaling evaluated on the face f and ell is the average of ell_coef evaluated in the two elements sharing the face f.

Parameters
[in]exsolPointer to a Coefficient object reproducing the anticipated values of the scalar field, u_ex.
[in]ell_coeffPointer to a Coefficient object used to compute the averaged value ell in the above integral.
[in]jump_scalingCan be configured to provide scaling by nu, nu/h, or nu*p^2/h
[in]irsOptional pointer to a custom integration rule e.g. higher order than the default rule.
Note
Quadratures with negative weights (as in some simplex integration rules in MFEM) can produce negative integrals even with non-negative integrands. To avoid returning negative errors this function uses the absolute values of the element-wise integrals. This may lead to results which are not entirely consistent with such integration rules.

Reimplemented from mfem::GridFunction.

Definition at line 773 of file pgridfunc.cpp.

◆ ComputeDivError()

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

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

See also
GridFunction::ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 402 of file pgridfunc.hpp.

◆ ComputeFlux()

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

Reimplemented from mfem::GridFunction.

Definition at line 1264 of file pgridfunc.cpp.

◆ ComputeGradError()

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

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

See also
GridFunction::ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 378 of file pgridfunc.hpp.

◆ ComputeH1Error() [1/2]

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

Returns either the H1-seminorm or the DG Face Jumps error or both depending on norm_type = 1, 2, 3

See also
GridFunction::ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t NU, int norm_type) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 447 of file pgridfunc.hpp.

◆ ComputeH1Error() [2/2]

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

Returns the error measured in H1-norm in parallel for H1 or L2 elements.

See also
GridFunction::ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 464 of file pgridfunc.hpp.

◆ ComputeHCurlError()

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

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

See also
GridFunction::ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 492 of file pgridfunc.hpp.

◆ ComputeHDivError()

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

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

See also
GridFunction::ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 477 of file pgridfunc.hpp.

◆ ComputeL1Error() [1/3]

real_t mfem::ParGridFunction::ComputeL1Error ( Coefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlineoverridevirtual

Returns ||u_ex - u_h||_L1 in parallel for scalar fields.

See also
GridFunction::ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 317 of file pgridfunc.hpp.

◆ ComputeL1Error() [2/3]

MFEM_DEPRECATED real_t mfem::ParGridFunction::ComputeL1Error ( Coefficient * exsol[],
const IntegrationRule * irs[] = NULL ) const
inlineoverridevirtual

Returns ||u_ex - u_h||_L1 in parallel for H1 or L2 elements.

See also
GridFunction::ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]) const for more detailed documentation.
Warning
While this function is nominally equivalent to ComputeLpError, with appropriate arguments, the returned errors may differ noticeably because ComputeLpError uses a higher order integration rule by default.
Deprecated
See ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]) const for the preferred implementation.

Reimplemented from mfem::GridFunction.

Definition at line 300 of file pgridfunc.hpp.

◆ ComputeL1Error() [3/3]

real_t mfem::ParGridFunction::ComputeL1Error ( VectorCoefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlineoverridevirtual

Returns ||u_ex - u_h||_L1 in parallel for vector fields.

See also
GridFunction::ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 326 of file pgridfunc.hpp.

◆ ComputeL2Error() [1/3]

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

Returns ||u_ex - u_h||_L2 in parallel for scalar fields.

See also
GridFunction::ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[], const Array<int> *elems) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 350 of file pgridfunc.hpp.

◆ ComputeL2Error() [2/3]

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

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

See also
GridFunction::ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[], const Array<int> *elems) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 336 of file pgridfunc.hpp.

◆ ComputeL2Error() [3/3]

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

Returns ||u_ex - u_h||_L2 in parallel for vector fields.

See also
GridFunction::ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[], const Array<int> *elems) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 364 of file pgridfunc.hpp.

◆ ComputeLpError() [1/2]

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

Returns ||u_ex - u_h||_Lp in parallel for scalar fields.

See also
GridFunction::ComputeLpError(const real_t p, Coefficient &exsol, Coefficient *weight, const IntegrationRule *irs[], const Array<int> *elems) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 549 of file pgridfunc.hpp.

◆ ComputeLpError() [2/2]

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

Returns ||u_ex - u_h||_Lp in parallel for vector fields.

See also
GridFunction::ComputeLpError(const real_t p, VectorCoefficient &exsol, Coefficient *weight, VectorCoefficient *v_weight, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 567 of file pgridfunc.hpp.

◆ ComputeMaxError() [1/3]

real_t mfem::ParGridFunction::ComputeMaxError ( Coefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlineoverridevirtual

Returns Max|u_ex - u_h| error in parallel for scalar fields.

See also
GridFunction::ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 524 of file pgridfunc.hpp.

◆ ComputeMaxError() [2/3]

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

Returns Max|u_ex - u_h| error in parallel for scalar or vector fields.

Note
This implementation of the max error of a vector field computes the max norm over vector components rather than the magnitude of the vector.
See also
GridFunction::ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 511 of file pgridfunc.hpp.

◆ ComputeMaxError() [3/3]

real_t mfem::ParGridFunction::ComputeMaxError ( VectorCoefficient & exsol,
const IntegrationRule * irs[] = NULL ) const
inlineoverridevirtual

Returns Max|u_ex - u_h| error in parallel for vector fields.

See also
GridFunction::ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]) const for more detailed documentation.

Reimplemented from mfem::GridFunction.

Definition at line 535 of file pgridfunc.hpp.

◆ CountElementsPerVDof()

void mfem::ParGridFunction::CountElementsPerVDof ( Array< int > & elem_per_vdof) const
overridevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 518 of file pgridfunc.cpp.

◆ Distribute() [1/2]

void mfem::ParGridFunction::Distribute ( const Vector & tv)
inline

Definition at line 181 of file pgridfunc.hpp.

◆ Distribute() [2/2]

void mfem::ParGridFunction::Distribute ( const Vector * tv)

Set the grid function on (all) dofs from a given vector on the true dofs, i.e. P tv.

Definition at line 142 of file pgridfunc.cpp.

◆ ExchangeFaceNbrData()

void mfem::ParGridFunction::ExchangeFaceNbrData ( )

Definition at line 215 of file pgridfunc.cpp.

◆ FaceNbrData() [1/2]

Vector & mfem::ParGridFunction::FaceNbrData ( )
inline

Definition at line 225 of file pgridfunc.hpp.

◆ FaceNbrData() [2/2]

const Vector & mfem::ParGridFunction::FaceNbrData ( ) const
inline

Definition at line 226 of file pgridfunc.hpp.

◆ GetDerivative()

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

Parallel version of GridFunction::GetDerivative(); see its documentation.

Definition at line 527 of file pgridfunc.cpp.

◆ GetElementDofValues()

void mfem::ParGridFunction::GetElementDofValues ( int el,
Vector & dof_vals ) const
overridevirtual

Sets the output vector dof_vals to the values of the degrees of freedom of element el. If el is greater than or equal to the number of local elements, it will be interpreted as a shifted index of a face neighbor element.

Reimplemented from mfem::GridFunction.

Definition at line 548 of file pgridfunc.cpp.

◆ GetSerialGridFunction() [1/2]

GridFunction mfem::ParGridFunction::GetSerialGridFunction ( int save_rank,
FiniteElementSpace & serial_fes ) const

Returns a GridFunction on MPI rank save_rank that does not have any duplication of vertices/nodes at processor boundaries.

The given serial_fes must be defined on the mesh returned by ParMesh::GetSerialMesh (with save_rank ranks), for example using the space belonging to the GridFunction obtained from ParGridFunction::GetSerialGridFunction(int,Mesh &) const.

Note
The returned GridFunction does not assume ownership of serial_fes.

Definition at line 992 of file pgridfunc.cpp.

◆ GetSerialGridFunction() [2/2]

GridFunction mfem::ParGridFunction::GetSerialGridFunction ( int save_rank,
Mesh & serial_mesh ) const

Returns a GridFunction on MPI rank save_rank that does not have any duplication of vertices/nodes at processor boundaries.

The serial_mesh is obtained using ParMesh::GetSerialMesh. Note that the save_rank must be the same as that used in ParMesh::GetSerialMesh.

Note
The returned GridFunction will own the newly created FiniteElementCollection and FiniteElementSpace objects.

Definition at line 1071 of file pgridfunc.cpp.

◆ GetTrueDofs() [1/2]

HypreParVector * mfem::ParGridFunction::GetTrueDofs ( ) const

Returns the true dofs in a new HypreParVector.

Definition at line 153 of file pgridfunc.cpp.

◆ GetTrueDofs() [2/2]

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

Extract the true-dofs from the GridFunction.

Definition at line 144 of file gridfunc.cpp.

◆ GetValue() [1/3]

real_t mfem::ParGridFunction::GetValue ( ElementTransformation & T)
inline

Definition at line 231 of file pgridfunc.hpp.

◆ GetValue() [2/3]

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

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

Reimplemented from mfem::GridFunction.

Definition at line 405 of file pgridfunc.cpp.

◆ GetValue() [3/3]

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

Return a scalar value from within the given element.

Reimplemented from mfem::GridFunction.

Definition at line 276 of file pgridfunc.cpp.

◆ GetVectorValue() [1/2]

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

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

Reimplemented from mfem::GridFunction.

Definition at line 454 of file pgridfunc.cpp.

◆ GetVectorValue() [2/2]

void mfem::ParGridFunction::GetVectorValue ( int i,
const IntegrationPoint & ip,
Vector & val ) const
overridevirtual

Return a vector value from within the given element.

Reimplemented from mfem::GridFunction.

Definition at line 348 of file pgridfunc.cpp.

◆ MakeRef() [1/6]

void mfem::ParGridFunction::MakeRef ( FiniteElementSpace * f,
real_t * v )
overridevirtual

Make the ParGridFunction reference external data on a new FiniteElementSpace.

This method changes the FiniteElementSpace associated with the ParGridFunction and sets the pointer v as external data in the ParGridFunction. The new space f is expected to be a ParFiniteElementSpace.

Reimplemented from mfem::GridFunction.

Definition at line 112 of file pgridfunc.cpp.

◆ MakeRef() [2/6]

void mfem::ParGridFunction::MakeRef ( FiniteElementSpace * f,
Vector & v,
int v_offset )
overridevirtual

Make the ParGridFunction reference external data on a new FiniteElementSpace.

This method changes the FiniteElementSpace associated with the ParGridFunction and sets the data of the Vector v (plus the v_offset) as external data in the ParGridFunction. The new space f is expected to be a ParFiniteElementSpace.

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

Reimplemented from mfem::GridFunction.

Definition at line 127 of file pgridfunc.cpp.

◆ MakeRef() [3/6]

void mfem::ParGridFunction::MakeRef ( ParFiniteElementSpace * f,
real_t * v )

Make the ParGridFunction reference external data on a new ParFiniteElementSpace.

This method changes the ParFiniteElementSpace associated with the ParGridFunction and sets the pointer v as external data in the ParGridFunction.

Definition at line 120 of file pgridfunc.cpp.

◆ MakeRef() [4/6]

void mfem::ParGridFunction::MakeRef ( ParFiniteElementSpace * f,
Vector & v,
int v_offset )

Make the ParGridFunction reference external data on a new ParFiniteElementSpace.

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

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

Definition at line 135 of file pgridfunc.cpp.

◆ MakeRef() [5/6]

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 1491 of file vector.hpp.

◆ MakeRef() [6/6]

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 1491 of file vector.hpp.

◆ operator=() [1/4]

ParGridFunction & mfem::ParGridFunction::operator= ( const HypreParVector & tv)
inline

Short semantic for Distribute()

Definition at line 189 of file pgridfunc.hpp.

◆ operator=() [2/4]

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

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

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

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

Definition at line 119 of file pgridfunc.hpp.

◆ operator=() [3/4]

ParGridFunction & mfem::ParGridFunction::operator= ( const Vector & v)
inline

Copy the data from a Vector to the ParGridFunction data.

Definition at line 127 of file pgridfunc.hpp.

◆ operator=() [4/4]

ParGridFunction & mfem::ParGridFunction::operator= ( real_t value)
inline

Assign constant values to the ParGridFunction data.

Definition at line 123 of file pgridfunc.hpp.

◆ ParallelAssemble() [1/3]

HypreParVector * mfem::ParGridFunction::ParallelAssemble ( ) const

Returns a new vector assembled on the true dofs.

Definition at line 208 of file pgridfunc.cpp.

◆ ParallelAssemble() [2/3]

void mfem::ParGridFunction::ParallelAssemble ( HypreParVector & tv) const

Returns the vector assembled on the true dofs.

Definition at line 203 of file pgridfunc.cpp.

◆ ParallelAssemble() [3/3]

void mfem::ParGridFunction::ParallelAssemble ( Vector & tv) const

Returns the vector assembled on the true dofs.

Definition at line 198 of file pgridfunc.cpp.

◆ ParallelAverage() [1/3]

HypreParVector * mfem::ParGridFunction::ParallelAverage ( ) const

Returns a new vector averaged on the true dofs.

Definition at line 174 of file pgridfunc.cpp.

◆ ParallelAverage() [2/3]

void mfem::ParGridFunction::ParallelAverage ( HypreParVector & tv) const

Returns the vector averaged on the true dofs.

Definition at line 167 of file pgridfunc.cpp.

◆ ParallelAverage() [3/3]

void mfem::ParGridFunction::ParallelAverage ( Vector & tv) const

Returns the vector averaged on the true dofs.

Definition at line 160 of file pgridfunc.cpp.

◆ ParallelProject() [1/3]

HypreParVector * mfem::ParGridFunction::ParallelProject ( ) const

Returns a new vector restricted to the true dofs.

Definition at line 191 of file pgridfunc.cpp.

◆ ParallelProject() [2/3]

void mfem::ParGridFunction::ParallelProject ( HypreParVector & tv) const

Returns the vector restricted to the true dofs.

Definition at line 186 of file pgridfunc.cpp.

◆ ParallelProject() [3/3]

void mfem::ParGridFunction::ParallelProject ( Vector & tv) const

Returns the vector restricted to the true dofs.

Definition at line 181 of file pgridfunc.cpp.

◆ ParFESpace()

ParFiniteElementSpace * mfem::ParGridFunction::ParFESpace ( ) const
inline

Definition at line 130 of file pgridfunc.hpp.

◆ ProjectBdrCoefficient() [1/4]

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

◆ ProjectBdrCoefficient() [2/4]

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

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 from mfem::GridFunction.

Definition at line 278 of file pgridfunc.hpp.

◆ ProjectBdrCoefficient() [3/4]

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

Definition at line 676 of file pgridfunc.cpp.

◆ ProjectBdrCoefficient() [4/4]

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

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

Reimplemented from mfem::GridFunction.

Definition at line 274 of file pgridfunc.hpp.

◆ ProjectBdrCoefficientTangent()

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

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 from mfem::GridFunction.

Definition at line 732 of file pgridfunc.cpp.

◆ ProjectCoefficient() [1/6]

void mfem::ParGridFunction::ProjectCoefficient ( Coefficient & coeff)
overridevirtual

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). For NURBS spaces these degrees of freedom are not available and L2 projection is resorted to as fallback.

Reimplemented from mfem::GridFunction.

Definition at line 567 of file pgridfunc.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 400 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 420 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). For NURBS spaces these degrees of freedom are not available and L2 projection is resorted to as fallback.

Definition at line 407 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 412 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 416 of file gridfunc.cpp.

◆ ProjectDiscCoefficient() [1/4]

void mfem::ParGridFunction::ProjectDiscCoefficient ( Coefficient & coeff,
AvgType type )
overridevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 633 of file pgridfunc.cpp.

◆ ProjectDiscCoefficient() [2/4]

void mfem::ParGridFunction::ProjectDiscCoefficient ( VectorCoefficient & coeff)
overridevirtual

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 from mfem::GridFunction.

Definition at line 589 of file pgridfunc.cpp.

◆ ProjectDiscCoefficient() [3/4]

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

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

◆ ProjectDiscCoefficient() [4/4]

void mfem::ParGridFunction::ProjectDiscCoefficient ( VectorCoefficient & coeff,
AvgType type )
overridevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 654 of file pgridfunc.cpp.

◆ ProlongateToMaxOrder()

std::unique_ptr< ParGridFunction > mfem::ParGridFunction::ProlongateToMaxOrder ( ) const

Return a GridFunction with the values of this, prolongated to the maximum order of all elements in the mesh.

Definition at line 1300 of file pgridfunc.cpp.

◆ Save() [1/3]

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

Save the local portion of the ParGridFunction. This differs from the serial GridFunction::Save in that it takes into account the signs of the local dofs.

Reimplemented from mfem::GridFunction.

Definition at line 1085 of file pgridfunc.cpp.

◆ Save() [2/3]

void mfem::ParGridFunction::Save ( const char * fname,
int precision = 16 ) const
overridevirtual

Save the ParGridFunction to files (one for each MPI rank). The files will be given suffixes according to the MPI rank. The given precision will be used for ASCII output.

Reimplemented from mfem::GridFunction.

Definition at line 956 of file pgridfunc.cpp.

◆ Save() [3/3]

void mfem::ParGridFunction::Save ( std::ostream & out) const
overridevirtual

Save the local portion of the ParGridFunction. This differs from the serial GridFunction::Save in that it takes into account the signs of the local dofs.

Reimplemented from mfem::GridFunction.

Definition at line 940 of file pgridfunc.cpp.

◆ SaveAsOne() [1/2]

void mfem::ParGridFunction::SaveAsOne ( const char * fname,
int precision = 16 ) const

Save the ParGridFunction to a single file (written using MPI rank 0). The given precision will be used for ASCII output.

Definition at line 966 of file pgridfunc.cpp.

◆ SaveAsOne() [2/2]

void mfem::ParGridFunction::SaveAsOne ( std::ostream & out = mfem::out) const

Merge the local grid functions.

Definition at line 1104 of file pgridfunc.cpp.

◆ SaveAsSerial()

void mfem::ParGridFunction::SaveAsSerial ( const char * fname,
int precision = 16,
int save_rank = 0 ) const

Write the serial GridFunction a single file (written using MPI rank 0). The given precision will be used for ASCII output.

Definition at line 978 of file pgridfunc.cpp.

◆ SetFromTrueDofs()

void mfem::ParGridFunction::SetFromTrueDofs ( const Vector & tv)
inlineoverridevirtual

Set the GridFunction from the given true-dof vector.

Reimplemented from mfem::GridFunction.

Definition at line 186 of file pgridfunc.hpp.

◆ SetSpace() [1/2]

void mfem::ParGridFunction::SetSpace ( FiniteElementSpace * f)
overridevirtual

Associate a new FiniteElementSpace with the ParGridFunction.

The ParGridFunction is resized using the SetSize() method. The new space f is expected to be a ParFiniteElementSpace.

Reimplemented from mfem::GridFunction.

Definition at line 97 of file pgridfunc.cpp.

◆ SetSpace() [2/2]

void mfem::ParGridFunction::SetSpace ( ParFiniteElementSpace * f)

Associate a new parallel space with the ParGridFunction.

Definition at line 105 of file pgridfunc.cpp.

◆ Update()

void mfem::ParGridFunction::Update ( )
overridevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 91 of file pgridfunc.cpp.

Member Data Documentation

◆ face_nbr_data

Vector mfem::ParGridFunction::face_nbr_data
protected

Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().

Definition at line 56 of file pgridfunc.hpp.

◆ pfes

ParFiniteElementSpace* mfem::ParGridFunction::pfes
protected

Points to the same object as fes.

Definition at line 52 of file pgridfunc.hpp.

◆ send_data

Vector mfem::ParGridFunction::send_data
protected

Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring processors.

Definition at line 61 of file pgridfunc.hpp.


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