15 #include "../config/config.hpp"
38 void SaveSTLTri(std::ostream &out,
double p1[],
double p2[],
double p3[]);
173 double _min = 0.0,
double _max = std::numeric_limits<double>::infinity());
230 int norm_type)
const;
246 exsol, NULL, NULL, irs);
275 int wcoef = 1,
int subdomain = -1);
315 virtual void Save(std::ostream &out)
const;
320 void SaveVTK(std::ostream &out,
const std::string &field_name,
int ref);
322 void SaveSTL(std::ostream &out,
int TimesToRefine = 1);
331 std::ostream &
operator<<(std::ostream &out,
const GridFunction &sol);
353 :
Vector(vdim_*qspace_->GetSize()),
361 :
Vector(qf_data, vdim_*qspace_->GetSize()),
433 void Save(std::ostream &out)
const;
437 std::ostream &
operator<<(std::ostream &out,
const QuadratureFunction &qf);
443 Vector &error_estimates,
444 Array<int> *aniso_flags = NULL,
445 int with_subdomains = 1);
449 GridFunction& gf1, GridFunction& gf2);
461 : n(_n), mesh_in(m), sol_in(s) { }
468 GridFunction *sol,
const int ny);
486 double *qf_data,
int vdim_)
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
int GetVDim() const
Get the vector dimension.
virtual ~GridFunction()
Destroys grid function.
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
QuadratureFunction(QuadratureSpace *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace.
Class for an integration rule - an Array of IntegrationPoint.
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Class for grid function - Vector with associated FE space.
double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Class used for extruding scalar GridFunctions.
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains)
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
void SetSize(int s)
Resize the vector to size s.
GridFunction * Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d, GridFunction *sol, const int ny)
Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
void MakeOwner(FiniteElementCollection *_fec)
Make the GridFunction the owner of 'fec' and 'fes'.
int GetSize()
Return the total number of quadrature points.
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh)
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
bool own_qspace
QuadratureSpace ownership flag.
void GetBdrValuesFrom(GridFunction &)
void GetGradient(ElementTransformation &tr, Vector &grad)
virtual ~QuadratureFunction()
Data type dense matrix using column-major storage.
Delta function coefficient.
void ImposeBounds(int i, const Vector &weights, const Vector &_lo, const Vector &_hi)
int vdim
Vector dimension.
void SetSpace(QuadratureSpace *qspace_, int vdim_=-1)
Change the QuadratureSpace and optionally the vector dimension.
bool OwnsSpace()
Get the QuadratureSpace ownership flag.
void GetNodalValues(int i, Array< double > &nval, int vdim=1) const
Returns the values in the vertices of i'th element for dimension vdim.
QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace, using the external data...
void GetVectorFieldNodalValues(Vector &val, int comp) const
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
const FiniteElementSpace * FESpace() const
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
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)
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
double GetDivergence(ElementTransformation &tr)
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
void SetOwnsSpace(bool own)
Set the QuadratureSpace ownership flag.
double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
GridFunction(const GridFunction &orig)
Copy constructor.
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
QuadratureFunction()
Create an empty QuadratureFunction.
void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
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...
void GetCurl(ElementTransformation &tr, Vector &curl)
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
FiniteElementCollection * OwnFEC()
virtual ~ExtrudeCoefficient()
FiniteElementSpace * FESpace()
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
Base class Coefficient that may optionally depend on time.
FiniteElementSpace * fes
FE space on which grid function lives.
void ProjectGridFunction(const GridFunction &src)
double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
QuadratureSpace * qspace
Associated QuadratureSpace.
FiniteElementCollection * fec
Used when the grid function is read from a file.
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.
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Class for integration point with weight.
void SaveSTL(std::ostream &out, int TimesToRefine=1)
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
void GetElementValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
void GetElementAverages(GridFunction &avgs)
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
void ProjectCoefficient(Coefficient &coeff)
void SetSpace(FiniteElementSpace *f)
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
ExtrudeCoefficient(Mesh *m, Coefficient &s, int _n)
GridFunction & operator=(double value)
Redefine '=' for GridFunction = constant.
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad)
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Class representing the storage layout of a QuadratureFunction.
long GetSequence() const
Return update counter (see Mesh::sequence)
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Class representing a function through its values (scalar or vector) at quadrature points...
void GetValuesFrom(GridFunction &)
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
void SaveVTK(std::ostream &out, const std::string &field_name, int ref)
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)