104 const int *partitioning = NULL);
230 int vdim = 1)
const override;
236 int comp = 0,
Vector *tr = NULL)
const override;
239 Vector &val)
const override;
303#pragma GCC diagnostic push
304#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
308#pragma GCC diagnostic pop
338 const Array<int> *elems = NULL)
const override
352 const Array<int> *elems = NULL)
const override
366 const Array<int> *elems = NULL)
const override
449 int norm_type)
const override
552 const Array<int> *elems = NULL)
const override
578 bool wcoef =
true,
int subdomain = -1)
override;
583 void Save(std::ostream &
out)
const override;
587 void SaveAsOne(
const char *fname,
int precision=16)
const;
592 void Save(
const char *fname,
int precision=16)
const override;
619 void SaveAsSerial(
const char *fname,
int precision=16,
int save_rank=0)
const;
621#ifdef MFEM_USE_ADIOS2
651 Vector &errors,
int norm_p = 2,
real_t solver_tol = 1e-12,
652 int solver_max_it = 200);
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Class for grid function - Vector with associated FE space.
virtual real_t ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(curl)-norm for ND elements.
virtual real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t Nu, int norm_type) const
virtual real_t ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Returns Max|u_ex - u_h| error for H1 or L2 elements.
virtual void MakeRef(FiniteElementSpace *f, real_t *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||exsol - u_h||_L2 for scalar or vector H1 or L2 elements.
virtual real_t ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Returns ||u_ex - u_h||_L1 for H1 or L2 elements.
virtual real_t ComputeLpError(const real_t p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||u_ex - u_h||_Lp for H1 or L2 elements.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
virtual real_t ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
virtual real_t ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(div)-norm for RT elements.
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
virtual real_t ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
virtual real_t ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
void ProjectBdrCoefficient(Coefficient &coeff, const Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Wrapper for hypre's parallel vector class.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
Abstract parallel finite element space.
Class for parallel grid function.
void CountElementsPerVDof(Array< int > &elem_per_vdof) const override
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
void ExchangeFaceNbrData()
void GetDerivative(int comp, int der_comp, ParGridFunction &der) const
Parallel version of GridFunction::GetDerivative(); see its documentation.
real_t ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns Max|u_ex - u_h| error in parallel for vector fields.
ParGridFunction(ParFiniteElementSpace *pf, Vector &base, int base_offset=0)
Construct a ParGridFunction using previously allocated Vector base starting at the given offset,...
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t Nu, int norm_type) const override
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.
void Save(std::ostream &out) const override
real_t ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns ||u_ex - u_h||_L1 in parallel for scalar fields.
ParGridFunction(ParFiniteElementSpace *pf, real_t *data)
Construct a ParGridFunction using previously allocated array data.
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.
void ProjectDiscCoefficient(VectorCoefficient &coeff) override
Project a discontinuous vector coefficient as a grid function on a continuous finite element space....
real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) 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.
void ProjectBdrCoefficient(VectorCoefficient &vcoeff, const Array< int > &attr) override
Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with ...
const Vector & FaceNbrData() const
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
ParGridFunction(const ParGridFunction &orig)
Copy constructor. The internal vector face_nbr_data is not copied.
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.
void AddDistribute(real_t a, const Vector &tv)
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 SetFromTrueDofs(const Vector &tv) override
Set the GridFunction from the given true-dof vector.
void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1) override
ParGridFunction & operator=(real_t value)
Assign constant values to the ParGridFunction data.
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.
ParGridFunction & operator=(const HypreParVector &tv)
Short semantic for Distribute()
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
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 ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const override
Returns the error measured H(div)-norm in parallel for RT elements.
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 bound...
ParFiniteElementSpace * pfes
Points to the same object as fes.
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
void Distribute(const Vector &tv)
real_t ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns ||u_ex - u_h||_L1 in parallel for vector fields.
ParFiniteElementSpace * ParFESpace() const
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.
virtual ~ParGridFunction()=default
void AddDistribute(real_t a, const Vector *tv)
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const override
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 ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const override
Returns ||div u_ex - div u_h||_L2 in parallel for RT elements.
void SaveAsSerial(const char *fname, int precision=16, int save_rank=0) const
ParGridFunction(ParFiniteElementSpace *pf)
Vector face_nbr_data
Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void SaveAsOne(const char *fname, int precision=16) const
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 ...
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
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.
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
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.
void Distribute(const Vector *tv)
std::unique_ptr< ParGridFunction > ProlongateToMaxOrder() const
Return a GridFunction with the values of this, prolongated to the maximum order of all elements in th...
void GetElementDofValues(int el, Vector &dof_vals) const override
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
real_t ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const override
Returns ||curl u_ex - curl u_h||_L2 in parallel for ND elements.
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, const Array< int > &bdr_attr) override
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
real_t GetValue(ElementTransformation &T)
ParGridFunction & operator=(const ParGridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Class for parallel meshes.
Base class for vector Coefficients that optionally depend on time and space.
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
real_t L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, real_t solver_tol, int solver_max_it)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
real_t GlobalLpNorm(const real_t p, real_t loc_norm, MPI_Comm comm)
Compute a global Lp norm from the local Lp norms computed by each processor.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
real_t p(const Vector &x, real_t t)