52 #ifndef AMESOS2_SOLVER_DECL_HPP
53 #define AMESOS2_SOLVER_DECL_HPP
55 #include <Teuchos_Describable.hpp>
56 #include <Teuchos_ParameterList.hpp>
57 #include <Teuchos_RCP.hpp>
58 #include <Teuchos_Comm.hpp>
77 template <
class Matrix,
class Vector>
78 class Solver :
public Teuchos::Describable {
132 virtual void solve(
void ) = 0;
149 virtual void solve(
const Teuchos::Ptr<Vector> X,
150 const Teuchos::Ptr<const Vector> B)
const = 0;
167 virtual void solve(Vector* X,
const Vector* B)
const = 0;
188 virtual type&
setParameters(
const Teuchos::RCP<Teuchos::ParameterList> & parameterList ) = 0;
195 virtual Teuchos::RCP<const Teuchos::ParameterList>
getValidParameters(
void )
const = 0;
223 virtual void setA(
const Teuchos::RCP<const Matrix> a,
EPhase keep_phase = CLEAN ) = 0;
244 virtual void setA(
const Matrix* a,
EPhase keep_phase = CLEAN ) = 0;
252 virtual void setX(
const Teuchos::RCP<Vector> x ) = 0;
256 virtual void setX( Vector* x ) = 0;
260 virtual const Teuchos::RCP<Vector>
getX(
void ) = 0;
264 virtual Vector*
getXRaw(
void ) = 0;
268 virtual void setB(
const Teuchos::RCP<const Vector> b ) = 0;
272 virtual void setB(
const Vector* b ) = 0;
276 virtual const Teuchos::RCP<const Vector>
getB(
void ) = 0;
280 virtual const Vector*
getBRaw(
void ) = 0;
284 virtual Teuchos::RCP<const Teuchos::Comm<int> >
getComm(
void )
const = 0;
292 virtual std::string
name(
void )
const = 0;
307 virtual void describe( Teuchos::FancyOStream &out,
308 const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default )
const = 0;
317 virtual void printTiming( Teuchos::FancyOStream &out,
319 const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default )
const = 0;
330 virtual void getTiming( Teuchos::ParameterList& timingParameterList )
const = 0;
339 #endif // AMESOS2_SOLVER_BASE_DECL_HPP
virtual void setX(const Teuchos::RCP< Vector > x)=0
Sets the LHS vector X.
Holds internal status data about the owning Amesos2 solver.
Definition: Amesos2_Status.hpp:73
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:65
virtual Status & getStatus() const =0
Returns a reference to this solver's internal status object.
virtual type & setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶meterList)=0
Set/update internal variables and solver options.
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters(void) const =0
Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.
virtual void printTiming(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0
Prints timing information about the current solver.
virtual std::string name(void) const =0
Return the name of this solver.
virtual std::string description(void) const =0
Returns a short description of this Solver.
virtual void getTiming(Teuchos::ParameterList &timingParameterList) const =0
Extracts timing information from the current solver.
virtual type & preOrdering(void)=0
Pre-orders the matrix.
virtual Vector * getXRaw(void)=0
Returns a raw pointer to the LHS of the linear system.
virtual type & numericFactorization(void)=0
Performs numeric factorization on the matrix.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0
virtual void setB(const Teuchos::RCP< const Vector > b)=0
Sets the RHS vector B.
virtual const Vector * getBRaw(void)=0
Returns a raw pointer to the RHS of the linear system.
virtual type & symbolicFactorization(void)=0
Performs symbolic factorization on the matrix.
virtual const Teuchos::RCP< const Vector > getB(void)=0
Returns the vector that is the RHS of the linear system.
Interface to Amesos2 solver objects.
Definition: Amesos2_Solver_decl.hpp:78
Enum and other types declarations for Amesos2.
virtual const Teuchos::RCP< Vector > getX(void)=0
Returns the vector that is the LHS of the linear system.
Container class for status variables.
virtual bool matrixShapeOK(void)=0
Returns true if the solver can handle the matrix shape.
virtual void solve(void)=0
Solves (or )
virtual void setA(const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN)=0
Sets the matrix A of this solver.
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm(void) const =0
Returns a pointer to the Teuchos::Comm communicator with this matrix.