/[escript]/trunk/paso/src/Solver_AMG.c
ViewVC logotype

Contents of /trunk/paso/src/Solver_AMG.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2509 - (show annotations)
Thu Jul 2 04:51:51 2009 UTC (10 years, 7 months ago) by artak
File MIME type: text/plain
File size: 15781 byte(s)
Now one can choose coarsening method from python for AMG. Some commented lines are removed.
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 /**************************************************************/
16
17 /* Paso: AMG preconditioner */
18
19 /**************************************************************/
20
21 /* Author: artak@access.edu.au */
22
23 /**************************************************************/
24
25 #include "Paso.h"
26 #include "Solver.h"
27 #include "Options.h"
28 #include "PasoUtil.h"
29 #include "UMFPACK.h"
30 #include "MKL.h"
31 #include "SystemMatrix.h"
32 #include "Pattern_coupling.h"
33
34 /**************************************************************/
35
36 /* free all memory used by AMG */
37
38 void Paso_Solver_AMG_free(Paso_Solver_AMG * in) {
39 if (in!=NULL) {
40 Paso_Solver_AMG_free(in->AMG_of_Schur);
41 Paso_Solver_Jacobi_free(in->GS);
42 MEMFREE(in->inv_A_FF);
43 MEMFREE(in->A_FF_pivot);
44 Paso_SparseMatrix_free(in->A_FC);
45 Paso_SparseMatrix_free(in->A_CF);
46 MEMFREE(in->rows_in_F);
47 MEMFREE(in->rows_in_C);
48 MEMFREE(in->mask_F);
49 MEMFREE(in->mask_C);
50 MEMFREE(in->x_F);
51 MEMFREE(in->b_F);
52 MEMFREE(in->x_C);
53 MEMFREE(in->b_C);
54 MEMFREE(in);
55 }
56 }
57
58 /**************************************************************/
59
60 /* constructs the block-block factorization of
61
62 [ A_FF A_FC ]
63 A_p=
64 [ A_CF A_FF ]
65
66 to
67
68 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ]
69 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ]
70
71
72 where S=A_FF-ACF*invA_FF*A_FC within the shape of S
73
74 then AMG is applied to S again until S becomes empty
75
76 */
77 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose,dim_t level, double coarsening_threshold, dim_t coarsening_method) {
78 Paso_Solver_AMG* out=NULL;
79 dim_t n=A_p->numRows;
80 dim_t n_block=A_p->row_block_size;
81 index_t* mis_marker=NULL;
82 index_t* counter=NULL;
83 index_t iPtr,*index, *where_p, iPtr_s;
84 dim_t i,k,j;
85 Paso_SparseMatrix * schur=NULL;
86 Paso_SparseMatrix * schur_withFillIn=NULL;
87 double S=0;
88
89 /*Make sure we have block sizes 1*/
90 A_p->pattern->input_block_size=A_p->col_block_size;
91 A_p->pattern->output_block_size=A_p->row_block_size;
92 A_p->pattern->block_size=A_p->block_size;
93
94 /* identify independend set of rows/columns */
95 mis_marker=TMPMEMALLOC(n,index_t);
96 counter=TMPMEMALLOC(n,index_t);
97 out=MEMALLOC(1,Paso_Solver_AMG);
98 out->AMG_of_Schur=NULL;
99 out->inv_A_FF=NULL;
100 out->A_FF_pivot=NULL;
101 out->A_FC=NULL;
102 out->A_CF=NULL;
103 out->rows_in_F=NULL;
104 out->rows_in_C=NULL;
105 out->mask_F=NULL;
106 out->mask_C=NULL;
107 out->x_F=NULL;
108 out->b_F=NULL;
109 out->x_C=NULL;
110 out->b_C=NULL;
111 out->GS=NULL;
112 out->A=Paso_SparseMatrix_getReference(A_p);
113 /*out->GS=Paso_Solver_getGS(A_p,verbose);*/
114 out->GS=Paso_Solver_getJacobi(A_p);
115 /*out->GS->sweeps=2;*/
116 out->level=level;
117
118 if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) || level==0 ) ) {
119 /* identify independend set of rows/columns */
120 #pragma omp parallel for private(i) schedule(static)
121 for (i=0;i<n;++i) mis_marker[i]=-1;
122
123 if (coarsening_method == PASO_YAIR_SHAPIRA_COARSENING) {
124 Paso_Pattern_coup(A_p,mis_marker,coarsening_threshold);
125 }
126 else if (coarsening_method == PASO_RUGE_STUEBEN_COARSENING) {
127 Paso_Pattern_RS(A_p,mis_marker,coarsening_threshold);
128 }
129 else if (coarsening_method == PASO_AGGREGATION_COARSENING) {
130 Paso_Pattern_Aggregiation(A_p,mis_marker,coarsening_threshold);
131 }
132 else {
133 /*Default coarseneing*/
134 Paso_Pattern_RS(A_p,mis_marker,coarsening_threshold);
135 }
136
137 if (Paso_noError()) {
138 #pragma omp parallel for private(i) schedule(static)
139 for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
140 out->n=n;
141 out->n_block=n_block;
142 out->n_F=Paso_Util_cumsum(n,counter);
143 out->mask_F=MEMALLOC(n,index_t);
144 out->rows_in_F=MEMALLOC(out->n_F,index_t);
145 out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
146 out->A_FF_pivot=NULL; /* later use for block size>3 */
147 if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
148 /* creates an index for F from mask */
149 #pragma omp parallel for private(i) schedule(static)
150 for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
151 #pragma omp parallel for private(i) schedule(static)
152 for (i = 0; i < n; ++i) {
153 if (mis_marker[i]) {
154 out->rows_in_F[counter[i]]=i;
155 out->mask_F[i]=counter[i];
156 } else {
157 out->mask_F[i]=-1;
158 }
159 }
160
161 /* Compute row-sum for getting rs(A_FF)^-1*/
162 #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
163 for (i = 0; i < out->n_F; ++i) {
164 S=0;
165 for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
166 j=A_p->pattern->index[iPtr];
167 if (mis_marker[j])
168 S+=A_p->val[iPtr];
169 }
170 out->inv_A_FF[i]=1./S;
171 }
172
173 if( Paso_noError()) {
174 /* if there are no nodes in the coarse level there is no more work to do */
175 out->n_C=n-out->n_F;
176 if (level>0 && out->n_C>500) {
177 /*if (out->n_F>500) {*/
178 out->rows_in_C=MEMALLOC(out->n_C,index_t);
179 out->mask_C=MEMALLOC(n,index_t);
180 if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
181 /* creates an index for C from mask */
182 #pragma omp parallel for private(i) schedule(static)
183 for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
184 Paso_Util_cumsum(n,counter);
185 #pragma omp parallel for private(i) schedule(static)
186 for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
187 #pragma omp parallel for private(i) schedule(static)
188 for (i = 0; i < n; ++i) {
189 if (! mis_marker[i]) {
190 out->rows_in_C[counter[i]]=i;
191 out->mask_C[i]=counter[i];
192 } else {
193 out->mask_C[i]=-1;
194 }
195 }
196 /* get A_CF block: */
197 out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
198 if (Paso_noError()) {
199 /* get A_FC block: */
200 out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
201 /* get A_CC block: */
202 if (Paso_noError()) {
203 schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
204 /*find the pattern of the schur complement with fill in*/
205
206 schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern)),1,1);
207 /* copy values over*/
208 #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
209 for (i = 0; i < schur_withFillIn->numRows; ++i) {
210 for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
211 j=schur_withFillIn->pattern->index[iPtr];
212 iPtr_s=schur->pattern->ptr[i];
213 schur_withFillIn->val[iPtr]=0.;
214 index=&(schur->pattern->index[iPtr_s]);
215 where_p=(index_t*)bsearch(&j,
216 index,
217 schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
218 sizeof(index_t),
219 Paso_comparIndex);
220 if (where_p!=NULL) {
221 schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
222 }
223 }
224 }
225 if (Paso_noError()) {
226 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
227 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level-1,coarsening_threshold,coarsening_method);
228 Paso_SparseMatrix_free(schur);
229 }
230 /* allocate work arrays for AMG application */
231 if (Paso_noError()) {
232 out->x_F=MEMALLOC(n_block*out->n_F,double);
233 out->b_F=MEMALLOC(n_block*out->n_F,double);
234 out->x_C=MEMALLOC(n_block*out->n_C,double);
235 out->b_C=MEMALLOC(n_block*out->n_C,double);
236
237 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
238 #pragma omp parallel for private(i,k) schedule(static)
239 for (i = 0; i < out->n_F; ++i) {
240 for (k=0; k<n_block;++k) {
241 out->x_F[i*n_block+k]=0.;
242 out->b_F[i*n_block+k]=0.;
243 }
244 }
245 #pragma omp parallel for private(i,k) schedule(static)
246 for (i = 0; i < out->n_C; ++i) {
247 for (k=0; k<n_block;++k) {
248 out->x_C[i*n_block+k]=0.;
249 out->b_C[i*n_block+k]=0.;
250 }
251 }
252 }
253 }
254 }
255 }
256 }
257 }
258 }
259 }
260 }
261 }
262 TMPMEMFREE(mis_marker);
263 TMPMEMFREE(counter);
264 if (Paso_noError()) {
265 if (verbose && level>0) {
266 printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,out->n_C);
267 }
268 return out;
269 } else {
270 Paso_Solver_AMG_free(out);
271 return NULL;
272 }
273 }
274
275 /**************************************************************/
276
277 /* apply AMG precondition b-> x
278
279 in fact it solves
280
281 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FC ] [ x_F ] = [b_F]
282 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
283
284 in the form
285
286 b->[b_F,b_C]
287 x_F=invA_FF*b_F
288 b_C=b_C-A_CF*x_F
289 x_C=AMG(b_C)
290 b_F=b_F-A_FC*x_C
291 x_F=invA_FF*b_F
292 x<-[x_F,x_C]
293
294 should be called within a parallel region
295 barrier synconization should be performed to make sure that the input vector available
296
297 */
298
299 void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
300 dim_t i;
301 double *r=MEMALLOC(amg->n,double);
302 double *x0=MEMALLOC(amg->n,double);
303 double time0=0;
304 bool_t verbose=0;
305
306 #ifdef MKL
307 Paso_SparseMatrix *temp=NULL;
308 #endif
309
310
311 if (amg->level==0 || amg->n_C<=500) {
312
313 time0=Paso_timer();
314 #ifdef MKL
315 temp=Paso_SparseMatrix_alloc(MATRIX_FORMAT_BLK1 + MATRIX_FORMAT_OFFSET1, amg->A->pattern,1,1);
316 #pragma omp parallel for private(i) schedule(static)
317 for (i=0;i<amg->A->len;++i) {
318 temp->val[i]=amg->A->val[i];
319 }
320 Paso_MKL1(temp,x,b,0);
321 Paso_SparseMatrix_free(temp);
322 #else
323 #ifdef UMFPACK
324 Paso_UMFPACK1(amg->A,x,b,0);
325 #else
326 Paso_Solver_solveJacobi(amg->GS,x,b);
327 #endif
328 #endif
329
330 time0=Paso_timer()-time0;
331 if (verbose) fprintf(stderr,"timing: DIRECT SOLVER: %e\n\n\n",time0);
332
333 } else {
334 /* presmoothing */
335 time0=Paso_timer();
336 Paso_Solver_solveJacobi(amg->GS,x,b);
337 time0=Paso_timer()-time0;
338 if (verbose) fprintf(stderr,"timing: Presmooting: %e\n",time0);
339 /* end of presmoothing */
340
341
342 time0=Paso_timer();
343 #pragma omp parallel for private(i) schedule(static)
344 for (i=0;i<amg->n;++i) r[i]=b[i];
345
346 /*r=b-Ax*/
347 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
348
349 /* b->[b_F,b_C] */
350 #pragma omp parallel for private(i) schedule(static)
351 for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];
352
353 #pragma omp parallel for private(i) schedule(static)
354 for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
355
356 /* x_F=invA_FF*b_F */
357 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
358
359 /* b_C=b_C-A_CF*x_F */
360 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
361
362 time0=Paso_timer()-time0;
363 if (verbose) fprintf(stderr,"timing: Before next level: %e\n",time0);
364
365 /* x_C=AMG(b_C) */
366 Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
367
368 time0=Paso_timer();
369 /* b_F=b_F-A_FC*x_C */
370 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
371 /* x_F=invA_FF*b_F */
372 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
373 /* x<-[x_F,x_C] */
374
375 #pragma omp parallel for private(i) schedule(static)
376 for (i=0;i<amg->n;++i) {
377 if (amg->mask_C[i]>-1) {
378 x[i]=amg->x_C[amg->mask_C[i]];
379 } else {
380 x[i]=amg->x_F[amg->mask_F[i]];
381 }
382 }
383
384 time0=Paso_timer()-time0;
385 if (verbose) fprintf(stderr,"timing: After next level: %e\n",time0);
386
387 /*postsmoothing*/
388 time0=Paso_timer();
389 #pragma omp parallel for private(i) schedule(static)
390 for (i=0;i<amg->n;++i) r[i]=b[i];
391
392 /*r=b-Ax */
393 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
394 Paso_Solver_solveJacobi(amg->GS,x0,r);
395
396 #pragma omp parallel for private(i) schedule(static)
397 for (i=0;i<amg->n;++i) x[i]+=x0[i];
398
399 time0=Paso_timer()-time0;
400 if (verbose) fprintf(stderr,"timing: Postsmoothing: %e\n",time0);
401
402 /*end of postsmoothing*/
403
404 }
405 MEMFREE(r);
406 MEMFREE(x0);
407 return;
408 }
409
410 /*
411 * $Log$
412 *
413 */

  ViewVC Help
Powered by ViewVC 1.1.26