17 #include "../config/config.hpp"
29 class ParFiniteElementSpace;
53 PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscInt *col = NULL);
60 PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscScalar *_data,
74 bool transpose =
false,
bool allocate =
true);
80 bool allocate =
true);
94 MPI_Comm
GetComm()
const {
return PetscObjectComm((PetscObject)
x); }
100 operator Vec()
const {
return x; }
128 void Print(
const char *fname = NULL,
bool binary =
false)
const;
167 void BlockDiagonalConstructor(MPI_Comm comm, PetscInt *row_starts,
169 bool assembled, Mat *
A);
208 PetscParMatrix(MPI_Comm comm, PetscInt glob_size, PetscInt *row_starts,
217 PetscInt global_num_cols, PetscInt *row_starts,
238 {
Mult(1.0, x, 0.0, y); }
244 MPI_Comm
GetComm()
const {
return PetscObjectComm((PetscObject)
A); }
247 operator Mat()
const {
return A; }
270 PetscInt
NNZ()
const;
286 void Print(
const char *fname = NULL,
bool binary =
false)
const;
316 PetscParMatrix *
RAP(PetscParMatrix *Rt, PetscParMatrix *A, PetscParMatrix *P);
319 PetscParMatrix *
RAP(PetscParMatrix *A, PetscParMatrix *P);
326 void EliminateBC(PetscParMatrix &A, PetscParMatrix &Ae,
327 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
341 bctype(_type), setup(false), eval_t(0.0),
342 eval_t_cached(std::numeric_limits<double>::min()) {}
357 {
mfem_error(
"PetscBCHandler::Eval method not overloaded"); }
369 void SetUp(PetscInt n);
382 double eval_t_cached;
396 const char*
GetName() {
return name.c_str(); }
402 class PetscSolverMonitor;
450 void Customize(
bool customize =
true)
const;
465 operator PetscObject()
const {
return obj; }
486 const std::string &prefix = std::string());
492 const std::string &prefix = std::string());
508 operator KSP()
const {
return (KSP)
obj; }
515 PetscPCGSolver(MPI_Comm comm,
const std::string &prefix = std::string());
518 const std::string &prefix = std::string());
527 const std::string &prefix = std::string());
529 const std::string &prefix = std::string());
531 const std::string &prefix = std::string());
540 operator PC()
const {
return (PC)
obj; }
588 const std::string &prefix = std::string());
591 const std::string &prefix = std::string());
599 const std::string &prefix = std::string());
608 const std::string &prefix = std::string());
610 const std::string &prefix = std::string());
624 operator SNES()
const {
return (SNES)
obj; }
640 PetscODESolver(MPI_Comm comm,
const std::string &prefix = std::string());
652 virtual void Step(
Vector &x,
double &t,
double &dt);
653 virtual void Run(
Vector &x,
double &t,
double &dt,
double t_final);
656 operator TS()
const {
return (TS)
obj; }
673 MFEM_ABORT(
"MonitorSolution() is not implemented!")
679 MFEM_ABORT(
"MonitorResidual() is not implemented!")
686 #endif // MFEM_USE_MPI
687 #endif // MFEM_USE_PETSC
PetscPreconditionerFactory(const std::string &_name="MFEM Factory")
void SetPreconditioner(Solver &precond)
Sets the solver to perform preconditioning.
void SetPrintLevel(int plev)
PetscParVector * B
Right-hand side and solution vector.
void CreatePrivateContext()
void Print(const char *fname=NULL, bool binary=false) const
Prints the matrix (to stdout if fname is NULL)
virtual void SetOperator(const Operator &op)
Specification of the nonlinear operator.
Abstract class for PETSc's preconditioners.
PetscParMatrix & operator+=(const PetscParMatrix &B)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
bool clcustom
Boolean to handle SetFromOptions calls.
void ConvertOperator(MPI_Comm comm, const Operator &op, Mat *B, Operator::Type tid)
Wrapper for PETSc's matrix class.
const Array< int > * nat_dof
void ApplyBC(const Vector &x, Vector &y)
y = x on ess_tdof_list_c and y = g (internally evaluated) on ess_tdof_list
Abstract class for PETSc's linear solvers.
Vec x
The actual PETSc object.
Abstract class for PETSc's solvers.
virtual void Run(Vector &x, double &t, double &dt, double t_final)
Perform time integration from time t [in] to time tf [in].
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
ParFiniteElementSpace * fespace
Type GetType() const
Returns the type of boundary conditions.
Base abstract class for time dependent operators.
PetscInt GetNumCols() const
Returns the local number of columns.
Pointer to an Operator of a specified type.
void Mult(double a, const Vector &x, double b, Vector &y) const
Matvec: y = a A x + b y.
PetscBDDCSolver(MPI_Comm comm, Operator &op, const PetscBDDCSolverParams &opts=PetscBDDCSolverParams(), const std::string &prefix=std::string())
virtual ~PetscSolver()
Destroy the PetscParVectors allocated (if any).
PetscNonlinearSolver(MPI_Comm comm, const std::string &prefix=std::string())
PetscInt GetGlobalNumRows() const
Returns the global number of rows.
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]...
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
void Init()
Initialize with defaults. Does not initialize inherited members.
PetscPreconditioner(MPI_Comm comm, const std::string &prefix=std::string())
PetscClassId cid
The class id of the actual PETSc object.
Abstract parallel finite element space.
void SetRelTol(double tol)
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 ...
PetscPCGSolver(MPI_Comm comm, const std::string &prefix=std::string())
void Destroy()
Delete all owned data. Does not perform re-initialization with defaults.
void SetMonitor(PetscSolverMonitor *ctx)
Sets user-defined monitoring routine.
void SetMaxIter(int max_iter)
Abstract class for PETSc's ODE solvers.
Wrapper for PETSc's vector class.
void Customize(bool customize=true) const
Customize object with options set.
PetscInt M() const
Returns the global number of rows.
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
virtual ~PetscLinearSolver()
PetscInt GetNumRows() const
Returns the local number of rows.
void Randomize(PetscInt seed)
Set random values.
void Print(const char *fname=NULL, bool binary=false) const
Prints the vector (to stdout if fname is NULL)
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
void ResetArray()
Reset the PETSc Vec object to use its default data. Call this method after the use of PlaceArray()...
virtual void Init(TimeDependentOperator &f_, enum PetscODESolver::Type type)
Initialize the ODE solver.
PetscParMatrix * Transpose(bool action=false)
Returns the transpose of the PetscParMatrix.
PetscObject obj
The actual PETSc object (KSP, PC, SNES or TS).
Abstract class for monitoring PETSc's solvers.
void SetJacobianType(Operator::Type type)
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Auxiliary class for BDDC customization.
virtual ~PetscPreconditioner()
MPI_Comm GetComm() const
Get the associated MPI communicator.
void FixResidualBC(const Vector &x, Vector &y)
y = x-g on ess_tdof_list, the rest of y is unchanged
PetscParMatrix & operator=(const PetscParMatrix &B)
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...
virtual ~PetscParVector()
Calls PETSc's destroy function.
PetscFieldSplitSolver(MPI_Comm comm, Operator &op, const std::string &prefix=std::string())
PetscBCHandler * bchandler
Handler for boundary conditions.
virtual Solver * NewPreconditioner(const OperatorHandle &oh)=0
PetscLinearSolver(MPI_Comm comm, const std::string &prefix=std::string(), bool wrap=false)
void SetTime(double t)
Sets the current time.
ID for class PetscParMatrix, MATSHELL format.
PetscParMatrix()
Create an empty matrix to be used as a reference to an existing matrix.
void SetType(enum Type _type)
Sets the type of boundary conditions.
Abstract class for PETSc's nonlinear solvers.
void SetBCHandler(PetscBCHandler *bch)
Sets the object to handle essential boundary conditions.
Wrapper for hypre's parallel vector class.
virtual void Eval(double t, Vector &g)
Boundary conditions evaluation.
Type
Enumeration defining IDs for some classes derived from Operator.
PetscBCHandler(Type _type=ZERO)
PetscInt NNZ() const
Returns the number of nonzeros.
virtual void Mult(const Vector &b, Vector &x) const
Application of the preconditioner.
Mat A
The actual PETSc object.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
PetscInt N() const
Returns the global number of columns.
PetscParVector * GetY() const
Returns the inner vector in the range of A (it creates it if needed)
void SetJacobianType(Operator::Type type)
void MakeWrapper(MPI_Comm comm, const Operator *op, Mat *B)
Creates a wrapper around a mfem::Operator op using PETSc's MATSHELL object and returns the Mat in B...
bool operatorset
Boolean to handle SetOperator calls.
PetscParVector * X
Auxiliary vectors for typecasting.
void MultTranspose(double a, const Vector &x, double b, Vector &y) const
Matvec transpose: y = a A^T x + b y.
virtual ~PetscPreconditionerFactory()
virtual ~PetscSolverMonitor()
virtual ~PetscNonlinearSolver()
virtual ~PetscODESolver()
virtual ~PetscBCHandler()
PetscSolver()
Construct an empty PetscSolver. Initialize protected objects to NULL.
void SetPreconditionerFactory(PetscPreconditionerFactory *factory)
Sets the object for the creation of the preconditioner.
Mat ReleaseMat(bool dereference)
Release the PETSc Mat object. If dereference is true, decrement the refcount of the Mat object...
PetscParVector(MPI_Comm comm, PetscInt glob_size, PetscInt *col=NULL)
Creates vector with given global size and partitioning of the columns.
virtual void MonitorSolution(PetscInt it, PetscReal norm, const Vector &x)
Monitor the solution vector x.
ID for class PetscParMatrix, MATAIJ format.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetSpace(ParFiniteElementSpace *fe)
MPI_Comm GetComm() const
Get the associated MPI communicator.
PetscInt GlobalSize() const
Returns the global number of rows.
const Array< int > * ess_dof
void SetAbsTol(double tol)
Helper class for handling essential boundary conditions.
PetscParVector * GetX() const
Returns the inner vector in the domain of A (it creates it if needed)
Array< int > & GetTDofs()
Gets essential dofs (local, per-process numbering)
PetscParVector & operator=(PetscScalar d)
Set constant values.
PetscSolverMonitor(bool monitor_sol=false, bool monitor_res=true)
virtual void MonitorResidual(PetscInt it, PetscReal norm, const Vector &r)
Monitor the residual vector r.
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array...
void FreePrivateContext()
Vector * GlobalVector() const
Returns the global vector in each processor.
PetscODESolver(MPI_Comm comm, const std::string &prefix=std::string())
PetscInt GetGlobalNumCols() const
Returns the global number of columns.
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Wrapper for hypre's ParCSR matrix class.
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
void SetTDofs(Array< int > &list)
Sets essential dofs (local, per-process numbering)
void * private_ctx
Private context for solver.
void SetUp(PetscInt n)
SetUp the helper object, where n is the size of the solution vector.
void MakeRef(const PetscParMatrix &master)
Makes this object a reference to another PetscParMatrix.
virtual ~PetscParMatrix()
Calls PETSc's destroy function.