MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex19p.cpp
Go to the documentation of this file.
1 // MFEM Example 19 - Parallel Version
2 //
3 // Compile with: make ex19p
4 //
5 // Sample runs:
6 // mpirun -np 2 ex19p -m ../data/beam-quad.mesh
7 // mpirun -np 2 ex19p -m ../data/beam-tri.mesh
8 // mpirun -np 2 ex19p -m ../data/beam-hex.mesh
9 // mpirun -np 2 ex19p -m ../data/beam-tet.mesh
10 //
11 // Description: This examples solves a quasi-static incompressible nonlinear
12 // elasticity problem of the form 0 = H(x), where H is an
13 // incompressible hyperelastic model and x is a block state vector
14 // containing displacement and pressure variables. The geometry of
15 // the domain is assumed to be as follows:
16 //
17 // +---------------------+
18 // boundary --->| |<--- boundary
19 // attribute 1 | | attribute 2
20 // (fixed) +---------------------+ (fixed, nonzero)
21 //
22 // The example demonstrates the use of block nonlinear operators
23 // (the class RubberOperator defining H(x)) as well as a nonlinear
24 // Newton solver for the quasi-static problem. Each Newton step
25 // requires the inversion of a Jacobian matrix, which is done
26 // through a (preconditioned) inner solver. The specialized block
27 // preconditioner is implemented as a user-defined solver.
28 //
29 // We recommend viewing examples 2, 5, and 10 before viewing this
30 // example.
31 
32 #include "mfem.hpp"
33 #include <memory>
34 #include <iostream>
35 #include <fstream>
36 
37 using namespace std;
38 using namespace mfem;
39 
40 // Custom block preconditioner for the Jacobian of the incompressible nonlinear
41 // elasticity operator. It has the form
42 //
43 // P^-1 = [ K^-1 0 ][ I -B^T ][ I 0 ]
44 // [ 0 I ][ 0 I ][ 0 -\gamma S^-1 ]
45 //
46 // where the original Jacobian has the form
47 //
48 // J = [ K B^T ]
49 // [ B 0 ]
50 //
51 // and K^-1 is an approximation of the inverse of the displacement part of the
52 // Jacobian and S^-1 is an approximation of the inverse of the Schur
53 // complement S = B K^-1 B^T. The Schur complement is approximated using
54 // a mass matrix of the pressure variables.
55 class JacobianPreconditioner : public Solver
56 {
57 protected:
58  // Finite element spaces for setting up preconditioner blocks
60 
61  // Offsets for extracting block vector segments
62  Array<int> &block_trueOffsets;
63 
64  // Jacobian for block access
65  BlockOperator *jacobian;
66 
67  // Scaling factor for the pressure mass matrix in the block preconditioner
68  double gamma;
69 
70  // Objects for the block preconditioner application
71  Operator *pressure_mass;
72  Solver *mass_pcg;
73  Solver *mass_prec;
74  Solver *stiff_pcg;
75  Solver *stiff_prec;
76 
77 public:
78  JacobianPreconditioner(Array<ParFiniteElementSpace *> &fes,
79  Operator &mass, Array<int> &offsets);
80 
81  virtual void Mult(const Vector &k, Vector &y) const;
82  virtual void SetOperator(const Operator &op);
83 
84  virtual ~JacobianPreconditioner();
85 };
86 
87 // After spatial discretization, the rubber model can be written as:
88 // 0 = H(x)
89 // where x is the block vector representing the deformation and pressure and
90 // H(x) is the nonlinear incompressible neo-Hookean operator.
91 class RubberOperator : public Operator
92 {
93 protected:
94  // Finite element spaces
96 
97  // Block nonlinear form
98  ParBlockNonlinearForm *Hform;
99 
100  // Pressure mass matrix for the preconditioner
101  Operator *pressure_mass;
102 
103  // Newton solver for the hyperelastic operator
104  NewtonSolver newton_solver;
105 
106  // Solver for the Jacobian solve in the Newton method
107  Solver *j_solver;
108 
109  // Preconditioner for the Jacobian
110  Solver *j_prec;
111 
112  // Shear modulus coefficient
113  Coefficient &mu;
114 
115  // Block offsets for variable access
116  Array<int> &block_trueOffsets;
117 
118 public:
119  RubberOperator(Array<ParFiniteElementSpace *> &fes, Array<Array<int> *>&ess_bdr,
120  Array<int> &block_trueOffsets, double rel_tol, double abs_tol,
121  int iter, Coefficient &mu);
122 
123  // Required to use the native newton solver
124  virtual Operator &GetGradient(const Vector &xp) const;
125  virtual void Mult(const Vector &k, Vector &y) const;
126 
127  // Driver for the newton solver
128  void Solve(Vector &xp) const;
129 
130  virtual ~RubberOperator();
131 };
132 
133 // Visualization driver
134 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
135  ParGridFunction *field, const char *field_name = NULL,
136  bool init_vis = false);
137 
138 // Configuration definition functions
139 void ReferenceConfiguration(const Vector &x, Vector &y);
140 void InitialDeformation(const Vector &x, Vector &y);
141 
142 
143 int main(int argc, char *argv[])
144 {
145  // 1. Initialize MPI
146  int num_procs, myid;
147  MPI_Init(&argc, &argv);
148  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
149  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
150 
151  // 2. Parse command-line options
152  const char *mesh_file = "../data/beam-hex.mesh";
153  int ser_ref_levels = 0;
154  int par_ref_levels = 0;
155  int order = 2;
156  bool visualization = true;
157  double newton_rel_tol = 1e-4;
158  double newton_abs_tol = 1e-6;
159  int newton_iter = 500;
160  double mu = 1.0;
161 
162  OptionsParser args(argc, argv);
163  args.AddOption(&mesh_file, "-m", "--mesh",
164  "Mesh file to use.");
165  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
166  "Number of times to refine the mesh uniformly in serial.");
167  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
168  "Number of times to refine the mesh uniformly in parallel.");
169  args.AddOption(&order, "-o", "--order",
170  "Order (degree) of the finite elements.");
171  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
172  "--no-visualization",
173  "Enable or disable GLVis visualization.");
174  args.AddOption(&newton_rel_tol, "-rel", "--relative-tolerance",
175  "Relative tolerance for the Newton solve.");
176  args.AddOption(&newton_abs_tol, "-abs", "--absolute-tolerance",
177  "Absolute tolerance for the Newton solve.");
178  args.AddOption(&newton_iter, "-it", "--newton-iterations",
179  "Maximum iterations for the Newton solve.");
180  args.AddOption(&mu, "-mu", "--shear-modulus",
181  "Shear modulus for the neo-Hookean material.");
182  args.Parse();
183  if (!args.Good())
184  {
185  if (myid == 0)
186  {
187  args.PrintUsage(cout);
188  }
189  MPI_Finalize();
190  return 1;
191  }
192  if (myid == 0)
193  {
194  args.PrintOptions(cout);
195  }
196 
197  // 3. Read the (serial) mesh from the given mesh file on all processors. We
198  // can handle triangular, quadrilateral, tetrahedral and hexahedral meshes
199  // with the same code.
200  Mesh *mesh = new Mesh(mesh_file, 1, 1);
201  int dim = mesh->Dimension();
202 
203  // 4. Refine the mesh in serial to increase the resolution. In this example
204  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
205  // a command-line parameter.
206  for (int lev = 0; lev < ser_ref_levels; lev++)
207  {
208  mesh->UniformRefinement();
209  }
210 
211  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
212  // this mesh further in parallel to increase the resolution. Once the
213  // parallel mesh is defined, the serial mesh can be deleted.
214  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
215  delete mesh;
216  for (int lev = 0; lev < par_ref_levels; lev++)
217  {
218  pmesh->UniformRefinement();
219  }
220 
221  // 6. Define the shear modulus for the incompressible Neo-Hookean material
222  ConstantCoefficient c_mu(mu);
223 
224  // 7. Define the finite element spaces for displacement and pressure
225  // (Taylor-Hood elements). By default, the displacement (u/x) is a second
226  // order vector field, while the pressure (p) is a linear scalar function.
227  H1_FECollection quad_coll(order, dim);
228  H1_FECollection lin_coll(order-1, dim);
229 
230  ParFiniteElementSpace R_space(pmesh, &quad_coll, dim, Ordering::byVDIM);
231  ParFiniteElementSpace W_space(pmesh, &lin_coll);
232 
234  spaces[0] = &R_space;
235  spaces[1] = &W_space;
236 
237  HYPRE_Int glob_R_size = R_space.GlobalTrueVSize();
238  HYPRE_Int glob_W_size = W_space.GlobalTrueVSize();
239 
240  // 8. Define the Dirichlet conditions (set to boundary attribute 1 and 2)
241  Array<Array<int> *> ess_bdr(2);
242 
243  Array<int> ess_bdr_u(R_space.GetMesh()->bdr_attributes.Max());
244  Array<int> ess_bdr_p(W_space.GetMesh()->bdr_attributes.Max());
245 
246  ess_bdr_p = 0;
247  ess_bdr_u = 0;
248  ess_bdr_u[0] = 1;
249  ess_bdr_u[1] = 1;
250 
251  ess_bdr[0] = &ess_bdr_u;
252  ess_bdr[1] = &ess_bdr_p;
253 
254  // 9. Print the mesh statistics
255  if (myid == 0)
256  {
257  std::cout << "***********************************************************\n";
258  std::cout << "dim(u) = " << glob_R_size << "\n";
259  std::cout << "dim(p) = " << glob_W_size << "\n";
260  std::cout << "dim(u+p) = " << glob_R_size + glob_W_size << "\n";
261  std::cout << "***********************************************************\n";
262  }
263 
264  // 10. Define the block structure of the solution vector (u then p)
265  Array<int> block_trueOffsets(3);
266  block_trueOffsets[0] = 0;
267  block_trueOffsets[1] = R_space.TrueVSize();
268  block_trueOffsets[2] = W_space.TrueVSize();
269  block_trueOffsets.PartialSum();
270 
271  BlockVector xp(block_trueOffsets);
272 
273  // 11. Define grid functions for the current configuration, reference
274  // configuration, final deformation, and pressure
275  ParGridFunction x_gf(&R_space);
276  ParGridFunction x_ref(&R_space);
277  ParGridFunction x_def(&R_space);
278  ParGridFunction p_gf(&W_space);
279 
282 
283  x_gf.ProjectCoefficient(deform);
284  x_ref.ProjectCoefficient(refconfig);
285  p_gf = 0.0;
286 
287  // 12. Set up the block solution vectors
288  x_gf.GetTrueDofs(xp.GetBlock(0));
289  p_gf.GetTrueDofs(xp.GetBlock(1));
290 
291  // 13. Initialize the incompressible neo-Hookean operator
292  RubberOperator oper(spaces, ess_bdr, block_trueOffsets,
293  newton_rel_tol, newton_abs_tol, newton_iter, c_mu);
294 
295  // 14. Solve the Newton system
296  oper.Solve(xp);
297 
298  // 15. Distribute the shared degrees of freedom
299  x_gf.Distribute(xp.GetBlock(0));
300  p_gf.Distribute(xp.GetBlock(1));
301 
302  // 16. Compute the final deformation
303  subtract(x_gf, x_ref, x_def);
304 
305  // 17. Visualize the results if requested
306  socketstream vis_u, vis_p;
307  if (visualization)
308  {
309  char vishost[] = "localhost";
310  int visport = 19916;
311  vis_u.open(vishost, visport);
312  vis_u.precision(8);
313  visualize(vis_u, pmesh, &x_gf, &x_def, "Deformation", true);
314  // Make sure all ranks have sent their 'u' solution before initiating
315  // another set of GLVis connections (one from each rank):
316  MPI_Barrier(pmesh->GetComm());
317  vis_p.open(vishost, visport);
318  vis_p.precision(8);
319  visualize(vis_p, pmesh, &x_gf, &p_gf, "Pressure", true);
320  }
321 
322  // 18. Save the displaced mesh, the final deformation, and the pressure
323  {
324  GridFunction *nodes = &x_gf;
325  int owns_nodes = 0;
326  pmesh->SwapNodes(nodes, owns_nodes);
327 
328  ostringstream mesh_name, pressure_name, deformation_name;
329  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
330  pressure_name << "pressure." << setfill('0') << setw(6) << myid;
331  deformation_name << "deformation." << setfill('0') << setw(6) << myid;
332 
333  ofstream mesh_ofs(mesh_name.str().c_str());
334  mesh_ofs.precision(8);
335  pmesh->Print(mesh_ofs);
336 
337  ofstream pressure_ofs(pressure_name.str().c_str());
338  pressure_ofs.precision(8);
339  p_gf.Save(pressure_ofs);
340 
341  ofstream deformation_ofs(deformation_name.str().c_str());
342  deformation_ofs.precision(8);
343  x_def.Save(deformation_ofs);
344  }
345 
346  // 19. Free the used memory
347  delete pmesh;
348 
349  MPI_Finalize();
350 
351  return 0;
352 }
353 
354 
355 JacobianPreconditioner::JacobianPreconditioner(Array<ParFiniteElementSpace *>
356  &fes,
357  Operator &mass,
358  Array<int> &offsets)
359  : Solver(offsets[2]), block_trueOffsets(offsets), pressure_mass(&mass)
360 {
361  fes.Copy(spaces);
362 
363  gamma = 0.00001;
364 
365  // The mass matrix and preconditioner do not change every Newton cycle, so
366  // we only need to define them once
367  HypreBoomerAMG *mass_prec_amg = new HypreBoomerAMG();
368  mass_prec_amg->SetPrintLevel(0);
369 
370  mass_prec = mass_prec_amg;
371 
372  CGSolver *mass_pcg_iter = new CGSolver(spaces[0]->GetComm());
373  mass_pcg_iter->SetRelTol(1e-12);
374  mass_pcg_iter->SetAbsTol(1e-12);
375  mass_pcg_iter->SetMaxIter(200);
376  mass_pcg_iter->SetPrintLevel(0);
377  mass_pcg_iter->SetPreconditioner(*mass_prec);
378  mass_pcg_iter->SetOperator(*pressure_mass);
379  mass_pcg_iter->iterative_mode = false;
380 
381  mass_pcg = mass_pcg_iter;
382 
383  // The stiffness matrix does change every Newton cycle, so we will define it
384  // during SetOperator
385  stiff_pcg = NULL;
386  stiff_prec = NULL;
387 }
388 
389 void JacobianPreconditioner::Mult(const Vector &k, Vector &y) const
390 {
391  // Extract the blocks from the input and output vectors
392  Vector disp_in(k.GetData() + block_trueOffsets[0],
393  block_trueOffsets[1]-block_trueOffsets[0]);
394  Vector pres_in(k.GetData() + block_trueOffsets[1],
395  block_trueOffsets[2]-block_trueOffsets[1]);
396 
397  Vector disp_out(y.GetData() + block_trueOffsets[0],
398  block_trueOffsets[1]-block_trueOffsets[0]);
399  Vector pres_out(y.GetData() + block_trueOffsets[1],
400  block_trueOffsets[2]-block_trueOffsets[1]);
401 
402  Vector temp(block_trueOffsets[1]-block_trueOffsets[0]);
403  Vector temp2(block_trueOffsets[1]-block_trueOffsets[0]);
404 
405  // Perform the block elimination for the preconditioner
406  mass_pcg->Mult(pres_in, pres_out);
407  pres_out *= -gamma;
408 
409  jacobian->GetBlock(0,1).Mult(pres_out, temp);
410  subtract(disp_in, temp, temp2);
411 
412  stiff_pcg->Mult(temp2, disp_out);
413 }
414 
415 void JacobianPreconditioner::SetOperator(const Operator &op)
416 {
417  jacobian = (BlockOperator *) &op;
418 
419  // Initialize the stiffness preconditioner and solver
420  if (stiff_prec == NULL)
421  {
422  HypreBoomerAMG *stiff_prec_amg = new HypreBoomerAMG();
423  stiff_prec_amg->SetPrintLevel(0);
424  stiff_prec_amg->SetElasticityOptions(spaces[0]);
425 
426  stiff_prec = stiff_prec_amg;
427 
428  GMRESSolver *stiff_pcg_iter = new GMRESSolver(spaces[0]->GetComm());
429  stiff_pcg_iter->SetRelTol(1e-8);
430  stiff_pcg_iter->SetAbsTol(1e-8);
431  stiff_pcg_iter->SetMaxIter(200);
432  stiff_pcg_iter->SetPrintLevel(0);
433  stiff_pcg_iter->SetPreconditioner(*stiff_prec);
434  stiff_pcg_iter->iterative_mode = false;
435 
436  stiff_pcg = stiff_pcg_iter;
437  }
438 
439  // At each Newton cycle, compute the new stiffness AMG preconditioner by
440  // updating the iterative solver which, in turn, updates its preconditioner
441  stiff_pcg->SetOperator(jacobian->GetBlock(0,0));
442 }
443 
444 JacobianPreconditioner::~JacobianPreconditioner()
445 {
446  delete mass_pcg;
447  delete mass_prec;
448  delete stiff_prec;
449  delete stiff_pcg;
450 }
451 
452 
453 RubberOperator::RubberOperator(Array<ParFiniteElementSpace *> &fes,
454  Array<Array<int> *> &ess_bdr,
455  Array<int> &trueOffsets,
456  double rel_tol,
457  double abs_tol,
458  int iter,
459  Coefficient &c_mu)
460  : Operator(fes[0]->TrueVSize() + fes[1]->TrueVSize()),
461  newton_solver(fes[0]->GetComm()), mu(c_mu), block_trueOffsets(trueOffsets)
462 {
463  Array<Vector *> rhs(2);
464  rhs = NULL; // Set all entries in the array
465 
466  fes.Copy(spaces);
467 
468  // Define the block nonlinear form
469  Hform = new ParBlockNonlinearForm(spaces);
470 
471  // Add the incompressible neo-Hookean integrator
472  Hform->AddDomainIntegrator(new IncompressibleNeoHookeanIntegrator(mu));
473 
474  // Set the essential boundary conditions
475  Hform->SetEssentialBC(ess_bdr, rhs);
476 
477  // Compute the pressure mass stiffness matrix
478  ParBilinearForm *a = new ParBilinearForm(spaces[1]);
479  ConstantCoefficient one(1.0);
480  OperatorHandle mass(Operator::Hypre_ParCSR);
481  a->AddDomainIntegrator(new MassIntegrator(one));
482  a->Assemble();
483  a->Finalize();
484  a->ParallelAssemble(mass);
485  delete a;
486 
487  mass.SetOperatorOwner(false);
488  pressure_mass = mass.Ptr();
489 
490  // Initialize the Jacobian preconditioner
491  JacobianPreconditioner *jac_prec =
492  new JacobianPreconditioner(fes, *pressure_mass, block_trueOffsets);
493  j_prec = jac_prec;
494 
495  // Set up the Jacobian solver
496  GMRESSolver *j_gmres = new GMRESSolver(spaces[0]->GetComm());
497  j_gmres->iterative_mode = false;
498  j_gmres->SetRelTol(1e-12);
499  j_gmres->SetAbsTol(1e-12);
500  j_gmres->SetMaxIter(300);
501  j_gmres->SetPrintLevel(0);
502  j_gmres->SetPreconditioner(*j_prec);
503  j_solver = j_gmres;
504 
505  // Set the newton solve parameters
506  newton_solver.iterative_mode = true;
507  newton_solver.SetSolver(*j_solver);
508  newton_solver.SetOperator(*this);
509  newton_solver.SetPrintLevel(1);
510  newton_solver.SetRelTol(rel_tol);
511  newton_solver.SetAbsTol(abs_tol);
512  newton_solver.SetMaxIter(iter);
513 }
514 
515 // Solve the Newton system
516 void RubberOperator::Solve(Vector &xp) const
517 {
518  Vector zero;
519  newton_solver.Mult(zero, xp);
520  MFEM_VERIFY(newton_solver.GetConverged(),
521  "Newton Solver did not converge.");
522 }
523 
524 // compute: y = H(x,p)
525 void RubberOperator::Mult(const Vector &k, Vector &y) const
526 {
527  Hform->Mult(k, y);
528 }
529 
530 // Compute the Jacobian from the nonlinear form
531 Operator &RubberOperator::GetGradient(const Vector &xp) const
532 {
533  return Hform->GetGradient(xp);
534 }
535 
536 RubberOperator::~RubberOperator()
537 {
538  delete Hform;
539  delete pressure_mass;
540  delete j_solver;
541  delete j_prec;
542 }
543 
544 
545 // Inline visualization
546 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
547  ParGridFunction *field, const char *field_name, bool init_vis)
548 {
549  if (!out)
550  {
551  return;
552  }
553 
554  GridFunction *nodes = deformed_nodes;
555  int owns_nodes = 0;
556 
557  mesh->SwapNodes(nodes, owns_nodes);
558 
559  out << "parallel " << mesh->GetNRanks() << " " << mesh->GetMyRank() << "\n";
560  out << "solution\n" << *mesh << *field;
561 
562  mesh->SwapNodes(nodes, owns_nodes);
563 
564  if (init_vis)
565  {
566  out << "window_size 800 800\n";
567  out << "window_title '" << field_name << "'\n";
568  if (mesh->SpaceDimension() == 2)
569  {
570  out << "view 0 0\n"; // view from top
571  out << "keys jlA\n"; // turn off perspective and light, +anti-aliasing
572  }
573  out << "keys cmA\n"; // show colorbar and mesh, +anti-aliasing
574  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
575  }
576  out << flush;
577 }
578 
579 void ReferenceConfiguration(const Vector &x, Vector &y)
580 {
581  // Set the reference, stress free, configuration
582  y = x;
583 }
584 
585 void InitialDeformation(const Vector &x, Vector &y)
586 {
587  // Set the initial configuration. Having this different from the reference
588  // configuration can help convergence
589  y = x;
590  y[1] = x[1] + 0.25*x[0];
591 }
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:376
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10.cpp:588
Conjugate gradient method.
Definition: solvers.hpp:111
void ReferenceConfiguration(const Vector &x, Vector &y)
Definition: ex19.cpp:525
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Subclass constant coefficient.
Definition: coefficient.hpp:57
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
Definition: operator.hpp:57
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:5460
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:478
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:190
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:398
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:290
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:272
int main(int argc, char *argv[])
Definition: ex1.cpp:45
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:129
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
int GetNRanks() const
Definition: pmesh.hpp:124
The BoomerAMG solver in hypre.
Definition: hypre.hpp:796
int dim
Definition: ex3.cpp:47
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:72
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:99
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:224
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
void SetPrintLevel(int print_level)
Definition: hypre.hpp:837
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:247
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:259
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3428
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2666
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int SpaceDimension() const
Definition: mesh.hpp:646
void SetAbsTol(double atol)
Definition: solvers.hpp:62
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
int GetMyRank() const
Definition: pmesh.hpp:125
void SetRelTol(double rtol)
Definition: solvers.hpp:61
MPI_Comm GetComm() const
Definition: pmesh.hpp:123
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:386
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:74
void PartialSum()
Partial Sum.
Definition: array.cpp:140
GMRES method.
Definition: solvers.hpp:143
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:131
A class representing a general parallel block nonlinear operator defined on the Cartesian product of ...
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:141
Class for parallel bilinear form.
int open(const char hostname[], int port)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:125
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:93
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
Base class for solvers.
Definition: operator.hpp:268
Class for parallel grid function.
Definition: pgridfunc.hpp:32
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
Abstract operator.
Definition: operator.hpp:21
A class to handle Block systems in a matrix-free implementation.
Class for parallel meshes.
Definition: pmesh.hpp:32
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:77
bool Good() const
Definition: optparser.hpp:120