1 |
/* $Id$ */ |
2 |
|
3 |
/* |
4 |
******************************************************************************** |
5 |
* Copyright 2006 by ACcESS MNRF * |
6 |
* * |
7 |
* http://www.access.edu.au * |
8 |
* Primary Business: Queensland, Australia * |
9 |
* Licensed under the Open Software License version 3.0 * |
10 |
* http://www.opensource.org/licenses/osl-3.0.php * |
11 |
******************************************************************************** |
12 |
*/ |
13 |
|
14 |
#ifndef INC_SOLVER |
15 |
#define INC_SOLVER |
16 |
|
17 |
#include "../SystemMatrix.h" |
18 |
#include "../performance.h" |
19 |
|
20 |
#define FINLEY_SOLVER_TRACE |
21 |
/* error codes used in the solver */ |
22 |
#define SOLVER_NO_ERROR 0 |
23 |
#define SOLVER_MAXITER_REACHED 1 |
24 |
#define SOLVER_INPUT_ERROR -1 |
25 |
#define SOLVER_MEMORY_ERROR -9 |
26 |
#define SOLVER_BREAKDOWN -10 |
27 |
|
28 |
static double ONE=1.; |
29 |
static double ZERO=0.; |
30 |
static double TOLERANCE_FOR_SCALARS=0.; |
31 |
|
32 |
/* ILU preconditioner */ |
33 |
struct Paso_Solver_ILU { |
34 |
dim_t n_block; |
35 |
dim_t n; |
36 |
index_t num_colors; |
37 |
index_t* colorOf; |
38 |
index_t* main_iptr; |
39 |
double* factors; |
40 |
Paso_SystemMatrixPattern* pattern; |
41 |
}; |
42 |
typedef struct Paso_Solver_ILU Paso_Solver_ILU; |
43 |
|
44 |
/* RILU preconditioner */ |
45 |
struct Paso_Solver_RILU { |
46 |
dim_t n; |
47 |
dim_t n_block; |
48 |
dim_t n_F; |
49 |
dim_t n_C; |
50 |
double* inv_A_FF; |
51 |
index_t* A_FF_pivot; |
52 |
Paso_SystemMatrix * A_FC; |
53 |
Paso_SystemMatrix * A_CF; |
54 |
index_t* rows_in_F; |
55 |
index_t* rows_in_C; |
56 |
index_t* mask_F; |
57 |
index_t* mask_C; |
58 |
double* x_F; |
59 |
double* b_F; |
60 |
double* x_C; |
61 |
double* b_C; |
62 |
struct Paso_Solver_RILU * RILU_of_Schur; |
63 |
}; |
64 |
typedef struct Paso_Solver_RILU Paso_Solver_RILU; |
65 |
|
66 |
|
67 |
/* jacobi preconditioner */ |
68 |
|
69 |
typedef struct Paso_Solver_Jacobi { |
70 |
dim_t n_block; |
71 |
dim_t n; |
72 |
double* values; |
73 |
index_t* pivot; |
74 |
} Paso_Solver_Jacobi; |
75 |
|
76 |
/* general preconditioner interface */ |
77 |
|
78 |
typedef struct Paso_Solver_Preconditioner { |
79 |
dim_t type; |
80 |
/* jacobi preconditioner */ |
81 |
Paso_Solver_Jacobi* jacobi; |
82 |
/* ilu preconditioner */ |
83 |
Paso_Solver_ILU* ilu; |
84 |
/* ilu preconditioner */ |
85 |
Paso_Solver_RILU* rilu; |
86 |
} Paso_Solver_Preconditioner; |
87 |
|
88 |
void Paso_Solver(Paso_SystemMatrix*,double*,double*,Paso_Options*,Paso_Performance* pp); |
89 |
void Paso_Solver_free(Paso_SystemMatrix*); |
90 |
err_t Paso_Solver_BiCGStab( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance,Paso_Performance* pp); |
91 |
err_t Paso_Solver_PCG( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance,Paso_Performance* pp); |
92 |
err_t Paso_Solver_GMRES(Paso_SystemMatrix * A, double * r, double * x, dim_t *num_iter, double * tolerance,dim_t length_of_recursion,dim_t restart,Paso_Performance* pp); |
93 |
void Paso_Preconditioner_free(Paso_Solver_Preconditioner*); |
94 |
void Paso_Solver_setPreconditioner(Paso_SystemMatrix* A,Paso_Options* options); |
95 |
void Paso_Solver_solvePreconditioner(Paso_SystemMatrix* A,double*,double*); |
96 |
void Paso_Solver_applyBlockDiagonalMatrix(dim_t n_block,dim_t n,double* D,index_t* pivot,double* x,double* b); |
97 |
|
98 |
void Paso_Solver_ILU_free(Paso_Solver_ILU * in); |
99 |
Paso_Solver_ILU* Paso_Solver_getILU(Paso_SystemMatrix * A_p,bool_t verbose); |
100 |
void Paso_Solver_solveILU(Paso_Solver_ILU * ilu, double * x, double * b); |
101 |
|
102 |
void Paso_Solver_RILU_free(Paso_Solver_RILU * in); |
103 |
Paso_Solver_RILU* Paso_Solver_getRILU(Paso_SystemMatrix * A_p,bool_t verbose); |
104 |
void Paso_Solver_solveRILU(Paso_Solver_RILU * rilu, double * x, double * b); |
105 |
|
106 |
void Paso_Solver_updateIncompleteSchurComplement(Paso_SystemMatrix* A_CC,Paso_SystemMatrix *A_CF,double* |
107 |
invA_FF,index_t* A_FF_pivot,Paso_SystemMatrix *A_FC); |
108 |
Paso_Solver_Jacobi* Paso_Solver_getJacobi(Paso_SystemMatrix * A_p); |
109 |
void Paso_Solver_solveJacobi(Paso_Solver_Jacobi * prec, double * x, double * b); |
110 |
void Paso_Solver_Jacobi_free(Paso_Solver_Jacobi * in); |
111 |
|
112 |
#endif /* #ifndef INC_SOLVER */ |