17 #include "../config/config.hpp" 27 #include "../general/mem_manager.hpp" 29 #include "petscconf.h" 30 #if !defined(PETSC_USE_REAL_DOUBLE) 31 #error "MFEM does not work with PETSc compiled without double precision" 33 #if defined(PETSC_USE_COMPLEX) 34 #error "MFEM does not work with PETSc compiled with complex numbers support" 36 #if defined(PETSC_USE_64BIT_INDICES) && !defined(HYPRE_BIGINT) 37 #error "Mismatch between HYPRE (32bit) and PETSc (64bit) integer types" 39 #if !defined(PETSC_USE_64BIT_INDICES) && defined(HYPRE_BIGINT) 40 #error "Mismatch between HYPRE (64bit) and PETSc (32bit) integer types" 43 #include "petscversion.h" 44 #if PETSC_VERSION_GE(3,12,0) 45 #include "petscsystypes.h" 69 typedef struct ::_p_Vec *
Vec;
70 typedef struct ::_p_Mat *
Mat;
71 typedef struct ::_p_KSP *
KSP;
72 typedef struct ::_p_PC *
PC;
73 typedef struct ::_p_SNES *
SNES;
74 typedef struct ::_p_TS *
TS;
100 int size_,
bool usedev_)
110 bool read_,
bool write_,
bool usedev_)
150 class ParFiniteElementSpace;
151 class PetscParMatrix;
205 bool transpose =
false,
bool allocate =
true);
211 bool allocate =
true);
301 void Print(
const char *fname = NULL,
bool binary =
false)
const;
303 const double *
Read(
bool=
true)
const override;
304 const double *
HostRead()
const override;
305 double *
Write(
bool=
true)
override;
349 void BlockDiagonalConstructor(MPI_Comm comm,
PetscInt *row_starts,
353 void SetUpForDevice();
445 {
Mult(1.0, x, 0.0, y); }
451 const double a = 1.0)
const override 452 {
Mult(
a, x, 1.0, y); }
455 const double a = 1.0)
const override 510 void Print(
const char *fname = NULL,
bool binary =
false)
const;
564 PetscParMatrix *
ParMult(
const PetscParMatrix *A,
const PetscParMatrix *B);
567 PetscParMatrix *
RAP(PetscParMatrix *Rt, PetscParMatrix *A, PetscParMatrix *P);
574 PetscParMatrix *
RAP(PetscParMatrix *A, PetscParMatrix *P);
577 PetscParMatrix *
RAP(HypreParMatrix *A, PetscParMatrix *P);
584 void EliminateBC(PetscParMatrix &A, PetscParMatrix &Ae,
585 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
599 bctype(type_), setup(false), eval_t(0.0),
600 eval_t_cached(
std::numeric_limits<double>::min()) {}
615 {
mfem_error(
"PetscBCHandler::Eval method not overloaded"); }
649 double eval_t_cached;
663 const char*
GetName() {
return name.c_str(); }
669 class PetscSolverMonitor;
717 void Customize(
bool customize =
true)
const;
755 bool wrap =
true,
bool iter_mode =
false);
757 const std::string &prefix = std::string(),
bool iter_mode =
false);
764 const std::string &prefix = std::string(),
bool iter_mode =
false);
791 PetscPCGSolver(MPI_Comm comm,
const std::string &prefix = std::string(),
792 bool iter_mode =
false);
794 bool iter_mode =
false);
796 const std::string &prefix = std::string(),
bool iter_mode =
false);
808 const std::string &prefix = std::string());
810 const std::string &prefix = std::string());
812 const std::string &prefix = std::string());
876 const std::string &prefix = std::string());
879 const std::string &prefix = std::string());
887 const std::string &prefix = std::string());
898 const std::string &prefix = std::string());
907 const std::string &prefix = std::string());
909 const std::string &prefix = std::string());
930 Vector &W,
bool &changed_y,
bool &changed_w));
957 PetscODESolver(MPI_Comm comm,
const std::string &prefix = std::string());
972 virtual void Step(
Vector &x,
double &
t,
double &dt);
973 virtual void Run(
Vector &x,
double &
t,
double &dt,
double t_final);
993 MFEM_ABORT(
"MonitorSolution() is not implemented!")
999 MFEM_ABORT(
"MonitorResidual() is not implemented!")
1009 #endif // MFEM_USE_MPI 1010 #endif // MFEM_USE_PETSC const double * GetDevicePointer() const
void SetPreconditioner(Solver &precond)
void SetPrintLevel(int plev)
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
PetscParVector * B
Right-hand side and solution vector.
void CreatePrivateContext()
void trans(const Vector &u, Vector &x)
void EliminateBC(const HypreParMatrix &A, const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s. B.
void SetUpdate(void(*update)(Operator *op, int it, const Vector &F, const Vector &X, const Vector &D, const Vector &P))
virtual void SetOperator(const Operator &op)
Specification of the nonlinear operator.
Abstract class for PETSc's preconditioners.
PetscParMatrix & operator+=(const PetscParMatrix &B)
petsc::Vec x
The actual PETSc object.
virtual void SetOperator(const Operator &op)
bool clcustom
Boolean to handle SetFromOptions calls.
void SetObjective(void(*obj)(Operator *op, const Vector &x, double *f))
Specification of an objective function to be used for line search.
PetscODESolver::Type GetType() const
void MakeAliasForSync(Memory< double > &base_, int offset_, int size_, bool read_, bool write_, bool usedev_)
petsc::Mat A
The actual PETSc object.
Wrapper for PETSc's matrix class.
PetscParVector & operator+=(const PetscParVector &y)
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.
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].
void SetPostCheck(void(*post)(Operator *op, const Vector &X, Vector &Y, Vector &W, bool &changed_y, bool &changed_w))
PetscBCHandler(Type type_=ZERO)
ParFiniteElementSpace * fespace
Base abstract class for first order time dependent operators.
void MakeAliasForSync(const Memory< double > &base_, int offset_, int size_, bool usedev_)
void SetType(PetscODESolver::Type)
void ConvertOperator(MPI_Comm comm, const Operator &op, petsc::Mat *B, Operator::Type tid)
virtual void MonitorSolver(PetscSolver *solver)
Generic monitor to take access to the solver.
Pointer to an Operator of a specified type.
PetscBDDCSolver(MPI_Comm comm, Operator &op, const PetscBDDCSolverParams &opts=PetscBDDCSolverParams(), const std::string &prefix=std::string())
PetscParVector & AddValues(const Array< PetscInt > &, const Array< PetscScalar > &)
Add values in a vector.
virtual ~PetscSolver()
Destroy the PetscParVectors allocated (if any).
PetscNonlinearSolver(MPI_Comm comm, const std::string &prefix=std::string())
PetscParVector * GetX() const
Returns the inner vector in the domain of A (it creates it if needed)
void SetBlockSize(PetscInt rbs, PetscInt cbs=-1)
Set row and column block sizes of a matrix.
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.
void MultTranspose(double a, const Vector &x, double b, Vector &y) const
Matvec transpose: y = a A^T x + b y.
PetscInt GetRowStart() const
Returns the global index of the first local row.
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Abstract parallel finite element space.
void SetRelTol(double tol)
PetscInt M() const
Returns the global number of rows.
Wrapper for syncing PETSc's vector memory.
void Randomize(PetscInt seed=0)
Set random values.
void Destroy()
Delete all owned data. Does not perform re-initialization with defaults.
void AddMult(const Vector &x, Vector &y, const double a=1.0) const override
Operator application: y+=A(x) (default) or y+=a*A(x).
void SetMonitor(PetscSolverMonitor *ctx)
Sets user-defined monitoring routine.
void SetMaxIter(int max_iter)
MPI_Comm GetComm() const
Get the associated MPI communicator.
Abstract class for PETSc's ODE solvers.
PetscParVector & operator*=(PetscScalar d)
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
void EliminateRowsCols(const Array< int > &rows_cols, const PetscParVector &X, PetscParVector &B, double diag=1.)
Eliminate rows and columns from the matrix, and rows from the vector B. Modify B with the BC values i...
const double * HostRead() const override
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
unsigned flags
Bit flags defined from the FlagMask enum.
void Print(const char *fname=NULL, bool binary=false) const
Prints the matrix (to stdout if fname is NULL)
virtual ~PetscLinearSolver()
std::function< double(const Vector &)> f(double mass_coeff)
void ResetArray()
Reset the PETSc Vec object to use its default data. Call this method after the use of PlaceArray()...
PetscParVector & SetValues(const Array< PetscInt > &, const Array< PetscScalar > &)
Set values in a vector.
MPI_Comm GetComm() const
Get the associated MPI communicator.
double * Write(bool=true) override
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
virtual void Init(TimeDependentOperator &f_, enum PetscODESolver::Type type)
Initialize the ODE solver.
PetscParMatrix * Transpose(bool action=false)
Returns the transpose of the PetscParMatrix.
void write(std::ostream &os, T value)
Write 'value' to stream.
PetscInt GetGlobalNumRows() const
Returns the global number of rows.
PetscObject obj
The actual PETSc object (KSP, PC, SNES or TS).
Abstract class for monitoring PETSc's solvers.
bool DeviceRequested() const
void Sync(const Memory &other) const
Copy the host/device pointer validity flags from other to *this.
PetscLinearSolver(MPI_Comm comm, const std::string &prefix=std::string(), bool wrap=true, bool iter_mode=false)
struct s_NavierContext ctx
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.
PetscInt GetNumRows() const
Returns the local number of rows.
virtual ~PetscPreconditioner()
void FixResidualBC(const Vector &x, Vector &y)
y = x-g on ess_tdof_list, the rest of y is unchanged
bool UseDevice() const override
Return the device flag of the Memory object used by the Vector.
void SetType(enum Type type_)
Sets the type of boundary conditions.
PetscParMatrix & operator=(const PetscParMatrix &B)
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.
PetscParMatrix & operator-=(const PetscParMatrix &B)
virtual Solver * NewPreconditioner(const OperatorHandle &oh)=0
double * ReadWrite(bool=true) override
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
void Zero(Vector &x)
Replace boundary dofs with 0.
T read(std::istream &is)
Read a value from the stream and return it.
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void Shift(double s)
Shift diagonal by a constant.
PetscInt GlobalSize() const
Returns the global number of rows.
PetscParVector * GetY() const
Returns the inner vector in the range of A (it creates it if needed)
double * HostReadWrite() override
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
void UpdateVecFromFlags()
Update PETSc's Vec after having accessed its data via GetMemory()
void SetTime(double t)
Sets the current time.
ID for class PetscParMatrix, MATSHELL format.
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
void SetMat(petsc::Mat newA)
Replace the inner Mat Object. The reference count of newA is increased.
PetscParMatrix()
Create an empty matrix to be used as a reference to an existing matrix.
void MakeWrapper(MPI_Comm comm, const Operator *op, petsc::Mat *B)
Creates a wrapper around a mfem::Operator op using PETSc's MATSHELL object and returns the Mat in B...
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
Abstract class for PETSc's nonlinear solvers.
void Mult(double a, const Vector &x, double b, Vector &y) const
Matvec: y = a A x + b y.
void SetBCHandler(PetscBCHandler *bch)
Sets the object to handle essential boundary conditions.
Wrapper for hypre's parallel vector class.
void SetBlockSize(PetscInt bs)
Set block size of a vector.
void SetDeviceInvalid() const
void SetHostValid() const
const double * Read(bool=true) const override
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void ZeroBC(const Vector &x, Vector &y)
y = x on ess_tdof_list_c and y = 0 on ess_tdof_list
virtual void Eval(double t, Vector &g)
Boundary conditions evaluation.
Type
Enumeration defining IDs for some classes derived from Operator.
PetscPreconditionerFactory(const std::string &name_="MFEM Factory")
PetscInt N() const
Returns the global number of columns.
void SetDeviceValid() const
void SetFlagsFromMask_() const
PetscInt GetNumCols() const
Returns the local number of columns.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
void SetJacobianType(Operator::Type type)
double * HostWrite() override
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
bool operatorset
Boolean to handle SetOperator calls.
PetscParVector * X
Auxiliary vectors for typecasting.
void MakeAlias(const Memory &base, int offset, int size)
Create a memory object that points inside the memory object base.
PetscInt NNZ() const
Returns the number of nonzeros.
PetscParVector & operator-=(const PetscParVector &y)
Type GetType() const
Returns the type of boundary conditions.
virtual ~PetscPreconditionerFactory()
virtual ~PetscSolverMonitor()
virtual void Mult(const Vector &b, Vector &x) const
Application of the preconditioner.
const double * GetHostPointer() const
virtual ~PetscNonlinearSolver()
void SetComputeNetFlux(bool net=true)
Setup BDDC with no-net-flux local solvers. Needs a ParFiniteElementSpace attached.
virtual ~PetscODESolver()
virtual ~PetscBCHandler()
PetscSolver()
Construct an empty PetscSolver. Initialize protected objects to NULL.
void SetHostInvalid() const
void SetPreconditionerFactory(PetscPreconditionerFactory *factory)
Sets the object for the creation of the preconditioner.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
petsc::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)
const Array< int > * ess_dof
struct _p_PetscObject * PetscObject
PetscH2Solver(Operator &op, ParFiniteElementSpace *fes, const std::string &prefix=std::string())
void ResetMemory()
Completes the operation started with PlaceMemory.
void SetAbsTol(double tol)
Helper class for handling essential boundary conditions.
Class used by MFEM to store pointers to host and/or device memory.
Array< int > & GetTDofs()
Gets essential dofs (local, per-process numbering)
void Customize(bool customize=true) const
Customize object with options set.
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.
bool WriteRequested() const
PetscParMatrix * TripleMatrixProduct(PetscParMatrix *R, PetscParMatrix *A, PetscParMatrix *P)
Returns the matrix R * A * P.
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array...
bool IsAliasForSync() const
void FreePrivateContext()
PetscODESolver(MPI_Comm comm, const std::string &prefix=std::string())
bool ReadRequested() const
void EliminateRows(const Array< int > &rows)
Eliminate only the rows from the matrix.
virtual void MultTranspose(const Vector &b, Vector &x) 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(), bool iter_mode=false)
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Wrapper for hypre's ParCSR matrix class.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
void Print(const char *fname=NULL, bool binary=false) const
Prints the vector (to stdout if fname is NULL)
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
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)
Vector * GlobalVector() const
Returns the global vector in each processor.
PetscInt GetGlobalNumCols() const
Returns the global number of columns.
PetscInt GetColStart() const
Returns the global index of the first local column.
void PlaceMemory(Memory< double > &, bool=false)
This requests write access from where the memory is valid and temporarily replaces the corresponding ...
MPI_Comm GetComm() const
Get the associated MPI communicator.
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.
void ScaleCols(const Vector &s)
Scale the local col i by s(i).
virtual ~PetscParMatrix()
Calls PETSc's destroy function.