MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex15p.cpp
Go to the documentation of this file.
1 // MFEM Example 15 - Parallel Version
2 //
3 // Compile with: make ex15p
4 //
5 // Sample runs: mpirun -np 4 ex15p
6 // mpirun -np 4 ex15p -o 1 -y 0.2
7 // mpirun -np 4 ex15p -o 4 -y 0.1
8 // mpirun -np 4 ex15p -n 5
9 // mpirun -np 4 ex15p -p 1 -n 3
10 //
11 // Other meshes:
12 //
13 // mpirun -np 4 ex15p -m ../data/square-disc-nurbs.mesh
14 // mpirun -np 4 ex15p -m ../data/disc-nurbs.mesh
15 // mpirun -np 4 ex15p -m ../data/fichera.mesh -tf 0.5
16 // mpirun -np 4 ex15p -m ../data/ball-nurbs.mesh -tf 0.5
17 // mpirun -np 4 ex15p -m ../data/mobius-strip.mesh
18 // mpirun -np 4 ex15p -m ../data/amr-quad.mesh
19 // mpirun -np 4 ex15p -m ../data/square-disc.mesh
20 // mpirun -np 4 ex15p -m ../data/escher.mesh -r 2 -tf 0.3
21 //
22 // Different estimators:
23 //
24 // mpirun -np 4 ex15p -est 0 -e 1e-4
25 // mpirun -np 4 ex15p -est 1 -e 1e-6
26 // mpirun -np 4 ex15p -est 1 -o 3 -tf 0.3
27 // mpirun -np 4 ex15p -est 2 -o 2
28 //
29 // Description: Building on Example 6, this example demonstrates dynamic AMR.
30 // The mesh is adapted to a time-dependent solution by refinement
31 // as well as by derefinement. For simplicity, the solution is
32 // prescribed and no time integration is done. However, the error
33 // estimation and refinement/derefinement decisions are realistic.
34 //
35 // At each outer iteration the right hand side function is changed
36 // to mimic a time dependent problem. Within each inner iteration
37 // the problem is solved on a sequence of meshes which are locally
38 // refined according to a chosen error estimator. Currently there
39 // are three error estimators supported: A L2 formulation of the
40 // Zienkiewicz-Zhu error estimator (0), a Kelly error indicator (1)
41 // and a traditional Zienkiewicz-Zhu error estimator (2). At the
42 // end of the inner iteration the error estimates are also used to
43 // identify any elements which may be over-refined and a single
44 // derefinement step is performed. After each refinement or
45 // derefinement step a rebalance operation is performed to keep
46 // the mesh evenly distributed among the available processors.
47 //
48 // The example demonstrates MFEM's capability to refine, derefine
49 // and load balance nonconforming meshes, in 2D and 3D, and on
50 // linear, curved and surface meshes. Interpolation of functions
51 // between coarse and fine meshes, persistent GLVis visualization,
52 // and saving of time-dependent fields for external visualization
53 // with VisIt (visit.llnl.gov) are also illustrated.
54 //
55 // We recommend viewing Examples 1, 6 and 9 before viewing this
56 // example.
57 
58 #include "mfem.hpp"
59 #include <fstream>
60 #include <iostream>
61 
62 using namespace std;
63 using namespace mfem;
64 
65 // Choices for the problem setup. Affect bdr_func and rhs_func.
66 int problem;
68 
69 // Prescribed time-dependent boundary and right-hand side functions.
70 double bdr_func(const Vector &pt, double t);
71 double rhs_func(const Vector &pt, double t);
72 
73 // Update the finite element space, interpolate the solution and perform
74 // parallel load balancing.
77  ParLinearForm &b);
78 
79 
80 int main(int argc, char *argv[])
81 {
82  // 1. Initialize MPI.
83  int num_procs, myid;
84  MPI_Init(&argc, &argv);
85  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
86  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
87 
88  // 2. Parse command-line options.
89  problem = 0;
90  nfeatures = 1;
91  const char *mesh_file = "../data/star-hilbert.mesh";
92  int order = 2;
93  double t_final = 1.0;
94  double max_elem_error = 1.0e-4;
95  double hysteresis = 0.25; // derefinement safety coefficient
96  int ref_levels = 0;
97  int nc_limit = 3; // maximum level of hanging nodes
98  bool visualization = true;
99  bool visit = false;
100  int which_estimator = 0;
101 
102  OptionsParser args(argc, argv);
103  args.AddOption(&mesh_file, "-m", "--mesh",
104  "Mesh file to use.");
105  args.AddOption(&problem, "-p", "--problem",
106  "Problem setup to use: 0 = spherical front, 1 = ball.");
107  args.AddOption(&nfeatures, "-n", "--nfeatures",
108  "Number of solution features (fronts/balls).");
109  args.AddOption(&order, "-o", "--order",
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.");
130  args.Parse();
131  if (!args.Good())
132  {
133  if (myid == 0)
134  {
135  args.PrintUsage(cout);
136  }
137  MPI_Finalize();
138  return 1;
139  }
140  if (myid == 0)
141  {
142  args.PrintOptions(cout);
143  }
144 
145  // 3. Read the (serial) mesh from the given mesh file on all processors. We
146  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
147  // and volume meshes with the same code.
148  Mesh *mesh = new Mesh(mesh_file, 1, 1);
149  int dim = mesh->Dimension();
150  int sdim = mesh->SpaceDimension();
151 
152  // 4. Project a NURBS mesh to a piecewise-quadratic curved mesh. Make sure
153  // that the mesh is non-conforming if it has quads or hexes and refine it.
154  if (mesh->NURBSext)
155  {
156  mesh->UniformRefinement();
157  if (ref_levels > 0) { ref_levels--; }
158  mesh->SetCurvature(2);
159  }
160  mesh->EnsureNCMesh(true);
161  for (int l = 0; l < ref_levels; l++)
162  {
163  mesh->UniformRefinement();
164  }
165  // Make sure tet-only meshes are marked for local refinement.
166  mesh->Finalize(true);
167 
168  // 5. Define a parallel mesh by partitioning the serial mesh. Once the
169  // parallel mesh is defined, the serial mesh can be deleted.
170  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
171  delete mesh;
172 
173  MFEM_VERIFY(pmesh.bdr_attributes.Size() > 0,
174  "Boundary attributes required in the mesh.");
175  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
176  ess_bdr = 1;
177 
178  // 6. Define a finite element space on the mesh. The polynomial order is one
179  // (linear) by default, but this can be changed on the command line.
180  H1_FECollection fec(order, dim);
181  ParFiniteElementSpace fespace(&pmesh, &fec);
182 
183  // 7. As in Example 1p, we set up bilinear and linear forms corresponding to
184  // the Laplace problem -\Delta u = 1. We don't assemble the discrete
185  // problem yet, this will be done in the inner loop.
186  ParBilinearForm a(&fespace);
187  ParLinearForm b(&fespace);
188 
189  ConstantCoefficient one(1.0);
192 
194  a.AddDomainIntegrator(integ);
196 
197  // 8. The solution vector x and the associated finite element grid function
198  // will be maintained over the AMR iterations.
199  ParGridFunction x(&fespace);
200 
201  // 9. Connect to GLVis. Prepare for VisIt output.
202  char vishost[] = "localhost";
203  int visport = 19916;
204 
205  socketstream sout;
206  if (visualization)
207  {
208  sout.open(vishost, visport);
209  if (!sout)
210  {
211  if (myid == 0)
212  {
213  cout << "Unable to connect to GLVis server at "
214  << vishost << ':' << visport << endl;
215  cout << "GLVis visualization disabled.\n";
216  }
217  visualization = false;
218  }
219  sout.precision(8);
220  }
221 
222  VisItDataCollection visit_dc("Example15-Parallel", &pmesh);
223  visit_dc.RegisterField("solution", &x);
224  int vis_cycle = 0;
225 
226  // 10. As in Example 6p, we set up an estimator that will be used to obtain
227  // element error indicators. The integrator needs to provide the method
228  // ComputeElementFlux. We supply an L2 space for the discontinuous flux
229  // and an H(div) space for the smoothed flux.
230  L2_FECollection flux_fec(order, dim);
231  RT_FECollection smooth_flux_fec(order-1, dim);
232  ErrorEstimator* estimator{nullptr};
233 
234  switch (which_estimator)
235  {
236  case 1:
237  {
238  auto flux_fes = new ParFiniteElementSpace(&pmesh, &flux_fec, sdim);
239  estimator = new KellyErrorEstimator(*integ, x, flux_fes);
240  break;
241  }
242  case 2:
243  {
244  auto flux_fes = new ParFiniteElementSpace(&pmesh, &fec, sdim);
245  estimator = new ZienkiewiczZhuEstimator(*integ, x, flux_fes);
246  break;
247  }
248 
249  default:
250  if (myid == 0)
251  {
252  std::cout << "Unknown estimator. Falling back to L2ZZ." << std::endl;
253  }
254  case 0:
255  {
256  auto flux_fes = new ParFiniteElementSpace(&pmesh, &flux_fec, sdim);
257  auto smooth_flux_fes = new ParFiniteElementSpace(&pmesh, &smooth_flux_fec);
258  estimator = new L2ZienkiewiczZhuEstimator(*integ, x, flux_fes, smooth_flux_fes);
259  break;
260  }
261  }
262 
263  // 11. As in Example 6p, we also need a refiner. This time the refinement
264  // strategy is based on a fixed threshold that is applied locally to each
265  // element. The global threshold is turned off by setting the total error
266  // fraction to zero. We also enforce a maximum refinement ratio between
267  // adjacent elements.
268  ThresholdRefiner refiner(*estimator);
269  refiner.SetTotalErrorFraction(0.0); // use purely local threshold
270  refiner.SetLocalErrorGoal(max_elem_error);
271  refiner.PreferConformingRefinement();
272  refiner.SetNCLimit(nc_limit);
273 
274  // 12. A derefiner selects groups of elements that can be coarsened to form
275  // a larger element. A conservative enough threshold needs to be set to
276  // prevent derefining elements that would immediately be refined again.
277  ThresholdDerefiner derefiner(*estimator);
278  derefiner.SetThreshold(hysteresis * max_elem_error);
279  derefiner.SetNCLimit(nc_limit);
280 
281  // 13. The outer time loop. In each iteration we update the right hand side,
282  // solve the problem on the current mesh, visualize the solution and
283  // refine the mesh as many times as necessary. Then we derefine any
284  // elements which have very small errors.
285  for (double time = 0.0; time < t_final + 1e-10; time += 0.01)
286  {
287  if (myid == 0)
288  {
289  cout << "\nTime " << time << "\n\nRefinement:" << endl;
290  }
291 
292  // Set the current time in the coefficients
293  bdr.SetTime(time);
294  rhs.SetTime(time);
295 
296  // Make sure errors will be recomputed in the following.
297  refiner.Reset();
298  derefiner.Reset();
299 
300  // 14. The inner refinement loop. At the end we want to have the current
301  // time step resolved to the prescribed tolerance in each element.
302  for (int ref_it = 1; ; ref_it++)
303  {
304  HYPRE_BigInt global_dofs = fespace.GlobalTrueVSize();
305  if (myid == 0)
306  {
307  cout << "Iteration: " << ref_it << ", number of unknowns: "
308  << global_dofs << flush;
309  }
310 
311  // 15. Recompute the field on the current mesh: assemble the stiffness
312  // matrix and the right-hand side.
313  a.Assemble();
314  b.Assemble();
315 
316  // 16. Project the exact solution to the essential DOFs.
317  x.ProjectBdrCoefficient(bdr, ess_bdr);
318 
319  // 17. Create and solve the parallel linear system.
320  Array<int> ess_tdof_list;
321  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
322 
323  HypreParMatrix A;
324  Vector B, X;
325  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
326 
327  HypreBoomerAMG amg(A);
328  amg.SetPrintLevel(0);
329  HyprePCG pcg(A);
330  pcg.SetTol(1e-12);
331  pcg.SetMaxIter(200);
332  pcg.SetPrintLevel(0);
333  pcg.SetPreconditioner(amg);
334  pcg.Mult(B, X);
335 
336  // 18. Extract the local solution on each processor.
337  a.RecoverFEMSolution(X, b, x);
338 
339  // 19. Send the solution by socket to a GLVis server and optionally
340  // save it in VisIt format.
341  if (visualization)
342  {
343  sout << "parallel " << num_procs << " " << myid << "\n";
344  sout << "solution\n" << pmesh << x << flush;
345  }
346  if (visit)
347  {
348  visit_dc.SetCycle(vis_cycle++);
349  visit_dc.SetTime(time);
350  visit_dc.Save();
351  }
352 
353  // 20. Apply the refiner on the mesh. The refiner calls the error
354  // estimator to obtain element errors, then it selects elements to
355  // be refined and finally it modifies the mesh. The Stop() method
356  // determines if all elements satisfy the local threshold.
357  refiner.Apply(pmesh);
358  if (myid == 0)
359  {
360  cout << ", total error: " << estimator->GetTotalError() << endl;
361  }
362 
363  // 21. Quit the AMR loop if the termination criterion has been met
364  if (refiner.Stop())
365  {
366  a.Update(); // Free the assembled data
367  break;
368  }
369 
370  // 22. Update the space, interpolate the solution, rebalance the mesh.
371  UpdateAndRebalance(pmesh, fespace, x, a, b);
372  }
373 
374  // 23. Use error estimates from the last inner iteration to check for
375  // possible derefinements. The derefiner works similarly as the
376  // refiner. The errors are not recomputed because the mesh did not
377  // change (and also the estimator was not Reset() at this time).
378  if (derefiner.Apply(pmesh))
379  {
380  if (myid == 0)
381  {
382  cout << "\nDerefined elements." << endl;
383  }
384 
385  // 24. Update the space and the solution, rebalance the mesh.
386  UpdateAndRebalance(pmesh, fespace, x, a, b);
387  }
388  }
389 
390  delete estimator;
391 
392  // 25. Exit
393  MPI_Finalize();
394  return 0;
395 }
396 
397 
400  ParLinearForm &b)
401 {
402  // Update the space: recalculate the number of DOFs and construct a matrix
403  // that will adjust any GridFunctions to the new mesh state.
404  fespace.Update();
405 
406  // Interpolate the solution on the new mesh by applying the transformation
407  // matrix computed in the finite element space. Multiple GridFunctions could
408  // be updated here.
409  x.Update();
410 
411  if (pmesh.Nonconforming())
412  {
413  // Load balance the mesh.
414  pmesh.Rebalance();
415 
416  // Update the space again, this time a GridFunction redistribution matrix
417  // is created. Apply it to the solution.
418  fespace.Update();
419  x.Update();
420  }
421 
422  // Inform the linear and bilinear forms that the space has changed.
423  a.Update();
424  b.Update();
425 
426  // Free any transformation matrices to save memory.
427  fespace.UpdatesFinished();
428 }
429 
430 
431 const double alpha = 0.02;
432 
433 // Spherical front with a Gaussian cross section and radius t
434 double front(double x, double y, double z, double t, int)
435 {
436  double r = sqrt(x*x + y*y + z*z);
437  return exp(-0.5*pow((r - t)/alpha, 2));
438 }
439 
440 double front_laplace(double x, double y, double z, double t, int dim)
441 {
442  double x2 = x*x, y2 = y*y, z2 = z*z, t2 = t*t;
443  double r = sqrt(x2 + y2 + z2);
444  double a2 = alpha*alpha, a4 = a2*a2;
445  return -exp(-0.5*pow((r - t)/alpha, 2)) / a4 *
446  (-2*t*(x2 + y2 + z2 - (dim-1)*a2/2)/r + x2 + y2 + z2 + t2 - dim*a2);
447 }
448 
449 // Smooth spherical step function with radius t
450 double ball(double x, double y, double z, double t, int)
451 {
452  double r = sqrt(x*x + y*y + z*z);
453  return -atan(2*(r - t)/alpha);
454 }
455 
456 double ball_laplace(double x, double y, double z, double t, int dim)
457 {
458  double x2 = x*x, y2 = y*y, z2 = z*z, t2 = 4*t*t;
459  double r = sqrt(x2 + y2 + z2);
460  double a2 = alpha*alpha;
461  double den = pow(-a2 - 4*(x2 + y2 + z2 - 2*r*t) - t2, 2.0);
462  return (dim == 2) ? 2*alpha*(a2 + t2 - 4*x2 - 4*y2)/r/den
463  /* */ : 4*alpha*(a2 + t2 - 4*r*t)/r/den;
464 }
465 
466 // Composes several features into one function
467 template<typename F0, typename F1>
468 double composite_func(const Vector &pt, double t, F0 f0, F1 f1)
469 {
470  int dim = pt.Size();
471  double x = pt(0), y = pt(1), z = 0.0;
472  if (dim == 3) { z = pt(2); }
473 
474  if (problem == 0)
475  {
476  if (nfeatures <= 1)
477  {
478  return f0(x, y, z, t, dim);
479  }
480  else
481  {
482  double sum = 0.0;
483  for (int i = 0; i < nfeatures; i++)
484  {
485  double x0 = 0.5*cos(2*M_PI * i / nfeatures);
486  double y0 = 0.5*sin(2*M_PI * i / nfeatures);
487  sum += f0(x - x0, y - y0, z, t, dim);
488  }
489  return sum;
490  }
491  }
492  else
493  {
494  double sum = 0.0;
495  for (int i = 0; i < nfeatures; i++)
496  {
497  double x0 = 0.5*cos(2*M_PI * i / nfeatures + M_PI*t);
498  double y0 = 0.5*sin(2*M_PI * i / nfeatures + M_PI*t);
499  sum += f1(x - x0, y - y0, z, 0.25, dim);
500  }
501  return sum;
502  }
503 }
504 
505 // Exact solution, used for the Dirichlet BC.
506 double bdr_func(const Vector &pt, double t)
507 {
508  return composite_func(pt, t, front, ball);
509 }
510 
511 // Laplace of the exact solution, used for the right hand side.
512 double rhs_func(const Vector &pt, double t)
513 {
514  return composite_func(pt, t, front_laplace, ball_laplace);
515 }
void SetTol(double tol)
Definition: hypre.cpp:3569
double front_laplace(double x, double y, double z, double t, int dim)
Definition: ex15.cpp:378
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
double bdr_func(const Vector &pt, double t)
Definition: ex15.cpp:444
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel pr...
Definition: estimators.hpp:443
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2988
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
void SetNCLimit(int nc_limit)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
virtual void Save()
Save the collection and a VisIt root file.
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
double rhs_func(const Vector &pt, double t)
Definition: ex15.cpp:450
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Definition: estimators.hpp:217
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3589
void Update(ParFiniteElementSpace *pf=NULL)
Update the object according to the given new FE space *pf.
Definition: plinearform.cpp:21
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
bool Apply(Mesh &mesh)
Perform the mesh operation.
void Rebalance()
Definition: pmesh.cpp:3602
bool Nonconforming() const
Definition: mesh.hpp:1367
Class for parallel linear form.
Definition: plinearform.hpp:26
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 ...
Definition: optparser.cpp:150
constexpr char vishost[]
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:8800
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:81
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
virtual void Reset()
Reset the associated estimator.
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1467
constexpr int visport
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)
Definition: mesh.cpp:4882
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
int nfeatures
Definition: ex15.cpp:63
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Base class for all element based error estimators.
Definition: estimators.hpp:41
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetTime(double t)
Set physical time (for time-dependent simulations)
void SetTime(double t)
Set the time for time dependent coefficients.
Definition: coefficient.hpp:48
void SetNCLimit(int nc_limit)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3579
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
void UpdateAndRebalance(ParMesh &pmesh, ParFiniteElementSpace &fespace, ParGridFunction &x, ParBilinearForm &a, ParLinearForm &b)
Definition: ex15p.cpp:398
int SpaceDimension() const
Definition: mesh.hpp:912
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
double ball_laplace(double x, double y, double z, double t, int dim)
Definition: ex15.cpp:394
PCG solver in hypre.
Definition: hypre.hpp:1060
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:601
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
int problem
Definition: ex15.cpp:62
HYPRE_Int HYPRE_BigInt
void PreferConformingRefinement()
Use conforming refinement, if possible (triangles, tetrahedra) – this is the default.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2600
int dim
Definition: ex24.cpp:53
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3594
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
Class for parallel bilinear form.
double ball(double x, double y, double z, double t, int)
Definition: ex15.cpp:388
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
void SetThreshold(double thresh)
Set the de-refinement threshold. The default value is zero.
RefCoord t[3]
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure...
Definition: estimators.hpp:88
void SetLocalErrorGoal(double err_goal)
Set the local stopping criterion: stop when local_err_i &lt;= local_err_goal. The default value is zero...
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetTotalErrorFraction(double fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2...
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
virtual void UpdatesFinished()
Free ParGridFunction transformation matrix (if any), to save memory.
Definition: pfespace.hpp:411
double front(double x, double y, double z, double t, int)
Definition: ex15.cpp:372
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3617
Class for parallel meshes.
Definition: pmesh.hpp:32
double composite_func(const Vector &pt, double t, F0 f0, F1 f1)
Definition: ex15.cpp:406
int main()
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
De-refinement operator using an error threshold.
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150