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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2761 - (show annotations)
Thu Nov 19 06:11:07 2009 UTC (11 years, 3 months ago) by artak
File MIME type: text/plain
File size: 19487 byte(s)
AMLI preconditioner (previously called  AMG)
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2009 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: AMLI preconditioner */
18
19 /**************************************************************/
20
21 /* Author: artak@uq.edu.au */
22
23 /**************************************************************/
24
25 #include "Paso.h"
26 #include "Solver.h"
27 #include "Options.h"
28 #include "PasoUtil.h"
29 #include "UMFPACK.h"
30 #include "MKL.h"
31 #include "SystemMatrix.h"
32 #include "Pattern_coupling.h"
33
34 /**************************************************************/
35
36 /* free all memory used by AMLI */
37
38 void Paso_Solver_AMLI_System_free(Paso_Solver_AMLI_System * in) {
39 dim_t i;
40 if (in!=NULL) {
41 for (i=0;i<in->block_size;++i) {
42 Paso_Solver_AMLI_free(in->amliblock[i]);
43 Paso_SparseMatrix_free(in->block[i]);
44 }
45 MEMFREE(in);
46 }
47 }
48
49
50 /* free all memory used by AMLI */
51
52 void Paso_Solver_AMLI_free(Paso_Solver_AMLI * in) {
53 if (in!=NULL) {
54 Paso_Solver_Jacobi_free(in->GS);
55 MEMFREE(in->inv_A_FF);
56 MEMFREE(in->A_FF_pivot);
57 Paso_SparseMatrix_free(in->A_FC);
58 Paso_SparseMatrix_free(in->A_CF);
59 Paso_SparseMatrix_free(in->A);
60 if(in->coarsest_level==TRUE) {
61 #ifdef MKL
62 Paso_MKL_free1(in->AOffset1);
63 Paso_SparseMatrix_free(in->AOffset1);
64 #else
65 #ifdef UMFPACK
66 Paso_UMFPACK1_free((Paso_UMFPACK_Handler*)(in->solver));
67 #endif
68 #endif
69 }
70 MEMFREE(in->rows_in_F);
71 MEMFREE(in->rows_in_C);
72 MEMFREE(in->mask_F);
73 MEMFREE(in->mask_C);
74 MEMFREE(in->x_F);
75 MEMFREE(in->b_F);
76 MEMFREE(in->x_C);
77 MEMFREE(in->b_C);
78 in->solver=NULL;
79 Paso_Solver_AMLI_free(in->AMLI_of_Schur);
80 MEMFREE(in->b_C);
81 MEMFREE(in);
82 }
83 }
84
85 /**************************************************************/
86
87 /* constructs the block-block factorization of
88
89 [ A_FF A_FC ]
90 A_p=
91 [ A_CF A_FF ]
92
93 to
94
95 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ]
96 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ]
97
98
99 where S=A_FF-ACF*invA_FF*A_FC within the shape of S
100
101 then AMLI is applied to S again until S becomes empty
102
103 */
104 Paso_Solver_AMLI* Paso_Solver_getAMLI(Paso_SparseMatrix *A_p,dim_t level,Paso_Options* options) {
105 Paso_Solver_AMLI* out=NULL;
106 Paso_Pattern* temp1=NULL;
107 Paso_Pattern* temp2=NULL;
108 bool_t verbose=options->verbose;
109 dim_t n=A_p->numRows;
110 dim_t n_block=A_p->row_block_size;
111 index_t* mis_marker=NULL;
112 index_t* counter=NULL;
113 index_t iPtr,*index, *where_p, iPtr_s;
114 dim_t i,j;
115 Paso_SparseMatrix * schur=NULL;
116 Paso_SparseMatrix * schur_withFillIn=NULL;
117 double S=0;
118
119
120 /* char filename[8];
121 sprintf(filename,"AMLILevel%d",level);
122
123 Paso_SparseMatrix_saveMM(A_p,filename);
124 */
125
126 /*Make sure we have block sizes 1*/
127 if (A_p->col_block_size>1) {
128 Paso_setError(TYPE_ERROR,"Paso_Solver_getAMLI: AMLI requires column block size 1.");
129 return NULL;
130 }
131 if (A_p->row_block_size>1) {
132 Paso_setError(TYPE_ERROR,"Paso_Solver_getAMLI: AMLI requires row block size 1.");
133 return NULL;
134 }
135 out=MEMALLOC(1,Paso_Solver_AMLI);
136 /* identify independend set of rows/columns */
137 mis_marker=TMPMEMALLOC(n,index_t);
138 counter=TMPMEMALLOC(n,index_t);
139 if ( !( Paso_checkPtr(mis_marker) || Paso_checkPtr(counter) || Paso_checkPtr(out)) ) {
140 out->AMLI_of_Schur=NULL;
141 out->inv_A_FF=NULL;
142 out->A_FF_pivot=NULL;
143 out->A_FC=NULL;
144 out->A_CF=NULL;
145 out->rows_in_F=NULL;
146 out->rows_in_C=NULL;
147 out->mask_F=NULL;
148 out->mask_C=NULL;
149 out->x_F=NULL;
150 out->b_F=NULL;
151 out->x_C=NULL;
152 out->b_C=NULL;
153 out->GS=NULL;
154 out->A=Paso_SparseMatrix_getReference(A_p);
155 out->GS=NULL;
156 out->solver=NULL;
157 /*out->GS=Paso_Solver_getGS(A_p,verbose);*/
158 out->level=level;
159 out->n=n;
160 out->n_F=n+1;
161 out->n_block=n_block;
162
163 if (level==0 || n<=options->min_coarse_matrix_size) {
164 out->coarsest_level=TRUE;
165 #ifdef MKL
166 out->AOffset1=Paso_SparseMatrix_alloc(MATRIX_FORMAT_BLK1 + MATRIX_FORMAT_OFFSET1, out->A->pattern,1,1, FALSE);
167 #pragma omp parallel for private(i) schedule(static)
168 for (i=0;i<out->A->len;++i) {
169 out->AOffset1->val[i]=out->A->val[i];
170 }
171 #else
172 #ifdef UMFPACK
173 #else
174 out->GS=Paso_Solver_getJacobi(A_p);
175 #endif
176 #endif
177 } else {
178 out->coarsest_level=FALSE;
179 out->GS=Paso_Solver_getJacobi(A_p);
180
181 /* identify independend set of rows/columns */
182 #pragma omp parallel for private(i) schedule(static)
183 for (i=0;i<n;++i) mis_marker[i]=-1;
184
185 if (options->coarsening_method == PASO_YAIR_SHAPIRA_COARSENING) {
186 Paso_Pattern_YS(A_p,mis_marker,options->coarsening_threshold);
187 }
188 else if (options->coarsening_method == PASO_RUGE_STUEBEN_COARSENING) {
189 Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
190 }
191 else if (options->coarsening_method == PASO_AGGREGATION_COARSENING) {
192 Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);
193 }
194 else {
195 /*Default coarseneing*/
196 Paso_Pattern_YS(A_p,mis_marker,options->coarsening_threshold);
197 /*Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);*/
198 /*Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);*/
199
200 }
201
202 #pragma omp parallel for private(i) schedule(static)
203 for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
204
205 out->n_F=Paso_Util_cumsum(n,counter);
206
207 if (out->n_F==0) {
208 out->coarsest_level=TRUE;
209 level=0;
210 if (verbose) {
211 /*printf("AMLI coarsening eliminates all unknowns, switching to Jacobi preconditioner.\n");*/
212 printf("AMLI coarsening does not eliminate any of the unknowns, switching to Jacobi preconditioner.\n");
213 }
214 }
215 else if (out->n_F==n) {
216 out->coarsest_level=TRUE;
217 level=0;
218 if (verbose) {
219 /*printf("AMLI coarsening eliminates all unknowns, switching to Jacobi preconditioner.\n");*/
220 printf("AMLI coarsening eliminates all of the unknowns, switching to Jacobi preconditioner.\n");
221
222 }
223 } else {
224
225 if (Paso_noError()) {
226
227 /*#pragma omp parallel for private(i) schedule(static)
228 for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
229 out->n_F=Paso_Util_cumsum(n,counter);
230 */
231
232 /*if(level==3) {
233 printf("##TOTAL: %d, ELIMINATED: %d\n",n,out->n_F);
234 for (i = 0; i < n; ++i) {
235 printf("##%d %d\n",i,mis_marker[i]);
236 }
237 }*/
238
239 out->mask_F=MEMALLOC(n,index_t);
240 out->rows_in_F=MEMALLOC(out->n_F,index_t);
241 out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
242 out->A_FF_pivot=NULL; /* later use for block size>3 */
243 if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
244 /* creates an index for F from mask */
245 #pragma omp parallel for private(i) schedule(static)
246 for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
247 #pragma omp parallel for private(i) schedule(static)
248 for (i = 0; i < n; ++i) {
249 if (mis_marker[i]) {
250 out->rows_in_F[counter[i]]=i;
251 out->mask_F[i]=counter[i];
252 } else {
253 out->mask_F[i]=-1;
254 }
255 }
256
257 /* Compute row-sum for getting rs(A_FF)^-1*/
258 #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
259 for (i = 0; i < out->n_F; ++i) {
260 S=0;
261 /*printf("[%d ]: [%d] -> ",i, out->rows_in_F[i]);*/
262 for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
263 j=A_p->pattern->index[iPtr];
264 /*if (j==out->rows_in_F[i]) printf("diagonal %e",A_p->val[iPtr]);*/
265 if (mis_marker[j])
266 S+=A_p->val[iPtr];
267 }
268 /*printf("-> %e \n",S);*/
269 out->inv_A_FF[i]=1./S;
270 }
271 }
272 }
273
274 /*check whether coarsening process actually makes sense to continue.
275 if coarse matrix at least smaller by 30% then continue, otherwise we stop.*/
276 if ((out->n_F*100/n)<30) {
277 level=1;
278 }
279
280 if ( Paso_noError()) {
281 /* if there are no nodes in the coarse level there is no more work to do */
282 out->n_C=n-out->n_F;
283
284 /*if (out->n_F>500) */
285 out->rows_in_C=MEMALLOC(out->n_C,index_t);
286 out->mask_C=MEMALLOC(n,index_t);
287 if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
288 /* creates an index for C from mask */
289 #pragma omp parallel for private(i) schedule(static)
290 for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
291 Paso_Util_cumsum(n,counter);
292 #pragma omp parallel for private(i) schedule(static)
293 for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
294 #pragma omp parallel for private(i) schedule(static)
295 for (i = 0; i < n; ++i) {
296 if (! mis_marker[i]) {
297 out->rows_in_C[counter[i]]=i;
298 out->mask_C[i]=counter[i];
299 } else {
300 out->mask_C[i]=-1;
301 }
302 }
303 }
304 }
305 if ( Paso_noError()) {
306 /* get A_CF block: */
307 out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
308 /* get A_FC block: */
309 out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
310 /* get A_CC block: */
311 schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
312
313 }
314 if ( Paso_noError()) {
315 /*find the pattern of the schur complement with fill in*/
316 temp1=Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern);
317 temp2=Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, temp1);
318 schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,temp2,1,1, TRUE);
319 Paso_Pattern_free(temp1);
320 Paso_Pattern_free(temp2);
321 }
322 if ( Paso_noError()) {
323 /* copy values over*/
324 #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
325 for (i = 0; i < schur_withFillIn->numRows; ++i) {
326 for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
327 j=schur_withFillIn->pattern->index[iPtr];
328 iPtr_s=schur->pattern->ptr[i];
329 index=&(schur->pattern->index[iPtr_s]);
330 where_p=(index_t*)bsearch(&j,
331 index,
332 schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
333 sizeof(index_t),
334 Paso_comparIndex);
335 if (where_p!=NULL) {
336 schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
337 }
338 }
339 }
340 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
341 out->AMLI_of_Schur=Paso_Solver_getAMLI(schur_withFillIn,level-1,options);
342 }
343 /* allocate work arrays for AMLI application */
344 if (Paso_noError()) {
345 out->x_F=MEMALLOC(n_block*out->n_F,double);
346 out->b_F=MEMALLOC(n_block*out->n_F,double);
347 out->x_C=MEMALLOC(n_block*out->n_C,double);
348 out->b_C=MEMALLOC(n_block*out->n_C,double);
349
350 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
351 #pragma omp parallel for private(i) schedule(static)
352 for (i = 0; i < out->n_F; ++i) {
353 out->x_F[i]=0.;
354 out->b_F[i]=0.;
355 }
356 #pragma omp parallel for private(i) schedule(static)
357 for (i = 0; i < out->n_C; ++i) {
358 out->x_C[i]=0.;
359 out->b_C[i]=0.;
360 }
361 }
362 }
363 Paso_SparseMatrix_free(schur);
364 Paso_SparseMatrix_free(schur_withFillIn);
365 }
366 }
367 }
368 TMPMEMFREE(mis_marker);
369 TMPMEMFREE(counter);
370
371 if (Paso_noError()) {
372 if (verbose && level>0 && !out->coarsest_level) {
373 printf("AMLI: level: %d: %d unknowns eliminated. %d left.\n",level, out->n_F,out->n_C);
374 }
375 return out;
376 } else {
377 Paso_Solver_AMLI_free(out);
378 return NULL;
379 }
380 }
381
382 /**************************************************************/
383
384 /* apply AMLI precondition b-> x
385
386 in fact it solves
387
388 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FC ] [ x_F ] = [b_F]
389 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
390
391 in the form
392
393 b->[b_F,b_C]
394 x_F=invA_FF*b_F
395 b_C=b_C-A_CF*x_F
396 x_C=AMLI(b_C)
397 b_F=b_F-A_FC*x_C
398 x_F=invA_FF*b_F
399 x<-[x_F,x_C]
400
401 should be called within a parallel region
402 barrier synconization should be performed to make sure that the input vector available
403
404 */
405
406 void Paso_Solver_solveAMLI(Paso_Solver_AMLI * amli, double * x, double * b) {
407 dim_t i;
408 double time0=0;
409 double *r=NULL, *x0=NULL;
410 bool_t verbose=0;
411 #ifdef UMFPACK
412 Paso_UMFPACK_Handler * ptr=NULL;
413 #endif
414
415 r=MEMALLOC(amli->n,double);
416 x0=MEMALLOC(amli->n,double);
417
418 if (amli->coarsest_level) {
419
420 time0=Paso_timer();
421 /*If all unknown are eliminated then Jacobi is the best preconditioner*/
422 if (amli->n_F==0 || amli->n_F==amli->n) {
423 Paso_Solver_solveJacobi(amli->GS,x,b);
424 }
425 else {
426 #ifdef MKL
427 Paso_MKL1(amli->AOffset1,x,b,verbose);
428 #else
429 #ifdef UMFPACK
430 ptr=(Paso_UMFPACK_Handler *)(amli->solver);
431 Paso_UMFPACK1(&ptr,amli->A,x,b,verbose);
432 amli->solver=(void*) ptr;
433 #else
434 Paso_Solver_solveJacobi(amli->GS,x,b);
435 #endif
436 #endif
437 }
438 time0=Paso_timer()-time0;
439 if (verbose) fprintf(stderr,"timing: DIRECT SOLVER: %e\n",time0);
440
441 } else {
442 /* presmoothing */
443 time0=Paso_timer();
444 Paso_Solver_solveJacobi(amli->GS,x,b);
445 time0=Paso_timer()-time0;
446 if (verbose) fprintf(stderr,"timing: Presmooting: %e\n",time0);
447 /* end of presmoothing */
448
449
450 time0=Paso_timer();
451 #pragma omp parallel for private(i) schedule(static)
452 for (i=0;i<amli->n;++i) r[i]=b[i];
453
454 /*r=b-Ax*/
455 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A,x,1.,r);
456
457 /* r->[b_F,b_C] */
458 #pragma omp parallel for private(i) schedule(static)
459 for (i=0;i<amli->n_F;++i) amli->b_F[i]=r[amli->rows_in_F[i]];
460
461 #pragma omp parallel for private(i) schedule(static)
462 for (i=0;i<amli->n_C;++i) amli->b_C[i]=r[amli->rows_in_C[i]];
463
464 /* x_F=invA_FF*b_F */
465 Paso_Solver_applyBlockDiagonalMatrix(1,amli->n_F,amli->inv_A_FF,amli->A_FF_pivot,amli->x_F,amli->b_F);
466
467 /* b_C=b_C-A_CF*x_F */
468 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A_CF,amli->x_F,1.,amli->b_C);
469
470 time0=Paso_timer()-time0;
471 if (verbose) fprintf(stderr,"timing: Before next level: %e\n",time0);
472
473 /* x_C=AMLI(b_C) */
474 Paso_Solver_solveAMLI(amli->AMLI_of_Schur,amli->x_C,amli->b_C);
475
476 time0=Paso_timer();
477
478 /* b_F=b_F-A_FC*x_C */
479 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A_FC,amli->x_C,1.,amli->b_F);
480 /* x_F=invA_FF*b_F */
481 Paso_Solver_applyBlockDiagonalMatrix(1,amli->n_F,amli->inv_A_FF,amli->A_FF_pivot,amli->x_F,amli->b_F);
482
483 /* x<-[x_F,x_C] */
484 #pragma omp parallel for private(i) schedule(static)
485 for (i=0;i<amli->n;++i) {
486 if (amli->mask_C[i]>-1) {
487 x[i]=amli->x_C[amli->mask_C[i]];
488 } else {
489 x[i]=amli->x_F[amli->mask_F[i]];
490 }
491 }
492
493 time0=Paso_timer()-time0;
494 if (verbose) fprintf(stderr,"timing: After next level: %e\n",time0);
495
496 /*postsmoothing*/
497 time0=Paso_timer();
498 #pragma omp parallel for private(i) schedule(static)
499 for (i=0;i<amli->n;++i) r[i]=b[i];
500
501 /*r=b-Ax */
502 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A,x,1.,r);
503 Paso_Solver_solveJacobi(amli->GS,x0,r);
504
505 #pragma omp parallel for private(i) schedule(static)
506 for (i=0;i<amli->n;++i) x[i]+=x0[i];
507
508 time0=Paso_timer()-time0;
509 if (verbose) fprintf(stderr,"timing: Postsmoothing: %e\n",time0);
510
511 /*end of postsmoothing*/
512
513 }
514 MEMFREE(r);
515 MEMFREE(x0);
516 return;
517 }

  ViewVC Help
Powered by ViewVC 1.1.26