MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
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. More...
 
 ParGridFunction (ParFiniteElementSpace *pf)
 
 ParGridFunction (ParFiniteElementSpace *pf, double *data)
 Construct a ParGridFunction using previously allocated array data. More...
 
 ParGridFunction (ParFiniteElementSpace *pf, GridFunction *gf)
 Construct a ParGridFunction using a GridFunction as external data. More...
 
 ParGridFunction (ParFiniteElementSpace *pf, HypreParVector *tv)
 Creates grid function on (all) dofs from a given vector on the true dofs, i.e. P tv. More...
 
 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. More...
 
 ParGridFunction (ParMesh *pmesh, std::istream &input)
 Construct a ParGridFunction on a given ParMesh, pmesh, reading from an std::istream. More...
 
ParGridFunctionoperator= (const ParGridFunction &rhs)
 Copy assignment. Only the data of the base class Vector is copied. More...
 
ParGridFunctionoperator= (double value)
 Assign constant values to the ParGridFunction data. More...
 
ParGridFunctionoperator= (const Vector &v)
 Copy the data from a Vector to the ParGridFunction data. More...
 
ParFiniteElementSpaceParFESpace () const
 
virtual void Update ()
 Transform by the Space UpdateMatrix (e.g., on Mesh change). More...
 
virtual void SetSpace (FiniteElementSpace *f)
 Associate a new FiniteElementSpace with the ParGridFunction. More...
 
void SetSpace (ParFiniteElementSpace *f)
 Associate a new parallel space with the ParGridFunction. More...
 
virtual void MakeRef (FiniteElementSpace *f, double *v)
 Make the ParGridFunction reference external data on a new FiniteElementSpace. More...
 
void MakeRef (ParFiniteElementSpace *f, double *v)
 Make the ParGridFunction reference external data on a new ParFiniteElementSpace. More...
 
virtual void MakeRef (FiniteElementSpace *f, Vector &v, int v_offset)
 Make the ParGridFunction reference external data on a new FiniteElementSpace. More...
 
void MakeRef (ParFiniteElementSpace *f, Vector &v, int v_offset)
 Make the ParGridFunction reference external data on a new ParFiniteElementSpace. More...
 
void Distribute (const Vector *tv)
 
void Distribute (const Vector &tv)
 
void AddDistribute (double a, const Vector *tv)
 
void AddDistribute (double a, const Vector &tv)
 
virtual void SetFromTrueDofs (const Vector &tv)
 Set the GridFunction from the given true-dof vector. More...
 
ParGridFunctionoperator= (const HypreParVector &tv)
 Short semantic for Distribute() More...
 
HypreParVectorGetTrueDofs () const
 Returns the true dofs in a new HypreParVector. More...
 
void ParallelAverage (Vector &tv) const
 Returns the vector averaged on the true dofs. More...
 
void ParallelAverage (HypreParVector &tv) const
 Returns the vector averaged on the true dofs. More...
 
HypreParVectorParallelAverage () const
 Returns a new vector averaged on the true dofs. More...
 
void ParallelProject (Vector &tv) const
 Returns the vector restricted to the true dofs. More...
 
void ParallelProject (HypreParVector &tv) const
 Returns the vector restricted to the true dofs. More...
 
HypreParVectorParallelProject () const
 Returns a new vector restricted to the true dofs. More...
 
void ParallelAssemble (Vector &tv) const
 Returns the vector assembled on the true dofs. More...
 
void ParallelAssemble (HypreParVector &tv) const
 Returns the vector assembled on the true dofs. More...
 
HypreParVectorParallelAssemble () const
 Returns a new vector assembled on the true dofs. More...
 
void ExchangeFaceNbrData ()
 
VectorFaceNbrData ()
 
const VectorFaceNbrData () const
 
virtual double GetValue (int i, const IntegrationPoint &ip, int vdim=1) const
 
double GetValue (ElementTransformation &T)
 
virtual double GetValue (ElementTransformation &T, const IntegrationPoint &ip, int comp=0, Vector *tr=NULL) const
 
virtual void GetVectorValue (int i, const IntegrationPoint &ip, Vector &val) const
 
virtual void GetVectorValue (ElementTransformation &T, const IntegrationPoint &ip, Vector &val, Vector *tr=NULL) const
 
virtual void ProjectCoefficient (Coefficient &coeff)
 
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 &vcoeff, 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...
 
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...
 
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 ComputeL1Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeL1Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
virtual double ComputeL1Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
 
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, double Nu, const IntegrationRule *irs[]=NULL) const
 Returns the Face Jumps error for L2 elements. 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 H(div)-norm for RT elements. More...
 
