/[escript]/branches/symbolic_from_3470/paso/src/Preconditioner.h
ViewVC logotype

Contents of /branches/symbolic_from_3470/paso/src/Preconditioner.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3868 - (show annotations)
Thu Mar 15 06:07:08 2012 UTC (7 years, 2 months ago) by caltinay
File MIME type: text/plain
File size: 12657 byte(s)
Update to latest trunk

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_PRECS
16 #define INC_PRECS
17
18 #include "SystemMatrix.h"
19 #include "performance.h"
20 #include "BOOMERAMG.h"
21
22 #define PRECONDITIONER_NO_ERROR 0
23 #define PRECONDITIONER_MAXITER_REACHED 1
24 #define PRECONDITIONER_INPUT_ERROR -1
25 #define PRECONDITIONER_MEMORY_ERROR -9
26 #define PRECONDITIONER_BREAKDOWN -10
27 #define PRECONDITIONER_NEGATIVE_NORM_ERROR -11
28 #define PRECONDITIONER_DIVERGENCE -12
29
30
31 #define PASO_AMG_UNDECIDED -1
32 #define PASO_AMG_IN_F 0
33 #define PASO_AMG_IN_C 1
34
35 /* GAUSS SEIDEL & Jacobi */
36 typedef struct Paso_Preconditioner_LocalSmoother {
37 bool_t Jacobi;
38 double* diag;
39 double* buffer;
40 index_t* pivot;
41 } Paso_Preconditioner_LocalSmoother;
42
43 typedef struct Paso_Preconditioner_Smoother {
44 Paso_Preconditioner_LocalSmoother* localSmoother;
45 bool_t is_local;
46 } Paso_Preconditioner_Smoother;
47
48 void Paso_Preconditioner_Smoother_free(Paso_Preconditioner_Smoother * in);
49 void Paso_Preconditioner_LocalSmoother_free(Paso_Preconditioner_LocalSmoother * in);
50
51
52 Paso_Preconditioner_Smoother* Paso_Preconditioner_Smoother_alloc(Paso_SystemMatrix * A_p, const bool_t jacobi, const bool_t is_local, const bool_t verbose);
53 Paso_Preconditioner_LocalSmoother* Paso_Preconditioner_LocalSmoother_alloc(Paso_SparseMatrix * A_p,const bool_t jacobi, const bool_t verbose);
54
55 void Paso_Preconditioner_Smoother_solve(Paso_SystemMatrix* A, Paso_Preconditioner_Smoother * gs, double * x, const double * b, const dim_t sweeps, const bool_t x_is_initial);
56 void Paso_Preconditioner_LocalSmoother_solve(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x, const double * b, const dim_t sweeps, const bool_t x_is_initial);
57 err_t Paso_Preconditioner_Smoother_solve_byTolerance(Paso_SystemMatrix* A, Paso_Preconditioner_Smoother * gs, double * x, const double * b, const double atol, dim_t *sweeps, const bool_t x_is_initial);
58
59 void Paso_Preconditioner_LocalSmoother_Sweep(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
60 void Paso_Preconditioner_LocalSmoother_Sweep_sequential(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
61 void Paso_Preconditioner_LocalSmoother_Sweep_tiled(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
62 void Paso_Preconditioner_LocalSmoother_Sweep_colored(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
63
64
65
66 /* Local preconditioner */
67 struct Paso_Preconditioner_AMG {
68 dim_t level;
69 dim_t n;
70 dim_t n_F;
71 dim_t n_block;
72 Paso_SystemMatrix * A_C; /* coarse level matrix */
73 Paso_SystemMatrix * P; /* prolongation n x n_C*/
74 Paso_SystemMatrix * R; /* restriction n_C x n */
75
76 Paso_Preconditioner_Smoother* Smoother;
77 dim_t post_sweeps;
78 dim_t pre_sweeps;
79 dim_t options_smoother; /* used in direct solver */
80 bool_t verbose; /* used in direct solver */
81 index_t reordering; /* applied reordering in direct solver */
82 dim_t refinements; /* number of refinements in direct solver (typically =0) */
83 double* r; /* buffer for residual */
84 double* x_C; /* solution of coarse level system */
85 double* b_C; /* right hand side of coarse level system */
86 struct Paso_Preconditioner_AMG * AMG_C;
87 };
88 typedef struct Paso_Preconditioner_AMG Paso_Preconditioner_AMG;
89
90 void Paso_Preconditioner_AMG_free(Paso_Preconditioner_AMG * in);
91 Paso_Preconditioner_AMG* Paso_Preconditioner_AMG_alloc(Paso_SystemMatrix * A_p,dim_t level,Paso_Options* options);
92 void Paso_Preconditioner_AMG_solve(Paso_SystemMatrix* A, Paso_Preconditioner_AMG * amg, double * x, double * b);
93 void Paso_Preconditioner_AMG_setStrongConnections(Paso_SystemMatrix* A, dim_t *degree_S, index_t* offset_S, index_t *S, const double theta, const double tau);
94 void Paso_Preconditioner_AMG_setStrongConnections_Block(Paso_SystemMatrix* A, dim_t *degree_S, index_t* offset_S, index_t *S, const double theta, const double tau);
95 Paso_SystemMatrix* Paso_Preconditioner_AMG_getProlongation(Paso_SystemMatrix* A_p, const index_t* offset_S, const dim_t* degree_S, const index_t* S, const dim_t n_C, index_t* counter_C, const index_t interpolation_method);
96 void Paso_Preconditioner_AMG_setClassicProlongation(Paso_SystemMatrix* P, Paso_SystemMatrix* A, const index_t* offset_S, const dim_t* degree_S, const index_t* S,const index_t *counter_C);
97 void Paso_Preconditioner_AMG_setClassicProlongation_Block(Paso_SystemMatrix* P, Paso_SystemMatrix* A, const index_t* offset_S, const dim_t* degree_S, const index_t* S,const index_t *counter_C);
98 void Paso_Preconditioner_AMG_setDirectProlongation(Paso_SystemMatrix* P, Paso_SystemMatrix* A, const index_t* offset_S, const dim_t* degree_S, const index_t* S,const index_t *counter_C);
99 void Paso_Preconditioner_AMG_setDirectProlongation_Block(Paso_SystemMatrix* P, Paso_SystemMatrix* A, const index_t* offset_S, const dim_t* degree_S, const index_t* S,const index_t *counter_C);
100 double Paso_Preconditioner_AMG_getCoarseLevelSparsity(const Paso_Preconditioner_AMG * in);
101 dim_t Paso_Preconditioner_AMG_getNumCoarseUnknwons(const Paso_Preconditioner_AMG * in);
102 index_t Paso_Preconditioner_AMG_getMaxLevel(const Paso_Preconditioner_AMG * in);
103 void Paso_Preconditioner_AMG_transposeStrongConnections(const dim_t n, const dim_t* degree_S, const index_t* offset_S, const index_t* S, const dim_t nT, dim_t *degree_ST, index_t* offset_ST,index_t* ST);
104 void Paso_Preconditioner_AMG_CIJPCoarsening(const dim_t n, const dim_t my_n, index_t*split_marker,
105 const dim_t* degree_S, const index_t* offset_S, const index_t* S,
106 const dim_t* degree_ST, const index_t* offset_ST, const index_t* ST,
107 Paso_Connector* col_connector, Paso_Distribution* col_dist);
108 Paso_SystemMatrix* Paso_Preconditioner_AMG_getRestriction(Paso_SystemMatrix* P);
109 Paso_SystemMatrix* Paso_Preconditioner_AMG_buildInterpolationOperator(Paso_SystemMatrix* A, Paso_SystemMatrix* P, Paso_SystemMatrix* R);
110 Paso_SystemMatrix* Paso_Preconditioner_AMG_buildInterpolationOperatorBlock(Paso_SystemMatrix* A, Paso_SystemMatrix* P, Paso_SystemMatrix* R);
111 Paso_SparseMatrix* Paso_Preconditioner_AMG_mergeSystemMatrix(Paso_SystemMatrix* A);
112 void Paso_Preconditioner_AMG_mergeSolve(Paso_Preconditioner_AMG* amg);
113
114 /* Local AMG preconditioner */
115 struct Paso_Preconditioner_LocalAMG {
116 dim_t level;
117 dim_t n;
118 dim_t n_F;
119 dim_t n_block;
120
121 Paso_SparseMatrix * A_C; /* coarse level matrix */
122 Paso_SparseMatrix * P; /* prolongation n x n_C*/
123 Paso_SparseMatrix * R; /* restriction n_C x n */
124
125 Paso_Preconditioner_LocalSmoother* Smoother;
126 dim_t post_sweeps;
127 dim_t pre_sweeps;
128 index_t reordering; /* applied reordering in direct solver */
129 dim_t refinements; /* number of refinements in direct solver (typically =0) */
130 double* r; /* buffer for residual */
131 double* x_C; /* solution of coarse level system */
132 double* b_C; /* right hand side of coarse level system */
133 struct Paso_Preconditioner_LocalAMG * AMG_C;
134 };
135 typedef struct Paso_Preconditioner_LocalAMG Paso_Preconditioner_LocalAMG;
136
137 void Paso_Preconditioner_LocalAMG_free(Paso_Preconditioner_LocalAMG * in);
138 Paso_Preconditioner_LocalAMG* Paso_Preconditioner_LocalAMG_alloc(Paso_SparseMatrix * A_p,dim_t level,Paso_Options* options);
139 void Paso_Preconditioner_LocalAMG_solve(Paso_SparseMatrix* A, Paso_Preconditioner_LocalAMG * amg, double * x, double * b);
140
141 void Paso_Preconditioner_LocalAMG_RungeStuebenSearch(const dim_t n, const index_t* offset, const dim_t* degree, const index_t* S, index_t*split_marker, const bool_t usePanel);
142 void Paso_Preconditioner_LocalAMG_setStrongConnections_Block(Paso_SparseMatrix* A, dim_t *degree, index_t *S, const double theta, const double tau);
143 void Paso_Preconditioner_LocalAMG_setStrongConnections(Paso_SparseMatrix* A, dim_t *degree, index_t *S, const double theta, const double tau);
144 Paso_SparseMatrix* Paso_Preconditioner_LocalAMG_getProlongation(Paso_SparseMatrix* A_p, const index_t* offset_S, const dim_t* degree_S, const index_t* S, const dim_t n_C, const index_t* counter_C, const index_t interpolation_method);
145 void Paso_Preconditioner_LocalAMG_setDirectProlongation_Block(Paso_SparseMatrix* P_p, const Paso_SparseMatrix* A_p, const index_t *counter_C);
146 void Paso_Preconditioner_LocalAMG_setDirectProlongation(Paso_SparseMatrix* P_p, const Paso_SparseMatrix* A_p, const index_t *counter_C);
147 void Paso_Preconditioner_LocalAMG_setClassicProlongation(Paso_SparseMatrix* P_p, Paso_SparseMatrix* A_p, const index_t* offset_S, const dim_t* degree_S, const index_t* S, const index_t *counter_C);
148 void Paso_Preconditioner_LocalAMG_setClassicProlongation_Block(Paso_SparseMatrix* P_p, Paso_SparseMatrix* A_p, const index_t* offset_S, const dim_t* degree_S, const index_t* S, const index_t *counter_C);
149 index_t Paso_Preconditioner_LocalAMG_getMaxLevel(const Paso_Preconditioner_LocalAMG * in);
150 double Paso_Preconditioner_LocalAMG_getCoarseLevelSparsity(const Paso_Preconditioner_LocalAMG * in);
151 dim_t Paso_Preconditioner_LocalAMG_getNumCoarseUnknwons(const Paso_Preconditioner_LocalAMG * in);
152 void Paso_Preconditioner_LocalAMG_enforceFFConnectivity(const dim_t n, const index_t* offset_S, const dim_t* degree_S, const index_t* S, index_t*split_marker);
153
154
155 struct Paso_Preconditioner_BoomerAMG
156 {
157 Paso_BOOMERAMG_Handler* pt;
158 };
159 typedef struct Paso_Preconditioner_BoomerAMG Paso_Preconditioner_BoomerAMG;
160 void Paso_Preconditioner_BoomerAMG_free(Paso_Preconditioner_BoomerAMG * in);
161 Paso_Preconditioner_BoomerAMG* Paso_Preconditioner_BoomerAMG_alloc(Paso_SystemMatrix * A_p,Paso_Options* options);
162 void Paso_Preconditioner_BoomerAMG_solve(Paso_SystemMatrix* A, Paso_Preconditioner_BoomerAMG * amg, double * x, double * b);
163
164
165 struct Paso_Preconditioner_AMG_Root
166 {
167 bool_t is_local;
168 Paso_Preconditioner_AMG* amg;
169 Paso_Preconditioner_LocalAMG* localamg;
170 Paso_Preconditioner_BoomerAMG* boomeramg;
171 dim_t sweeps;
172 Paso_Preconditioner_Smoother* amgsubstitute;
173 };
174 typedef struct Paso_Preconditioner_AMG_Root Paso_Preconditioner_AMG_Root;
175
176 Paso_Preconditioner_AMG_Root* Paso_Preconditioner_AMG_Root_alloc(Paso_SystemMatrix *A, Paso_Options* options);
177 void Paso_Preconditioner_AMG_Root_free(Paso_Preconditioner_AMG_Root * in);
178 void Paso_Preconditioner_AMG_Root_solve(Paso_SystemMatrix* A, Paso_Preconditioner_AMG_Root * amg, double * x, double * b);
179
180 /*===============================================*/
181 /* ILU preconditioner */
182 struct Paso_Solver_ILU {
183 double* factors;
184 };
185 typedef struct Paso_Solver_ILU Paso_Solver_ILU;
186
187
188
189 /* RILU preconditioner */
190 struct Paso_Solver_RILU {
191 dim_t n;
192 dim_t n_block;
193 dim_t n_F;
194 dim_t n_C;
195 double* inv_A_FF;
196 index_t* A_FF_pivot;
197 Paso_SparseMatrix * A_FC;
198 Paso_SparseMatrix * A_CF;
199 index_t* rows_in_F;
200 index_t* rows_in_C;
201 index_t* mask_F;
202 index_t* mask_C;
203 double* x_F;
204 double* b_F;
205 double* x_C;
206 double* b_C;
207 struct Paso_Solver_RILU * RILU_of_Schur;
208 };
209 typedef struct Paso_Solver_RILU Paso_Solver_RILU;
210
211
212
213 /* general preconditioner interface */
214
215 typedef struct Paso_Preconditioner {
216 dim_t type;
217 dim_t sweeps;
218 /* jacobi preconditioner */
219 Paso_Preconditioner_Smoother* jacobi;
220 /* Gauss-Seidel preconditioner */
221 Paso_Preconditioner_Smoother* gs;
222 /* amg preconditioner */
223 Paso_Preconditioner_AMG_Root *amg;
224
225 /* ilu preconditioner */
226 Paso_Solver_ILU* ilu;
227 /* rilu preconditioner */
228 Paso_Solver_RILU* rilu;
229
230 } Paso_Preconditioner;
231
232 void Paso_Preconditioner_free(Paso_Preconditioner*);
233 Paso_Preconditioner* Paso_Preconditioner_alloc(Paso_SystemMatrix* A,Paso_Options* options);
234 void Paso_Preconditioner_solve(Paso_Preconditioner* prec, Paso_SystemMatrix* A,double*,double*);
235
236
237 /*******************************************/
238 void Paso_Solver_ILU_free(Paso_Solver_ILU * in);
239 Paso_Solver_ILU* Paso_Solver_getILU(Paso_SparseMatrix * A_p,bool_t verbose);
240 void Paso_Solver_solveILU(Paso_SparseMatrix * A, Paso_Solver_ILU * ilu, double * x, const double * b);
241
242 void Paso_Solver_RILU_free(Paso_Solver_RILU * in);
243 Paso_Solver_RILU* Paso_Solver_getRILU(Paso_SparseMatrix * A_p,bool_t verbose);
244 void Paso_Solver_solveRILU(Paso_Solver_RILU * rilu, double * x, double * b);
245
246 void Paso_Solver_updateIncompleteSchurComplement(Paso_SparseMatrix* A_CC, Paso_SparseMatrix *A_CF,double* invA_FF,index_t* A_FF_pivot, Paso_SparseMatrix *A_FC);
247
248
249 #endif /* #ifndef INC_PRECS */

  ViewVC Help
Powered by ViewVC 1.1.26