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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2107 - (show annotations)
Fri Nov 28 04:39:07 2008 UTC (10 years, 10 months ago) by artak
File MIME type: text/plain
File size: 5999 byte(s)
Current version of AMG uses ILU for relaxation, however it is not stable when schur matrix becames more denser. For example it is not stable for -e paramenter is more 400 in 2D for convection problem 
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2008 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
31 #define TOLERANCE_FOR_SCALARS (double)(0.)
32 #define PASO_ONE (double)(1.0)
33 #define PASO_ZERO (double)(0.0)
34
35 /* static double ONE=1.; */
36 /* static double ZERO=0.;*/
37 /*static double TOLERANCE_FOR_SCALARS=0.;*/
38
39 /* jacobi preconditioner */
40
41 typedef struct Paso_Solver_Jacobi {
42 dim_t n_block;
43 dim_t n;
44 double* values;
45 index_t* pivot;
46 } Paso_Solver_Jacobi;
47
48
49 /* ILU preconditioner */
50 struct Paso_Solver_ILU {
51 dim_t n_block;
52 dim_t n;
53 index_t num_colors;
54 index_t* colorOf;
55 index_t* main_iptr;
56 double* factors;
57 Paso_Pattern* pattern;
58 };
59 typedef struct Paso_Solver_ILU Paso_Solver_ILU;
60
61 /* GS preconditioner */
62 struct Paso_Solver_GS {
63 dim_t n_block;
64 dim_t n;
65 index_t num_colors;
66 index_t* colorOf;
67 index_t* main_iptr;
68 double* diag;
69 Paso_SparseMatrix * factors;
70 Paso_Pattern* pattern;
71 dim_t sweeps;
72 double* x_old;
73 };
74 typedef struct Paso_Solver_GS Paso_Solver_GS;
75
76 /* RILU preconditioner */
77 struct Paso_Solver_RILU {
78 dim_t n;
79 dim_t n_block;
80 dim_t n_F;
81 dim_t n_C;
82 double* inv_A_FF;
83 index_t* A_FF_pivot;
84 Paso_SparseMatrix * A_FC;
85 Paso_SparseMatrix * A_CF;
86 index_t* rows_in_F;
87 index_t* rows_in_C;
88 index_t* mask_F;
89 index_t* mask_C;
90 double* x_F;
91 double* b_F;
92 double* x_C;
93 double* b_C;
94 struct Paso_Solver_RILU * RILU_of_Schur;
95 };
96 typedef struct Paso_Solver_RILU Paso_Solver_RILU;
97
98 /* AMG preconditioner */
99 struct Paso_Solver_AMG {
100 dim_t n;
101 dim_t level;
102 dim_t n_block;
103 dim_t n_F;
104 dim_t n_C;
105 double* inv_A_FF;
106 index_t* A_FF_pivot;
107 Paso_SparseMatrix * A_FC;
108 Paso_SparseMatrix * A_CF;
109 index_t* rows_in_F;
110 index_t* rows_in_C;
111 index_t* mask_F;
112 index_t* mask_C;
113 double* x_F;
114 double* b_F;
115 double* x_C;
116 double* b_C;
117 Paso_SparseMatrix * A;
118 Paso_Solver_ILU* GS;
119 struct Paso_Solver_AMG * AMG_of_Schur;
120 };
121 typedef struct Paso_Solver_AMG Paso_Solver_AMG;
122
123
124 /* general preconditioner interface */
125
126 typedef struct Paso_Solver_Preconditioner {
127 dim_t type;
128 /* jacobi preconditioner */
129 Paso_Solver_Jacobi* jacobi;
130 /* ilu preconditioner */
131 Paso_Solver_ILU* ilu;
132 /* rilu preconditioner */
133 Paso_Solver_RILU* rilu;
134 /* Gauss-Seidel preconditioner */
135 Paso_Solver_GS* gs;
136 /* amg preconditioner */
137 Paso_Solver_AMG* amg;
138
139 } Paso_Solver_Preconditioner;
140
141 void Paso_Solver(Paso_SystemMatrix*,double*,double*,Paso_Options*,Paso_Performance* pp);
142 void Paso_Solver_free(Paso_SystemMatrix*);
143 err_t Paso_Solver_BiCGStab( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
144 err_t Paso_Solver_PCG( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
145 err_t Paso_Solver_TFQMR( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
146 err_t Paso_Solver_MINRES( Paso_SystemMatrix * A, double* B, double * X, dim_t *iter, double * tolerance, Paso_Performance* pp);
147 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);
148 void Paso_Preconditioner_free(Paso_Solver_Preconditioner*);
149 void Paso_Solver_setPreconditioner(Paso_SystemMatrix* A,Paso_Options* options);
150 void Paso_Solver_solvePreconditioner(Paso_SystemMatrix* A,double*,double*);
151 void Paso_Solver_applyBlockDiagonalMatrix(dim_t n_block,dim_t n,double* D,index_t* pivot,double* x,double* b);
152
153 void Paso_Solver_ILU_free(Paso_Solver_ILU * in);
154 Paso_Solver_ILU* Paso_Solver_getILU(Paso_SparseMatrix * A_p,bool_t verbose);
155 void Paso_Solver_solveILU(Paso_Solver_ILU * ilu, double * x, double * b);
156
157 void Paso_Solver_GS_free(Paso_Solver_GS * in);
158 Paso_Solver_GS* Paso_Solver_getGS(Paso_SparseMatrix * A_p,bool_t verbose);
159 void Paso_Solver_solveGS(Paso_Solver_GS * gs, double * x, double * b);
160
161 void Paso_Solver_RILU_free(Paso_Solver_RILU * in);
162 Paso_Solver_RILU* Paso_Solver_getRILU(Paso_SparseMatrix * A_p,bool_t verbose);
163 void Paso_Solver_solveRILU(Paso_Solver_RILU * rilu, double * x, double * b);
164
165 void Paso_Solver_AMG_free(Paso_Solver_AMG * in);
166 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix * A_p,bool_t verbose,dim_t level);
167 void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b);
168
169 void Paso_Solver_updateIncompleteSchurComplement(Paso_SparseMatrix* A_CC, Paso_SparseMatrix *A_CF,double* invA_FF,index_t* A_FF_pivot, Paso_SparseMatrix *A_FC);
170 Paso_Solver_Jacobi* Paso_Solver_getJacobi(Paso_SparseMatrix * A_p);
171 void Paso_Solver_solveJacobi(Paso_Solver_Jacobi * prec, double * x, double * b);
172 void Paso_Solver_Jacobi_free(Paso_Solver_Jacobi * in);
173
174 err_t Paso_Solver_GMRES2(Paso_Function * F, const double* f0, const double* x0, double * x, dim_t *iter, double* tolerance, Paso_Performance* pp);
175 err_t Paso_Solver_NewtonGMRES(Paso_Function *F, double *x, Paso_Options* options, Paso_Performance* pp);
176
177 Paso_Function * Paso_Function_LinearSystem_alloc(Paso_SystemMatrix* A, double* b, Paso_Options* options);
178 err_t Paso_Function_LinearSystem_call(Paso_Function * F,double* value, const double* arg, Paso_Performance *pp);
179 void Paso_Function_LinearSystem_free(Paso_Function * F);
180 err_t Paso_Function_LinearSystem_setInitialGuess(Paso_SystemMatrix* A, double* x, Paso_Performance *pp);
181
182 #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