38 S->
Mult(tmp, fieldtrue);
50 double dx, loc_area = 0.0;
51 for (
int i = 0; i < pmesh.
GetNE(); i++)
56 MPI_Allreduce(&loc_area, &glob_area, 1, MPI_DOUBLE,
62 dx = glob_area / glob_zones;
break;
64 dx = sqrt(glob_area / glob_zones);
break;
66 dx = sqrt(2.0 * glob_area / glob_zones);
break;
68 dx = pow(glob_area / glob_zones, 1.0/3.0);
break;
70 dx = pow(6.0 * glob_area / glob_zones, 1.0/3.0);
break;
71 default: MFEM_ABORT(
"Unknown zone type!"); dx = 0.0;
82 "Only Ordering::byNODES is supported.");
85 const int size = dist_s.
Size();
90 for (
int d = 0; d <
dim; d++)
93 for (
int i = 0; i < size; i++)
95 magn(i) += der(i) * der(i);
97 dist_v(i + d*size) = (dist_s(i) > 0.0) ? -der(i) : der(i);
101 for (
int i = 0; i < size; i++)
103 const double vec_magn = std::sqrt(magn(i) + 1e-12);
104 for (
int d = 0; d <
dim; d++)
106 dist_v(i + d*size) *= fabs(dist_s(i)) / vec_magn;
116 "This function expects a vector ParGridFunction!");
131 MFEM_VERIFY(check_h1 && pfes.
GetVDim() == 1,
132 "This solver supports only scalar H1 spaces.");
139 source.ProjectCoefficient(zero_level_set);
146 for (
int i = 0; i <
source.Size(); i++)
148 const double x =
source(i);
149 source(i) = ((x < -1.0) || (x > 1.0)) ? 0.0 : (1.0 - x) * (1.0 + x);
153 int amg_print_level = 0;
182 if (pmesh.bdr_attributes.Size())
184 Array<int> ess_bdr(pmesh.bdr_attributes.Max());
217 for (
int ii = 0; ii < diffused_source.
Size(); ii++)
223 diffused_source(ii) = 0.5 * (u_neumann(ii) + u_dirichlet(ii));
258 double d_min_loc = distance.
Min();
260 MPI_Allreduce(&d_min_loc, &d_min_glob, 1, MPI_DOUBLE,
262 distance -= d_min_glob;
273 x.ProjectCoefficient(grad_u);
276 sol_sock_x <<
"parallel " << pfes.
GetNRanks() <<
" " 278 sol_sock_x.precision(8);
279 sol_sock_x <<
"solution\n" << pmesh << x;
280 sol_sock_x <<
"window_geometry " << 0 <<
" " << 0 <<
" " 281 << 500 <<
" " << 500 <<
"\n" 282 <<
"window_title '" <<
"Heat Directions" <<
"'\n" 283 <<
"keys evvRj*******A\n" << std::flush;
287 double NormalizationDistanceSolver::NormalizationCoeff::
292 u.GetGradient(T, u_grad);
293 const double u_value =
u.GetValue(T, ip);
295 return u_value / sqrt(u_value * u_value + u_grad * u_grad + 1e-12);
306 NormalizationCoeff rv_coeff(u_gf);
317 MFEM_VERIFY((check_h1 || check_l2) && fesd->
GetVDim() == 1,
318 "This solver supports only scalar H1 or L2 spaces.");
321 const int dim=mesh->Dimension();
323 MPI_Comm lcomm=fesd->
GetComm();
325 MPI_Comm_rank(lcomm,&myrank);
327 const int order = fesd->
GetOrder(0);
371 for (
int pp=3; pp<maxp; pp++)
375 std::cout <<
"pp = " << pp << std::endl;
418 trans.SetIntPoint(&ip);
430 for (
int jj=0; jj<ndim; jj++)
432 ngrad2 = ngrad2 + qval(jj)*qval(jj);
435 energy = energy + w * ngrad2 *
diffcoef * 0.5;
440 energy = energy + w * pval * pval * 0.5;
444 energy = energy - w*pval;
448 energy = energy + w*pval;
486 trans.SetIntPoint(&ip);
502 elvect.
Add(w * pval, shapef);
510 elvect.
Add( -w, shapef);
514 elvect.
Add( w, shapef);
543 trans.SetIntPoint(&ip);
585 trans.SetIntPoint(&ip);
601 for (
int jj=0; jj<ndim; jj++)
605 tmpv.
Add(vgrad[jj],shapef);
611 for (
int jj=0; jj<ndim; jj++)
613 ngrad2 = ngrad2 + qval(jj)*qval(jj);
616 energy = energy + w * std::pow(ngrad2+
ee*
ee,
pp/2.0)/
pp;
622 energy = energy - w * pval * tval;
626 energy = energy + w * pval * tval;
669 trans.SetIntPoint(&ip);
685 for (
int jj=0; jj<ndim; jj++)
689 tmpv.
Add(vgrad[jj],shapef);
696 for (
int jj=0; jj<ndim; jj++)
698 ngrad2 = ngrad2 + qval(jj)*qval(jj);
702 aa = std::pow(aa, (
pp - 2.0) / 2.0);
704 elvect.
Add(w * aa,lvec);
710 elvect.
Add( -w*fval, shapef);
714 elvect.
Add( w*fval, shapef);
752 trans.SetIntPoint(&ip);
767 for (
int jj=0; jj<ndim; jj++)
771 tmpv.
Add(vgrad[jj],shapef);
779 for (
int jj=0; jj<ndim; jj++)
781 ngrad2 = ngrad2 + qval(jj)*qval(jj);
784 aa = ngrad2 +
ee *
ee;
785 aa1 = std::pow(aa, (
pp - 2.0) / 2.0);
786 aa0 = (
pp-2.0) * std::pow(aa, (
pp - 4.0) / 2.0);
810 gmres->
Mult(rhs, *sv);
Abstract class for all finite elements.
Class for domain integration L(v) := (f, v)
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Conjugate gradient method.
void trans(const Vector &u, Vector &x)
int GetNPoints() const
Returns the number of the points in the integration rule.
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun) override
Compute the local energy.
Class for an integration rule - an Array of IntegrationPoint.
void GetDerivative(int comp, int der_comp, ParGridFunction &der)
Parallel version of GridFunction::GetDerivative(); see its documentation.
void SetCol(int c, const double *col)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
A coefficient that is constant across space and time.
Geometry::Type GetElementBaseGeometry(int i) const
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 ...
int Dimension() const
Dimension of the reference space used within the elements.
void SetSize(int s)
Resize the vector to size s.
Class for domain integrator L(v) := (f, grad v)
void ComputeScalarDistance(Coefficient &u_coeff, ParGridFunction &dist)
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
void Filter(ParGridFunction &func, ParGridFunction &ffield)
Pointer to an Operator of a specified type.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Data type dense matrix using column-major storage.
void ScalarDistToVector(ParGridFunction &dist_s, ParGridFunction &dist_v)
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
virtual void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)=0
void DeleteAll()
Delete the whole array.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
The BoomerAMG solver in hypre.
void source(const Vector &x, Vector &f)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
Assemble the local gradient matrix.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
long long GetGlobalNE() const
Return the total (global) number of elements.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
const FiniteElementCollection * FEColl() const
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
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.
ParMesh * GetParMesh() const
void SetMaxIter(int max_it)
Vector coefficient defined as the Gradient of a scalar GridFunction.
Newton's method for solving F(x)=b for a given operator F.
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Parallel smoothers in hypre.
void GetColumn(int c, Vector &col) const
VectorCoefficient * fgrad
void Mult(const double *x, double *y) const
Matrix vector multiplication.
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
IterativeSolver::PrintLevel print_level
int GetVDim() const
Returns vector dimension.
Wrapper for hypre's parallel vector class.
Mesh * GetMesh() const
Returns the mesh.
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
void SetAbsTol(double atol)
void SetRelTol(double rtol)
int GetDim() const
Returns the reference space dimension for the finite element.
void SetPower(double pp_)
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...
double Min() const
Returns the minimal element of the vector.
ParFiniteElementSpace * ParFESpace() const
int GetOrder(int i) const
Returns the polynomial degree of the i'th finite element.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
int GetNE() const
Returns number of elements.
Solver wrapper which orthogonalizes the input and output vector.
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Class for integration point with weight.
void SetInput(Coefficient &nfunc)
void ComputeScalarDistance(Coefficient &func, ParGridFunction &fdist)
void DiffuseField(ParGridFunction &field, int smooth_steps)
Ordering::Type GetOrdering() const
Return the ordering method.
virtual void ComputeVectorDistance(Coefficient &zero_level_set, ParGridFunction &distance)
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
bool iterations
Detailed information about each iteration will be reported to mfem::out.
double u(const Vector &xvec)
Class for parallel grid function.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
Perform the local action of the NonlinearFormIntegrator.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
void ComputeScalarDistance(Coefficient &zero_level_set, ParGridFunction &distance)
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
double GetElementVolume(int i)
double AvgElementSize(ParMesh &pmesh)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Arbitrary order "L2-conforming" discontinuous finite elements.
void Neg()
(*this) = -(*this)