/[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 1902 - (hide annotations)
Wed Oct 22 03:54:14 2008 UTC (10 years, 10 months ago) by artak
File MIME type: text/plain
File size: 19991 byte(s)
some bugs fixed in coresening methods
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     index_t iPtr,*index, *where_p, iPtr_s;
78     dim_t i,k,j,j0;
79 jfenwick 1851 Paso_SparseMatrix * schur=NULL;
80 artak 1881 Paso_SparseMatrix * schur_withFillIn=NULL;
81 artak 1902 Paso_Pattern* mult=NULL;
82     double A11,A12,A13,A21,A22,A23,A31,A32,A33,D,time0,time1,time2,S;
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     out=MEMALLOC(1,Paso_Solver_AMG);
88     out->AMG_of_Schur=NULL;
89     out->inv_A_FF=NULL;
90     out->A_FF_pivot=NULL;
91     out->A_FC=NULL;
92     out->A_CF=NULL;
93     out->rows_in_F=NULL;
94     out->rows_in_C=NULL;
95     out->mask_F=NULL;
96     out->mask_C=NULL;
97     out->x_F=NULL;
98     out->b_F=NULL;
99     out->x_C=NULL;
100     out->b_C=NULL;
101 artak 1890 out->A=Paso_SparseMatrix_getReference(A_p);
102     out->level=level;
103 artak 1881
104 jfenwick 1851 if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) {
105     /* identify independend set of rows/columns */
106     time0=Paso_timer();
107     #pragma omp parallel for private(i) schedule(static)
108     for (i=0;i<n;++i) mis_marker[i]=-1;
109 artak 1902 /*Paso_Pattern_RS(A_p,mis_marker,0.25);*/
110     Paso_Pattern_coup(A_p,mis_marker);
111 jfenwick 1851 time2=Paso_timer()-time0;
112     if (Paso_noError()) {
113     #pragma omp parallel for private(i) schedule(static)
114     for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
115     out->n=n;
116     out->n_block=n_block;
117     out->n_F=Paso_Util_cumsum(n,counter);
118     out->mask_F=MEMALLOC(n,index_t);
119     out->rows_in_F=MEMALLOC(out->n_F,index_t);
120     out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
121     out->A_FF_pivot=NULL; /* later use for block size>3 */
122     if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
123     #pragma omp parallel
124     {
125     /* creates an index for F from mask */
126     #pragma omp for private(i) schedule(static)
127     for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
128     #pragma omp for private(i) schedule(static)
129     for (i = 0; i < n; ++i) {
130     if (mis_marker[i]) {
131     out->rows_in_F[counter[i]]=i;
132     out->mask_F[i]=counter[i];
133     } else {
134     out->mask_F[i]=-1;
135     }
136     }
137 artak 1881 /* Compute row-sum for getting rs(A_FF)*/
138     #pragma omp for private(i,iPtr) schedule(static)
139     for (i = 0; i < out->n_F; ++i) {
140 artak 1902 out->inv_A_FF[i]=0;
141 artak 1881 for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
142 artak 1902 out->inv_A_FF[i]+=A_p->val[iPtr];
143 artak 1881 }
144     }
145 artak 1902
146     /* fprintf(stderr,"\n MATRIX START\n");
147     j=0;
148     for (i = 0; i < n; ++i) {
149     for (iPtr=A_p->pattern->ptr[i];iPtr<A_p->pattern->ptr[i + 1]; ++iPtr) {
150     fprintf(stderr,"A[%d,%d]=%.2f ",i,A_p->pattern->index[iPtr],A_p->val[iPtr]);
151     }
152     if (mis_marker[i]) {
153     fprintf(stderr,"ROW SUM %.2f\n ",out->inv_A_FF[j]);
154     j++;
155     }
156     else {
157     fprintf(stderr,"\n ");
158     }
159     }
160     fprintf(stderr,"\n MATRIX END\n");
161     */
162 jfenwick 1851 #pragma omp for private(i, where_p,iPtr,A11,A12,A13,A21,A22,A23,A31,A32,A33,D,index) schedule(static)
163     for (i = 0; i < out->n_F; i++) {
164     /* find main diagonal */
165     iPtr=A_p->pattern->ptr[out->rows_in_F[i]];
166     index=&(A_p->pattern->index[iPtr]);
167     where_p=(index_t*)bsearch(&out->rows_in_F[i],
168     index,
169     A_p->pattern->ptr[out->rows_in_F[i] + 1]-A_p->pattern->ptr[out->rows_in_F[i]],
170     sizeof(index_t),
171     Paso_comparIndex);
172     if (where_p==NULL) {
173     Paso_setError(VALUE_ERROR, "Paso_Solver_getAMG: main diagonal element missing.");
174     } else {
175     iPtr+=(index_t)(where_p-index);
176     /* get inverse of A_FF block: */
177 artak 1902 S=out->inv_A_FF[i];
178     if (ABS(A_p->val[iPtr])>0.) {
179     if(ABS(S)>0.)
180     out->inv_A_FF[i]=1./S;
181 artak 1881 } else {
182 artak 1902 Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getAMG: Break-down in AMG decomposition: non-regular main diagonal block.");
183 artak 1890 }
184 artak 1881 }
185 jfenwick 1851 }
186     } /* end parallel region */
187    
188     if( Paso_noError()) {
189     /* if there are no nodes in the coarse level there is no more work to do */
190     out->n_C=n-out->n_F;
191 artak 1890 /*if (out->n_C>11) {*/
192 artak 1902 /*fprintf(stderr,"SPARSITY in level %d \n", level);*/
193 artak 1890 if (level>0) {
194 jfenwick 1851 out->rows_in_C=MEMALLOC(out->n_C,index_t);
195     out->mask_C=MEMALLOC(n,index_t);
196     if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
197     /* creates an index for C from mask */
198     #pragma omp parallel for private(i) schedule(static)
199     for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
200     Paso_Util_cumsum(n,counter);
201     #pragma omp parallel
202     {
203     #pragma omp for private(i) schedule(static)
204     for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
205     #pragma omp for private(i) schedule(static)
206     for (i = 0; i < n; ++i) {
207     if (! mis_marker[i]) {
208     out->rows_in_C[counter[i]]=i;
209     out->mask_C[i]=counter[i];
210     } else {
211     out->mask_C[i]=-1;
212     }
213     }
214     } /* end parallel region */
215     /* get A_CF block: */
216     out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
217     if (Paso_noError()) {
218     /* get A_FC block: */
219     out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
220 artak 1881 /* get A_CC block: */
221 jfenwick 1851 if (Paso_noError()) {
222     schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
223 artak 1881
224     /*find the pattern of the schur complement with fill in*/
225     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);
226 artak 1902 /* fprintf(stderr,"\nSHURE MATRIX START\n");
227 artak 1881
228 artak 1902 for (i = 0; i < schur->numRows; ++i) {
229     for (iPtr=schur->pattern->ptr[i];iPtr<schur->pattern->ptr[i + 1]; ++iPtr) {
230     fprintf(stderr,"A[%d,%d]=%.2f ",i,schur->pattern->index[iPtr],schur->val[iPtr]);
231     }
232     fprintf(stderr,"\n");
233     }
234     fprintf(stderr,"MATRIX END\n");
235    
236     fprintf(stderr,"\nA_CF MATRIX START\n");
237     for (i = 0; i < out->n_C; ++i) {
238     for (iPtr=out->A_CF->pattern->ptr[i];iPtr<out->A_CF->pattern->ptr[i + 1]; ++iPtr) {
239     fprintf(stderr,"A_CF[%d,%d]=%.2f ",i,out->A_CF->pattern->index[iPtr],out->A_CF->val[iPtr]);
240     }
241     fprintf(stderr,"\n");
242     }
243     fprintf(stderr,"MATRIX END\n");
244    
245     fprintf(stderr,"\nA_FC MATRIX START\n");
246     for (i = 0; i < out->n_F; ++i) {
247     for (iPtr=out->A_FC->pattern->ptr[i];iPtr<out->A_FC->pattern->ptr[i + 1]; ++iPtr) {
248     fprintf(stderr,"A_FC[%d,%d]=%.2f ",i,out->A_FC->pattern->index[iPtr],out->A_FC->val[iPtr]);
249     }
250     fprintf(stderr,"\n");
251     }
252     fprintf(stderr,"MATRIX END\n");
253    
254     mult=Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern);
255     fprintf(stderr,"\nPATTERN MATRIX START\n");
256     for (i = 0; i < mult->numOutput; ++i) {
257     for (iPtr=mult->ptr[i];iPtr<mult->ptr[i + 1]; ++iPtr) {
258     fprintf(stderr," P[%d,%d]=X ",i,mult->index[iPtr]);
259     }
260     fprintf(stderr,"\n");
261     }
262     fprintf(stderr,"MATRIX END\n");
263     */
264 artak 1881 /* copy values over*/
265     #pragma omp for private(i,iPtr,iPtr_s,j,j0) schedule(static)
266     for (i = 0; i < schur_withFillIn->numRows; ++i) {
267     for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
268     j=schur_withFillIn->pattern->index[iPtr];
269     schur_withFillIn->val[iPtr]=0.;
270     for (iPtr_s=schur->pattern->ptr[i];iPtr_s<schur->pattern->ptr[i + 1]; ++iPtr_s){
271     j0=schur->pattern->index[iPtr_s];
272     if (j==j0) {
273     schur_withFillIn->val[iPtr]=schur->val[iPtr_s];
274     break;
275     }
276     }
277     }
278     }
279 artak 1902 time0=Paso_timer()-time0;
280 jfenwick 1851 if (Paso_noError()) {
281     time1=Paso_timer();
282     /* update A_CC block to get Schur complement and then apply AMG to it */
283 artak 1881 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
284 artak 1902
285     /*
286     fprintf(stderr,"\n SHURE WITH FILL IN MATRIX START\n");
287     for (i = 0; i < schur_withFillIn->numRows; ++i) {
288     for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
289     fprintf(stderr,"A[%d,%d]=%.2f ",i,schur_withFillIn->pattern->index[iPtr],schur_withFillIn->val[iPtr]);
290     }
291     fprintf(stderr,"\n");
292     }
293     fprintf(stderr,"MATRIX END\n");
294     */
295 jfenwick 1851 time1=Paso_timer()-time1;
296 artak 1890 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level-1);
297 artak 1881
298 jfenwick 1851 Paso_SparseMatrix_free(schur);
299 artak 1890 /* Paso_SparseMatrix_free(schur_withFillIn);*/
300 jfenwick 1851 }
301     /* allocate work arrays for AMG application */
302     if (Paso_noError()) {
303     out->x_F=MEMALLOC(n_block*out->n_F,double);
304     out->b_F=MEMALLOC(n_block*out->n_F,double);
305     out->x_C=MEMALLOC(n_block*out->n_C,double);
306     out->b_C=MEMALLOC(n_block*out->n_C,double);
307 artak 1890
308 jfenwick 1851 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
309     #pragma omp parallel
310     {
311     #pragma omp for private(i,k) schedule(static)
312     for (i = 0; i < out->n_F; ++i) {
313     for (k=0; k<n_block;++k) {
314     out->x_F[i*n_block+k]=0.;
315     out->b_F[i*n_block+k]=0.;
316     }
317     }
318     #pragma omp for private(i,k) schedule(static)
319     for (i = 0; i < out->n_C; ++i) {
320     for (k=0; k<n_block;++k) {
321     out->x_C[i*n_block+k]=0.;
322     out->b_C[i*n_block+k]=0.;
323     }
324     }
325     } /* end parallel region */
326     }
327     }
328     }
329     }
330     }
331     }
332     }
333     }
334     }
335     }
336     TMPMEMFREE(mis_marker);
337     TMPMEMFREE(counter);
338     if (Paso_noError()) {
339     if (verbose) {
340     printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F);
341     if (out->n_C>0) {
342     printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1);
343     } else {
344     printf("timing: AMG: MIS: %e\n",time2);
345     }
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_FF ] [ 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,k;
380     dim_t n_block=amg->n_block;
381 artak 1890 double *r=MEMALLOC(amg->n,double);
382 artak 1892 Paso_Solver_GS* GS=NULL;
383 artak 1890
384     if (amg->level==0) {
385 jfenwick 1851 /* x=invA_FF*b */
386 artak 1902 /*Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,x,b);*/
387     GS=Paso_Solver_getGS(amg->A,-1);
388     Paso_Solver_solveGS(GS,x,b);
389     Paso_Solver_GS_free(GS);
390 jfenwick 1851 } else {
391 artak 1902
392     /* presmoothing */
393 artak 1890 GS=Paso_Solver_getGS(amg->A,-1);
394     Paso_Solver_solveGS(GS,x,b);
395    
396     #pragma omp parallel for private(i) schedule(static)
397     for (i=0;i<amg->n;++i) r[i]=b[i];
398    
399     /*r=b-Ax*/
400     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
401     /****************/
402 jfenwick 1851 /* b->[b_F,b_C] */
403     if (n_block==1) {
404     #pragma omp parallel for private(i) schedule(static)
405 artak 1902 for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];
406 artak 1890 #pragma omp parallel for private(i) schedule(static)
407     for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
408 jfenwick 1851 } else {
409     #pragma omp parallel for private(i,k) schedule(static)
410     for (i=0;i<amg->n_F;++i)
411 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];
412 jfenwick 1851 #pragma omp parallel for private(i,k) schedule(static)
413     for (i=0;i<amg->n_C;++i)
414 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];
415 jfenwick 1851 }
416 artak 1890
417 jfenwick 1851 /* x_F=invA_FF*b_F */
418     Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
419 artak 1890
420 jfenwick 1851 /* b_C=b_C-A_CF*x_F */
421     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
422 artak 1890
423 jfenwick 1851 /* x_C=AMG(b_C) */
424     Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
425 artak 1890
426 jfenwick 1851 /* b_F=b_F-A_FC*x_C */
427     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
428     /* x_F=invA_FF*b_F */
429     Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
430     /* x<-[x_F,x_C] */
431 artak 1890
432 jfenwick 1851 if (n_block==1) {
433     #pragma omp parallel for private(i) schedule(static)
434     for (i=0;i<amg->n;++i) {
435     if (amg->mask_C[i]>-1) {
436 artak 1890 x[i]+=amg->x_C[amg->mask_C[i]];
437 jfenwick 1851 } else {
438 artak 1890 x[i]+=amg->x_F[amg->mask_F[i]];
439 jfenwick 1851 }
440     }
441     } else {
442     #pragma omp parallel for private(i,k) schedule(static)
443     for (i=0;i<amg->n;++i) {
444     if (amg->mask_C[i]>-1) {
445 artak 1890 for (k=0;k<n_block;k++) x[n_block*i+k]+=amg->x_C[n_block*amg->mask_C[i]+k];
446 jfenwick 1851 } else {
447 artak 1890 for (k=0;k<n_block;k++) x[n_block*i+k]+=amg->x_F[n_block*amg->mask_F[i]+k];
448 jfenwick 1851 }
449     }
450     }
451     /* all done */
452 artak 1890 /*post smoothing*/
453 artak 1902 Paso_Solver_solveGS1(GS,x,b);
454    
455 artak 1890 Paso_Solver_GS_free(GS);
456 jfenwick 1851 }
457 artak 1890 MEMFREE(r);
458 jfenwick 1851 return;
459     }
460    
461     /*
462     * $Log$
463     *
464     */

  ViewVC Help
Powered by ViewVC 1.1.26