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

  ViewVC Help
Powered by ViewVC 1.1.26