INMOST
A toolkit for distributed mathematical modeling
inmost_solver.h
1 #ifndef INMOST_SOLVER_INCLUDED
2 #define INMOST_SOLVER_INCLUDED
3 
4 #include "inmost_common.h"
5 #include "inmost_sparse.h"
6 
7 #if defined(USE_SOLVER)
8 namespace INMOST {
9 
10  class SolverInterface;
11 
12  class SolverParameters;
13 
14  enum SolverVerbosityLevel {
15  SolverVerbosityLevel0 = 0,
16  SolverVerbosityLevel1 = 1,
17  SolverVerbosityLevel2 = 2
18  };
19 
29  class Solver {
30  public:
31  //Backward-compatibility
32  typedef std::string Type;
33  static const Type INNER_ILU2;
34  static const Type INNER_DDPQILUC;
35  static const Type INNER_MPTILUC;
36  static const Type INNER_MLMPTILUC;
37  static const Type INNER_MPTILU2;
38  static const Type Trilinos_Aztec;
39  static const Type Trilinos_Belos;
40  static const Type Trilinos_ML;
41  static const Type Trilinos_Ifpack;
42  static const Type PETSc;
43  static const Type ANI;
44  static const Type FCBIILU2;
45  static const Type K3BIILU2;
46  static const Type SUPERLU;
52  void SetParameterEnum(std::string name, INMOST_DATA_ENUM_TYPE value);
53 
59  void SetParameterReal(std::string name, INMOST_DATA_REAL_TYPE value);
60 
62  std::string GetReason() { return ReturnReason(); }
63  //Backward-compatibility
64  private:
65  static std::vector<SolverParameters> parameters;
66  static int *argc;
67  static char ***argv;
68  static bool is_initialized;
69  static bool is_finalized;
70  SolverInterface *solver;
71  std::string prefix;
72  SolverVerbosityLevel versbosity;
73  double preconditioner_time;
74  double iterations_time;
75  bool is_solved;
76 
79  static void parseXMLDatabase(const char *xml_database);
80 
81  public:
87  class OrderInfo {
88  private:
89  typedef std::vector<INMOST_DATA_ENUM_TYPE> storage_type;
90  storage_type global_to_proc;
91  storage_type global_overlap;
92  std::vector<INMOST_DATA_ENUM_TYPE> vector_exchange_recv;
93  std::vector<INMOST_DATA_ENUM_TYPE> vector_exchange_send;
94  std::vector<INMOST_DATA_REAL_TYPE> send_storage;
95  std::vector<INMOST_DATA_REAL_TYPE> recv_storage;
96  std::vector<INMOST_MPI_Request> send_requests;
97  std::vector<INMOST_MPI_Request> recv_requests;
98  std::vector<INMOST_DATA_ENUM_TYPE> extended_indexes;
99  INMOST_DATA_ENUM_TYPE local_vector_begin;
100  INMOST_DATA_ENUM_TYPE local_vector_end;
101  INMOST_DATA_ENUM_TYPE initial_matrix_begin;
102  INMOST_DATA_ENUM_TYPE initial_matrix_end;
103  INMOST_DATA_ENUM_TYPE local_matrix_begin;
104  INMOST_DATA_ENUM_TYPE local_matrix_end;
105  bool have_matrix;
106  INMOST_MPI_Comm comm;
107  int rank;
108  int size;
109  public:
111  void Clear();
112 
114  bool &HaveMatrix() { return have_matrix; }
115 
118 
121  OrderInfo(const OrderInfo &other);
122 
125  OrderInfo &operator=(OrderInfo const &other);
126 
129 
135  void PrepareMatrix(Sparse::Matrix &m, INMOST_DATA_ENUM_TYPE overlap);
136 
140 
144 
151  INMOST_DATA_ENUM_TYPE GetProcessor(INMOST_DATA_ENUM_TYPE gind) const;
152 
157  void GetOverlapRegion(INMOST_DATA_ENUM_TYPE proc, INMOST_DATA_ENUM_TYPE &mbeg, INMOST_DATA_ENUM_TYPE &mend) const;
158 
163  void GetLocalRegion(INMOST_DATA_ENUM_TYPE proc, INMOST_DATA_ENUM_TYPE &mbeg, INMOST_DATA_ENUM_TYPE &mend) const;
164 
168  void GetVectorRegion(INMOST_DATA_ENUM_TYPE &mbeg, INMOST_DATA_ENUM_TYPE &mend) const {
169  mbeg = local_vector_begin;
170  mend = local_vector_end;
171  }
173  INMOST_DATA_ENUM_TYPE GetRank() const { return rank; }
175  INMOST_DATA_ENUM_TYPE GetSize() const { return size; }
176 
180 
184 
188  void Integrate(INMOST_DATA_REAL_TYPE *inout, INMOST_DATA_ENUM_TYPE num) const;
190  INMOST_MPI_Comm GetComm() const { return comm; }
193  INMOST_MPI_Request *GetSendRequests() {
194  assert(!send_requests.empty());
195  return &send_requests[0];
196  }
199  INMOST_MPI_Request *GetRecvRequests() {
200  assert(!recv_requests.empty());
201  return &recv_requests[0];
202  }
204  INMOST_DATA_ENUM_TYPE GetSendRequestsSize() { return static_cast<INMOST_DATA_ENUM_TYPE>(send_requests.size()); }
206  INMOST_DATA_ENUM_TYPE GetRecvRequestsSize() { return static_cast<INMOST_DATA_ENUM_TYPE>(recv_requests.size()); }
208  INMOST_DATA_ENUM_TYPE *GetSendExchangeArray() {
209  assert(!vector_exchange_send.empty());
210  return &vector_exchange_send[0];
211  }
213  INMOST_DATA_ENUM_TYPE GetSendExchangeSize() { return static_cast<INMOST_DATA_ENUM_TYPE>(send_storage.size()); }
215  INMOST_DATA_ENUM_TYPE *GetRecvExchangeArray() {
216  assert(!vector_exchange_recv.empty());
217  return &vector_exchange_recv[0];
218  }
220  INMOST_DATA_ENUM_TYPE GetRecvExchangeSize() { return static_cast<INMOST_DATA_ENUM_TYPE>(recv_storage.size()); }
223  //void BeginSequentialCode() {for(int i = 0; i < rank; i++) MPI_Barrier(comm);}
226  //void EndSequentialCode() {for(int i = rank; i < size; i++) MPI_Barrier(comm);}
227  };
228 
237  Solver(std::string solverName, std::string prefix = "", INMOST_MPI_Comm _comm = INMOST_MPI_COMM_WORLD);
238 
241  Solver(const Solver &other);
242 
245  Solver &operator=(const Solver &other);
246 
249  std::string SolverName() const;
250 
253  std::string SolverPrefix() const;
254 
257  SolverVerbosityLevel VerbosityLevel() const;
258 
260  void SetVerbosityLevel(SolverVerbosityLevel level);
261 
283  static void Initialize(int *argc, char ***argv, const char *database = NULL);
284 
290  static void Finalize();
291 
293  static bool isInitialized();
294 
296  static bool isFinalized();
297 
312  void SetMatrix(Sparse::Matrix &A, bool ModifiedPattern = true, bool OldPreconditioner = false);
313 
325 
327  bool Clear();
328 
332  std::string GetParameter(std::string name) const;
333 
404  void SetParameter(std::string name, std::string value);
407  INMOST_DATA_ENUM_TYPE Iterations() const;
410  INMOST_DATA_REAL_TYPE Residual() const;
411 
414  const std::string ReturnReason() const;
415 
418  double PreconditionerTime() const;
419 
422  double IterationsTime() const;
423 
426  bool IsSolved() const;
427 
430  const std::string SolutionMetadataLine(const std::string &delimiter) const;
443  INMOST_DATA_REAL_TYPE Condest(INMOST_DATA_REAL_TYPE tol, INMOST_DATA_ENUM_TYPE maxits = 100);
444 
448  static bool isSolverAvailable(std::string name);
449 
452  static std::vector<std::string> getAvailableSolvers();
453 
456  };
457 
458  typedef std::vector<std::string>::iterator solvers_names_iterator_t;
459  typedef std::vector<SolverParameters>::iterator solver_parameters_iterator_t;
460 
462  namespace TTSP {
463 
464  void Initialize(const std::string &path);
465 
466  void Deinitialize();
467 
468  void Enable(const std::string &solver_name, const std::string &solver_prefix);
469 
470  void Disable(const std::string &solver_name, const std::string &solver_prefix);
471 
472  bool isEnabled(const std::string &solver_name, const std::string &solver_prefix);
473 
474  bool isDisabled(const std::string &solver_name, const std::string &solver_prefix);
475 
476  void SolverOptimize(Solver &solver, bool use_last_suggestion = false);
477 
478  void SolverOptimizeSaveResult(Solver &solver, double metrics, bool is_good);
479 
480  void DestroySavedOptimizer(Solver &solver);
481 
482  }
483 
484 }
485 
486 #endif
487 
488 #endif //INMOST_SOLVER_INCLUDED
Base class for low level parallel operations with objects of Solver class.
Definition: inmost_solver.h:87
INMOST_MPI_Comm GetComm() const
Get the communicator which the solver is associated with.
INMOST_DATA_ENUM_TYPE GetProcessor(INMOST_DATA_ENUM_TYPE gind) const
Retrieve the processor number by binary search for the specified global index.
INMOST_MPI_Request * GetSendRequests()
MPI structures that hold information on sent data.
void Clear()
Clear all structures and data.
INMOST_DATA_ENUM_TYPE GetRank() const
Get the rank of the current communicator, i.e. the current process index.
OrderInfo & operator=(OrderInfo const &other)
Assign all structures.
INMOST_DATA_ENUM_TYPE GetSendRequestsSize()
Total number of send requests.
void Update(Sparse::Vector &x)
Update the shared data in parallel vector.
~OrderInfo()
Clears all data.
void RestoreVector(Sparse::Vector &v) const
Restore initial nonparallel state of the Vector.
void GetLocalRegion(INMOST_DATA_ENUM_TYPE proc, INMOST_DATA_ENUM_TYPE &mbeg, INMOST_DATA_ENUM_TYPE &mend) const
Get the local index region for the specified processor.
void RestoreMatrix(Sparse::Matrix &m)
Restore initial nonparallel state of the Matrix with no overlap.
INMOST_DATA_ENUM_TYPE GetRecvExchangeSize()
Size of the array used to receive data.
bool & HaveMatrix()
Return true if Matrix data have already been specified.
void PrepareVector(Sparse::Vector &v) const
Prepare parallel state of the Vector.
OrderInfo(const OrderInfo &other)
Copy all structures.
INMOST_DATA_ENUM_TYPE * GetRecvExchangeArray()
Request raw acces to array used to receive data.
INMOST_MPI_Request * GetRecvRequests()
MPI structures that hold information on posted receive requests.
void Integrate(INMOST_DATA_REAL_TYPE *inout, INMOST_DATA_ENUM_TYPE num) const
Get the sum of num elements of real array on all processes.
void PrepareMatrix(Sparse::Matrix &m, INMOST_DATA_ENUM_TYPE overlap)
Prepare parallel state of the Matrix with specified overlap size.
INMOST_DATA_ENUM_TYPE GetSendExchangeSize()
Size of the array used to send data.
void Accumulate(Sparse::Vector &x)
Sum shared values in parallel vector.
INMOST_DATA_ENUM_TYPE * GetSendExchangeArray()
Request raw access to array used to send data.
void GetOverlapRegion(INMOST_DATA_ENUM_TYPE proc, INMOST_DATA_ENUM_TYPE &mbeg, INMOST_DATA_ENUM_TYPE &mend) const
Get the overlap index region for the specified processor.
INMOST_DATA_ENUM_TYPE GetRecvRequestsSize()
Total number of posted receive requests.
OrderInfo()
Default initialize all structures.
void GetVectorRegion(INMOST_DATA_ENUM_TYPE &mbeg, INMOST_DATA_ENUM_TYPE &mend) const
Get the local index region for the current processor.
INMOST_DATA_ENUM_TYPE GetSize() const
Get the size of the current communicator, i.e. the total number of processes used.
Main class to set and solve linear system.
Definition: inmost_solver.h:29
void SetParameterEnum(std::string name, INMOST_DATA_ENUM_TYPE value)
Backward-compatibility access to integer-typed parameters.
static const Type INNER_DDPQILUC
inner Solver based on BiCGStab(L) solver with second order Crout-ILU with inversed-based condition es...
Definition: inmost_solver.h:34
INMOST_DATA_ENUM_TYPE Iterations() const
Return the number of iterations performed by the last solution.
SolverVerbosityLevel VerbosityLevel() const
Return the solver verbosity specified level of the current solver.
static bool isFinalized()
Checks the stage of parallel solution is finalized.
const std::string SolutionMetadataLine(const std::string &delimiter) const
Get the useful info of the last solution.
static void Finalize()
Finalize the stage of parallel solution.
static const Type PETSc
external Solver PETSc,
Definition: inmost_solver.h:42
bool Solve(Sparse::Vector &RHS, Sparse::Vector &SOL)
Solver the linear system: A*x = b.
const std::string ReturnReason() const
Get the reason of convergence or divergence of the last solution.
static const Type ANI
external Solver from ANI3D based on ILU2 (sequential Fortran version),
Definition: inmost_solver.h:43
void SetParameter(std::string name, std::string value)
static const Type Trilinos_Ifpack
external Solver AztecOO with Ifpack preconditioner.
Definition: inmost_solver.h:41
static const Type Trilinos_Belos
external Solver Belos from Trilinos package, currently without preconditioner.
Definition: inmost_solver.h:39
INMOST_DATA_REAL_TYPE Residual() const
Return the final precondioned residual achieved by the last solution.
static bool isSolverAvailable(std::string name)
Checks if solver available.
static bool isInitialized()
Checks the stage of parallel solution is initialized.
bool IsSolved() const
Return the last solution successfullness.
~Solver()
Delete solver and associated data.
void SetParameterReal(std::string name, INMOST_DATA_REAL_TYPE value)
Backward-compatibility access to real-typed parameters.
static const Type Trilinos_ML
external Solver AztecOO with ML preconditioner.
Definition: inmost_solver.h:40
double PreconditionerTime() const
Return the time of preconditioner evaluation by the last solution.
static const Type K3BIILU2
inner K3BIILU2 Solver (BIILU2 parallel version).
Definition: inmost_solver.h:45
Solver(std::string solverName, std::string prefix="", INMOST_MPI_Comm _comm=INMOST_MPI_COMM_WORLD)
Main constructor of the solver.
std::string SolverPrefix() const
Return the solver user specified name of the current solver.
static const Type INNER_MPTILU2
inner Solver based on BiCGStab(L) solver with second order ILU and maximum product transversal reorde...
Definition: inmost_solver.h:37
Solver & operator=(const Solver &other)
Assign a solver.
static const Type SUPERLU
external Solver SuperLU
Definition: inmost_solver.h:46
static const Type Trilinos_Aztec
external Solver AztecOO from Trilinos package.
Definition: inmost_solver.h:38
std::string SolverName() const
Return the solver name.
static const Type INNER_ILU2
inner Solver based on BiCGStab(L) solver with second order ILU factorization as preconditioner.
Definition: inmost_solver.h:33
Solver(const Solver &other)
Copy a solver.
static const Type INNER_MLMPTILUC
inner Solver based on BiCGStab(L) solver with second order Crout-ILU with inversed-based condition es...
Definition: inmost_solver.h:36
INMOST_DATA_REAL_TYPE Condest(INMOST_DATA_REAL_TYPE tol, INMOST_DATA_ENUM_TYPE maxits=100)
Computes the smallest and the largest eigenvalue with the power method.
static const Type FCBIILU2
external FCBIILU2 Solver (BIILU2 parallel F2C version).
Definition: inmost_solver.h:44
static void Initialize(int *argc, char ***argv, const char *database=NULL)
Initialize the stage of parallel solution.
static std::vector< std::string > getAvailableSolvers()
Return the list of all available solvers.
static const Type INNER_MPTILUC
inner Solver based on BiCGStab(L) solver with second order Crout-ILU with inversed-based condition es...
Definition: inmost_solver.h:35
std::string GetParameter(std::string name) const
Get the solver output parameter.
void SetMatrix(Sparse::Matrix &A, bool ModifiedPattern=true, bool OldPreconditioner=false)
Set the matrix and construct the preconditioner.
std::string GetReason()
Backward-compatibility acces to return reason.
Definition: inmost_solver.h:62
double IterationsTime() const
Return the time of iteration stage by the last solution.
bool Clear()
Clear all internal data of the current solver including matrix, preconditioner etc.
Class to store the distributed sparse matrix by compressed rows.
Distributed vector class.
Definition: inmost_sparse.h:75