12 #ifndef MFEM_DIST_SOLVER_HPP
13 #define MFEM_DIST_SOLVER_HPP
104 double rtol = 1e-7,
double atol = 1e-12)
105 : maxp(maxp_), newton_iter(newton_iter_),
106 newton_rel_tol(rtol), newton_abs_tol(atol) { }
115 const int newton_iter;
116 const double newton_rel_tol, newton_abs_tol;
135 const double norm = V.
Norml2() + 1e-12;
149 : basef(basec_), corrf(corrc_) { }
154 double u = basef.
Eval(T,ip);
155 double c = corrf.
Eval(T,ip);
156 if (u<0.0) { u*=-1.0; }
176 double rd=rh/(2*std::sqrt(3.0));
186 const Vector &elfun)
override;
215 bool ownership_=
true)
232 const Vector &elfun)
override;
252 int maxiter = 100,
double rtol = 1e-12,
253 double atol = 1e-15,
int print_lv = 0)
255 fecp(order, mesh.Dimension()),
256 fesp(&mesh, &fecp, 1),
Abstract class for all finite elements.
void trans(const Vector &u, Vector &x)
HypreParVector * NewTrueDofVector()
Class for grid function - Vector with associated FE space.
void ComputeScalarDistance(Coefficient &func, ParGridFunction &fdist)
Base class for vector Coefficients that optionally depend on time and space.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
double Norml2() const
Returns the l2 norm of the vector.
virtual ~DistanceSolver()
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Data type dense matrix using column-major storage.
Abstract parallel finite element space.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
HeatDistanceSolver(double diff_coeff)
void ScalarDistToVector(ParGridFunction &dist_s, ParGridFunction &dist_v)
VectorCoefficient * fgrad
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
PUMPLaplacian(Coefficient *nfunc, VectorCoefficient *nfgrad, bool ownership_=true)
The BoomerAMG solver in hypre.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetPrintLevel(int print_level)
void SetInput(Coefficient &nfunc)
void SetMaxIter(int max_it)
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
ScreenedPoisson(Coefficient &nfunc, double rh)
NormalizationDistanceSolver()
Wrapper for hypre's parallel vector class.
PProductCoefficient(Coefficient &basec_, Coefficient &corrc_)
IterativeSolver::PrintLevel print_level
void SetAbsTol(double atol)
void SetRelTol(double rtol)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void ComputeVectorDistance(Coefficient &zero_level_set, ParGridFunction &distance)
NormalizedGradCoefficient(const GridFunction &u_gf, int dim)
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)
Class for integration point with weight.
PLapDistanceSolver(int maxp_=30, int newton_iter_=10, double rtol=1e-7, double atol=1e-12)
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
double u(const Vector &xvec)
void SetPower(double pp_)
Class for parallel grid function.
void SetMaxPower(int new_pp)
Class for parallel meshes.
void ComputeScalarDistance(Coefficient &u_coeff, ParGridFunction &dist)
Settings for the output behavior of the IterativeSolver.
void Filter(ParGridFunction &func, ParGridFunction &ffield)
PDEFilter(ParMesh &mesh, double rh, int order=2, int maxiter=100, double rtol=1e-12, double atol=1e-15, int print_lv=0)
void GetGradient(ElementTransformation &tr, Vector &grad) const
double AvgElementSize(ParMesh &pmesh)
virtual void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)=0