15 #include "../config/config.hpp" 17 #include "coefficient.hpp" 19 #ifdef MFEM_USE_ADIOS2 20 #include "../general/adios2stream.hpp" 49 void SaveSTLTri(std::ostream &
out,
double p1[],
double p2[],
double p3[]);
99 :
Vector(base, base_offset,
f->GetVSize())
132 MFEM_VERIFY(
t_vec.
Size() > 0,
"SetTrueVector() before GetTrueVector()");
221 int comp = 0,
Vector *tr = NULL)
const;
374 double min_ = 0.0,
double max_ =
infinity());
557 int norm_type)
const;
594 int norm_type,
const Array<int> *elems = NULL,
674 bool wcoef =
true,
int subdomain = -1);
736 virtual void Save(std::ostream &
out)
const;
740 virtual void Save(
const char *fname,
int precision=16)
const;
742 #ifdef MFEM_USE_ADIOS2 752 void SaveVTK(std::ostream &
out,
const std::string &field_name,
int ref);
756 void SaveSTL(std::ostream &
out,
int TimesToRefine = 1);
765 std::ostream &
operator<<(std::ostream &
out,
const GridFunction &sol);
783 : nu(nu_), type(type_) { }
794 std::ostream &
operator<<(std::ostream &
out,
const QuadratureFunction &qf);
800 Vector &error_estimates,
801 Array<int> *aniso_flags = NULL,
802 int with_subdomains = 1,
803 bool with_coeff =
false);
816 const Vector *midpoint=NULL);
834 FiniteElementSpace *ufes,
852 Vector &error_estimates,
853 bool subdomain_reconstruction =
true,
854 bool with_coeff =
false,
855 double tichonov_coeff = 0.0);
859 GridFunction& gf1, GridFunction& gf2);
871 : n(n_), mesh_in(m), sol_in(
s) { }
878 GridFunction *sol,
const int ny);
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual ~GridFunction()
Destroys grid function.
virtual double ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
virtual void CountElementsPerVDof(Array< int > &elem_per_vdof) const
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
Class used for extruding scalar GridFunctions.
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Base class for vector Coefficients that optionally depend on time and space.
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
void RestrictConforming()
GridFunction * Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d, GridFunction *sol, const int ny)
Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
virtual void ComputeElementL1Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
virtual void ComputeElementLpErrors(const double p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
void GetGradient(ElementTransformation &tr, Vector &grad) const
Gradient of a scalar function at a quadrature point.
int Size() const
Returns the size of the vector.
void GetElementAverages(GridFunction &avgs) const
Data type dense matrix using column-major storage.
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
void AccumulateAndCountBdrTangentValues(VectorCoefficient &vcoeff, Array< int > &bdr_attr, Array< int > &values_counter)
Vector & GetTrueVector()
Read and write access to the (optional) internal true-dof Vector.
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
double Eval(double h, int p) const
ExtrudeCoefficient(Mesh *m, Coefficient &s, int n_)
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
virtual void ComputeElementMaxErrors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
void GetHessians(int i, const IntegrationRule &ir, DenseMatrix &hess, int vdim=1) const
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
void ImposeBounds(int i, const Vector &weights, const Vector &lo_, const Vector &hi_)
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
std::function< double(const Vector &)> f(double mass_coeff)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
GridFunction(FiniteElementSpace *f, double *data)
Construct a GridFunction using previously allocated array data.
void GetDerivative(int comp, int der_comp, GridFunction &der)
Compute a certain derivative of a function's component. Derivatives of the function are computed at t...
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
void GetCurl(ElementTransformation &tr, Vector &curl) const
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
GridFunction(const GridFunction &orig)
Copy constructor. The internal true-dof vector t_vec is not copied.
JumpScaling(double nu_=1.0, JumpScalingType type_=CONSTANT)
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
void BoundingBox(const Array< int > &patch, FiniteElementSpace *ufes, int order, Vector &xmin, Vector &xmax, double &angle, Vector &midpoint, int iface)
Defines the bounding box for the face patches used by NewZZErorrEstimator.
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.
double LSZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, Vector &error_estimates, bool subdomain_reconstruction, bool with_coeff, double tichonov_coeff)
A ‘‘true’’ ZZ error estimator that uses face-based patches for flux reconstruction.
void AccumulateAndCountDerivativeValues(int comp, int der_comp, GridFunction &der, Array< int > &zones_per_dof)
Used for the serial and parallel implementations of the GetDerivative() method; see its documentation...
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 ...
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
double GetDivergence(ElementTransformation &tr) const
FiniteElementCollection * OwnFEC()
virtual ~ExtrudeCoefficient()
FiniteElementSpace * FESpace()
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
void GetNodalValues(int i, Array< double > &nval, int vdim=1) const
Returns the values in the vertices of i'th element for dimension vdim.
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
virtual void GetElementDofValues(int el, Vector &dof_vals) const
double p(const Vector &x, double t)
virtual double 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.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual double ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, class JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
FiniteElementSpace * fes
FE space on which the grid function lives. Owned if fec is not NULL.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
void GetVectorFieldNodalValues(Vector &val, int comp) const
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.
FiniteElementCollection * fec
Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner()...
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh) const
Compute the vector gradient with respect to the reference element variable.
void TensorProductLegendre(int dim, int order, const Vector &x_in, const Vector &xmax, const Vector &xmin, Vector &poly, double angle, const Vector *midpoint)
Defines the global tensor product polynomial space used by NewZZErorrEstimator.
virtual void ComputeElementL2Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
double ComputeElementLpDistance(double p, int i, GridFunction &gf1, GridFunction &gf2)
Compute the Lp distance between two grid functions on the given element.
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
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
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains, bool with_coeff)
Class for integration point with weight.
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 w...
virtual void ComputeElementL1Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
void LegacyNCReorder()
Loading helper.
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
GridFunction(FiniteElementSpace *f, Vector &base, int base_offset=0)
Construct a GridFunction using previously allocated Vector base starting at the given offset...
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
void GetLaplacians(int i, const IntegrationRule &ir, Vector &laps, int vdim=1) const
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
double u(const Vector &xvec)
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
virtual void ComputeElementMaxErrors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
const FiniteElementSpace * FESpace() const
void GetValuesFrom(const GridFunction &orig_func)
void GetBdrValuesFrom(const GridFunction &orig_func)
void AccumulateAndCountBdrValues(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr, Array< int > &values_counter)
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 ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
virtual void ComputeElementL2Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.