virtual double ComputeHCurlError (VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
 Returns the error measured 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 ComputeLpError (const double p, Coefficient &exsol, Coefficient *weight=NULL, 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 ComputeFlux (BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
 
virtual void Save (std::ostream &out) const
 
virtual void Save (adios2stream &out, const std::string &variable_name, const adios2stream::data_type type=adios2stream::data_type::point_data) const
 
void SaveAsOne (std::ostream &out=mfem::out)
 Merge the local grid functions. More...
 
virtual ~ParGridFunction ()
 
- Public Member Functions inherited from mfem::GridFunction
 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. If all dofs are true, then tv will be set to point to the data of *this. More...
 
void SetTrueVector ()
 Shortcut for calling GetTrueDofs() with GetTrueVector() as argument. 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
 
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...
 
void ProjectCoefficient (Coefficient &coeff, Array< int > &dofs, int vd=0)
 
void ProjectCoefficient (VectorCoefficient &vcoeff)
 
void ProjectCoefficient (VectorCoefficient &vcoeff, Array< int > &dofs)
 
void ProjectCoefficient (Coefficient *coeff[])
 
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...
 
void ProjectBdrCoefficientNormal (VectorCoefficient &vcoeff, Array< int > &bdr_attr)
 
virtual double ComputeW11Error (Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=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 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
 
GridFunctionoperator= (double value)
 Redefine '=' for GridFunction = constant. More...
 
GridFunctionoperator= (const Vector &v)
 Copy the data from v. More...
 
FiniteElementSpaceFESpace ()
 
const FiniteElementSpaceFESpace () const
 
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...
 
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...
 
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 ()
 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...
 
void UseDevice (bool use_dev) const
 Enable execution of Vector operations using the mfem::Device. More...
 
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...
 
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...
 
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)
 Update the memory location of the vector to match v. More...
 
void SyncAliasMemory (const Vector &v)
 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/= (double c)
 
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 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...
 
const double * Read (bool on_dev=true) const
 Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev). More...
 
const double * HostRead () const
 Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false). More...
 
double * Write (bool on_dev=true)
 Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev). More...
 
double * HostWrite ()
 Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false). More...
 
