12 #ifndef MFEM_PGRIDFUNC 13 #define MFEM_PGRIDFUNC 15 #include "../config/config.hpp" 19 #include "../general/globals.hpp" 29 double GlobalLpNorm(
const double p,
double loc_norm, MPI_Comm comm);
87 const int *partitioning = NULL);
215 {
return GetValue(T.ElementNo, T.GetIntPoint()); }
219 int comp = 0,
Vector *tr = NULL)
const;
420 bool wcoef =
true,
int subdomain = -1);
425 virtual void Save(std::ostream &
out)
const;
429 void SaveAsOne(
const char *fname,
int precision=16)
const;
434 virtual void Save(
const char *fname,
int precision=16)
const;
436 #ifdef MFEM_USE_ADIOS2 458 const ParGridFunction &x,
459 ParFiniteElementSpace &smooth_flux_fes,
460 ParFiniteElementSpace &flux_fes,
461 Vector &errors,
int norm_p = 2,
double solver_tol = 1e-12,
462 int solver_max_it = 200);
466 #endif // MFEM_USE_MPI virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
ParGridFunction(ParFiniteElementSpace *pf, Vector &base, int base_offset=0)
Construct a ParGridFunction using previously allocated Vector base starting at the given offset...
virtual void GetElementDofValues(int el, Vector &dof_vals) const
virtual double ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
Class for an integration rule - an Array of IntegrationPoint.
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Class for grid function - Vector with associated FE space.
void GetDerivative(int comp, int der_comp, ParGridFunction &der)
Parallel version of GridFunction::GetDerivative(); see its documentation.
void ExchangeFaceNbrData()
void AddDistribute(double a, const Vector &tv)
virtual double ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
ParGridFunction(const ParGridFunction &orig)
Copy constructor. The internal vector face_nbr_data is not copied.
Base class for vector Coefficients that optionally depend on time and space.
void SaveAsOne(const char *fname, int precision=16) const
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
virtual ~ParGridFunction()
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
ParFiniteElementSpace * pfes
Points to the same object as fes.
ParGridFunction & operator=(const HypreParVector &tv)
Short semantic for Distribute()
ParGridFunction(ParFiniteElementSpace *pf)
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
void Distribute(const Vector &tv)
virtual double ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
virtual void ProjectBdrCoefficient(VectorCoefficient &vcoeff, Array< int > &attr)
Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with ...
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual double ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
double GetValue(ElementTransformation &T)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
const Vector & FaceNbrData() const
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
double f(const Vector &xvec)
ParGridFunction & operator=(double value)
Assign constant values to the ParGridFunction data.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
virtual double ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(div)-norm for RT elements.
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured H(curl)-norm for ND elements.
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, const Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
Wrapper for hypre's parallel vector class.
double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm)
Compute a global Lp norm from the local Lp norms computed by each processor.
virtual double ComputeL1Error(Coefficient &exsol, 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.
double p(const Vector &x, double t)
ParGridFunction & operator=(const ParGridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
ParFiniteElementSpace * ParFESpace() const
virtual double ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
Returns the Face Jumps error for L2 elements.
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, double solver_tol, int solver_max_it)
void Distribute(const Vector *tv)
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
virtual double ComputeLpError(const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
virtual void Save(std::ostream &out) const
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(curl)-norm for ND elements.
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array< int > &attr)
Project a set of Coefficients on the components of the GridFunction, modifying only DOFs on the bound...
Vector face_nbr_data
Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
Class for integration point with weight.
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
void AddDistribute(double a, const Vector *tv)
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
virtual double ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||u_ex - u_h||_L2 for H1 or L2 elements.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Class for parallel grid function.
Class for parallel meshes.
ParGridFunction(ParFiniteElementSpace *pf, double *data)
Construct a ParGridFunction using previously allocated array data.
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const