MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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/fichera-mixed.mesh -tf 0.5
17// mpirun -np 4 ex15p -m ../data/ball-nurbs.mesh -tf 0.5
18// mpirun -np 4 ex15p -m ../data/mobius-strip.mesh
19// mpirun -np 4 ex15p -m ../data/amr-quad.mesh
20// mpirun -np 4 ex15p -m ../data/square-disc.mesh
21// mpirun -np 4 ex15p -m ../data/escher.mesh -r 2 -tf 0.3
22//
23// Different estimators:
24//
25// mpirun -np 4 ex15p -est 0 -e 1e-4
26// mpirun -np 4 ex15p -est 1 -e 1e-6
27// mpirun -np 4 ex15p -est 1 -o 3 -tf 0.3
28// mpirun -np 4 ex15p -est 2 -o 2
29//
30// Description: Building on Example 6, this example demonstrates dynamic AMR.
31// The mesh is adapted to a time-dependent solution by refinement
32// as well as by derefinement. For simplicity, the solution is
33// prescribed and no time integration is done. However, the error
34// estimation and refinement/derefinement decisions are realistic.
35//
36// At each outer iteration the right hand side function is changed
37// to mimic a time dependent problem. Within each inner iteration
38// the problem is solved on a sequence of meshes which are locally
39// refined according to a chosen error estimator. Currently there
40// are three error estimators supported: A L2 formulation of the
41// Zienkiewicz-Zhu error estimator (0), a Kelly error indicator (1)
42// and a traditional Zienkiewicz-Zhu error estimator (2). At the
43// end of the inner iteration the error estimates are also used to
44// identify any elements which may be over-refined and a single
45// derefinement step is performed. After each refinement or
46// derefinement step a rebalance operation is performed to keep
47// the mesh evenly distributed among the available processors.
48//
49// The example demonstrates MFEM's capability to refine, derefine
50// and load balance nonconforming meshes, in 2D and 3D, and on
51// linear, curved and surface meshes. Interpolation of functions
52// between coarse and fine meshes, persistent GLVis visualization,
53// and saving of time-dependent fields for external visualization
54// with VisIt (visit.llnl.gov) are also illustrated.
55//
56// We recommend viewing Examples 1, 6 and 9 before viewing this
57// example.
58
59#include "mfem.hpp"
60#include <fstream>
61#include <iostream>
62
63using namespace std;
64using namespace mfem;
65
66// Choices for the problem setup. Affect bdr_func and rhs_func.
69
70// Prescribed time-dependent boundary and right-hand side functions.
71real_t bdr_func(const Vector &pt, real_t t);
72real_t rhs_func(const Vector &pt, real_t t);
73
74// Update the finite element space, interpolate the solution and perform
75// parallel load balancing.
79
80
81int main(int argc, char *argv[])
82{
83 // 1. Initialize MPI and HYPRE.
84 Mpi::Init(argc, argv);
85 int num_procs = Mpi::WorldSize();
86 int myid = Mpi::WorldRank();
88
89 // 2. Parse command-line options.
90 problem = 0;
91 nfeatures = 1;
92 const char *mesh_file = "../data/star-hilbert.mesh";
93 int order = 2;
94 real_t t_final = 1.0;
95 real_t max_elem_error = 1.0e-4;
96 real_t hysteresis = 0.25; // derefinement safety coefficient
97 int ref_levels = 0;
98 int nc_limit = 3; // maximum level of hanging nodes
99 bool visualization = true;
100 bool visit = false;
101 int which_estimator = 0;
102
103 OptionsParser args(argc, argv);
104 args.AddOption(&mesh_file, "-m", "--mesh",
105 "Mesh file to use.");
106 args.AddOption(&problem, "-p", "--problem",
107 "Problem setup to use: 0 = spherical front, 1 = ball.");
108 args.AddOption(&nfeatures, "-n", "--nfeatures",
109 "Number of solution features (fronts/balls).");
110 args.AddOption(&order, "-o", "--order",
111 "Finite element order (polynomial degree).");
112 args.AddOption(&max_elem_error, "-e", "--max-err",
113 "Maximum element error");
114 args.AddOption(&hysteresis, "-y", "--hysteresis",
115 "Derefinement safety coefficient.");
116 args.AddOption(&ref_levels, "-r", "--ref-levels",
117 "Number of initial uniform refinement levels.");
118 args.AddOption(&nc_limit, "-l", "--nc-limit",
119 "Maximum level of hanging nodes.");
120 args.AddOption(&t_final, "-tf", "--t-final",
121 "Final time; start time is 0.");
122 args.AddOption(&which_estimator, "-est", "--estimator",
123 "Which estimator to use: "
124 "0 = L2ZZ, 1 = Kelly, 2 = ZZ. Defaults to L2ZZ.");
125 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
126 "--no-visualization",
127 "Enable or disable GLVis visualization.");
128 args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
129 "--no-visit-datafiles",
130 "Save data files for VisIt (visit.llnl.gov) visualization.");
131 args.Parse();
132 if (!args.Good())
133 {
134 if (myid == 0)
135 {
136 args.PrintUsage(cout);
137 }
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);
195 b.AddDomainIntegrator(new DomainLFIntegrator(rhs));
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);
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 (real_t 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
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 return 0;
394}
395
396
400{
401 // Update the space: recalculate the number of DOFs and construct a matrix
402 // that will adjust any GridFunctions to the new mesh state.
403 fespace.Update();
404
405 // Interpolate the solution on the new mesh by applying the transformation
406 // matrix computed in the finite element space. Multiple GridFunctions could
407 // be updated here.
408 x.Update();
409
410 if (pmesh.Nonconforming())
411 {
412 // Load balance the mesh.
413 pmesh.Rebalance();
414
415 // Update the space again, this time a GridFunction redistribution matrix
416 // is created. Apply it to the solution.
417 fespace.Update();
418 x.Update();
419 }
420
421 // Inform the linear and bilinear forms that the space has changed.
422 a.Update();
423 b.Update();
424
425 // Free any transformation matrices to save memory.
426 fespace.UpdatesFinished();
427}
428
429
430const real_t alpha = 0.02;
431
432// Spherical front with a Gaussian cross section and radius t
434{
435 real_t r = sqrt(x*x + y*y + z*z);
436 return exp(-0.5*pow((r - t)/alpha, 2));
437}
438
440{
441 real_t x2 = x*x, y2 = y*y, z2 = z*z, t2 = t*t;
442 real_t r = sqrt(x2 + y2 + z2);
443 real_t a2 = alpha*alpha, a4 = a2*a2;
444 return -exp(-0.5*pow((r - t)/alpha, 2)) / a4 *
445 (-2*t*(x2 + y2 + z2 - (dim-1)*a2/2)/r + x2 + y2 + z2 + t2 - dim*a2);
446}
447
448// Smooth spherical step function with radius t
450{
451 real_t r = sqrt(x*x + y*y + z*z);
452 return -atan(2*(r - t)/alpha);
453}
454
456{
457 real_t x2 = x*x, y2 = y*y, z2 = z*z, t2 = 4*t*t;
458 real_t r = sqrt(x2 + y2 + z2);
459 real_t a2 = alpha*alpha;
460 real_t den = pow(-a2 - 4*(x2 + y2 + z2 - 2*r*t) - t2, 2.0);
461 return (dim == 2) ? 2*alpha*(a2 + t2 - 4*x2 - 4*y2)/r/den
462 /* */ : 4*alpha*(a2 + t2 - 4*r*t)/r/den;
463}
464
465// Composes several features into one function
466template<typename F0, typename F1>
467real_t composite_func(const Vector &pt, real_t t, F0 f0, F1 f1)
468{
469 int dim = pt.Size();
470 real_t x = pt(0), y = pt(1), z = 0.0;
471 if (dim == 3) { z = pt(2); }
472
473 if (problem == 0)
474 {
475 if (nfeatures <= 1)
476 {
477 return f0(x, y, z, t, dim);
478 }
479 else
480 {
481 real_t sum = 0.0;
482 for (int i = 0; i < nfeatures; i++)
483 {
484 real_t x0 = 0.5*cos(2*M_PI * i / nfeatures);
485 real_t y0 = 0.5*sin(2*M_PI * i / nfeatures);
486 sum += f0(x - x0, y - y0, z, t, dim);
487 }
488 return sum;
489 }
490 }
491 else
492 {
493 real_t sum = 0.0;
494 for (int i = 0; i < nfeatures; i++)
495 {
496 real_t x0 = 0.5*cos(2*M_PI * i / nfeatures + M_PI*t);
497 real_t y0 = 0.5*sin(2*M_PI * i / nfeatures + M_PI*t);
498 sum += f1(x - x0, y - y0, z, 0.25, dim);
499 }
500 return sum;
501 }
502}
503
504// Exact solution, used for the Dirichlet BC.
506{
507 return composite_func(pt, t, front, ball);
508}
509
510// Laplace of the exact solution, used for the right hand side.
512{
514}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
Abstract base class BilinearFormIntegrator.
virtual void SetTime(real_t t)
Set the time for time dependent coefficients.
A coefficient that is constant across space and time.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
Class for domain integration .
Definition lininteg.hpp:109
Base class for all element based error estimators.
A general function coefficient.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
void SetPrintLevel(int print_level)
Definition hypre.hpp:1771
PCG solver in hypre.
Definition hypre.hpp:1275
void SetPrintLevel(int print_lvl)
Definition hypre.cpp:4156
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition hypre.cpp:4161
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Definition hypre.cpp:4184
void SetMaxIter(int max_iter)
Definition hypre.cpp:4146
void SetTol(real_t tol)
Definition hypre.cpp:4136
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel pr...
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
bool Apply(Mesh &mesh)
Perform the mesh operation.
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
bool Nonconforming() const
Definition mesh.hpp:2229
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition mesh.cpp:3241
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition mesh.cpp:10626
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.
Definition mesh.cpp:6211
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
void UpdatesFinished() override
Free ParGridFunction transformation matrix (if any), to save memory.
Definition pfespace.hpp:425
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
void Update(bool want_transform=true) override
Class for parallel grid function.
Definition pgridfunc.hpp:33
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition pgridfunc.cpp:90
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
void Rebalance()
Definition pmesh.cpp:4009
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
De-refinement operator using an error threshold.
void SetThreshold(real_t thresh)
Set the de-refinement threshold. The default value is zero.
void SetNCLimit(int nc_limit_)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
virtual void Reset()
Reset the associated estimator.
Mesh refinement operator using an error threshold.
void SetNCLimit(int nc_limit_)
Set the maximum ratio of refinement levels of adjacent elements (0 = unlimited).
void SetTotalErrorFraction(real_t fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2.
virtual void Reset()
Reset the associated estimator.
void PreferConformingRefinement()
Use conforming refinement, if possible (triangles, tetrahedra) – this is the default.
void SetLocalErrorGoal(real_t err_goal)
Set the local stopping criterion: stop when local_err_i <= local_err_goal. The default value is zero.
Vector data type.
Definition vector.hpp:80
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
Data collection with VisIt I/O routines.
virtual void Save() override
Save the collection and a VisIt root file.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
real_t front_laplace(real_t x, real_t y, real_t z, real_t t, int dim)
Definition ex15p.cpp:439
real_t ball(real_t x, real_t y, real_t z, real_t t, int)
Definition ex15p.cpp:449
void UpdateAndRebalance(ParMesh &pmesh, ParFiniteElementSpace &fespace, ParGridFunction &x, ParBilinearForm &a, ParLinearForm &b)
Definition ex15p.cpp:397
real_t rhs_func(const Vector &pt, real_t t)
Definition ex15p.cpp:511
real_t front(real_t x, real_t y, real_t z, real_t t, int)
Definition ex15p.cpp:433
real_t composite_func(const Vector &pt, real_t t, F0 f0, F1 f1)
Definition ex15p.cpp:467
int problem
Definition ex15p.cpp:67
real_t bdr_func(const Vector &pt, real_t t)
Definition ex15p.cpp:505
const real_t alpha
Definition ex15p.cpp:430
real_t ball_laplace(real_t x, real_t y, real_t z, real_t t, int dim)
Definition ex15p.cpp:455
int nfeatures
Definition ex15p.cpp:68
int dim
Definition ex24.cpp:53
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const int visport
float real_t
Definition config.hpp:43
const char vishost[]
RefCoord t[3]