Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_Solver_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Amesos2: Templated Direct Sparse Solver Package
4//
5// Copyright 2011 NTESS and the Amesos2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
18#ifndef AMESOS2_SOLVER_DECL_HPP
19#define AMESOS2_SOLVER_DECL_HPP
20
21#include <Teuchos_Describable.hpp>
22#include <Teuchos_ParameterList.hpp>
23#include <Teuchos_RCP.hpp>
24#include <Teuchos_Comm.hpp>
25
26#include "Amesos2_TypeDecl.hpp"
27#include "Amesos2_Status.hpp"
28
29namespace Amesos2 {
30
31
43 template <class Matrix, class Vector>
44 class Solver : public Teuchos::Describable {
45
46 public:
47
48 typedef Solver<Matrix,Vector> type;
49
51
52
60 virtual type& preOrdering( void ) = 0;
61
62
68 virtual type& symbolicFactorization( void ) = 0;
69
70
83 virtual type& numericFactorization( void ) = 0;
84
85
98 virtual void solve( void ) = 0;
99
100
115 virtual void solve(const Teuchos::Ptr<Vector> X,
116 const Teuchos::Ptr<const Vector> B) const = 0;
117
118
133 virtual void solve(Vector* X, const Vector* B) const = 0;
134
136
153 virtual type& setParameters( const Teuchos::RCP<Teuchos::ParameterList> & parameterList ) = 0;
154
155
160 virtual Teuchos::RCP<const Teuchos::ParameterList> getValidParameters( void ) const = 0;
161
163
164
188 virtual void setA( const Teuchos::RCP<const Matrix> a, EPhase keep_phase = CLEAN ) = 0;
189
209 virtual void setA( const Matrix* a, EPhase keep_phase = CLEAN ) = 0;
210
211
213 virtual bool matrixShapeOK( void ) = 0;
214
215
217 virtual void setX( const Teuchos::RCP<Vector> x ) = 0;
218
219
221 virtual void setX( Vector* x ) = 0;
222
223
225 virtual const Teuchos::RCP<Vector> getX( void ) = 0;
226
227
229 virtual Vector* getXRaw( void ) = 0;
230
231
233 virtual void setB( const Teuchos::RCP<const Vector> b ) = 0;
234
235
237 virtual void setB( const Vector* b ) = 0;
238
239
241 virtual const Teuchos::RCP<const Vector> getB( void ) = 0;
242
243
245 virtual const Vector* getBRaw( void ) = 0;
246
247
249 virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm( void ) const = 0;
250
251
253 virtual Status& getStatus() const = 0;
254
255
257 virtual std::string name( void ) const = 0;
258
260
261
267 virtual std::string description( void ) const = 0;
268
269
272 virtual void describe( Teuchos::FancyOStream &out,
273 const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default ) const = 0;
274
276
277
283 virtual void printTiming( Teuchos::FancyOStream &out,
284 const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default ) const = 0;
285
286
295 virtual void getTiming( Teuchos::ParameterList& timingParameterList ) const = 0;
296
298
299 }; // End class Solver
300
301
302} // end namespace Amesos2
303
304#endif // AMESOS2_SOLVER_BASE_DECL_HPP
Container class for status variables.
Enum and other types declarations for Amesos2.
Interface to Amesos2 solver objects.
Definition Amesos2_Solver_decl.hpp:44
virtual std::string description(void) const =0
Returns a short description of this Solver.
virtual type & numericFactorization(void)=0
Performs numeric factorization on the matrix.
virtual void setB(const Teuchos::RCP< const Vector > b)=0
Sets the RHS vector B.
virtual void setX(Vector *x)=0
Sets the LHS vector X using a raw pointer.
virtual type & symbolicFactorization(void)=0
Performs symbolic factorization on the matrix.
virtual Vector * getXRaw(void)=0
Returns a raw pointer to the LHS of the linear system.
virtual void printTiming(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0
Prints timing information about the current solver.
virtual Status & getStatus() const =0
Returns a reference to this solver's internal status object.
virtual void setB(const Vector *b)=0
Sets the RHS vector B using a raw pointer.
virtual type & preOrdering(void)=0
Pre-orders the matrix.
virtual const Vector * getBRaw(void)=0
Returns a raw pointer to the RHS of the linear system.
virtual void setA(const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN)=0
Sets the matrix A of this solver.
virtual bool matrixShapeOK(void)=0
Returns true if the solver can handle the matrix shape.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0
virtual type & setParameters(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)=0
Set/update internal variables and solver options.
virtual void solve(Vector *X, const Vector *B) const =0
Solve using the given X and B vectors.
virtual void getTiming(Teuchos::ParameterList &timingParameterList) const =0
Extracts timing information from the current solver.
virtual const Teuchos::RCP< Vector > getX(void)=0
Returns the vector that is the LHS of the linear system.
virtual void solve(const Teuchos::Ptr< Vector > X, const Teuchos::Ptr< const Vector > B) const =0
Solve using the given X and B vectors.
virtual void setX(const Teuchos::RCP< Vector > x)=0
Sets the LHS vector X.
virtual void solve(void)=0
Solves (or )
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters(void) const =0
Return a const parameter list of all of the valid parameters that this->setParameterList(....
virtual void setA(const Matrix *a, EPhase keep_phase=CLEAN)=0
Sets the matrix A of this solver.
virtual std::string name(void) const =0
Return the name 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.
virtual const Teuchos::RCP< const Vector > getB(void)=0
Returns the vector that is the RHS of the linear system.
Holds internal status data about the owning Amesos2 solver.
Definition Amesos2_Status.hpp:39
EPhase
Used to indicate a phase in the direct solution.
Definition Amesos2_TypeDecl.hpp:31