Escript  Revision_4320
Preconditioner.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2013 by University of Queensland
5 * http://www.uq.edu.au
6 *
7 * Primary Business: Queensland, Australia
8 * Licensed under the Open Software License version 3.0
9 * http://www.opensource.org/licenses/osl-3.0.php
10 *
11 * Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 * Development since 2012 by School of Earth Sciences
13 *
14 *****************************************************************************/
15 
16 
17 #ifndef INC_PRECS
18 #define INC_PRECS
19 
20 #include "SystemMatrix.h"
21 #include "performance.h"
22 #include "BOOMERAMG.h"
23 #include "MergedSolver.h"
24 
25 #define PRECONDITIONER_NO_ERROR 0
26 #define PRECONDITIONER_MAXITER_REACHED 1
27 #define PRECONDITIONER_INPUT_ERROR -1
28 #define PRECONDITIONER_MEMORY_ERROR -9
29 #define PRECONDITIONER_BREAKDOWN -10
30 #define PRECONDITIONER_NEGATIVE_NORM_ERROR -11
31 #define PRECONDITIONER_DIVERGENCE -12
32 
33 
34 #define PASO_AMG_UNDECIDED -1
35 #define PASO_AMG_IN_F 0
36 #define PASO_AMG_IN_C 1
37 
38 /* GAUSS SEIDEL & Jacobi */
41  double* diag;
42  double* buffer;
45 
50 
53 
54 
57 
58 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);
59 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);
60 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);
61 
66 
67 
68 
69 /* Local preconditioner */
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 
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  Paso_MergedSolver* merged_solver; /* used on the coarsest level */
88 };
90 
93 void Paso_Preconditioner_AMG_solve(Paso_SystemMatrix* A, Paso_Preconditioner_AMG * amg, double * x, double * b);
94 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);
95 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);
96 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);
97 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);
98 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);
99 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);
100 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);
104 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);
105 void Paso_Preconditioner_AMG_CIJPCoarsening(const dim_t n, const dim_t my_n, index_t*split_marker,
106  const dim_t* degree_S, const index_t* offset_S, const index_t* S,
107  const dim_t* degree_ST, const index_t* offset_ST, const index_t* ST,
108  Paso_Connector* col_connector, Paso_Distribution* col_dist);
114 
115 /* Local AMG preconditioner */
118  Paso_SparseMatrix * A_C; /* coarse level matrix */
119  Paso_SparseMatrix * P; /* prolongation n x n_C*/
120  Paso_SparseMatrix * R; /* restriction n_C x n */
121 
125  index_t reordering; /* applied reordering in direct solver */
126  dim_t refinements; /* number of refinements in direct solver (typically =0) */
127  double* r; /* buffer for residual */
128  double* x_C; /* solution of coarse level system */
129  double* b_C; /* right hand side of coarse level system */
131 };
133 
137 
138 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);
139 void Paso_Preconditioner_LocalAMG_setStrongConnections_Block(Paso_SparseMatrix* A, dim_t *degree, index_t *S, const double theta, const double tau);
140 void Paso_Preconditioner_LocalAMG_setStrongConnections(Paso_SparseMatrix* A, dim_t *degree, index_t *S, const double theta, const double tau);
141 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);
144 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);
145 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 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);
150 
151 
153 {
155 };
160 
161 
163 {
170 };
172 
176 
177 /*===============================================*/
178 /* ILU preconditioner */
180  double* factors;
181 };
183 
184 
185 
186 /* RILU preconditioner */
192  double* inv_A_FF;
200  double* x_F;
201  double* b_F;
202  double* x_C;
203  double* b_C;
205 };
207 
208 
209 
210 /* general preconditioner interface */
211 
212 typedef struct Paso_Preconditioner {
215  /* jacobi preconditioner */
217  /* Gauss-Seidel preconditioner */
219  /* amg preconditioner */
221 
222  /* ilu preconditioner */
224  /* rilu preconditioner */
226 
228 
231 void Paso_Preconditioner_solve(Paso_Preconditioner* prec, Paso_SystemMatrix* A,double*,double*);
232 
233 
234 /*******************************************/
237 void Paso_Solver_solveILU(Paso_SparseMatrix * A, Paso_Solver_ILU * ilu, double * x, const double * b);
238 
241 void Paso_Solver_solveRILU(Paso_Solver_RILU * rilu, double * x, double * b);
242 
244 
245 
246 #endif /* #ifndef INC_PRECS */