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