escript  Revision_
Solver.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2016 by The University of Queensland
5 * http://www.uq.edu.au
6 *
7 * Primary Business: Queensland, Australia
8 * Licensed under the Apache License, version 2.0
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 * Development 2012-2013 by School of Earth Sciences
13 * Development from 2014 by Centre for Geoscience Computing (GeoComp)
14 *
15 *****************************************************************************/
16 
17 
18 #ifndef __PASO_SOLVER_H__
19 #define __PASO_SOLVER_H__
20 
21 #include "Paso.h"
22 #include "Functions.h"
23 #include "performance.h"
24 #include "SystemMatrix.h"
25 
26 namespace paso {
27 
28 #define TOLERANCE_FOR_SCALARS (double)(0.)
29 
30 void solve_free(SystemMatrix* A);
31 
32 SolverResult Solver(SystemMatrix_ptr, double*, double*, Options*, Performance*);
33 
34 void Solver_free(SystemMatrix*);
35 
36 SolverResult Solver_BiCGStab(SystemMatrix_ptr A, double* B, double* X,
37  dim_t* iter, double* tolerance, Performance* pp);
38 
39 SolverResult Solver_PCG(SystemMatrix_ptr A, double* B, double* X, dim_t* iter,
40  double* tolerance, Performance* pp);
41 
42 SolverResult Solver_TFQMR(SystemMatrix_ptr A, double* B, double* X, dim_t* iter,
43  double* tolerance, Performance* pp);
44 
45 SolverResult Solver_MINRES(SystemMatrix_ptr A, double* B, double* X,
46  dim_t* iter, double* tolerance, Performance* pp);
47 
48 SolverResult Solver_GMRES(SystemMatrix_ptr A, double* r, double* x,
49  dim_t* num_iter, double* tolerance,
50  dim_t length_of_recursion, dim_t restart,
51  Performance* pp);
52 
53 SolverResult Solver_GMRES2(Function* F, const double* f0, const double* x0,
54  double* x, dim_t* iter, double* tolerance,
55  Performance* pp);
56 
57 SolverResult Solver_NewtonGMRES(Function* F, double* x, Options* options,
58  Performance* pp);
59 
60 } // namespace paso
61 
62 #endif // __PASO_SOLVER_H__
63 
SolverResult Solver_TFQMR(SystemMatrix_ptr A, double *B, double *X, dim_t *iter, double *tolerance, Performance *pp)
Definition: TFQMR.cpp:62
void Solver_free(SystemMatrix *A)
Definition: Solver.cpp:40
SolverResult Solver_NewtonGMRES(Function *F, double *x, Options *options, Performance *pp)
Definition: NewtonGMRES.cpp:43
boost::shared_ptr< SystemMatrix > SystemMatrix_ptr
Definition: SystemMatrix.h:40
Definition: AMG.cpp:45
SolverResult Solver_GMRES(SystemMatrix_ptr A, double *r, double *x, dim_t *iter, double *tolerance, dim_t Length_of_recursion, dim_t restart, Performance *pp)
Definition: GMRES.cpp:68
SolverResult Solver_PCG(SystemMatrix_ptr A, double *r, double *x, dim_t *iter, double *tolerance, Performance *pp)
Definition: PCG.cpp:62
void solve_free(SystemMatrix *in)
Definition: solve.cpp:128
SolverResult
Definition: Paso.h:42
SolverResult Solver_BiCGStab(SystemMatrix_ptr A, double *r, double *x, dim_t *iter, double *tolerance, Performance *pp)
Definition: BiCGStab.cpp:77
SolverResult Solver_MINRES(SystemMatrix_ptr A, double *R, double *X, dim_t *iter, double *tolerance, Performance *pp)
Definition: MINRES.cpp:59
SolverResult Solver(SystemMatrix_ptr A, double *x, double *b, Options *options, Performance *pp)
calls the iterative solver
Definition: Solver.cpp:46
SolverResult Solver_GMRES2(Function *F, const double *f0, const double *x0, double *dx, dim_t *iter, double *tolerance, Performance *pp)
Definition: GMRES2.cpp:24
index_t dim_t
Definition: DataTypes.h:64