MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
petsc.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 // Author: Stefano Zampini <stefano.zampini@gmail.com>
13 
14 #ifndef MFEM_PETSC
15 #define MFEM_PETSC
16 
17 #include "../config/config.hpp"
18 #include "linalg.hpp"
19 
20 #ifdef MFEM_USE_PETSC
21 #ifdef MFEM_USE_MPI
22 
23 #include <petsc.h>
24 #include <limits>
25 
26 namespace mfem
27 {
28 
29 class ParFiniteElementSpace;
30 class PetscParMatrix;
31 
32 /// Wrapper for PETSc's vector class
33 class PetscParVector : public Vector
34 {
35 protected:
36  /// The actual PETSc object
37  Vec x;
38 
39  friend class PetscParMatrix;
40  friend class PetscODESolver;
41  friend class PetscLinearSolver;
42  friend class PetscPreconditioner;
43  friend class PetscNonlinearSolver;
44  friend class PetscBDDCSolver;
45 
46  // Set Vector::data and Vector::size from x
47  void _SetDataAndSize_();
48 
49 public:
50  /// Creates vector with given global size and partitioning of the columns.
51  /** If @a col is provided, processor P owns columns [col[P],col[P+1]).
52  Otherwise, PETSc decides the partitioning */
53  PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscInt *col = NULL);
54 
55  /** @brief Creates vector with given global size, partitioning of the
56  columns, and data.
57 
58  The data must be allocated and destroyed outside. If @a _data is NULL, a
59  dummy vector without a valid data array will be created. */
60  PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscScalar *_data,
61  PetscInt *col);
62 
63  /// Creates vector compatible with @a y
65 
66  /// Creates a PetscParVector from a Vector (data is not copied)
67  PetscParVector(MPI_Comm comm, const Vector &_x);
68 
69  /** @brief Creates vector compatible with the Operator (i.e. in the domain
70  of) @a op or its adjoint. */
71  /** The argument @a allocate determines if the memory is actually allocated
72  to store the data. */
73  explicit PetscParVector(MPI_Comm comm, const Operator &op,
74  bool transpose = false, bool allocate = true);
75 
76  /// Creates vector compatible with (i.e. in the domain of) @a A or @a A^T
77  /** The argument @a allocate determines if the memory is actually allocated
78  to store the data. */
79  explicit PetscParVector(const PetscParMatrix &A, bool transpose = false,
80  bool allocate = true);
81 
82  /// Creates PetscParVector out of PETSc Vec object.
83  /** @param[in] y The PETSc Vec object.
84  @param[in] ref If true, we increase the reference count of @a y. */
85  explicit PetscParVector(Vec y, bool ref=false);
86 
87  /// Create a true dof parallel vector on a given ParFiniteElementSpace
88  explicit PetscParVector(ParFiniteElementSpace *pfes);
89 
90  /// Calls PETSc's destroy function
91  virtual ~PetscParVector();
92 
93  /// Get the associated MPI communicator
94  MPI_Comm GetComm() const { return PetscObjectComm((PetscObject)x); }
95 
96  /// Returns the global number of rows
97  PetscInt GlobalSize() const;
98 
99  /// Conversion function to PETSc's Vec type
100  operator Vec() const { return x; }
101 
102  /// Returns the global vector in each processor
103  Vector* GlobalVector() const;
104 
105  /// Set constant values
106  PetscParVector& operator= (PetscScalar d);
107 
108  /// Define '=' for PETSc vectors.
110 
111  /** @brief Temporarily replace the data of the PETSc Vec object. To return to
112  the original data array, call ResetArray().
113 
114  @note This method calls PETSc's VecPlaceArray() function.
115  @note The inherited Vector::data pointer is not affected by this call. */
116  void PlaceArray(PetscScalar *temp_data);
117 
118  /** @brief Reset the PETSc Vec object to use its default data. Call this
119  method after the use of PlaceArray().
120 
121  @note This method calls PETSc's VecResetArray() function. */
122  void ResetArray();
123 
124  /// Set random values
125  void Randomize(PetscInt seed);
126 
127  /// Prints the vector (to stdout if @a fname is NULL)
128  void Print(const char *fname = NULL, bool binary = false) const;
129 };
130 
131 
132 /// Wrapper for PETSc's matrix class
133 class PetscParMatrix : public Operator
134 {
135 protected:
136  /// The actual PETSc object
137  Mat A;
138 
139  /// Auxiliary vectors for typecasting
140  mutable PetscParVector *X, *Y;
141 
142  /// Initialize with defaults. Does not initialize inherited members.
143  void Init();
144 
145  /// Delete all owned data. Does not perform re-initialization with defaults.
146  void Destroy();
147 
148  /** @brief Creates a wrapper around a mfem::Operator @a op using PETSc's
149  MATSHELL object and returns the Mat in @a B.
150 
151  This does not take any reference to @a op, that should not be destroyed
152  until @a B is needed. */
153  void MakeWrapper(MPI_Comm comm, const Operator* op, Mat *B);
154 
155  /// Convert an mfem::Operator into a Mat @a B; @a op can be destroyed unless
156  /// tid == PETSC_MATSHELL
157  /// if op is a BlockOperator, the operator type is relevant to the individual
158  /// blocks
159  void ConvertOperator(MPI_Comm comm, const Operator& op, Mat *B,
160  Operator::Type tid);
161 
162  friend class PetscLinearSolver;
163  friend class PetscPreconditioner;
164 
165 private:
166  /// Constructs a block-diagonal Mat object
167  void BlockDiagonalConstructor(MPI_Comm comm, PetscInt *row_starts,
168  PetscInt *col_starts, SparseMatrix *diag,
169  bool assembled, Mat *A);
170 
171 public:
172  /// Create an empty matrix to be used as a reference to an existing matrix.
173  PetscParMatrix();
174 
175  /// Creates PetscParMatrix out of PETSc's Mat.
176  /** @param[in] a The PETSc Mat object.
177  @param[in] ref If true, we increase the reference count of @a a. */
178  PetscParMatrix(Mat a, bool ref=false);
179 
180  /** @brief Convert a HypreParMatrix @a ha to a PetscParMatrix in the given
181  PETSc format @a tid. */
182  /** The supported type ids are: Operator::PETSC_MATAIJ,
183  Operator::PETSC_MATIS, and Operator::PETSC_MATSHELL. */
184  explicit PetscParMatrix(const HypreParMatrix *ha,
186 
187  /** @brief Convert an mfem::Operator into a PetscParMatrix in the given PETSc
188  format @a tid. */
189  /** If @a tid is Operator::PETSC_MATSHELL and @a op is not a PetscParMatrix,
190  it converts any mfem::Operator @a op implementing Operator::Mult() and
191  Operator::MultTranspose() into a PetscParMatrix. The Operator @a op
192  should not be deleted while the constructed PetscParMatrix is used.
193 
194  Otherwise, it tries to convert the operator in PETSc's classes.
195 
196  In particular, if @a op is a BlockOperator, then a MATNEST Mat object is
197  created using @a tid as the type for the blocks.
198  Note that if @a op is already a PetscParMatrix of the same type as
199  @a tid, the resulting PetscParMatrix will share the same Mat object */
200  PetscParMatrix(MPI_Comm comm, const Operator *op,
202 
203  /// Creates block-diagonal square parallel matrix.
204  /** The block-diagonal is given by @a diag which must be in CSR format
205  (finalized). The new PetscParMatrix does not take ownership of any of the
206  input arrays. The type id @a tid can be either PETSC_MATAIJ (parallel
207  distributed CSR) or PETSC_MATIS. */
208  PetscParMatrix(MPI_Comm comm, PetscInt glob_size, PetscInt *row_starts,
209  SparseMatrix *diag, Operator::Type tid);
210 
211  /// Creates block-diagonal rectangular parallel matrix.
212  /** The block-diagonal is given by @a diag which must be in CSR format
213  (finalized). The new PetscParMatrix does not take ownership of any of the
214  input arrays. The type id @a tid can be either PETSC_MATAIJ (parallel
215  distributed CSR) or PETSC_MATIS. */
216  PetscParMatrix(MPI_Comm comm, PetscInt global_num_rows,
217  PetscInt global_num_cols, PetscInt *row_starts,
218  PetscInt *col_starts, SparseMatrix *diag,
219  Operator::Type tid);
220 
221  /// Calls PETSc's destroy function.
222  virtual ~PetscParMatrix() { Destroy(); }
223 
224  /// @name Assignment operators
225  ///@{
229  ///@}
230 
231  /// Matvec: @a y = @a a A @a x + @a b @a y.
232  void Mult(double a, const Vector &x, double b, Vector &y) const;
233 
234  /// Matvec transpose: @a y = @a a A^T @a x + @a b @a y.
235  void MultTranspose(double a, const Vector &x, double b, Vector &y) const;
236 
237  virtual void Mult(const Vector &x, Vector &y) const
238  { Mult(1.0, x, 0.0, y); }
239 
240  virtual void MultTranspose(const Vector &x, Vector &y) const
241  { MultTranspose(1.0, x, 0.0, y); }
242 
243  /// Get the associated MPI communicator
244  MPI_Comm GetComm() const { return PetscObjectComm((PetscObject)A); }
245 
246  /// Typecasting to PETSc's Mat type
247  operator Mat() const { return A; }
248 
249  /// Returns the local number of rows
250  PetscInt GetNumRows() const;
251 
252  /// Returns the local number of columns
253  PetscInt GetNumCols() const;
254 
255  /// Returns the global number of rows
256  PetscInt M() const;
257 
258  /// Returns the global number of columns
259  PetscInt N() const;
260 
261  /// Returns the global number of rows
262  PetscInt GetGlobalNumRows() const { return M(); }
263 
264  /// Returns the global number of columns
265  PetscInt GetGlobalNumCols() const { return N(); }
266 
267  /// Returns the number of nonzeros.
268  /** Differently from HYPRE, this call is collective on the communicator,
269  as this number is not stored inside PETSc, but needs to be computed. */
270  PetscInt NNZ() const;
271 
272  /// Returns the inner vector in the domain of A (it creates it if needed)
273  PetscParVector* GetX() const;
274 
275  /// Returns the inner vector in the range of A (it creates it if needed)
276  PetscParVector* GetY() const;
277 
278  /// Returns the transpose of the PetscParMatrix.
279  /** If @a action is false, the new matrix is constructed with the PETSc
280  function MatTranspose().
281  If @a action is true, then the matrix is not actually transposed.
282  Instead, an object that behaves like the transpose is returned. */
283  PetscParMatrix* Transpose(bool action = false);
284 
285  /// Prints the matrix (to stdout if fname is NULL)
286  void Print(const char *fname = NULL, bool binary = false) const;
287 
288  /// Scale all entries by s: A_scaled = s*A.
289  void operator*=(double s);
290 
291  /** @brief Eliminate rows and columns from the matrix, and rows from the
292  vector @a B. Modify @a B with the BC values in @a X. */
293  void EliminateRowsCols(const Array<int> &rows_cols, const PetscParVector &X,
294  PetscParVector &B);
295  void EliminateRowsCols(const Array<int> &rows_cols, const HypreParVector &X,
296  HypreParVector &B);
297 
298  /** @brief Eliminate rows and columns from the matrix and store the
299  eliminated elements in a new matrix Ae (returned).
300 
301  The sum of the modified matrix and the returned matrix, Ae, is equal to
302  the original matrix. */
303  PetscParMatrix* EliminateRowsCols(const Array<int> &rows_cols);
304 
305  /// Makes this object a reference to another PetscParMatrix
306  void MakeRef(const PetscParMatrix &master);
307 
308  /** @brief Release the PETSc Mat object. If @a dereference is true, decrement
309  the refcount of the Mat object. */
310  Mat ReleaseMat(bool dereference);
311 
312  Type GetType() const;
313 };
314 
315 /// Returns the matrix Rt^t * A * P
316 PetscParMatrix * RAP(PetscParMatrix *Rt, PetscParMatrix *A, PetscParMatrix *P);
317 
318 /// Returns the matrix P^t * A * P
319 PetscParMatrix * RAP(PetscParMatrix *A, PetscParMatrix *P);
320 
321 /** @brief Eliminate essential BC specified by @a ess_dof_list from the solution
322  @a X to the r.h.s. @a B.
323 
324  Here, @a A is a matrix with eliminated BC, while @a Ae is such that
325  (@a A + @a Ae) is the original (Neumann) matrix before elimination. */
326 void EliminateBC(PetscParMatrix &A, PetscParMatrix &Ae,
327  const Array<int> &ess_dof_list, const Vector &X, Vector &B);
328 
329 /// Helper class for handling essential boundary conditions.
331 {
332 public:
333  enum Type
334  {
336  CONSTANT, ///< Constant in time b.c.
338  };
339 
341  bctype(_type), setup(false), eval_t(0.0),
342  eval_t_cached(std::numeric_limits<double>::min()) {}
343  PetscBCHandler(Array<int>& ess_tdof_list, Type _type = ZERO);
344 
345  virtual ~PetscBCHandler() {}
346 
347  /// Returns the type of boundary conditions
348  Type GetType() const { return bctype; }
349 
350  /// Sets the type of boundary conditions
351  void SetType(enum Type _type) { bctype = _type; setup = false; }
352 
353  /// Boundary conditions evaluation
354  /** In the result vector, @a g, only values at the essential dofs need to be
355  set. */
356  virtual void Eval(double t, Vector &g)
357  { mfem_error("PetscBCHandler::Eval method not overloaded"); }
358 
359  /// Sets essential dofs (local, per-process numbering)
360  void SetTDofs(Array<int>& list);
361 
362  /// Gets essential dofs (local, per-process numbering)
363  Array<int>& GetTDofs() { return ess_tdof_list; }
364 
365  /// Sets the current time
366  void SetTime(double t) { eval_t = t; }
367 
368  /// SetUp the helper object, where @a n is the size of the solution vector
369  void SetUp(PetscInt n);
370 
371  /// y = x on ess_tdof_list_c and y = g (internally evaluated) on ess_tdof_list
372  void ApplyBC(const Vector &x, Vector &y);
373 
374  /// y = x-g on ess_tdof_list, the rest of y is unchanged
375  void FixResidualBC(const Vector& x, Vector& y);
376 
377 private:
378  enum Type bctype;
379  bool setup;
380 
381  double eval_t;
382  double eval_t_cached;
383  Vector eval_g;
384 
385  Array<int> ess_tdof_list; //Essential true dofs
386 };
387 
388 // Helper class for user-defined preconditioners that needs to be setup
390 {
391 private:
392  std::string name;
393 public:
394  PetscPreconditionerFactory(const std::string &_name = "MFEM Factory")
395  : name(_name) { }
396  const char* GetName() { return name.c_str(); }
397  virtual Solver *NewPreconditioner(const OperatorHandle& oh) = 0;
399 };
400 
401 // Forward declarations of helper classes
402 class PetscSolverMonitor;
403 
404 /// Abstract class for PETSc's solvers.
406 {
407 protected:
408  /// Boolean to handle SetFromOptions calls.
409  mutable bool clcustom;
410 
411  /// The actual PETSc object (KSP, PC, SNES or TS).
412  PetscObject obj;
413 
414  /// The class id of the actual PETSc object
415  PetscClassId cid;
416 
417  /// Right-hand side and solution vector
418  mutable PetscParVector *B, *X;
419 
420  /// Handler for boundary conditions
422 
423  /// Private context for solver
424  void *private_ctx;
425 
426  /// Boolean to handle SetOperator calls.
427  mutable bool operatorset;
428 
429 public:
430  /// Construct an empty PetscSolver. Initialize protected objects to NULL.
431  PetscSolver();
432 
433  /// Destroy the PetscParVectors allocated (if any).
434  virtual ~PetscSolver();
435 
436  /** @name Update of PETSc options.
437  The following Set methods can be used to update the internal PETSc
438  options.
439  @note They will be overwritten by the options in the input PETSc file. */
440  ///@{
441  void SetTol(double tol);
442  void SetRelTol(double tol);
443  void SetAbsTol(double tol);
444  void SetMaxIter(int max_iter);
445  void SetPrintLevel(int plev);
446  ///@}
447 
448  /// Customize object with options set
449  /** If @a customize is false, it disables any options customization. */
450  void Customize(bool customize = true) const;
451  int GetConverged();
452  int GetNumIterations();
453  double GetFinalNorm();
454 
455  /// Sets user-defined monitoring routine.
456  void SetMonitor(PetscSolverMonitor *ctx);
457 
458  /// Sets the object to handle essential boundary conditions
459  void SetBCHandler(PetscBCHandler *bch);
460 
461  /// Sets the object for the creation of the preconditioner
463 
464  /// Conversion function to PetscObject.
465  operator PetscObject() const { return obj; }
466 
467 protected:
468  /// These two methods handle creation and destructions of
469  /// private data for the Solver objects
470  void CreatePrivateContext();
471  void FreePrivateContext();
472 };
473 
474 
475 /// Abstract class for PETSc's linear solvers.
476 class PetscLinearSolver : public PetscSolver, public Solver
477 {
478 private:
479  /// Internal flag to handle HypreParMatrix conversion or not.
480  bool wrap;
481 
482 public:
483  PetscLinearSolver(MPI_Comm comm, const std::string &prefix = std::string(),
484  bool wrap = false);
486  const std::string &prefix = std::string());
487  /// Constructs a solver using a HypreParMatrix.
488  /** If @a wrap is true, then the MatMult ops of HypreParMatrix are wrapped.
489  No preconditioner can be automatically constructed from PETSc. If
490  @a wrap is false, the HypreParMatrix is converted into PETSc format. */
491  PetscLinearSolver(const HypreParMatrix &A, bool wrap = true,
492  const std::string &prefix = std::string());
493  virtual ~PetscLinearSolver();
494 
495  virtual void SetOperator(const Operator &op);
496 
497  /// Allows to prescribe a different operator (@a pop) to construct
498  /// the preconditioner
499  void SetOperator(const Operator &op, const Operator &pop);
500 
501  /// Sets the solver to perform preconditioning
502  void SetPreconditioner(Solver &precond);
503 
504  /// Application of the solver.
505  virtual void Mult(const Vector &b, Vector &x) const;
506 
507  /// Conversion function to PETSc's KSP type.
508  operator KSP() const { return (KSP)obj; }
509 };
510 
511 
513 {
514 public:
515  PetscPCGSolver(MPI_Comm comm, const std::string &prefix = std::string());
516  PetscPCGSolver(PetscParMatrix &A, const std::string &prefix = std::string());
517  PetscPCGSolver(HypreParMatrix &A,bool wrap=true,
518  const std::string &prefix = std::string());
519 };
520 
521 
522 /// Abstract class for PETSc's preconditioners.
523 class PetscPreconditioner : public PetscSolver, public Solver
524 {
525 public:
526  PetscPreconditioner(MPI_Comm comm,
527  const std::string &prefix = std::string());
529  const std::string &prefix = std::string());
530  PetscPreconditioner(MPI_Comm comm, Operator &op,
531  const std::string &prefix = std::string());
532  virtual ~PetscPreconditioner();
533 
534  virtual void SetOperator(const Operator &op);
535 
536  /// Application of the preconditioner.
537  virtual void Mult(const Vector &b, Vector &x) const;
538 
539  /// Conversion function to PETSc's PC type.
540  operator PC() const { return (PC)obj; }
541 };
542 
543 
544 /// Auxiliary class for BDDC customization.
546 {
547 protected:
553  friend class PetscBDDCSolver;
554 
555 public:
557  nat_dof(NULL), nat_dof_local(false)
558  {}
560 
561  /// Specify dofs on the essential boundary.
562  /** If @a loc is false, it is a list of true dofs in local ordering.
563  If @a loc is true, it is a marker for Vdofs in local ordering. */
564  void SetEssBdrDofs(const Array<int> *essdofs, bool loc = false)
565  {
566  ess_dof = essdofs;
567  ess_dof_local = loc;
568  }
569  /// Specify dofs on the natural boundary.
570  /** If @a loc is false, it is a list of true dofs in local ordering.
571  If @a loc is true, it is a marker for Vdofs in local ordering. */
572  void SetNatBdrDofs(const Array<int> *natdofs, bool loc = false)
573  {
574  nat_dof = natdofs;
575  nat_dof_local = loc;
576  }
577 };
578 
579 
581 {
582 private:
583  void BDDCSolverConstructor(const PetscBDDCSolverParams &opts);
584 
585 public:
586  PetscBDDCSolver(MPI_Comm comm, Operator &op,
588  const std::string &prefix = std::string());
591  const std::string &prefix = std::string());
592 };
593 
594 
596 {
597 public:
598  PetscFieldSplitSolver(MPI_Comm comm, Operator &op,
599  const std::string &prefix = std::string());
600 };
601 
602 
603 /// Abstract class for PETSc's nonlinear solvers.
605 {
606 public:
607  PetscNonlinearSolver(MPI_Comm comm,
608  const std::string &prefix = std::string());
609  PetscNonlinearSolver(MPI_Comm comm, Operator &op,
610  const std::string &prefix = std::string());
611  virtual ~PetscNonlinearSolver();
612 
613  /// Specification of the nonlinear operator.
614  virtual void SetOperator(const Operator &op);
615 
616  /// Specifies the desired format of the Jacobian in case a PetscParMatrix
617  /// is not returned by the GetGradient method
618  void SetJacobianType(Operator::Type type);
619 
620  /// Application of the solver.
621  virtual void Mult(const Vector &b, Vector &x) const;
622 
623  /// Conversion function to PETSc's SNES type.
624  operator SNES() const { return (SNES)obj; }
625 };
626 
627 
628 /// Abstract class for PETSc's ODE solvers.
629 class PetscODESolver : public PetscSolver, public ODESolver
630 {
631 public:
632  /// The type of the ODE. Use ODE_SOLVER_LINEAR if the jacobians
633  /// are linear and independent of time.
634  enum Type
635  {
638  };
639 
640  PetscODESolver(MPI_Comm comm, const std::string &prefix = std::string());
641  virtual ~PetscODESolver();
642 
643  /// Initialize the ODE solver.
644  virtual void Init(TimeDependentOperator &f_,
645  enum PetscODESolver::Type type);
647 
648  /// Specifies the desired format of the Jacobian in case a PetscParMatrix
649  /// is not returned by the GetGradient methods
650  void SetJacobianType(Operator::Type type);
651 
652  virtual void Step(Vector &x, double &t, double &dt);
653  virtual void Run(Vector &x, double &t, double &dt, double t_final);
654 
655  /// Conversion function to PETSc's TS type.
656  operator TS() const { return (TS)obj; }
657 };
658 
659 
660 /// Abstract class for monitoring PETSc's solvers.
662 {
663 public:
664  bool mon_sol;
665  bool mon_res;
666  PetscSolverMonitor(bool monitor_sol = false, bool monitor_res = true)
667  : mon_sol(monitor_sol), mon_res(monitor_res) {}
668  virtual ~PetscSolverMonitor() {}
669 
670  /// Monitor the solution vector x
671  virtual void MonitorSolution(PetscInt it, PetscReal norm, const Vector &x)
672  {
673  MFEM_ABORT("MonitorSolution() is not implemented!")
674  }
675 
676  /// Monitor the residual vector r
677  virtual void MonitorResidual(PetscInt it, PetscReal norm, const Vector &r)
678  {
679  MFEM_ABORT("MonitorResidual() is not implemented!")
680  }
681 };
682 
683 
684 } // namespace mfem
685 
686 #endif // MFEM_USE_MPI
687 #endif // MFEM_USE_PETSC
688 
689 #endif
PetscPreconditionerFactory(const std::string &_name="MFEM Factory")
Definition: petsc.hpp:394
void SetPreconditioner(Solver &precond)
Sets the solver to perform preconditioning.
Definition: petsc.cpp:2003
void SetPrintLevel(int plev)
Definition: petsc.cpp:1407
PetscParVector * B
Right-hand side and solution vector.
Definition: petsc.hpp:418
void CreatePrivateContext()
Definition: petsc.cpp:1679
void Print(const char *fname=NULL, bool binary=false) const
Prints the matrix (to stdout if fname is NULL)
Definition: petsc.cpp:1072
virtual void SetOperator(const Operator &op)
Specification of the nonlinear operator.
Definition: petsc.cpp:2616
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:523
PetscParMatrix & operator+=(const PetscParMatrix &B)
Definition: petsc.cpp:484
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: petsc.cpp:1839
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1558
bool clcustom
Boolean to handle SetFromOptions calls.
Definition: petsc.hpp:409
void ConvertOperator(MPI_Comm comm, const Operator &op, Mat *B, Operator::Type tid)
Definition: petsc.cpp:656
double GetFinalNorm()
Definition: petsc.cpp:1654
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:133
const Array< int > * nat_dof
Definition: petsc.hpp:551
void ApplyBC(const Vector &x, Vector &y)
y = x on ess_tdof_list_c and y = g (internally evaluated) on ess_tdof_list
Definition: petsc.cpp:1758
Abstract class for PETSc&#39;s linear solvers.
Definition: petsc.hpp:476
Vec x
The actual PETSc object.
Definition: petsc.hpp:37
Abstract class for PETSc&#39;s solvers.
Definition: petsc.hpp:405
virtual void Run(Vector &x, double &t, double &dt, double t_final)
Perform time integration from time t [in] to time tf [in].
Definition: petsc.cpp:2849
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: petsc.hpp:237
ParFiniteElementSpace * fespace
Definition: petsc.hpp:548
Type GetType() const
Returns the type of boundary conditions.
Definition: petsc.hpp:348
Base abstract class for time dependent operators.
Definition: operator.hpp:151
PetscInt GetNumCols() const
Returns the local number of columns.
Definition: petsc.cpp:366
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void Mult(double a, const Vector &x, double b, Vector &y) const
Matvec: y = a A x + b y.
Definition: petsc.cpp:1039
PetscBDDCSolver(MPI_Comm comm, Operator &op, const PetscBDDCSolverParams &opts=PetscBDDCSolverParams(), const std::string &prefix=std::string())
Definition: petsc.cpp:2531
virtual ~PetscSolver()
Destroy the PetscParVectors allocated (if any).
Definition: petsc.cpp:1318
PetscNonlinearSolver(MPI_Comm comm, const std::string &prefix=std::string())
Definition: petsc.cpp:2576
PetscInt GetGlobalNumRows() const
Returns the global number of rows.
Definition: petsc.hpp:262
Type GetType() const
Definition: petsc.cpp:1260
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Definition: petsc.cpp:2822
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
void Init()
Initialize with defaults. Does not initialize inherited members.
Definition: petsc.cpp:394
PetscPreconditioner(MPI_Comm comm, const std::string &prefix=std::string())
Definition: petsc.cpp:2116
PetscClassId cid
The class id of the actual PETSc object.
Definition: petsc.hpp:415
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Constant in time b.c.
Definition: petsc.hpp:336
void SetRelTol(double tol)
Definition: petsc.cpp:1330
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: petsc.hpp:240
PetscPCGSolver(MPI_Comm comm, const std::string &prefix=std::string())
Definition: petsc.cpp:2086
void Destroy()
Delete all owned data. Does not perform re-initialization with defaults.
Definition: petsc.cpp:907
void SetMonitor(PetscSolverMonitor *ctx)
Sets user-defined monitoring routine.
Definition: petsc.cpp:1479
void SetMaxIter(int max_iter)
Definition: petsc.cpp:1380
Abstract class for PETSc&#39;s ODE solvers.
Definition: petsc.hpp:629
Wrapper for PETSc&#39;s vector class.
Definition: petsc.hpp:33
void Customize(bool customize=true) const
Customize object with options set.
Definition: petsc.cpp:1555
PetscInt M() const
Returns the global number of rows.
Definition: petsc.cpp:373
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Definition: petsc.cpp:2667
PetscParVector * Y
Definition: petsc.hpp:140
virtual ~PetscLinearSolver()
Definition: petsc.cpp:2076
PetscInt GetNumRows() const
Returns the local number of rows.
Definition: petsc.cpp:359
void Randomize(PetscInt seed)
Set random values.
Definition: petsc.cpp:320
void Print(const char *fname=NULL, bool binary=false) const
Prints the vector (to stdout if fname is NULL)
Definition: petsc.cpp:332
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Definition: hypre.cpp:1612
void ResetArray()
Reset the PETSc Vec object to use its default data. Call this method after the use of PlaceArray()...
Definition: petsc.cpp:315
virtual void Init(TimeDependentOperator &f_, enum PetscODESolver::Type type)
Initialize the ODE solver.
Definition: petsc.cpp:2725
PetscParMatrix * Transpose(bool action=false)
Returns the transpose of the PetscParMatrix.
Definition: petsc.cpp:1020
PetscObject obj
The actual PETSc object (KSP, PC, SNES or TS).
Definition: petsc.hpp:412
Abstract class for monitoring PETSc&#39;s solvers.
Definition: petsc.hpp:661
Data type sparse matrix.
Definition: sparsemat.hpp:38
void SetJacobianType(Operator::Type type)
Definition: petsc.cpp:2661
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:146
Auxiliary class for BDDC customization.
Definition: petsc.hpp:545
virtual ~PetscPreconditioner()
Definition: petsc.cpp:2227
MPI_Comm GetComm() const
Get the associated MPI communicator.
Definition: petsc.hpp:244
void FixResidualBC(const Vector &x, Vector &y)
y = x-g on ess_tdof_list, the rest of y is unchanged
Definition: petsc.cpp:1783
PetscParMatrix & operator=(const PetscParMatrix &B)
Definition: petsc.cpp:468
void EliminateRowsCols(const Array< int > &rows_cols, const PetscParVector &X, PetscParVector &B)
Eliminate rows and columns from the matrix, and rows from the vector B. Modify B with the BC values i...
Definition: petsc.cpp:1219
virtual ~PetscParVector()
Calls PETSc&#39;s destroy function.
Definition: petsc.cpp:180
PetscFieldSplitSolver(MPI_Comm comm, Operator &op, const std::string &prefix=std::string())
Definition: petsc.cpp:2540
PetscBCHandler * bchandler
Handler for boundary conditions.
Definition: petsc.hpp:421
virtual Solver * NewPreconditioner(const OperatorHandle &oh)=0
PetscLinearSolver(MPI_Comm comm, const std::string &prefix=std::string(), bool wrap=false)
Definition: petsc.cpp:1804
void SetTime(double t)
Sets the current time.
Definition: petsc.hpp:366
ID for class PetscParMatrix, MATSHELL format.
Definition: operator.hpp:131
PetscParMatrix()
Create an empty matrix to be used as a reference to an existing matrix.
Definition: petsc.cpp:401
void SetType(enum Type _type)
Sets the type of boundary conditions.
Definition: petsc.hpp:351
Abstract class for PETSc&#39;s nonlinear solvers.
Definition: petsc.hpp:604
int GetNumIterations()
Definition: petsc.cpp:1621
void SetBCHandler(PetscBCHandler *bch)
Sets the object to handle essential boundary conditions.
Definition: petsc.cpp:1502
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
virtual void Eval(double t, Vector &g)
Boundary conditions evaluation.
Definition: petsc.hpp:356
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:124
PetscBCHandler(Type _type=ZERO)
Definition: petsc.hpp:340
PetscInt NNZ() const
Returns the number of nonzeros.
Definition: petsc.cpp:387
virtual void Mult(const Vector &b, Vector &x) const
Application of the preconditioner.
Definition: petsc.cpp:2197
Mat A
The actual PETSc object.
Definition: petsc.hpp:137
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: petsc.hpp:646
PetscInt N() const
Returns the global number of columns.
Definition: petsc.cpp:380
PetscParVector * GetY() const
Returns the inner vector in the range of A (it creates it if needed)
Definition: petsc.cpp:1010
void SetJacobianType(Operator::Type type)
Definition: petsc.cpp:2816
void MakeWrapper(MPI_Comm comm, const Operator *op, Mat *B)
Creates a wrapper around a mfem::Operator op using PETSc&#39;s MATSHELL object and returns the Mat in B...
Definition: petsc.cpp:635
bool operatorset
Boolean to handle SetOperator calls.
Definition: petsc.hpp:427
PetscParVector * X
Auxiliary vectors for typecasting.
Definition: petsc.hpp:140
void MultTranspose(double a, const Vector &x, double b, Vector &y) const
Matvec transpose: y = a A^T x + b y.
Definition: petsc.cpp:1055
virtual ~PetscSolverMonitor()
Definition: petsc.hpp:668
virtual ~PetscNonlinearSolver()
Definition: petsc.cpp:2608
virtual ~PetscODESolver()
Definition: petsc.cpp:2717
virtual ~PetscBCHandler()
Definition: petsc.hpp:345
int GetConverged()
Definition: petsc.cpp:1588
PetscSolver()
Construct an empty PetscSolver. Initialize protected objects to NULL.
Definition: petsc.cpp:1308
void SetPreconditionerFactory(PetscPreconditionerFactory *factory)
Sets the object for the creation of the preconditioner.
Definition: petsc.cpp:1521
Mat ReleaseMat(bool dereference)
Release the PETSc Mat object. If dereference is true, decrement the refcount of the Mat object...
Definition: petsc.cpp:1247
PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscInt *col=NULL)
Creates vector with given global size and partitioning of the columns.
Definition: petsc.cpp:162
virtual void MonitorSolution(PetscInt it, PetscReal norm, const Vector &x)
Monitor the solution vector x.
Definition: petsc.hpp:671
ID for class PetscParMatrix, MATAIJ format.
Definition: operator.hpp:129
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: petsc.cpp:2151
void SetSpace(ParFiniteElementSpace *fe)
Definition: petsc.hpp:559
MPI_Comm GetComm() const
Get the associated MPI communicator.
Definition: petsc.hpp:94
PetscInt GlobalSize() const
Returns the global number of rows.
Definition: petsc.cpp:147
PetscParVector * X
Definition: petsc.hpp:418
const Array< int > * ess_dof
Definition: petsc.hpp:549
void _SetDataAndSize_()
Definition: petsc.cpp:136
void SetAbsTol(double tol)
Definition: petsc.cpp:1355
Helper class for handling essential boundary conditions.
Definition: petsc.hpp:330
PetscParVector * GetX() const
Returns the inner vector in the domain of A (it creates it if needed)
Definition: petsc.cpp:1000
Array< int > & GetTDofs()
Gets essential dofs (local, per-process numbering)
Definition: petsc.hpp:363
void SetTol(double tol)
Definition: petsc.cpp:1325
PetscParVector & operator=(PetscScalar d)
Set constant values.
Definition: petsc.cpp:298
PetscSolverMonitor(bool monitor_sol=false, bool monitor_res=true)
Definition: petsc.hpp:666
virtual void MonitorResidual(PetscInt it, PetscReal norm, const Vector &r)
Monitor the residual vector r.
Definition: petsc.hpp:677
Vector data type.
Definition: vector.hpp:48
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array...
Definition: petsc.cpp:310
void FreePrivateContext()
Definition: petsc.cpp:1708
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: petsc.cpp:273
PetscODESolver(MPI_Comm comm, const std::string &prefix=std::string())
Definition: petsc.cpp:2692
Base class for solvers.
Definition: operator.hpp:268
PetscInt GetGlobalNumCols() const
Returns the global number of columns.
Definition: petsc.hpp:265
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Definition: petsc.hpp:572
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Definition: petsc.cpp:2043
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: petsc.cpp:1034
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
Definition: petsc.hpp:564
void SetTDofs(Array< int > &list)
Sets essential dofs (local, per-process numbering)
Definition: petsc.cpp:1735
void * private_ctx
Private context for solver.
Definition: petsc.hpp:424
void SetUp(PetscInt n)
SetUp the helper object, where n is the size of the solution vector.
Definition: petsc.cpp:1742
void MakeRef(const PetscParMatrix &master)
Makes this object a reference to another PetscParMatrix.
Definition: petsc.cpp:990
virtual ~PetscParMatrix()
Calls PETSc&#39;s destroy function.
Definition: petsc.hpp:222