21 for (
int i =
sequence.Size()-1; i >= 0; i--)
36 case NONE:
if (last) {
return NONE; }
goto next_step;
46 for (
int i = 0; i <
sequence.Size(); i++)
94 const int NE = mesh.
GetNE();
96 MFEM_ASSERT(local_err.
Size() == NE,
"invalid size of local_err");
98 const double total_err =
GetNorm(local_err, mesh);
112 for (
int el = 0; el < NE; el++)
123 if (aniso_flags.
Size() > 0)
169 MFEM_VERIFY(max_it > 0,
"max_it must be strictly positive")
197 int order_quad = 2*
order + 3;
205 for (
int i = 0; i < max_it; i++)
208 int NE = mesh.
GetNE();
210 double norm_of_coeff = 0.0;
225 double av_norm_of_coeff = norm_of_coeff / sqrt(globalNE);
228 Vector element_norms_of_fine_scale(NE);
240 for (
int j = 0; j < NE; j++)
243 double element_osc = h * element_norms_of_fine_scale(j);
244 if ( element_osc >
threshold * av_norm_of_coeff )
254 MPI_Comm comm = pmesh->
GetComm();
255 MPI_Allreduce(MPI_IN_PLACE, &
global_osc, 1, MPI_DOUBLE, MPI_SUM, comm);
256 MPI_Comm_rank(comm, &rank);
267 MFEM_WARNING(
"Reached maximum number of elements " 268 "before resolving data to tolerance.");
const IntegrationRule * ir_default[Geometry::NumGeom]
char ref_type
refinement XYZ bit mask (7 = full isotropic)
Class for grid function - Vector with associated FE space.
Array< int > mesh_refinements
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
virtual int ApplyImpl(Mesh &mesh)
Apply the operator to the mesh.
int Dimension() const
Dimension of the reference space used within the elements.
void SetSize(int s)
Resize the vector to size s.
bool Nonconforming() const
int Size() const
Returns the size of the vector.
virtual int PreprocessMesh(Mesh &mesh, int max_it)
Apply the operator to the mesh max_it times or until tolerance achieved.
Abstract parallel finite element space.
virtual int ApplyImpl(Mesh &mesh)
Apply the operator to the mesh.
virtual void Reset()
Reset all MeshOperators in the sequence.
int index
Mesh element number.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
virtual int ApplyImpl(Mesh &mesh)
Apply the MeshOperatorSequence.
bit mask for all "action" bits
long long GetGlobalNE() const
Return the total (global) number of elements.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
double GetElementSize(ElementTransformation *T, int type=0)
virtual void Reset()
Reset the associated estimator.
virtual const Vector & GetLocalErrors()=0
Get a Vector with all element errors.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
const IntegrationRule ** irs
virtual int ApplyImpl(Mesh &mesh)
Rebalance a parallel mesh (only non-conforming parallel meshes are supported).
bit mask for all "info" bits
Base class for all element based error estimators.
Array< Refinement > marked_elements
double Normlp(double p) const
Returns the l_p norm of the vector.
virtual const Array< int > & GetAnisotropicFlags()=0
Get an Array<int> with anisotropic flags for all mesh elements.
ErrorEstimator & estimator
virtual void Reset()=0
Force recomputation of the estimates on the next call to GetLocalErrors.
virtual int ApplyImpl(Mesh &mesh)
Apply the operator to the mesh once.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
virtual ~MeshOperatorSequence()
Delete all operators from the sequence.
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
a stopping criterion was satisfied
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).
AnisotropicErrorEstimator * aniso_estimator
ErrorEstimator & estimator
int GetNE() const
Returns number of elements.
virtual long long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
bool DerefineByError(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
void Destroy()
Destroy a vector.
The AnisotropicErrorEstimator class is the base class for all error estimators that compute one non-n...
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
double ParNormlp(const Vector &vec, double p, MPI_Comm comm)
Compute the l_p norm of the Vector which is split without overlap across the given communicator...
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.
int Size() const
Return the logical size of the array.
long long num_marked_elements
ThresholdRefiner(ErrorEstimator &est)
Construct a ThresholdRefiner using the given ErrorEstimator.
double GetNorm(const Vector &local_err, Mesh &mesh) const
virtual void Reset()
Reset.
Class for parallel grid function.
Array< MeshOperator * > sequence
MeshOperators sequence, owned by us.
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Class for parallel meshes.
Arbitrary order "L2-conforming" discontinuous finite elements.