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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1892 - (hide annotations)
Fri Oct 17 00:55:59 2008 UTC (10 years, 11 months ago) by artak
File MIME type: text/plain
File size: 18246 byte(s)
memory leak is fixed
1 jfenwick 1851
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 with reordering */
18    
19     /**************************************************************/
20    
21     /* Author: artak@access.edu.au */
22    
23     /**************************************************************/
24    
25     #include "Paso.h"
26     #include "Solver.h"
27     #include "PasoUtil.h"
28    
29     /**************************************************************/
30    
31     /* free all memory used by AMG */
32    
33     void Paso_Solver_AMG_free(Paso_Solver_AMG * in) {
34     if (in!=NULL) {
35     Paso_Solver_AMG_free(in->AMG_of_Schur);
36     MEMFREE(in->inv_A_FF);
37     MEMFREE(in->A_FF_pivot);
38     Paso_SparseMatrix_free(in->A_FC);
39     Paso_SparseMatrix_free(in->A_CF);
40     MEMFREE(in->rows_in_F);
41     MEMFREE(in->rows_in_C);
42     MEMFREE(in->mask_F);
43     MEMFREE(in->mask_C);
44     MEMFREE(in->x_F);
45     MEMFREE(in->b_F);
46     MEMFREE(in->x_C);
47     MEMFREE(in->b_C);
48     MEMFREE(in);
49     }
50     }
51    
52     /**************************************************************/
53    
54     /* constructs the block-block factorization of
55    
56     [ A_FF A_FC ]
57     A_p=
58     [ A_CF A_FF ]
59    
60     to
61    
62     [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ]
63     [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ]
64    
65    
66     where S=A_FF-ACF*invA_FF*A_FC within the shape of S
67    
68     then AMG is applied to S again until S becomes empty
69    
70     */
71 artak 1890 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose,dim_t level) {
72 jfenwick 1851 Paso_Solver_AMG* out=NULL;
73     dim_t n=A_p->numRows;
74     dim_t n_block=A_p->row_block_size;
75     index_t* mis_marker=NULL;
76 artak 1881 index_t* counter=NULL;
77     double *rs=NULL;
78     index_t iPtr,*index, *where_p, iPtr_s;
79     dim_t i,k,j,j0;
80 jfenwick 1851 Paso_SparseMatrix * schur=NULL;
81 artak 1881 Paso_SparseMatrix * schur_withFillIn=NULL;
82 artak 1891 double A11,A12,A13,A21,A22,A23,A31,A32,A33,D,time0,time1,time2;
83 artak 1890
84 jfenwick 1851 /* identify independend set of rows/columns */
85     mis_marker=TMPMEMALLOC(n,index_t);
86     counter=TMPMEMALLOC(n,index_t);
87 artak 1881 rs=TMPMEMALLOC(n,double);
88 jfenwick 1851 out=MEMALLOC(1,Paso_Solver_AMG);
89     out->AMG_of_Schur=NULL;
90     out->inv_A_FF=NULL;
91     out->A_FF_pivot=NULL;
92     out->A_FC=NULL;
93     out->A_CF=NULL;
94     out->rows_in_F=NULL;
95     out->rows_in_C=NULL;
96     out->mask_F=NULL;
97     out->mask_C=NULL;
98     out->x_F=NULL;
99     out->b_F=NULL;
100     out->x_C=NULL;
101     out->b_C=NULL;
102 artak 1890 out->A=Paso_SparseMatrix_getReference(A_p);
103     out->level=level;
104 artak 1881
105     /* fprintf(stderr,"START OF MATRIX \n\n");
106     for (i = 0; i < A_p->numRows; ++i) {
107     for (iPtr=A_p->pattern->ptr[i];iPtr<A_p->pattern->ptr[i + 1]; ++iPtr) {
108     j=A_p->pattern->index[iPtr];
109     fprintf(stderr,"A[%d,%d]=%.2f ",i,j,A_p->val[iPtr]);
110     }
111     fprintf(stderr,"\n");
112     }
113     fprintf(stderr,"END OF MATRIX \n\n");
114     */
115 jfenwick 1851 if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) {
116     /* identify independend set of rows/columns */
117     time0=Paso_timer();
118     #pragma omp parallel for private(i) schedule(static)
119     for (i=0;i<n;++i) mis_marker[i]=-1;
120 artak 1890 Paso_Pattern_RS(A_p,mis_marker,0.25);
121     /*
122     for (i=0;i<n;++i) fprintf(stderr," i=%d mis[i]=%d \n",i,mis_marker[i]);
123     */
124 jfenwick 1851 time2=Paso_timer()-time0;
125     if (Paso_noError()) {
126     #pragma omp parallel for private(i) schedule(static)
127     for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
128     out->n=n;
129     out->n_block=n_block;
130     out->n_F=Paso_Util_cumsum(n,counter);
131     out->mask_F=MEMALLOC(n,index_t);
132     out->rows_in_F=MEMALLOC(out->n_F,index_t);
133     out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
134     out->A_FF_pivot=NULL; /* later use for block size>3 */
135     if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
136     #pragma omp parallel
137     {
138     /* creates an index for F from mask */
139     #pragma omp for private(i) schedule(static)
140     for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
141     #pragma omp for private(i) schedule(static)
142     for (i = 0; i < n; ++i) {
143     if (mis_marker[i]) {
144     out->rows_in_F[counter[i]]=i;
145     out->mask_F[i]=counter[i];
146     } else {
147     out->mask_F[i]=-1;
148     }
149     }
150 artak 1881 /* Compute row-sum for getting rs(A_FF)*/
151     #pragma omp for private(i,iPtr) schedule(static)
152     for (i = 0; i < out->n_F; ++i) {
153     rs[i]=0;
154     for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
155     rs[i]+=A_p->val[iPtr];
156     }
157     }
158    
159 jfenwick 1851 #pragma omp for private(i, where_p,iPtr,A11,A12,A13,A21,A22,A23,A31,A32,A33,D,index) schedule(static)
160     for (i = 0; i < out->n_F; i++) {
161     /* find main diagonal */
162     iPtr=A_p->pattern->ptr[out->rows_in_F[i]];
163     index=&(A_p->pattern->index[iPtr]);
164     where_p=(index_t*)bsearch(&out->rows_in_F[i],
165     index,
166     A_p->pattern->ptr[out->rows_in_F[i] + 1]-A_p->pattern->ptr[out->rows_in_F[i]],
167     sizeof(index_t),
168     Paso_comparIndex);
169     if (where_p==NULL) {
170     Paso_setError(VALUE_ERROR, "Paso_Solver_getAMG: main diagonal element missing.");
171     } else {
172     iPtr+=(index_t)(where_p-index);
173     /* get inverse of A_FF block: */
174 artak 1881 if (ABS(rs[i])>0.) {
175     out->inv_A_FF[i]=1./rs[i];
176     } else {
177 artak 1890 out->inv_A_FF[i]=0;
178     }
179    
180     /* } else {
181 jfenwick 1851 Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getAMG: Break-down in AMG decomposition: non-regular main diagonal block.");
182 artak 1890 }*/
183 artak 1881 }
184 jfenwick 1851 }
185     } /* end parallel region */
186    
187     if( Paso_noError()) {
188     /* if there are no nodes in the coarse level there is no more work to do */
189     out->n_C=n-out->n_F;
190 artak 1890 /*if (out->n_C>11) {*/
191     if (level>0) {
192 jfenwick 1851 out->rows_in_C=MEMALLOC(out->n_C,index_t);
193     out->mask_C=MEMALLOC(n,index_t);
194     if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
195     /* creates an index for C from mask */
196     #pragma omp parallel for private(i) schedule(static)
197     for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
198     Paso_Util_cumsum(n,counter);
199     #pragma omp parallel
200     {
201     #pragma omp for private(i) schedule(static)
202     for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
203     #pragma omp for private(i) schedule(static)
204     for (i = 0; i < n; ++i) {
205     if (! mis_marker[i]) {
206     out->rows_in_C[counter[i]]=i;
207     out->mask_C[i]=counter[i];
208     } else {
209     out->mask_C[i]=-1;
210     }
211     }
212     } /* end parallel region */
213     /* get A_CF block: */
214     out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
215     if (Paso_noError()) {
216     /* get A_FC block: */
217     out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
218 artak 1881 /* get A_CC block: */
219 jfenwick 1851 if (Paso_noError()) {
220     schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
221 artak 1881
222     /*find the pattern of the schur complement with fill in*/
223     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);
224    
225     /* copy values over*/
226     #pragma omp for private(i,iPtr,iPtr_s,j,j0) schedule(static)
227     for (i = 0; i < schur_withFillIn->numRows; ++i) {
228     for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
229     j=schur_withFillIn->pattern->index[iPtr];
230     schur_withFillIn->val[iPtr]=0.;
231     for (iPtr_s=schur->pattern->ptr[i];iPtr_s<schur->pattern->ptr[i + 1]; ++iPtr_s){
232     j0=schur->pattern->index[iPtr_s];
233     if (j==j0) {
234     schur_withFillIn->val[iPtr]=schur->val[iPtr_s];
235     break;
236     }
237     }
238     }
239     }
240    
241     /* for (i = 0; i < schur_withFillIn->numRows; ++i) {
242     for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
243     j=schur_withFillIn->pattern->index[iPtr];
244     fprintf(stderr,"A_CC[%d,%d]=%.2f ",i,j,schur_withFillIn->val[iPtr]);
245     }
246     fprintf(stderr,"\n");
247     }*/
248 jfenwick 1851 time0=Paso_timer()-time0;
249     if (Paso_noError()) {
250     time1=Paso_timer();
251     /* update A_CC block to get Schur complement and then apply AMG to it */
252 artak 1881 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
253 jfenwick 1851 time1=Paso_timer()-time1;
254 artak 1890 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level-1);
255 artak 1881
256 jfenwick 1851 Paso_SparseMatrix_free(schur);
257 artak 1890 /* Paso_SparseMatrix_free(schur_withFillIn);*/
258 jfenwick 1851 }
259     /* allocate work arrays for AMG application */
260     if (Paso_noError()) {
261     out->x_F=MEMALLOC(n_block*out->n_F,double);
262     out->b_F=MEMALLOC(n_block*out->n_F,double);
263     out->x_C=MEMALLOC(n_block*out->n_C,double);
264     out->b_C=MEMALLOC(n_block*out->n_C,double);
265 artak 1890
266 jfenwick 1851 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
267     #pragma omp parallel
268     {
269     #pragma omp for private(i,k) schedule(static)
270     for (i = 0; i < out->n_F; ++i) {
271     for (k=0; k<n_block;++k) {
272     out->x_F[i*n_block+k]=0.;
273     out->b_F[i*n_block+k]=0.;
274     }
275     }
276     #pragma omp for private(i,k) schedule(static)
277     for (i = 0; i < out->n_C; ++i) {
278     for (k=0; k<n_block;++k) {
279     out->x_C[i*n_block+k]=0.;
280     out->b_C[i*n_block+k]=0.;
281     }
282     }
283     } /* end parallel region */
284     }
285     }
286     }
287     }
288     }
289     }
290     }
291     }
292     }
293     }
294     TMPMEMFREE(mis_marker);
295     TMPMEMFREE(counter);
296 artak 1881 TMPMEMFREE(rs);
297 jfenwick 1851 if (Paso_noError()) {
298     if (verbose) {
299     printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F);
300     if (out->n_C>0) {
301     printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1);
302     } else {
303     printf("timing: AMG: MIS: %e\n",time2);
304     }
305     }
306     return out;
307     } else {
308     Paso_Solver_AMG_free(out);
309     return NULL;
310     }
311     }
312    
313     /**************************************************************/
314    
315     /* apply AMG precondition b-> x
316    
317     in fact it solves
318    
319     [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] [ x_F ] = [b_F]
320     [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
321    
322     in the form
323    
324     b->[b_F,b_C]
325     x_F=invA_FF*b_F
326     b_C=b_C-A_CF*x_F
327     x_C=AMG(b_C)
328     b_F=b_F-A_FC*x_C
329     x_F=invA_FF*b_F
330     x<-[x_F,x_C]
331    
332     should be called within a parallel region
333     barrier synconization should be performed to make sure that the input vector available
334    
335     */
336    
337     void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
338     dim_t i,k;
339     dim_t n_block=amg->n_block;
340 artak 1890 double *r=MEMALLOC(amg->n,double);
341     double *x0=MEMALLOC(amg->n,double);
342 artak 1892 Paso_Solver_GS* GS=NULL;
343 artak 1890
344     if (amg->level==0) {
345 jfenwick 1851 /* x=invA_FF*b */
346     Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,x,b);
347     } else {
348 artak 1890 fprintf(stdout,"LEVEL %d \n",amg->level);
349 artak 1881 /* presmoothing on (Shure, x, b, r) */
350 artak 1890 /****************/
351     GS=Paso_Solver_getGS(amg->A,-1);
352     Paso_Solver_solveGS(GS,x,b);
353    
354     #pragma omp parallel for private(i) schedule(static)
355     for (i=0;i<amg->n;++i) r[i]=b[i];
356    
357     /*r=b-Ax*/
358     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
359     /****************/
360 jfenwick 1851 /* b->[b_F,b_C] */
361     if (n_block==1) {
362     #pragma omp parallel for private(i) schedule(static)
363 artak 1890 for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]]; /* was b istead of r */
364     #pragma omp parallel for private(i) schedule(static)
365     for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
366 jfenwick 1851 } else {
367     #pragma omp parallel for private(i,k) schedule(static)
368     for (i=0;i<amg->n_F;++i)
369 artak 1890 for (k=0;k<n_block;k++) amg->b_F[amg->n_block*i+k]=r[n_block*amg->rows_in_F[i]+k];
370 jfenwick 1851 #pragma omp parallel for private(i,k) schedule(static)
371     for (i=0;i<amg->n_C;++i)
372 artak 1890 for (k=0;k<n_block;k++) amg->b_C[amg->n_block*i+k]=r[n_block*amg->rows_in_C[i]+k];
373 jfenwick 1851 }
374 artak 1890
375     /****************/
376 artak 1892 /* Corsening part */
377 artak 1890
378     /* r_F=A_FF^-1*r_F */
379     /*Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->b_F,amg->b_F);
380     fprintf(stderr,"2\n");
381     */
382     /* r_C=r_C-A_CF*r_F */
383     /*Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->b_F,1.,amg->b_C);
384     fprintf(stderr,"3\n");
385     */
386     /****************/
387    
388    
389 jfenwick 1851 /* x_F=invA_FF*b_F */
390     Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
391 artak 1890
392 jfenwick 1851 /* b_C=b_C-A_CF*x_F */
393     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
394 artak 1890
395 jfenwick 1851 /* x_C=AMG(b_C) */
396     Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
397 artak 1890
398 jfenwick 1851 /* b_F=b_F-A_FC*x_C */
399     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
400     /* x_F=invA_FF*b_F */
401     Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
402     /* x<-[x_F,x_C] */
403 artak 1890
404     /* */
405    
406    
407 jfenwick 1851 if (n_block==1) {
408     #pragma omp parallel for private(i) schedule(static)
409     for (i=0;i<amg->n;++i) {
410     if (amg->mask_C[i]>-1) {
411 artak 1890 x[i]+=amg->x_C[amg->mask_C[i]];
412 jfenwick 1851 } else {
413 artak 1890 x[i]+=amg->x_F[amg->mask_F[i]];
414 jfenwick 1851 }
415     }
416     } else {
417     #pragma omp parallel for private(i,k) schedule(static)
418     for (i=0;i<amg->n;++i) {
419     if (amg->mask_C[i]>-1) {
420 artak 1890 for (k=0;k<n_block;k++) x[n_block*i+k]+=amg->x_C[n_block*amg->mask_C[i]+k];
421 jfenwick 1851 } else {
422 artak 1890 for (k=0;k<n_block;k++) x[n_block*i+k]+=amg->x_F[n_block*amg->mask_F[i]+k];
423 jfenwick 1851 }
424     }
425     }
426     /* all done */
427 artak 1890 /*post smoothing*/
428     /* Paso_Solver_solveGS(GS,x0,b);
429     if (n_block==1) {
430     #pragma omp parallel for private(i) schedule(static)
431     for (i=0;i<amg->n;++i) x[i]+=x0[i];
432     } else {
433     #pragma omp parallel for private(i,k) schedule(static)
434     for (i=0;i<amg->n;++i) {
435     for (k=0;k<n_block;k++) x[n_block*i+k]+=x0[n_block*i+k];
436     }
437     }
438     */
439     Paso_Solver_GS_free(GS);
440 jfenwick 1851 }
441 artak 1890 MEMFREE(x0);
442     MEMFREE(r);
443 jfenwick 1851 return;
444     }
445    
446     /*
447     * $Log$
448     *
449     */

  ViewVC Help
Powered by ViewVC 1.1.26