/[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 2113 - (hide annotations)
Mon Dec 1 06:27:14 2008 UTC (11 years, 2 months ago) by artak
File MIME type: text/plain
File size: 16031 byte(s)
Temporary solution for nusty rows: a condition added 
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 1890 Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose,dim_t level) {
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 1954 /*Paso_Pattern_RS(A_p,mis_marker,0.25);*/
116     /*Paso_Pattern_Aggregiation(A_p,mis_marker,0.5);*/
117     Paso_Pattern_coup(A_p,mis_marker,0.05);
118 jfenwick 1851 if (Paso_noError()) {
119     #pragma omp parallel for private(i) schedule(static)
120     for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
121     out->n=n;
122     out->n_block=n_block;
123     out->n_F=Paso_Util_cumsum(n,counter);
124     out->mask_F=MEMALLOC(n,index_t);
125     out->rows_in_F=MEMALLOC(out->n_F,index_t);
126     out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
127     out->A_FF_pivot=NULL; /* later use for block size>3 */
128     if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
129     /* creates an index for F from mask */
130 artak 2107 #pragma omp parallel for private(i) schedule(static)
131 jfenwick 1851 for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
132 artak 2107 #pragma omp parallel for private(i) schedule(static)
133 jfenwick 1851 for (i = 0; i < n; ++i) {
134     if (mis_marker[i]) {
135     out->rows_in_F[counter[i]]=i;
136     out->mask_F[i]=counter[i];
137     } else {
138     out->mask_F[i]=-1;
139     }
140     }
141 artak 2107
142     /* Compute row-sum for getting rs(A_FF)^-1*/
143     #pragma omp parallel for private(i,iPtr,j) schedule(static)
144 artak 1881 for (i = 0; i < out->n_F; ++i) {
145 artak 2107 S=0;
146 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) {
147 artak 1975 j=A_p->pattern->index[iPtr];
148     if (mis_marker[j])
149 artak 2107 S+=A_p->val[iPtr];
150 artak 1881 }
151 artak 2113 /* This check is to make sure we dont get some nusty rows which were not removed durring coarsing process.*/
152     /* TODO: we have to mechanism that this does not happend at all, and get rid of this 'If'. */
153     if (ABS(S)>1.e-10)
154 artak 2107 out->inv_A_FF[i]=1./S;
155 artak 2113 else
156     out->inv_A_FF[i]=1.;
157 artak 1881 }
158 jfenwick 1851
159     if( Paso_noError()) {
160     /* if there are no nodes in the coarse level there is no more work to do */
161     out->n_C=n-out->n_F;
162 artak 2112 if (level<3) {
163     /*if (out->n_F>500) {*/
164 jfenwick 1851 out->rows_in_C=MEMALLOC(out->n_C,index_t);
165     out->mask_C=MEMALLOC(n,index_t);
166     if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
167     /* creates an index for C from mask */
168     #pragma omp parallel for private(i) schedule(static)
169     for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
170     Paso_Util_cumsum(n,counter);
171 artak 2107 #pragma omp parallel for private(i) schedule(static)
172 jfenwick 1851 for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
173 artak 2107 #pragma omp parallel for private(i) schedule(static)
174 jfenwick 1851 for (i = 0; i < n; ++i) {
175     if (! mis_marker[i]) {
176     out->rows_in_C[counter[i]]=i;
177     out->mask_C[i]=counter[i];
178     } else {
179     out->mask_C[i]=-1;
180     }
181     }
182     /* get A_CF block: */
183     out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
184     if (Paso_noError()) {
185     /* get A_FC block: */
186     out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
187 artak 1881 /* get A_CC block: */
188 jfenwick 1851 if (Paso_noError()) {
189     schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
190 artak 2107 /*find the pattern of the schur complement with fill in*/
191 artak 2112
192 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);
193 artak 2107
194     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));
195    
196 artak 1881 /* copy values over*/
197 artak 2107 #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
198 artak 1881 for (i = 0; i < schur_withFillIn->numRows; ++i) {
199     for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
200     j=schur_withFillIn->pattern->index[iPtr];
201 artak 2107 iPtr_s=schur->pattern->ptr[i];
202 artak 1881 schur_withFillIn->val[iPtr]=0.;
203 artak 2107 index=&(schur->pattern->index[iPtr_s]);
204     where_p=(index_t*)bsearch(&j,
205     index,
206     schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
207     sizeof(index_t),
208     Paso_comparIndex);
209     if (where_p!=NULL) {
210     schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
211 artak 1881 }
212     }
213     }
214 artak 2107
215 jfenwick 1851 if (Paso_noError()) {
216 artak 1881 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
217 artak 1939 out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level+1);
218 artak 2107 Paso_SparseMatrix_free(schur);
219 jfenwick 1851 }
220     /* allocate work arrays for AMG application */
221     if (Paso_noError()) {
222     out->x_F=MEMALLOC(n_block*out->n_F,double);
223     out->b_F=MEMALLOC(n_block*out->n_F,double);
224     out->x_C=MEMALLOC(n_block*out->n_C,double);
225     out->b_C=MEMALLOC(n_block*out->n_C,double);
226 artak 1890
227 jfenwick 1851 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
228 artak 2107 #pragma omp parallel for private(i,k) schedule(static)
229     for (i = 0; i < out->n_F; ++i) {
230 jfenwick 1851 for (k=0; k<n_block;++k) {
231     out->x_F[i*n_block+k]=0.;
232     out->b_F[i*n_block+k]=0.;
233     }
234     }
235 artak 2107 #pragma omp parallel for private(i,k) schedule(static)
236 jfenwick 1851 for (i = 0; i < out->n_C; ++i) {
237     for (k=0; k<n_block;++k) {
238     out->x_C[i*n_block+k]=0.;
239     out->b_C[i*n_block+k]=0.;
240     }
241     }
242     }
243     }
244     }
245     }
246     }
247     }
248     }
249     }
250     }
251     }
252     TMPMEMFREE(mis_marker);
253     TMPMEMFREE(counter);
254     if (Paso_noError()) {
255     if (verbose) {
256     printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F);
257 artak 2112 if (level<3) {
258     /*if (out->n_F<500) {*/
259 jfenwick 1851 printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1);
260     } else {
261     printf("timing: AMG: MIS: %e\n",time2);
262     }
263     }
264     return out;
265     } else {
266     Paso_Solver_AMG_free(out);
267     return NULL;
268     }
269     }
270    
271     /**************************************************************/
272    
273     /* apply AMG precondition b-> x
274    
275     in fact it solves
276    
277     [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] [ x_F ] = [b_F]
278     [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
279    
280     in the form
281    
282     b->[b_F,b_C]
283     x_F=invA_FF*b_F
284     b_C=b_C-A_CF*x_F
285     x_C=AMG(b_C)
286     b_F=b_F-A_FC*x_C
287     x_F=invA_FF*b_F
288     x<-[x_F,x_C]
289    
290     should be called within a parallel region
291     barrier synconization should be performed to make sure that the input vector available
292    
293     */
294    
295     void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) {
296 artak 2107 dim_t i;
297 artak 1890 double *r=MEMALLOC(amg->n,double);
298 artak 1939 /*Paso_Solver_GS* GS=NULL;*/
299 artak 1954 double *x0=MEMALLOC(amg->n,double);
300 artak 2107 double time0=0;
301    
302 artak 2112 if (amg->level==3) {
303     /*if (amg->n_F<=500) {*/
304 artak 2107 time0=Paso_timer();
305    
306 artak 2112 Paso_Solver_solveJacobi(amg->GS,x,b);
307 artak 2107
308     /* #pragma omp parallel for private(i) schedule(static)
309     for (i=0;i<amg->n;++i) r[i]=b[i];
310     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
311     Paso_Solver_solveGS(amg->GS,x0,r);
312     #pragma omp parallel for private(i) schedule(static)
313     for (i=0;i<amg->n;++i) {
314     x[i]+=x0[i];
315     }
316     */
317 artak 2112 /*Paso_UMFPACK1(amg->A,x,b,0);*/
318 artak 1939
319 artak 2107 time0=Paso_timer()-time0;
320     /*fprintf(stderr,"timing: DIRECT SOLVER: %e/\n",time0);*/
321 jfenwick 1851 } else {
322 artak 1954
323 artak 1902 /* presmoothing */
324 artak 2112 Paso_Solver_solveJacobi(amg->GS,x,b);
325 artak 2107 /*
326     #pragma omp parallel for private(i) schedule(static)
327     for (i=0;i<amg->n;++i) r[i]=b[i];
328    
329     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
330     Paso_Solver_solveGS(amg->GS,x0,r);
331    
332     #pragma omp parallel for private(i) schedule(static)
333     for (i=0;i<amg->n;++i) {
334     x[i]+=x0[i];
335     }
336     */
337 artak 1931 /* end of presmoothing */
338 artak 1939
339 artak 1890 #pragma omp parallel for private(i) schedule(static)
340     for (i=0;i<amg->n;++i) r[i]=b[i];
341    
342     /*r=b-Ax*/
343     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
344 artak 1954
345 jfenwick 1851 /* b->[b_F,b_C] */
346 artak 1954 #pragma omp parallel for private(i) schedule(static)
347     for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]];
348    
349     #pragma omp parallel for private(i) schedule(static)
350     for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]];
351 artak 1890
352 jfenwick 1851 /* x_F=invA_FF*b_F */
353 artak 1954 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
354 artak 1890
355 jfenwick 1851 /* b_C=b_C-A_CF*x_F */
356     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C);
357 artak 1890
358 jfenwick 1851 /* x_C=AMG(b_C) */
359     Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C);
360 artak 1890
361 jfenwick 1851 /* b_F=b_F-A_FC*x_C */
362     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F);
363     /* x_F=invA_FF*b_F */
364 artak 1954 Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F);
365 jfenwick 1851 /* x<-[x_F,x_C] */
366 artak 2113
367 artak 1954 #pragma omp parallel for private(i) schedule(static)
368     for (i=0;i<amg->n;++i) {
369     if (amg->mask_C[i]>-1) {
370 artak 2107 x[i]=amg->x_C[amg->mask_C[i]];
371 artak 1954 } else {
372 artak 2107 x[i]=amg->x_F[amg->mask_F[i]];
373 artak 1954 }
374 jfenwick 1851 }
375 artak 1954
376 artak 1931 /*postsmoothing*/
377 artak 2107
378 artak 1954 #pragma omp parallel for private(i) schedule(static)
379     for (i=0;i<amg->n;++i) r[i]=b[i];
380 artak 2107
381 artak 1954 /*r=b-Ax*/
382     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
383 artak 2112 Paso_Solver_solveJacobi(amg->GS,x0,r);
384 artak 2107
385     #pragma omp parallel for private(i) schedule(static)
386     for (i=0;i<amg->n;++i) x[i]+=x0[i];
387    
388    
389     /*#pragma omp parallel for private(i) schedule(static)
390     for (i=0;i<amg->n;++i) r[i]=b[i];
391    
392     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r);
393 artak 1954 Paso_Solver_solveGS(amg->GS,x0,r);
394    
395     #pragma omp parallel for private(i) schedule(static)
396     for (i=0;i<amg->n;++i) {
397     x[i]+=x0[i];
398     }
399 artak 2107 */
400 artak 1931 /*end of postsmoothing*/
401    
402 jfenwick 1851 }
403 artak 1890 MEMFREE(r);
404 artak 1954 MEMFREE(x0);
405 artak 1931 return;
406 jfenwick 1851 }
407    
408     /*
409     * $Log$
410     *
411     */

  ViewVC Help
Powered by ViewVC 1.1.26