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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.26