/[escript]/trunk/paso/src/Solver.h
ViewVC logotype

Contents of /trunk/paso/src/Solver.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3097 - (show annotations)
Fri Aug 20 04:59:12 2010 UTC (9 years, 10 months ago) by gross
File MIME type: text/plain
File size: 9013 byte(s)
some modifications to the GaussSeidel
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2010 by University of Queensland
5 * Earth Systems Science Computational Center (ESSCC)
6 * http://www.uq.edu.au/esscc
7 *
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
15 #ifndef INC_SOLVER
16 #define INC_SOLVER
17
18 #include "SystemMatrix.h"
19 #include "performance.h"
20 #include "Functions.h"
21
22 #define PASO_TRACE
23 /* error codes used in the solver */
24 #define SOLVER_NO_ERROR 0
25 #define SOLVER_MAXITER_REACHED 1
26 #define SOLVER_INPUT_ERROR -1
27 #define SOLVER_MEMORY_ERROR -9
28 #define SOLVER_BREAKDOWN -10
29 #define SOLVER_NEGATIVE_NORM_ERROR -11
30 #define SOLVER_DIVERGENCE -12
31
32 #define TOLERANCE_FOR_SCALARS (double)(0.)
33 #define PASO_ONE (double)(1.0)
34 #define PASO_ZERO (double)(0.0)
35
36 #define MAX_BLOCK_SIZE 3
37
38 /* static double ONE=1.; */
39 /* static double ZERO=0.;*/
40 /*static double TOLERANCE_FOR_SCALARS=0.;*/
41
42 /* jacobi preconditioner */
43
44 typedef struct Paso_Solver_Jacobi {
45 double* values;
46 index_t* pivot;
47 } Paso_Solver_Jacobi;
48
49 /* GS preconditioner */
50 typedef struct Paso_Solver_LocalGS {
51 double* diag;
52 double* buffer;
53 index_t* pivot;
54 dim_t sweeps;
55 } Paso_Solver_LocalGS;
56
57 typedef struct Paso_Solver_GS {
58 Paso_Solver_LocalGS* localGS;
59 bool_t is_local;
60 } Paso_Solver_GS;
61
62
63
64 /*===============================================*/
65 /* ILU preconditioner */
66 struct Paso_Solver_ILU {
67 double* factors;
68 };
69 typedef struct Paso_Solver_ILU Paso_Solver_ILU;
70
71
72
73
74 /* RILU preconditioner */
75 struct Paso_Solver_RILU {
76 dim_t n;
77 dim_t n_block;
78 dim_t n_F;
79 dim_t n_C;
80 double* inv_A_FF;
81 index_t* A_FF_pivot;
82 Paso_SparseMatrix * A_FC;
83 Paso_SparseMatrix * A_CF;
84 index_t* rows_in_F;
85 index_t* rows_in_C;
86 index_t* mask_F;
87 index_t* mask_C;
88 double* x_F;
89 double* b_F;
90 double* x_C;
91 double* b_C;
92 struct Paso_Solver_RILU * RILU_of_Schur;
93 };
94 typedef struct Paso_Solver_RILU Paso_Solver_RILU;
95
96 struct Paso_Solver_Smoother {
97 dim_t ID;
98 Paso_Solver_Jacobi* Jacobi;
99 Paso_Solver_LocalGS* GS;
100 };
101 typedef struct Paso_Solver_Smoother Paso_Solver_Smoother;
102
103 /* AMG preconditioner */
104 struct Paso_Solver_AMG {
105 dim_t n;
106 dim_t level;
107 bool_t coarsest_level;
108 dim_t n_block;
109 dim_t n_F;
110 dim_t n_C;
111
112 Paso_SparseMatrix * A_FF;
113 Paso_SparseMatrix * A_FC;
114 Paso_SparseMatrix * A_CF;
115 Paso_SparseMatrix * W_FC;
116
117 Paso_SparseMatrix * P;
118 Paso_SparseMatrix * R;
119
120 index_t* rows_in_F;
121 index_t* rows_in_C;
122 index_t* mask_F;
123 index_t* mask_C;
124 double* x_F;
125 double* b_F;
126 double* x_C;
127 double* b_C;
128
129 dim_t post_sweeps;
130 dim_t pre_sweeps;
131
132 Paso_SparseMatrix * A;
133 Paso_SparseMatrix * AOffset1;
134 Paso_SparseMatrix * AUnrolled;
135 void* solver;
136 Paso_Solver_Smoother* Smoother;
137 struct Paso_Solver_AMG * AMG_of_Coarse;
138 };
139 typedef struct Paso_Solver_AMG Paso_Solver_AMG;
140
141
142 /* AMLI preconditioner */
143 struct Paso_Solver_AMLI {
144 dim_t n;
145 dim_t level;
146 bool_t coarsest_level;
147 dim_t n_block;
148 dim_t n_F;
149 dim_t n_C;
150 double* inv_A_FF;
151 index_t* A_FF_pivot;
152 Paso_SparseMatrix * A_FC;
153 Paso_SparseMatrix * A_CF;
154 index_t* rows_in_F;
155 index_t* rows_in_C;
156 index_t* mask_F;
157 index_t* mask_C;
158 double* x_F;
159 double* b_F;
160 double* x_C;
161 double* b_C;
162
163 dim_t post_sweeps;
164 dim_t pre_sweeps;
165
166 Paso_SparseMatrix * A;
167 Paso_SparseMatrix * AOffset1;
168 void* solver;
169 Paso_Solver_Jacobi* GS;
170 struct Paso_Solver_AMLI * AMLI_of_Schur;
171 };
172 typedef struct Paso_Solver_AMLI Paso_Solver_AMLI;
173
174
175 /* AMLI preconditioner on blocks*/
176 struct Paso_Solver_AMLI_System {
177 dim_t block_size;
178 Paso_SparseMatrix *block[MAX_BLOCK_SIZE];
179 Paso_Solver_AMLI *amliblock[MAX_BLOCK_SIZE];
180 };
181 typedef struct Paso_Solver_AMLI_System Paso_Solver_AMLI_System;
182
183
184 /* AMG preconditioner on blocks*/
185 struct Paso_Solver_AMG_System {
186 dim_t block_size;
187 Paso_SparseMatrix *block[MAX_BLOCK_SIZE];
188 Paso_Solver_AMG *amgblock[MAX_BLOCK_SIZE];
189 };
190 typedef struct Paso_Solver_AMG_System Paso_Solver_AMG_System;
191
192 /* general preconditioner interface */
193
194 typedef struct Paso_Solver_Preconditioner {
195 dim_t type;
196
197 /* jacobi preconditioner */
198 Paso_Solver_Jacobi* jacobi;
199 /* Gauss-Seidel preconditioner */
200 Paso_Solver_GS* gs;
201
202 /* ilu preconditioner */
203 Paso_Solver_ILU* ilu;
204 /* rilu preconditioner */
205 Paso_Solver_RILU* rilu;
206 /* amg preconditioner */
207 Paso_Solver_AMG* amg;
208 /* amg on System */
209 Paso_Solver_AMG_System* amgSystem;
210 /* amg preconditioner */
211 Paso_Solver_AMLI* amli;
212 /* amg on System */
213 Paso_Solver_AMLI_System* amliSystem;
214
215 } Paso_Solver_Preconditioner;
216
217 void Paso_Solver(Paso_SystemMatrix*,double*,double*,Paso_Options*,Paso_Performance* pp);
218 void Paso_Solver_free(Paso_SystemMatrix*);
219 err_t Paso_Solver_BiCGStab( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
220 err_t Paso_Solver_PCG( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
221 err_t Paso_Solver_TFQMR( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
222 err_t Paso_Solver_MINRES( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
223 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);
224 err_t Paso_Solver_GMRES2(Paso_Function * F, const double* f0, const double* x0, double * x, dim_t *iter, double* tolerance, Paso_Performance* pp);
225
226 err_t Paso_Solver_NewtonGMRES(Paso_Function *F, double *x, Paso_Options* options, Paso_Performance* pp);
227 Paso_Function * Paso_Function_LinearSystem_alloc(Paso_SystemMatrix* A, double* b, Paso_Options* options);
228 err_t Paso_Function_LinearSystem_call(Paso_Function * F,double* value, const double* arg, Paso_Performance *pp);
229 void Paso_Function_LinearSystem_free(Paso_Function * F);
230 err_t Paso_Function_LinearSystem_setInitialGuess(Paso_SystemMatrix* A, double* x, Paso_Performance *pp);
231
232 void Paso_Preconditioner_free(Paso_Solver_Preconditioner*);
233 void Paso_Solver_setPreconditioner(Paso_SystemMatrix* A,Paso_Options* options);
234 void Paso_Solver_solvePreconditioner(Paso_SystemMatrix* A,double*,double*);
235 void Paso_Solver_applyBlockDiagonalMatrix(dim_t n_block,dim_t n,double* D,index_t* pivot,double* x,double* b);
236
237 Paso_Solver_Jacobi* Paso_Solver_getJacobi(Paso_SystemMatrix * A_p);
238 void Paso_Solver_solveJacobi(Paso_SystemMatrix * A_p, Paso_Solver_Jacobi * prec, double * x, double * b);
239 void Paso_Solver_solveLocalJacobi(Paso_SparseMatrix * A_p, Paso_Solver_Jacobi * prec, double * x, double * b);
240 void Paso_Solver_Jacobi_free(Paso_Solver_Jacobi * in);
241 Paso_Solver_Jacobi* Paso_Solver_getLocalJacobi(Paso_SparseMatrix * A_p);
242
243 void Paso_Solver_GS_free(Paso_Solver_GS * in);
244 void Paso_Solver_LocalGS_free(Paso_Solver_LocalGS * in);
245 Paso_Solver_GS* Paso_Solver_getGS(Paso_SystemMatrix * A_p, dim_t sweeps, bool_t is_local, bool_t verbose);
246 Paso_Solver_LocalGS* Paso_Solver_getLocalGS(Paso_SparseMatrix * A_p, dim_t sweeps, bool_t verbose);
247 void Paso_Solver_solveGS(Paso_SystemMatrix* A, Paso_Solver_GS * gs, double * x, const double * b);
248 void Paso_Solver_solveLocalGS(Paso_SparseMatrix* A, Paso_Solver_LocalGS * gs, double * x, const double * b);
249
250 void Paso_Solver_localGSSweep(Paso_SparseMatrix* A, Paso_Solver_LocalGS * gs, double * x);
251 void Paso_Solver_localGSSweep_sequential(Paso_SparseMatrix* A, Paso_Solver_LocalGS * gs, double * x);
252 void Paso_Solver_localGSSweep_tiled(Paso_SparseMatrix* A, Paso_Solver_LocalGS * gs, double * x);
253 void Paso_Solver_localGSSweep_colored(Paso_SparseMatrix* A, Paso_Solver_LocalGS * gs, double * x);
254
255 /*******************************************/
256 void Paso_Solver_ILU_free(Paso_Solver_ILU * in);
257 Paso_Solver_ILU* Paso_Solver_getILU(Paso_SparseMatrix * A_p,bool_t verbose);
258 void Paso_Solver_solveILU(Paso_SparseMatrix * A, Paso_Solver_ILU * ilu, double * x, const double * b);
259
260 void Paso_Solver_RILU_free(Paso_Solver_RILU * in);
261 Paso_Solver_RILU* Paso_Solver_getRILU(Paso_SparseMatrix * A_p,bool_t verbose);
262 void Paso_Solver_solveRILU(Paso_Solver_RILU * rilu, double * x, double * b);
263
264 void Paso_Solver_AMG_System_free(Paso_Solver_AMG_System * in);
265 void Paso_Solver_AMG_free(Paso_Solver_AMG * in);
266 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix * A_p,dim_t level,Paso_Options* options);
267 void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b);
268
269 void Paso_Solver_AMLI_System_free(Paso_Solver_AMLI_System * in);
270 void Paso_Solver_AMLI_free(Paso_Solver_AMLI * in);
271 Paso_Solver_AMLI* Paso_Solver_getAMLI(Paso_SparseMatrix * A_p,dim_t level,Paso_Options* options);
272 void Paso_Solver_solveAMLI(Paso_Solver_AMLI * amli, double * x, double * b);
273
274 void Paso_Solver_updateIncompleteSchurComplement(Paso_SparseMatrix* A_CC, Paso_SparseMatrix *A_CF,double* invA_FF,index_t* A_FF_pivot, Paso_SparseMatrix *A_FC);
275
276
277
278
279
280
281 #endif /* #ifndef INC_SOLVER */

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26