/[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 2661 - (show annotations)
Fri Sep 11 00:59:59 2009 UTC (10 years, 5 months ago) by artak
File MIME type: text/plain
File size: 17003 byte(s)
Finially unit tests for AMG preconditioner for single equations as well as for systems. Minor bug in AMG is also fixed
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2009 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@uq.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_System_free(Paso_Solver_AMG_System * in) {
39 if (in!=NULL) {
40 Paso_Solver_AMG_free(in->amgblock1);
41 Paso_Solver_AMG_free(in->amgblock2);
42 Paso_Solver_AMG_free(in->amgblock3);
43 Paso_SparseMatrix_free(in->block1);
44 Paso_SparseMatrix_free(in->block2);
45 Paso_SparseMatrix_free(in->block3);
46 MEMFREE(in);
47 }
48 }
49
50
51 /* free all memory used by AMG */
52
53 void Paso_Solver_AMG_free(Paso_Solver_AMG * in) {
54 if (in!=NULL) {
55 Paso_Solver_AMG_free(in->AMG_of_Schur);
56 Paso_Solver_Jacobi_free(in->GS);
57 MEMFREE(in->inv_A_FF);
58 MEMFREE(in->A_FF_pivot);
59 Paso_SparseMatrix_free(in->A_FC);
60 Paso_SparseMatrix_free(in->A_CF);
61 Paso_SparseMatrix_free(in->A);
62 MEMFREE(in->rows_in_F);
63 MEMFREE(in->rows_in_C);
64 MEMFREE(in->mask_F);
65 MEMFREE(in->mask_C);
66 MEMFREE(in->x_F);
67 MEMFREE(in->b_F);
68 MEMFREE(in->x_C);
69 MEMFREE(in->b_C);
70 #ifdef MKL
71 #else
72 #ifdef UMFPACK
73 Paso_UMFPACK1_free((Paso_UMFPACK_Handler*)(in->solver));
74 #endif
75 #endif
76 in->solver=NULL;
77 MEMFREE(in->b_C);
78 MEMFREE(in);
79 }
80 }
81
82 /**************************************************************/
83
84 /* constructs the block-block factorization of
85
86 [ A_FF A_FC ]
87 A_p=
88 [ A_CF A_FF ]
89
90 to
91
92 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ]
93 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ]
94
95
96 where S=A_FF-ACF*invA_FF*A_FC within the shape of S
97
98 then AMG is applied to S again until S becomes empty
99
100 */
101 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,dim_t level,Paso_Options* options) {
102 Paso_Solver_AMG* out=NULL;
103 Paso_Pattern* temp1=NULL;
104 Paso_Pattern* temp2=NULL;
105 bool_t verbose=options->verbose;
106 dim_t n=A_p->numRows;
107 dim_t n_block=A_p->row_block_size;
108 index_t* mis_marker=NULL;
109 index_t* counter=NULL;
110 index_t iPtr,*index, *where_p, iPtr_s;
111 dim_t i,k,j;
112 Paso_SparseMatrix * schur=NULL;
113 Paso_SparseMatrix * schur_withFillIn=NULL;
114 double S=0;
115
116 /*Make sure we have block sizes 1*/
117 if (A_p->col_block_size>1) {
118 Paso_setError(TYPE_ERROR,"Paso_Solver_getAMG: AMG requires column block size 1.");
119 return NULL;
120 }
121 if (A_p->row_block_size>1) {
122 Paso_setError(TYPE_ERROR,"Paso_Solver_getAMG: AMG requires row block size 1.");
123 return NULL;
124 }
125 out=MEMALLOC(1,Paso_Solver_AMG);
126 /* identify independend set of rows/columns */
127 mis_marker=TMPMEMALLOC(n,index_t);
128 counter=TMPMEMALLOC(n,index_t);
129 if ( !( Paso_checkPtr(mis_marker) || Paso_checkPtr(counter) || Paso_checkPtr(out)) ) {
130 out->AMG_of_Schur=NULL;
131 out->inv_A_FF=NULL;
132 out->A_FF_pivot=NULL;
133 out->A_FC=NULL;
134 out->A_CF=NULL;
135 out->rows_in_F=NULL;
136 out->rows_in_C=NULL;
137 out->mask_F=NULL;
138 out->mask_C=NULL;
139 out->x_F=NULL;
140 out->b_F=NULL;
141 out->x_C=NULL;
142 out->b_C=NULL;
143 out->GS=NULL;
144 out->A=Paso_SparseMatrix_getReference(A_p);
145 out->GS=NULL;
146 out->solver=NULL;
147 /*out->GS=Paso_Solver_getGS(A_p,verbose);*/
148 out->level=level;
149 out->n=n;
150 out->n_block=n_block;
151
152 if (level==0 || n<=options->min_coarse_matrix_size) {
153 out->coarsest_level=TRUE;
154 #ifdef MKL
155 #else
156 #ifdef UMFPACK
157 #else
158 out->GS=Paso_Solver_getJacobi(A_p);
159 #endif
160 #endif
161 } else {
162 out->coarsest_level=FALSE;
163 out->GS=Paso_Solver_getJacobi(A_p);
164
165 /* identify independend set of rows/columns */
166 #pragma omp parallel for private(i) schedule(static)
167 for (i=0;i<n;++i) mis_marker[i]=-1;
168
169 if (options->coarsening_method == PASO_YAIR_SHAPIRA_COARSENING) {
170 Paso_Pattern_YS(A_p,mis_marker,options->coarsening_threshold);
171 }
172 else if (options->coarsening_method == PASO_RUGE_STUEBEN_COARSENING) {
173 Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
174 }
175 else if (options->coarsening_method == PASO_AGGREGATION_COARSENING) {
176 Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);
177 }
178 else {
179 /*Default coarseneing*/
180 Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
181 /*Paso_Pattern_YS(A_p,mis_marker,options->coarsening_threshold);*/
182 /*Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);*/
183 }
184
185 if (Paso_noError()) {
186 #pragma omp parallel for private(i) schedule(static)
187 for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
188
189 out->n_F=Paso_Util_cumsum(n,counter);
190 out->mask_F=MEMALLOC(n,index_t);
191 out->rows_in_F=MEMALLOC(out->n_F,index_t);
192 out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
193 out->A_FF_pivot=NULL; /* later use for block size>3 */
194 if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
195 /* creates an index for F from mask */
196 #pragma omp parallel for private(i) schedule(static)
197 for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
198 #pragma omp parallel for private(i) schedule(static)
199 for (i = 0; i < n; ++i) {
200 if (mis_marker[i]) {
201 out->rows_in_F[counter[i]]=i;
202 out->mask_F[i]=counter[i];
203 } else {
204 out->mask_F[i]=-1;
205 }
206 }
207
208 /* Compute row-sum for getting rs(A_FF)^-1*/
209 #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
210 for (i = 0; i < out->n_F; ++i) {
211 S=0;
212 /*printf("[%d ]: [%d] -> ",i, out->rows_in_F[i]);*/
213 for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
214 j=A_p->pattern->index[iPtr];
215 /*if (j==out->rows_in_F[i]) printf("diagonal %e",A_p->val[iPtr]);*/
216 if (mis_marker[j])
217 S+=A_p->val[iPtr];
218 }
219 /*printf("-> %e \n",S);*/
220 out->inv_A_FF[i]=1./S;
221 }
222 }
223 }
224
225 if ( Paso_noError()) {
226 /* if there are no nodes in the coarse level there is no more work to do */
227 out->n_C=n-out->n_F;
228
229 /*if (out->n_F>500) */
230 out->rows_in_C=MEMALLOC(out->n_C,index_t);
231 out->mask_C=MEMALLOC(n,index_t);
232 if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
233 /* creates an index for C from mask */
234 #pragma omp parallel for private(i) schedule(static)
235 for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
236 Paso_Util_cumsum(n,counter);
237 #pragma omp parallel for private(i) schedule(static)
238 for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
239 #pragma omp parallel for private(i) schedule(static)
240 for (i = 0; i < n; ++i) {
241 if (! mis_marker[i]) {
242 out->rows_in_C[counter[i]]=i;
243 out->mask_C[i]=counter[i];
244 } else {
245 out->mask_C[i]=-1;
246 }
247 }
248 }
249 }
250 if ( Paso_noError()) {
251 /* get A_CF block: */
252 out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
253 /* get A_FC block: */
254 out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
255 /* get A_CC block: */
256 schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
257
258 }
259 if ( Paso_noError()) {
260 /*find the pattern of the schur complement with fill in*/
261 temp1=Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern);
262 temp2=Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, temp1);
263 schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,temp2,1,1, TRUE);
264 Paso_Pattern_free(temp1);
265 Paso_Pattern_free(temp2);
266 }
267 if ( Paso_noError()) {
268 /* copy values over*/
269 #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
270 for (i = 0; i < schur_withFillIn->numRows; ++i) {
271 for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
272 j=schur_withFillIn->pattern->index[iPtr];
273 iPtr_s=schur->pattern->ptr[i];
274 schur_withFillIn->val[iPtr]=0.;
275 index=&(schur->pattern->index[iPtr_s]);
276 where_p=(index_t*)bsearch(&j,
277 index,
278 schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
279 sizeof(index_t),
280 Paso_comparIndex);
281 if (where_p!=NULL) {
282 schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
283 }
284 }
285 }
286 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
287 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,level-1,options);
288 }
289 /* allocate work arrays for AMG application */
290 if (Paso_noError()) {
291 out->x_F=MEMALLOC(n_block*out->n_F,double);
292 out->b_F=MEMALLOC(n_block*out->n_F,double);
293 out->x_C=MEMALLOC(n_block*out->n_C,double);
294 out->b_C=MEMALLOC(n_block*out->n_C,double);
295
296 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
297 #pragma omp parallel for private(i,k) schedule(static)
298 for (i = 0; i < out->n_F; ++i) {
299 for (k=0; k<n_block;++k) {
300 out->x_F[i*n_block+k]=0.;
301 out->b_F[i*n_block+k]=0.;
302 }
303 }
304 #pragma omp parallel for private(i,k) schedule(static)
305 for (i = 0; i < out->n_C; ++i) {
306 for (k=0; k<n_block;++k) {
307 out->x_C[i*n_block+k]=0.;
308 out->b_C[i*n_block+k]=0.;
309 }
310 }
311 }
312 }
313 }
314 }
315 TMPMEMFREE(mis_marker);
316 TMPMEMFREE(counter);
317 Paso_SparseMatrix_free(schur);
318 Paso_SparseMatrix_free(schur_withFillIn);
319 if (Paso_noError()) {
320 if (verbose && level>0 && !out->coarsest_level) {
321 printf("AMG: level: %d: %d unknowns eliminated. %d left.\n",level, out->n_F,out->n_C);
322 }
323 return out;
324 } else {
325 Paso_Solver_AMG_free(out);
326 return NULL;
327 }
328 }
329
330 /**************************************************************/
331
332 /* apply AMG precondition b-> x
333
334 in fact it solves
335
336 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FC ] [ x_F ] = [b_F]
337 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
338
339 in the form
340
341 b->[b_F,b_C]
342 x_F=invA_FF*b_F
343 b_C=b_C-A_CF*x_F
344 x_C=AMG(b_C)
345 b_F=b_F-A_FC*x_C
346 x_F=invA_FF*b_F
347 x<-[x_F,x_C]
348
349 should be called within a parallel region
350 barrier synconization should be performed to make sure that the input vector available
351
352 */
353
354 void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
355 dim_t i;
356 double time0=0;
357 double *r=NULL, *x0=NULL;
358 bool_t verbose=0;
359 #ifdef MKL
360 Paso_SparseMatrix *temp=NULL;
361 #else
362 #ifdef UMFPACK
363 Paso_UMFPACK_Handler * ptr=NULL;
364 #endif
365 #endif
366 r=MEMALLOC(amg->n,double);
367 x0=MEMALLOC(amg->n,double);
368
369 if (amg->coarsest_level) {
370
371 time0=Paso_timer();
372 #ifdef MKL
373 temp=Paso_SparseMatrix_alloc(MATRIX_FORMAT_BLK1 + MATRIX_FORMAT_OFFSET1, amg->A->pattern,1,1, FALSE);
374 #pragma omp parallel for private(i) schedule(static)
375 for (i=0;i<amg->A->len;++i) {
376 temp->val[i]=amg->A->val[i];
377 }
378 Paso_MKL1(temp,x,b,0);
379 Paso_SparseMatrix_free(temp);
380 #else
381 #ifdef UMFPACK
382 ptr=(Paso_UMFPACK_Handler *)(amg->solver);
383 Paso_UMFPACK1(&ptr,amg->A,x,b,verbose);
384 amg->solver=(void*) ptr;
385 #else
386 Paso_Solver_solveJacobi(amg->GS,x,b);
387 #endif
388 #endif
389 time0=Paso_timer()-time0;
390 if (verbose) fprintf(stderr,"timing: DIRECT SOLVER: %e\n",time0);
391
392 } else {
393 /* presmoothing */
394 time0=Paso_timer();
395 Paso_Solver_solveJacobi(amg->GS,x,b);
396 time0=Paso_timer()-time0;
397 if (verbose) fprintf(stderr,"timing: Presmooting: %e\n",time0);
398 /* end of presmoothing */
399
400
401 time0=Paso_timer();
402 #pragma omp parallel for private(i) schedule(static)
403 for (i=0;i<amg->n;++i) r[i]=b[i];
404
405 /*r=b-Ax*/
406 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
407
408 /* b->[b_F,b_C] */
409 #pragma omp parallel for private(i) schedule(static)
410 for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];
411
412 #pragma omp parallel for private(i) schedule(static)
413 for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
414
415 /* x_F=invA_FF*b_F */
416 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
417
418 /* b_C=b_C-A_CF*x_F */
419 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
420
421 time0=Paso_timer()-time0;
422 if (verbose) fprintf(stderr,"timing: Before next level: %e\n",time0);
423
424 /* x_C=AMG(b_C) */
425 Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
426
427 time0=Paso_timer();
428 /* b_F=b_F-A_FC*x_C */
429 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
430 /* x_F=invA_FF*b_F */
431 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
432 /* x<-[x_F,x_C] */
433
434 #pragma omp parallel for private(i) schedule(static)
435 for (i=0;i<amg->n;++i) {
436 if (amg->mask_C[i]>-1) {
437 x[i]=amg->x_C[amg->mask_C[i]];
438 } else {
439 x[i]=amg->x_F[amg->mask_F[i]];
440 }
441 }
442
443 time0=Paso_timer()-time0;
444 if (verbose) fprintf(stderr,"timing: After next level: %e\n",time0);
445
446 /*postsmoothing*/
447 time0=Paso_timer();
448 #pragma omp parallel for private(i) schedule(static)
449 for (i=0;i<amg->n;++i) r[i]=b[i];
450
451 /*r=b-Ax */
452 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
453 Paso_Solver_solveJacobi(amg->GS,x0,r);
454
455 #pragma omp parallel for private(i) schedule(static)
456 for (i=0;i<amg->n;++i) x[i]+=x0[i];
457
458 time0=Paso_timer()-time0;
459 if (verbose) fprintf(stderr,"timing: Postsmoothing: %e\n",time0);
460
461 /*end of postsmoothing*/
462
463 }
464 MEMFREE(r);
465 MEMFREE(x0);
466 return;
467 }

  ViewVC Help
Powered by ViewVC 1.1.26