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

  ViewVC Help
Powered by ViewVC 1.1.26