/[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 2307 - (hide annotations)
Mon Mar 16 00:38:09 2009 UTC (10 years, 11 months ago) by artak
File MIME type: text/plain
File size: 15647 byte(s)
some clean up. RS and Aggregiation are removed.
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 artak 2112 Paso_Solver_Jacobi_free(in->GS);
39 jfenwick 1851 MEMFREE(in->inv_A_FF);
40     MEMFREE(in->A_FF_pivot);
41     Paso_SparseMatrix_free(in->A_FC);
42     Paso_SparseMatrix_free(in->A_CF);
43     MEMFREE(in->rows_in_F);
44     MEMFREE(in->rows_in_C);
45     MEMFREE(in->mask_F);
46     MEMFREE(in->mask_C);
47     MEMFREE(in->x_F);
48     MEMFREE(in->b_F);
49     MEMFREE(in->x_C);
50     MEMFREE(in->b_C);
51     MEMFREE(in);
52     }
53     }
54    
55     /**************************************************************/
56    
57     /* constructs the block-block factorization of
58    
59     [ A_FF A_FC ]
60     A_p=
61     [ A_CF A_FF ]
62    
63     to
64    
65     [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ]
66     [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ]
67    
68    
69     where S=A_FF-ACF*invA_FF*A_FC within the shape of S
70    
71     then AMG is applied to S again until S becomes empty
72    
73     */
74 artak 2178 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose,dim_t level, double couplingParam) {
75 jfenwick 1851 Paso_Solver_AMG* out=NULL;
76     dim_t n=A_p->numRows;
77     dim_t n_block=A_p->row_block_size;
78     index_t* mis_marker=NULL;
79 artak 1881 index_t* counter=NULL;
80     index_t iPtr,*index, *where_p, iPtr_s;
81 artak 2107 dim_t i,k,j;
82 jfenwick 1851 Paso_SparseMatrix * schur=NULL;
83 artak 1881 Paso_SparseMatrix * schur_withFillIn=NULL;
84 artak 2107 double time0=0,time1=0,time2=0,S=0;
85     /*Paso_Pattern* test;*/
86 artak 1890
87 jfenwick 1851 /* identify independend set of rows/columns */
88     mis_marker=TMPMEMALLOC(n,index_t);
89     counter=TMPMEMALLOC(n,index_t);
90     out=MEMALLOC(1,Paso_Solver_AMG);
91     out->AMG_of_Schur=NULL;
92     out->inv_A_FF=NULL;
93     out->A_FF_pivot=NULL;
94     out->A_FC=NULL;
95     out->A_CF=NULL;
96     out->rows_in_F=NULL;
97     out->rows_in_C=NULL;
98     out->mask_F=NULL;
99     out->mask_C=NULL;
100     out->x_F=NULL;
101     out->b_F=NULL;
102     out->x_C=NULL;
103     out->b_C=NULL;
104 artak 2107 out->GS=NULL;
105 artak 1890 out->A=Paso_SparseMatrix_getReference(A_p);
106 artak 2112 /*out->GS=Paso_Solver_getGS(A_p,verbose);*/
107     out->GS=Paso_Solver_getJacobi(A_p);
108 artak 2107 /*out->GS->sweeps=2;*/
109 artak 1890 out->level=level;
110 artak 2107
111 jfenwick 1851 if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) {
112     /* identify independend set of rows/columns */
113     #pragma omp parallel for private(i) schedule(static)
114     for (i=0;i<n;++i) mis_marker[i]=-1;
115 artak 2178 Paso_Pattern_coup(A_p,mis_marker,couplingParam);
116 jfenwick 1851 if (Paso_noError()) {
117     #pragma omp parallel for private(i) schedule(static)
118     for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
119     out->n=n;
120     out->n_block=n_block;
121     out->n_F=Paso_Util_cumsum(n,counter);
122     out->mask_F=MEMALLOC(n,index_t);
123     out->rows_in_F=MEMALLOC(out->n_F,index_t);
124     out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
125     out->A_FF_pivot=NULL; /* later use for block size>3 */
126     if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
127     /* creates an index for F from mask */
128 artak 2107 #pragma omp parallel for private(i) schedule(static)
129 jfenwick 1851 for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
130 artak 2107 #pragma omp parallel for private(i) schedule(static)
131 jfenwick 1851 for (i = 0; i < n; ++i) {
132     if (mis_marker[i]) {
133     out->rows_in_F[counter[i]]=i;
134     out->mask_F[i]=counter[i];
135     } else {
136     out->mask_F[i]=-1;
137     }
138     }
139 artak 2107
140     /* Compute row-sum for getting rs(A_FF)^-1*/
141 artak 2131 #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
142 artak 1881 for (i = 0; i < out->n_F; ++i) {
143 artak 2107 S=0;
144 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) {
145 artak 1975 j=A_p->pattern->index[iPtr];
146     if (mis_marker[j])
147 artak 2107 S+=A_p->val[iPtr];
148 artak 1881 }
149 artak 2107 out->inv_A_FF[i]=1./S;
150 artak 1881 }
151 jfenwick 1851
152     if( Paso_noError()) {
153     /* if there are no nodes in the coarse level there is no more work to do */
154     out->n_C=n-out->n_F;
155 artak 2280 if (level<3) {
156     /*if (out->n_F>500) {*/
157 jfenwick 1851 out->rows_in_C=MEMALLOC(out->n_C,index_t);
158     out->mask_C=MEMALLOC(n,index_t);
159     if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
160     /* creates an index for C from mask */
161     #pragma omp parallel for private(i) schedule(static)
162     for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
163     Paso_Util_cumsum(n,counter);
164 artak 2107 #pragma omp parallel for private(i) schedule(static)
165 jfenwick 1851 for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
166 artak 2107 #pragma omp parallel for private(i) schedule(static)
167 jfenwick 1851 for (i = 0; i < n; ++i) {
168     if (! mis_marker[i]) {
169     out->rows_in_C[counter[i]]=i;
170     out->mask_C[i]=counter[i];
171     } else {
172     out->mask_C[i]=-1;
173     }
174     }
175     /* get A_CF block: */
176     out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
177     if (Paso_noError()) {
178     /* get A_FC block: */
179     out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
180 artak 1881 /* get A_CC block: */
181 jfenwick 1851 if (Paso_noError()) {
182     schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
183 artak 2107 /*find the pattern of the schur complement with fill in*/
184 artak 2112
185 artak 1881 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);
186 artak 2107
187 artak 2280 /*fprintf(stderr,"Sparsity of Schure: %dx%d LEN %d Percentage %f\n",schur_withFillIn->pattern->numOutput,schur_withFillIn->pattern->numInput,schur_withFillIn->len,schur_withFillIn->len/(1.*schur_withFillIn->pattern->numOutput*schur_withFillIn->pattern->numInput));*/
188 artak 2107
189 artak 1881 /* copy values over*/
190 artak 2107 #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
191 artak 1881 for (i = 0; i < schur_withFillIn->numRows; ++i) {
192     for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
193     j=schur_withFillIn->pattern->index[iPtr];
194 artak 2107 iPtr_s=schur->pattern->ptr[i];
195 artak 1881 schur_withFillIn->val[iPtr]=0.;
196 artak 2107 index=&(schur->pattern->index[iPtr_s]);
197     where_p=(index_t*)bsearch(&j,
198     index,
199     schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
200     sizeof(index_t),
201     Paso_comparIndex);
202     if (where_p!=NULL) {
203     schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
204 artak 1881 }
205     }
206     }
207 artak 2107
208 jfenwick 1851 if (Paso_noError()) {
209 artak 1881 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
210 artak 2178 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level+1,couplingParam);
211 artak 2107 Paso_SparseMatrix_free(schur);
212 jfenwick 1851 }
213     /* allocate work arrays for AMG application */
214     if (Paso_noError()) {
215     out->x_F=MEMALLOC(n_block*out->n_F,double);
216     out->b_F=MEMALLOC(n_block*out->n_F,double);
217     out->x_C=MEMALLOC(n_block*out->n_C,double);
218     out->b_C=MEMALLOC(n_block*out->n_C,double);
219 artak 1890
220 jfenwick 1851 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
221 artak 2107 #pragma omp parallel for private(i,k) schedule(static)
222     for (i = 0; i < out->n_F; ++i) {
223 jfenwick 1851 for (k=0; k<n_block;++k) {
224     out->x_F[i*n_block+k]=0.;
225     out->b_F[i*n_block+k]=0.;
226     }
227     }
228 artak 2107 #pragma omp parallel for private(i,k) schedule(static)
229 jfenwick 1851 for (i = 0; i < out->n_C; ++i) {
230     for (k=0; k<n_block;++k) {
231     out->x_C[i*n_block+k]=0.;
232     out->b_C[i*n_block+k]=0.;
233     }
234     }
235     }
236     }
237     }
238     }
239     }
240     }
241     }
242     }
243     }
244     }
245     TMPMEMFREE(mis_marker);
246     TMPMEMFREE(counter);
247     if (Paso_noError()) {
248     if (verbose) {
249     printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F);
250 artak 2280 if (level<3) {
251     /* if (out->n_F<500) {*/
252 jfenwick 1851 printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1);
253     } else {
254     printf("timing: AMG: MIS: %e\n",time2);
255     }
256     }
257     return out;
258     } else {
259     Paso_Solver_AMG_free(out);
260     return NULL;
261     }
262     }
263    
264     /**************************************************************/
265    
266     /* apply AMG precondition b-> x
267    
268     in fact it solves
269    
270     [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] [ x_F ] = [b_F]
271     [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
272    
273     in the form
274    
275     b->[b_F,b_C]
276     x_F=invA_FF*b_F
277     b_C=b_C-A_CF*x_F
278     x_C=AMG(b_C)
279     b_F=b_F-A_FC*x_C
280     x_F=invA_FF*b_F
281     x<-[x_F,x_C]
282    
283     should be called within a parallel region
284     barrier synconization should be performed to make sure that the input vector available
285    
286     */
287    
288     void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
289 artak 2107 dim_t i;
290 artak 1890 double *r=MEMALLOC(amg->n,double);
291 artak 1939 /*Paso_Solver_GS* GS=NULL;*/
292 artak 1954 double *x0=MEMALLOC(amg->n,double);
293 artak 2107 double time0=0;
294    
295 artak 2280 if (amg->level==3) {
296     /*if (amg->n_F<=500) {*/
297 artak 2107 time0=Paso_timer();
298    
299 artak 2112 Paso_Solver_solveJacobi(amg->GS,x,b);
300 artak 2107
301     /* #pragma omp parallel for private(i) schedule(static)
302     for (i=0;i<amg->n;++i) r[i]=b[i];
303     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
304     Paso_Solver_solveGS(amg->GS,x0,r);
305     #pragma omp parallel for private(i) schedule(static)
306     for (i=0;i<amg->n;++i) {
307     x[i]+=x0[i];
308     }
309     */
310 artak 2112 /*Paso_UMFPACK1(amg->A,x,b,0);*/
311 artak 1939
312 artak 2107 time0=Paso_timer()-time0;
313     /*fprintf(stderr,"timing: DIRECT SOLVER: %e/\n",time0);*/
314 jfenwick 1851 } else {
315 artak 1954
316 artak 1902 /* presmoothing */
317 artak 2112 Paso_Solver_solveJacobi(amg->GS,x,b);
318 artak 2107 /*
319     #pragma omp parallel for private(i) schedule(static)
320     for (i=0;i<amg->n;++i) r[i]=b[i];
321    
322     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
323     Paso_Solver_solveGS(amg->GS,x0,r);
324    
325     #pragma omp parallel for private(i) schedule(static)
326     for (i=0;i<amg->n;++i) {
327     x[i]+=x0[i];
328     }
329     */
330 artak 1931 /* end of presmoothing */
331 artak 1939
332 artak 1890 #pragma omp parallel for private(i) schedule(static)
333     for (i=0;i<amg->n;++i) r[i]=b[i];
334    
335     /*r=b-Ax*/
336     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
337 artak 1954
338 jfenwick 1851 /* b->[b_F,b_C] */
339 artak 1954 #pragma omp parallel for private(i) schedule(static)
340     for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];
341    
342     #pragma omp parallel for private(i) schedule(static)
343     for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
344 artak 1890
345 jfenwick 1851 /* x_F=invA_FF*b_F */
346 artak 1954 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
347 artak 1890
348 jfenwick 1851 /* b_C=b_C-A_CF*x_F */
349     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
350 artak 1890
351 jfenwick 1851 /* x_C=AMG(b_C) */
352     Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
353 artak 1890
354 jfenwick 1851 /* b_F=b_F-A_FC*x_C */
355     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
356     /* x_F=invA_FF*b_F */
357 artak 1954 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
358 jfenwick 1851 /* x<-[x_F,x_C] */
359 artak 2113
360 artak 1954 #pragma omp parallel for private(i) schedule(static)
361     for (i=0;i<amg->n;++i) {
362     if (amg->mask_C[i]>-1) {
363 artak 2107 x[i]=amg->x_C[amg->mask_C[i]];
364 artak 1954 } else {
365 artak 2107 x[i]=amg->x_F[amg->mask_F[i]];
366 artak 1954 }
367 jfenwick 1851 }
368 artak 1954
369 artak 1931 /*postsmoothing*/
370 artak 2107
371 artak 1954 #pragma omp parallel for private(i) schedule(static)
372     for (i=0;i<amg->n;++i) r[i]=b[i];
373 artak 2107
374 artak 1954 /*r=b-Ax*/
375     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
376 artak 2112 Paso_Solver_solveJacobi(amg->GS,x0,r);
377 artak 2107
378     #pragma omp parallel for private(i) schedule(static)
379     for (i=0;i<amg->n;++i) x[i]+=x0[i];
380    
381    
382     /*#pragma omp parallel for private(i) schedule(static)
383     for (i=0;i<amg->n;++i) r[i]=b[i];
384    
385     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
386 artak 1954 Paso_Solver_solveGS(amg->GS,x0,r);
387    
388     #pragma omp parallel for private(i) schedule(static)
389     for (i=0;i<amg->n;++i) {
390     x[i]+=x0[i];
391     }
392 artak 2107 */
393 artak 1931 /*end of postsmoothing*/
394    
395 jfenwick 1851 }
396 artak 1890 MEMFREE(r);
397 artak 1954 MEMFREE(x0);
398 artak 1931 return;
399 jfenwick 1851 }
400    
401     /*
402     * $Log$
403     *
404     */

  ViewVC Help
Powered by ViewVC 1.1.26