/[escript]/branches/arrexp_2137_win/paso/src/Solver_AMG.c
ViewVC logotype

Contents of /branches/arrexp_2137_win/paso/src/Solver_AMG.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2202 - (show annotations)
Fri Jan 9 01:28:32 2009 UTC (10 years, 4 months ago) by jfenwick
File MIME type: text/plain
File size: 17649 byte(s)
Branching the array experiments from version 2137.
The idea is to make the changes required for the c++ changes to compile 
on windows without bringing in the later python changes.


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

  ViewVC Help
Powered by ViewVC 1.1.26