double * ReadWrite (bool on_dev=true)
 Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev). More...
 
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 ProjectBdrCoefficient (Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
 
- Protected Member Functions inherited from mfem::GridFunction
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 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

ParFiniteElementSpacepfes
 Points to the same object as fes. More...
 
Vector face_nbr_data
 Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData(). More...
 
Vector send_data
 Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring processors. More...
 
- Protected Attributes inherited from mfem::GridFunction
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 sequence
 
Vector t_vec
 
- Protected Attributes inherited from mfem::Vector
Memory< double > data
 
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 32 of file pgridfunc.hpp.

Constructor & Destructor Documentation

mfem::ParGridFunction::ParGridFunction ( )
inline

Definition at line 50 of file pgridfunc.hpp.

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

Copy constructor. The internal vector face_nbr_data is not copied.

Definition at line 53 of file pgridfunc.hpp.

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

Definition at line 56 of file pgridfunc.hpp.

mfem::ParGridFunction::ParGridFunction ( ParFiniteElementSpace pf,
double *  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 65 of file pgridfunc.hpp.

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.

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.

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.

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 71 of file pgridfunc.cpp.

virtual mfem::ParGridFunction::~ParGridFunction ( )
inlinevirtual

Definition at line 418 of file pgridfunc.hpp.

Member Function Documentation

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

Definition at line 138 of file pgridfunc.cpp.

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

Definition at line 161 of file pgridfunc.hpp.

virtual double mfem::ParGridFunction::ComputeCurlError ( VectorCoefficient excurl,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 295 of file pgridfunc.hpp.

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

Returns the Face Jumps error for L2 elements.

Reimplemented from mfem::GridFunction.

Definition at line 658 of file pgridfunc.cpp.

virtual double mfem::ParGridFunction::ComputeDivError ( Coefficient exdiv,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 303 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 1024 of file pgridfunc.cpp.

virtual double mfem::ParGridFunction::ComputeGradError ( VectorCoefficient exgrad,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 287 of file pgridfunc.hpp.

virtual double mfem::ParGridFunction::ComputeH1Error ( Coefficient exsol,
VectorCoefficient exgrad,
Coefficient ell_coef,
double  Nu,
int  norm_type 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 319 of file pgridfunc.hpp.

virtual double mfem::ParGridFunction::ComputeH1Error ( Coefficient exsol,
VectorCoefficient exgrad,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 331 of file pgridfunc.hpp.

virtual double mfem::ParGridFunction::ComputeHCurlError ( VectorCoefficient exsol,
VectorCoefficient excurl,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 348 of file pgridfunc.hpp.

virtual double mfem::ParGridFunction::ComputeHDivError ( VectorCoefficient exsol,
Coefficient exdiv,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Reimplemented from mfem::GridFunction.

Definition at line 339 of file pgridfunc.hpp.

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

Definition at line 252 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 259 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 263 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 267 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 274 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 278 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 377 of file pgridfunc.hpp.

virtual double mfem::ParGridFunction::ComputeLpError ( const double  p,
VectorCoefficient exsol,
Coefficient weight = NULL,
VectorCoefficient v_weight = NULL,
const IntegrationRule irs[] = NULL 
) const
inlinevirtual

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

Definition at line 388 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 357 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 365 of file pgridfunc.hpp.

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

Reimplemented from mfem::GridFunction.

Definition at line 371 of file pgridfunc.hpp.

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 132 of file pgridfunc.cpp.

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

Definition at line 159 of file pgridfunc.hpp.

void mfem::ParGridFunction::ExchangeFaceNbrData ( )

Definition at line 205 of file pgridfunc.cpp.

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

Definition at line 203 of file pgridfunc.hpp.

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

Definition at line 204 of file pgridfunc.hpp.

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

Returns the true dofs in a new HypreParVector.

Definition at line 143 of file pgridfunc.cpp.

double mfem::ParGridFunction::GetValue ( int  i,
const IntegrationPoint ip,
int  vdim = 1 
) const
virtual

Return a scalar value from within the given element.

Reimplemented from mfem::GridFunction.

Definition at line 264 of file pgridfunc.cpp.

double mfem::ParGridFunction::GetValue ( ElementTransformation T)
inline

Definition at line 209 of file pgridfunc.hpp.

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

Definition at line 372 of file pgridfunc.cpp.

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

Return a vector value from within the given element.

Reimplemented from mfem::GridFunction.

Definition at line 321 of file pgridfunc.cpp.

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

Definition at line 415 of file pgridfunc.cpp.

void mfem::ParGridFunction::MakeRef ( FiniteElementSpace f,
double *  v 
)
virtual

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 102 of file pgridfunc.cpp.

void mfem::ParGridFunction::MakeRef ( ParFiniteElementSpace f,
double *  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 110 of file pgridfunc.cpp.

void mfem::ParGridFunction::MakeRef ( FiniteElementSpace f,
Vector v,
int  v_offset 
)
virtual

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 117 of file pgridfunc.cpp.

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 125 of file pgridfunc.cpp.

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 97 of file pgridfunc.hpp.

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

Assign constant values to the ParGridFunction data.

Definition at line 101 of file pgridfunc.hpp.

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

Copy the data from a Vector to the ParGridFunction data.

Definition at line 105 of file pgridfunc.hpp.

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

Short semantic for Distribute()

Definition at line 167 of file pgridfunc.hpp.

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

Returns the vector assembled on the true dofs.

Definition at line 188 of file pgridfunc.cpp.

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

Returns the vector assembled on the true dofs.

Definition at line 193 of file pgridfunc.cpp.

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

Returns a new vector assembled on the true dofs.

Definition at line 198 of file pgridfunc.cpp.

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

Returns the vector averaged on the true dofs.

Definition at line 150 of file pgridfunc.cpp.

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

Returns the vector averaged on the true dofs.

Definition at line 157 of file pgridfunc.cpp.

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

Returns a new vector averaged on the true dofs.

Definition at line 164 of file pgridfunc.cpp.

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

Returns the vector restricted to the true dofs.

Definition at line 171 of file pgridfunc.cpp.

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

Returns the vector restricted to the true dofs.

Definition at line 176 of file pgridfunc.cpp.

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

Returns a new vector restricted to the true dofs.

Definition at line 181 of file pgridfunc.cpp.

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

Definition at line 108 of file pgridfunc.hpp.

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

Definition at line 582 of file pgridfunc.cpp.

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

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 240 of file pgridfunc.hpp.

virtual void mfem::ParGridFunction::ProjectBdrCoefficient ( Coefficient coeff[],
Array< int > &  attr 
)
inlinevirtual

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 245 of file pgridfunc.hpp.

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

Definition at line 620 of file pgridfunc.cpp.

void mfem::ParGridFunction::ProjectCoefficient ( Coefficient coeff)
virtual

Reimplemented from mfem::GridFunction.

Definition at line 474 of file pgridfunc.cpp.

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

Definition at line 495 of file pgridfunc.cpp.

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

Definition at line 539 of file pgridfunc.cpp.

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

Definition at line 560 of file pgridfunc.cpp.

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

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 819 of file pgridfunc.cpp.

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

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 836 of file pgridfunc.cpp.

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

Merge the local grid functions.

Definition at line 855 of file pgridfunc.cpp.

virtual void mfem::ParGridFunction::SetFromTrueDofs ( const Vector tv)
inlinevirtual

Set the GridFunction from the given true-dof vector.

Reimplemented from mfem::GridFunction.

Definition at line 164 of file pgridfunc.hpp.

void mfem::ParGridFunction::SetSpace ( FiniteElementSpace f)
virtual

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 87 of file pgridfunc.cpp.

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

Associate a new parallel space with the ParGridFunction.

Definition at line 95 of file pgridfunc.cpp.

void mfem::ParGridFunction::Update ( )
virtual

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

Reimplemented from mfem::GridFunction.

Definition at line 81 of file pgridfunc.cpp.

Member Data Documentation

Vector mfem::ParGridFunction::face_nbr_data
protected

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

Definition at line 39 of file pgridfunc.hpp.

ParFiniteElementSpace* mfem::ParGridFunction::pfes
protected

Points to the same object as fes.

Definition at line 35 of file pgridfunc.hpp.

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 44 of file pgridfunc.hpp.


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