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

  ViewVC Help
Powered by ViewVC 1.1.26