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

  ViewVC Help
Powered by ViewVC 1.1.26