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);
226 return GlobalLpNorm(std::numeric_limits<double>::infinity(),
242 exsol, NULL, NULL, irs);
262 p, exsol, weight, v_weight, irs),
pfes->
GetComm());
267 int wcoef = 1,
int subdomain = -1);
272 virtual void Save(std::ostream &
out)
const;
287 const ParGridFunction &x,
288 ParFiniteElementSpace &smooth_flux_fes,
289 ParFiniteElementSpace &flux_fes,
290 Vector &errors,
int norm_p = 2,
double solver_tol = 1e-12,
291 int solver_max_it = 200);
295 #endif // MFEM_USE_MPI
virtual double ComputeLpError(const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
void ExchangeFaceNbrData()
void AddDistribute(double a, const Vector &tv)
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
virtual ~ParGridFunction()
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
ParFiniteElementSpace * pfes
ParGridFunction & operator=(const HypreParVector &tv)
Short semantic for Distribute.
void SaveAsOne(std::ostream &out=mfem::out)
Merge the local grid functions.
ParGridFunction(ParFiniteElementSpace *pf)
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
virtual void Save(std::ostream &out) const
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)
void Distribute(const Vector &tv)
double GetValue(ElementTransformation &T)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
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.
const Vector & FaceNbrData() const
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
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 th...
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
Wrapper for hypre's parallel vector class.
virtual double ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
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 ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Base class Coefficient that may optionally depend on time.
virtual double ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
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 ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Class for integration point with weight.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void AddDistribute(double a, const Vector *tv)
virtual void ProjectCoefficient(Coefficient &coeff)
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
virtual double ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
GridFunction & operator=(double value)
Redefine '=' for GridFunction = constant.
Class for parallel grid function.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Class for parallel meshes.
ParGridFunction(ParFiniteElementSpace *pf, double *data)
Construct a ParGridFunction using previously allocated array data.
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
ParFiniteElementSpace * ParFESpace() const