80 int main(
int argc,
char *argv[])
83 Mpi::Init(argc, argv);
84 int num_procs = Mpi::WorldSize();
85 int myid = Mpi::WorldRank();
91 const char *mesh_file =
"../data/star-hilbert.mesh";
94 double max_elem_error = 1.0e-4;
95 double hysteresis = 0.25;
98 bool visualization =
true;
100 int which_estimator = 0;
103 args.
AddOption(&mesh_file,
"-m",
"--mesh",
104 "Mesh file to use.");
106 "Problem setup to use: 0 = spherical front, 1 = ball.");
108 "Number of solution features (fronts/balls).");
110 "Finite element order (polynomial degree).");
111 args.
AddOption(&max_elem_error,
"-e",
"--max-err",
112 "Maximum element error");
113 args.
AddOption(&hysteresis,
"-y",
"--hysteresis",
114 "Derefinement safety coefficient.");
115 args.
AddOption(&ref_levels,
"-r",
"--ref-levels",
116 "Number of initial uniform refinement levels.");
117 args.
AddOption(&nc_limit,
"-l",
"--nc-limit",
118 "Maximum level of hanging nodes.");
119 args.
AddOption(&t_final,
"-tf",
"--t-final",
120 "Final time; start time is 0.");
121 args.
AddOption(&which_estimator,
"-est",
"--estimator",
122 "Which estimator to use: " 123 "0 = L2ZZ, 1 = Kelly, 2 = ZZ. Defaults to L2ZZ.");
124 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
125 "--no-visualization",
126 "Enable or disable GLVis visualization.");
127 args.
AddOption(&visit,
"-visit",
"--visit-datafiles",
"-no-visit",
128 "--no-visit-datafiles",
129 "Save data files for VisIt (visit.llnl.gov) visualization.");
147 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
156 if (ref_levels > 0) { ref_levels--; }
160 for (
int l = 0; l < ref_levels; l++)
169 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
173 "Boundary attributes required in the mesh.");
193 a.AddDomainIntegrator(integ);
212 cout <<
"Unable to connect to GLVis server at " 214 cout <<
"GLVis visualization disabled.\n";
216 visualization =
false;
233 switch (which_estimator)
251 std::cout <<
"Unknown estimator. Falling back to L2ZZ." << std::endl;
284 for (
double time = 0.0; time < t_final + 1e-10; time += 0.01)
288 cout <<
"\nTime " << time <<
"\n\nRefinement:" << endl;
301 for (
int ref_it = 1; ; ref_it++)
306 cout <<
"Iteration: " << ref_it <<
", number of unknowns: " 307 << global_dofs << flush;
324 a.FormLinearSystem(ess_tdof_list, x,
b, A, X, B);
336 a.RecoverFEMSolution(X,
b, x);
342 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
343 sout <<
"solution\n" << pmesh << x << flush;
356 refiner.
Apply(pmesh);
359 cout <<
", total error: " << estimator->GetTotalError() << endl;
377 if (derefiner.
Apply(pmesh))
381 cout <<
"\nDerefined elements." << endl;
432 double front(
double x,
double y,
double z,
double t,
int)
434 double r = sqrt(x*x + y*y + z*z);
435 return exp(-0.5*pow((r -
t)/
alpha, 2));
440 double x2 = x*x, y2 = y*y, z2 = z*z, t2 =
t*
t;
441 double r = sqrt(x2 + y2 + z2);
443 return -exp(-0.5*pow((r -
t)/
alpha, 2)) / a4 *
444 (-2*
t*(x2 + y2 + z2 - (
dim-1)*a2/2)/r + x2 + y2 + z2 + t2 -
dim*a2);
448 double ball(
double x,
double y,
double z,
double t,
int)
450 double r = sqrt(x*x + y*y + z*z);
451 return -atan(2*(r -
t)/
alpha);
456 double x2 = x*x, y2 = y*y, z2 = z*z, t2 = 4*
t*
t;
457 double r = sqrt(x2 + y2 + z2);
459 double den = pow(-a2 - 4*(x2 + y2 + z2 - 2*r*
t) - t2, 2.0);
460 return (
dim == 2) ? 2*
alpha*(a2 + t2 - 4*x2 - 4*y2)/r/den
461 : 4*
alpha*(a2 + t2 - 4*r*
t)/r/den;
465 template<
typename F0,
typename F1>
469 double x = pt(0), y = pt(1), z = 0.0;
470 if (
dim == 3) { z = pt(2); }
476 return f0(x, y, z,
t,
dim);
483 double x0 = 0.5*cos(2*M_PI * i /
nfeatures);
484 double y0 = 0.5*sin(2*M_PI * i /
nfeatures);
485 sum += f0(x - x0, y - y0, z,
t,
dim);
495 double x0 = 0.5*cos(2*M_PI * i /
nfeatures + M_PI*
t);
496 double y0 = 0.5*sin(2*M_PI * i /
nfeatures + M_PI*
t);
497 sum += f1(x - x0, y - y0, z, 0.25,
dim);
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)
double bdr_func(const Vector &pt, double t)
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel pr...
void SetNCLimit(int nc_limit_)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
int Dimension() const
Dimension of the reference space used within the elements.
void PrintUsage(std::ostream &out) const
Print the usage message.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
virtual void Update(bool want_transform=true)
bool Nonconforming() const
int Size() const
Returns the size of the vector.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
virtual void SetTime(double t)
Set the time for time dependent coefficients.
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
double ball_laplace(double x, double y, double z, double t, int dim)
void SetPrintLevel(int print_lvl)
The BoomerAMG solver in hypre.
bool Apply(Mesh &mesh)
Perform the mesh operation.
double rhs_func(const Vector &pt, double t)
virtual void Reset()
Reset the associated estimator.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
double front_laplace(double x, double y, double z, double t, int dim)
void EnsureNCMesh(bool simplices_nonconforming=false)
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
virtual void Reset()
Reset the associated estimator.
void SetPrintLevel(int print_level)
Mesh refinement operator using an error threshold.
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Base class for all element based error estimators.
double composite_func(const Vector &pt, double t, F0 f0, F1 f1)
void SetTime(double t)
Set physical time (for time-dependent simulations)
HYPRE_BigInt GlobalTrueVSize() const
int main(int argc, char *argv[])
void SetMaxIter(int max_iter)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
void UpdateAndRebalance(ParMesh &pmesh, ParFiniteElementSpace &fespace, ParGridFunction &x, ParBilinearForm &a, ParLinearForm &b)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
virtual void Save() override
Save the collection and a VisIt root file.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
double ball(double x, double y, double z, double t, int)
void PreferConformingRefinement()
Use conforming refinement, if possible (triangles, tetrahedra) – this is the default.
int SpaceDimension() const
Dimension of the physical space containing the mesh.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
void SetThreshold(double thresh)
Set the de-refinement threshold. The default value is zero.
int Size() const
Return the logical size of the array.
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
void SetLocalErrorGoal(double err_goal)
Set the local stopping criterion: stop when local_err_i <= local_err_goal. The default value is zero...
A general function coefficient.
Arbitrary order H1-conforming (continuous) finite elements.
Class for parallel grid function.
double front(double x, double y, double z, double t, int)
void SetTotalErrorFraction(double fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2...
Wrapper for hypre's ParCSR matrix class.
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
Class for parallel meshes.
void SetNCLimit(int nc_limit_)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
Arbitrary order "L2-conforming" discontinuous finite elements.
De-refinement operator using an error threshold.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.