15 #include "../config/config.hpp"
43 void SaveSTLTri(std::ostream &
out,
double p1[],
double p2[],
double p3[]);
193 double _min = 0.0,
double _max =
infinity());
276 int norm_type)
const;
319 int wcoef = 1,
int subdomain = -1);
380 virtual void Save(std::ostream &
out)
const;
385 void SaveVTK(std::ostream &
out,
const std::string &field_name,
int ref);
387 void SaveSTL(std::ostream &
out,
int TimesToRefine = 1);
396 std::ostream &
operator<<(std::ostream &
out,
const GridFunction &sol);
418 :
Vector(vdim_*qspace_->GetSize()),
426 :
Vector(qf_data, vdim_*qspace_->GetSize()),
529 void Save(std::ostream &
out)
const;
533 std::ostream &
operator<<(std::ostream &
out,
const QuadratureFunction &qf);
539 Vector &error_estimates,
540 Array<int> *aniso_flags = NULL,
541 int with_subdomains = 1);
545 GridFunction& gf1, GridFunction& gf2);
557 : n(_n), mesh_in(m), sol_in(s) { }
564 GridFunction *sol,
const int ny);
582 double *qf_data,
int vdim_)
607 for (
int i = 0; i<values.
Size(); i++)
627 for (
int j = 0; j<sl_size; j++)
628 for (
int i = 0; i<
vdim; i++)
630 values(i,j) = *(q++);
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 GetElementAverages(GridFunction &avgs) const
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Class for grid function - Vector with associated FE space.
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Class used for extruding scalar GridFunctions.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
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'.
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.
int GetSize() const
Return the total number of quadrature points.
virtual ~QuadratureFunction()
Data type dense matrix using column-major storage.
Delta function coefficient.
int Size() const
Returns the size of the vector.
void ImposeBounds(int i, const Vector &weights, const Vector &_lo, const Vector &_hi)
int vdim
Vector dimension.
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh) const
Vector & GetTrueVector()
Read and write access to the (optional) internal true-dof Vector.
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
void SetSpace(QuadratureSpace *qspace_, int vdim_=-1)
Change the QuadratureSpace and optionally the vector dimension.
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
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)
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 SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad) 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.
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
GridFunction(const GridFunction &orig)
Copy constructor.
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
QuadratureFunction()
Create an empty QuadratureFunction.
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.
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 GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
FiniteElementCollection * OwnFEC()
virtual ~ExtrudeCoefficient()
FiniteElementSpace * FESpace()
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
QuadratureFunction & operator=(double value)
Redefine '=' for QuadratureFunction = constant.
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Base class Coefficient that may optionally depend on time.
FiniteElementSpace * fes
FE space on which grid function lives.
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
QuadratureSpace * qspace
Associated QuadratureSpace.
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
FiniteElementCollection * fec
Used when the grid function is read from a file.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
double GetDivergence(ElementTransformation &tr) 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.
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 ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
virtual void ProjectCoefficient(Coefficient &coeff)
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.
virtual 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 GetCurl(ElementTransformation &tr, Vector &curl) const
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Class representing the storage layout of a QuadratureFunction.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
long GetSequence() const
Return update counter (see Mesh::sequence)
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(const GridFunction &orig_func)
void GetBdrValuesFrom(const GridFunction &orig_func)
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 GetGradient(ElementTransformation &tr, Vector &grad) const
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)