12 #ifndef MFEM_DIST_SOLVER_HPP 13 #define MFEM_DIST_SOLVER_HPP 109 double rtol = 1e-7,
double atol = 1e-12)
110 : maxp(maxp_), newton_iter(newton_iter_),
111 newton_rel_tol(rtol), newton_abs_tol(atol) { }
120 const int newton_iter;
121 const double newton_rel_tol, newton_abs_tol;
140 const double norm = V.
Norml2() + 1e-12;
154 : basef(basec_), corrf(corrc_) { }
159 double u = basef.
Eval(T,ip);
160 double c = corrf.
Eval(T,ip);
161 if (
u<0.0) {
u*=-1.0; }
181 double rd=rh/(2*std::sqrt(3.0));
191 const Vector &elfun)
override;
220 bool ownership_=
true)
237 const Vector &elfun)
override;
257 int maxiter = 100,
double rtol = 1e-12,
258 double atol = 1e-15,
int print_lv = 0)
260 fecp(order, mesh.Dimension()),
261 fesp(&mesh, &fecp, 1),
315 #endif // MFEM_USE_MPI Abstract class for all finite elements.
void trans(const Vector &u, Vector &x)
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
HypreParVector * NewTrueDofVector()
Class for grid function - Vector with associated FE space.
virtual ~DistanceSolver()
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 ComputeScalarDistance(Coefficient &u_coeff, ParGridFunction &dist)
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)
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
HeatDistanceSolver(double diff_coeff)
NormalizedGradCoefficient(const GridFunction &u_gf, int dim)
Data type dense matrix using column-major storage.
void ScalarDistToVector(ParGridFunction &dist_s, ParGridFunction &dist_v)
Abstract parallel finite element space.
ScreenedPoisson(Coefficient &nfunc, double rh)
virtual void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)=0
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 ...
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
NormalizationDistanceSolver()
PLapDistanceSolver(int maxp_=30, int newton_iter_=10, double rtol=1e-7, double atol=1e-12)
The BoomerAMG solver in hypre.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetPrintLevel(int print_level)
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
void SetMaxIter(int max_it)
VectorCoefficient * fgrad
IterativeSolver::PrintLevel print_level
Wrapper for hypre's parallel vector class.
void SetAbsTol(double atol)
void SetRelTol(double rtol)
void SetPower(double pp_)
PProductCoefficient(Coefficient &basec_, Coefficient &corrc_)
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Class for integration point with weight.
void SetInput(Coefficient &nfunc)
void ComputeScalarDistance(Coefficient &func, ParGridFunction &fdist)
virtual void ComputeVectorDistance(Coefficient &zero_level_set, ParGridFunction &distance)
double Norml2() const
Returns the l2 norm of the vector.
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)
Class for parallel grid function.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
Class for parallel meshes.
void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)
PUMPLaplacian(Coefficient *nfunc, VectorCoefficient *nfgrad, bool ownership_=true)
Settings for the output behavior of the IterativeSolver.
double AvgElementSize(ParMesh &pmesh)
void SetMaxPower(int new_pp)