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