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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3402 - (show annotations)
Tue Dec 7 07:36:12 2010 UTC (8 years, 9 months ago) by gross
File MIME type: text/plain
File size: 6646 byte(s)
AMG support now block size > 3 (requires clapack or MKL).

additional diagnostics for AMG 


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
21
22 #define PASO_AMG_UNDECIDED -1
23 #define PASO_AMG_IN_F 0
24 #define PASO_AMG_IN_C 1
25
26 /* GAUSS SEIDEL & Jacobi */
27 typedef struct Paso_Preconditioner_LocalSmoother {
28 bool_t Jacobi;
29 double* diag;
30 double* buffer;
31 index_t* pivot;
32 } Paso_Preconditioner_LocalSmoother;
33
34 typedef struct Paso_Preconditioner_Smoother {
35 Paso_Preconditioner_LocalSmoother* localSmoother;
36 bool_t is_local;
37 } Paso_Preconditioner_Smoother;
38
39 void Paso_Preconditioner_Smoother_free(Paso_Preconditioner_Smoother * in);
40 void Paso_Preconditioner_LocalSmoother_free(Paso_Preconditioner_LocalSmoother * in);
41
42 Paso_Preconditioner_Smoother* Paso_Preconditioner_Smoother_alloc(Paso_SystemMatrix * A_p, const bool_t jacobi, const bool_t is_local, const bool_t verbose);
43 Paso_Preconditioner_LocalSmoother* Paso_Preconditioner_LocalSmoother_alloc(Paso_SparseMatrix * A_p,const bool_t jacobi, const bool_t verbose);
44
45 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);
46 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);
47
48 void Paso_Preconditioner_LocalSmoother_Sweep(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
49 void Paso_Preconditioner_LocalSmoother_Sweep_sequential(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
50 void Paso_Preconditioner_LocalSmoother_Sweep_tiled(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
51 void Paso_Preconditioner_LocalSmoother_Sweep_colored(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * gs, double * x);
52
53
54 /* AMG preconditioner */
55 struct Paso_Preconditioner_LocalAMG {
56 dim_t level;
57 dim_t n;
58 dim_t n_F;
59 dim_t n_block;
60 Paso_SparseMatrix * A_C; /* coarse level matrix */
61 Paso_SparseMatrix * P; /* prolongation n x n_C*/
62 Paso_SparseMatrix * R; /* restriction n_C x n */
63
64 Paso_Preconditioner_LocalSmoother* Smoother;
65 dim_t post_sweeps;
66 dim_t pre_sweeps;
67 index_t reordering; /* applied reordering in direct solver */
68 dim_t refinements; /* number of refinements in direct solver (typically =0) */
69 double* r; /* buffer for residual */
70 double* x_C; /* solution of coarse level system */
71 double* b_C; /* right hand side of coarse level system */
72 struct Paso_Preconditioner_LocalAMG * AMG_C;
73 };
74 typedef struct Paso_Preconditioner_LocalAMG Paso_Preconditioner_LocalAMG;
75
76 void Paso_Preconditioner_LocalAMG_free(Paso_Preconditioner_LocalAMG * in);
77 Paso_Preconditioner_LocalAMG* Paso_Preconditioner_LocalAMG_alloc(Paso_SparseMatrix * A_p,dim_t level,Paso_Options* options);
78 void Paso_Preconditioner_LocalAMG_solve(Paso_SparseMatrix* A, Paso_Preconditioner_LocalAMG * amg, double * x, double * b);
79 void Paso_Preconditioner_AMG_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);
80 void Paso_Preconditioner_AMG_setStrongConnections_Block(Paso_SparseMatrix* A, dim_t *degree, index_t *S, const double theta, const double tau);
81 void Paso_Preconditioner_AMG_setStrongConnections(Paso_SparseMatrix* A, dim_t *degree, index_t *S, const double theta, const double tau);
82 Paso_SparseMatrix* Paso_Preconditioner_AMG_getDirectProlongation(const Paso_SparseMatrix* A_p, const dim_t* degree, const index_t* S, const dim_t n_C, const index_t* counter_C);
83 void Paso_Preconditioner_AMG_setDirectProlongation_Block(Paso_SparseMatrix* P_p, const Paso_SparseMatrix* A_p, const index_t *counter_C);
84 void Paso_Preconditioner_AMG_setDirectProlongation(Paso_SparseMatrix* P_p, const Paso_SparseMatrix* A_p, const index_t *counter_C);
85 index_t Paso_Preconditioner_LocalAMG_getMaxLevel(const Paso_Preconditioner_LocalAMG * in);
86 double Paso_Preconditioner_LocalAMG_getCoarseLevelSparsity(const Paso_Preconditioner_LocalAMG * in);
87 dim_t Paso_Preconditioner_LocalAMG_getNumCoarseUnknwons(const Paso_Preconditioner_LocalAMG * in);
88 /*===============================================*/
89 /* ILU preconditioner */
90 struct Paso_Solver_ILU {
91 double* factors;
92 };
93 typedef struct Paso_Solver_ILU Paso_Solver_ILU;
94
95
96
97 /* RILU preconditioner */
98 struct Paso_Solver_RILU {
99 dim_t n;
100 dim_t n_block;
101 dim_t n_F;
102 dim_t n_C;
103 double* inv_A_FF;
104 index_t* A_FF_pivot;
105 Paso_SparseMatrix * A_FC;
106 Paso_SparseMatrix * A_CF;
107 index_t* rows_in_F;
108 index_t* rows_in_C;
109 index_t* mask_F;
110 index_t* mask_C;
111 double* x_F;
112 double* b_F;
113 double* x_C;
114 double* b_C;
115 struct Paso_Solver_RILU * RILU_of_Schur;
116 };
117 typedef struct Paso_Solver_RILU Paso_Solver_RILU;
118
119
120
121 /* general preconditioner interface */
122
123 typedef struct Paso_Preconditioner {
124 dim_t type;
125 dim_t sweeps;
126 /* jacobi preconditioner */
127 Paso_Preconditioner_Smoother* jacobi;
128 /* Gauss-Seidel preconditioner */
129 Paso_Preconditioner_Smoother* gs;
130 /* amg preconditioner */
131 Paso_Preconditioner_LocalAMG* localamg;
132 Paso_Preconditioner_LocalSmoother* localamgsubstitute;
133
134 /* ilu preconditioner */
135 Paso_Solver_ILU* ilu;
136 /* rilu preconditioner */
137 Paso_Solver_RILU* rilu;
138
139 } Paso_Preconditioner;
140
141 void Paso_Preconditioner_free(Paso_Preconditioner*);
142 Paso_Preconditioner* Paso_Preconditioner_alloc(Paso_SystemMatrix* A,Paso_Options* options);
143 void Paso_Preconditioner_solve(Paso_Preconditioner* prec, Paso_SystemMatrix* A,double*,double*);
144
145
146 /*******************************************/
147 void Paso_Solver_ILU_free(Paso_Solver_ILU * in);
148 Paso_Solver_ILU* Paso_Solver_getILU(Paso_SparseMatrix * A_p,bool_t verbose);
149 void Paso_Solver_solveILU(Paso_SparseMatrix * A, Paso_Solver_ILU * ilu, double * x, const double * b);
150
151 void Paso_Solver_RILU_free(Paso_Solver_RILU * in);
152 Paso_Solver_RILU* Paso_Solver_getRILU(Paso_SparseMatrix * A_p,bool_t verbose);
153 void Paso_Solver_solveRILU(Paso_Solver_RILU * rilu, double * x, double * b);
154
155 void Paso_Solver_updateIncompleteSchurComplement(Paso_SparseMatrix* A_CC, Paso_SparseMatrix *A_CF,double* invA_FF,index_t* A_FF_pivot, Paso_SparseMatrix *A_FC);
156
157
158 #endif /* #ifndef INC_PRECS */

  ViewVC Help
Powered by ViewVC 1.1.26