/[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 3094 - (show annotations)
Fri Aug 13 08:38:06 2010 UTC (8 years, 11 months ago) by gross
File MIME type: text/plain
File size: 21090 byte(s)
The MPI and sequational GAUSS_SEIDEL have been merged.
The couring and main diagonal pointer is now manged by the patternm which means that they are calculated once only even if the preconditioner is deleted.



1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2010 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 out->post_sweeps=options->post_sweeps;
163 out->pre_sweeps=options->pre_sweeps;
164
165 if (level==0 || n<=options->min_coarse_matrix_size) {
166 out->coarsest_level=TRUE;
167 #ifdef MKL
168 out->AOffset1=Paso_SparseMatrix_alloc(MATRIX_FORMAT_BLK1 + MATRIX_FORMAT_OFFSET1, out->A->pattern,1,1, FALSE);
169 #pragma omp parallel for private(i) schedule(static)
170 for (i=0;i<out->A->len;++i) {
171 out->AOffset1->val[i]=out->A->val[i];
172 }
173 #else
174 #ifdef UMFPACK
175 #else
176 out->GS=Paso_Solver_getLocalJacobi(A_p);
177 #endif
178 #endif
179 } else {
180 out->coarsest_level=FALSE;
181 out->GS=Paso_Solver_getLocalJacobi(A_p);
182
183 /* identify independend set of rows/columns */
184 #pragma omp parallel for private(i) schedule(static)
185 for (i=0;i<n;++i) mis_marker[i]=-1;
186
187 if (options->coarsening_method == PASO_YAIR_SHAPIRA_COARSENING) {
188 Paso_Pattern_YS(A_p,mis_marker,options->coarsening_threshold);
189 }
190 else if (options->coarsening_method == PASO_RUGE_STUEBEN_COARSENING) {
191 Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
192 }
193 else if (options->coarsening_method == PASO_AGGREGATION_COARSENING) {
194 Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);
195 }
196 else {
197 /*Default coarseneing*/
198 /*Paso_Pattern_YS(A_p,mis_marker,options->coarsening_threshold);*/
199 Paso_Pattern_RS(A_p,mis_marker,options->coarsening_threshold);
200 /*Paso_Pattern_Aggregiation(A_p,mis_marker,options->coarsening_threshold);*/
201
202 }
203
204 #pragma omp parallel for private(i) schedule(static)
205 for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
206
207 out->n_F=Paso_Util_cumsum(n,counter);
208
209 if (out->n_F==0) {
210 out->coarsest_level=TRUE;
211 level=0;
212 if (verbose) {
213 /*printf("AMLI coarsening eliminates all unknowns, switching to Jacobi preconditioner.\n");*/
214 printf("AMLI coarsening does not eliminate any of the unknowns, switching to Jacobi preconditioner.\n");
215 }
216 }
217 else if (out->n_F==n) {
218 out->coarsest_level=TRUE;
219 level=0;
220 if (verbose) {
221 /*printf("AMLI coarsening eliminates all unknowns, switching to Jacobi preconditioner.\n");*/
222 printf("AMLI coarsening eliminates all of the unknowns, switching to Jacobi preconditioner.\n");
223
224 }
225 } else {
226
227 if (Paso_noError()) {
228
229 /*#pragma omp parallel for private(i) schedule(static)
230 for (i = 0; i < n; ++i) counter[i]=mis_marker[i];
231 out->n_F=Paso_Util_cumsum(n,counter);
232 */
233
234 /*if(level==3) {
235 printf("##TOTAL: %d, ELIMINATED: %d\n",n,out->n_F);
236 for (i = 0; i < n; ++i) {
237 printf("##%d %d\n",i,mis_marker[i]);
238 }
239 }*/
240
241 out->mask_F=MEMALLOC(n,index_t);
242 out->rows_in_F=MEMALLOC(out->n_F,index_t);
243 out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double);
244 out->A_FF_pivot=NULL; /* later use for block size>3 */
245 if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) {
246 /* creates an index for F from mask */
247 #pragma omp parallel for private(i) schedule(static)
248 for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1;
249 #pragma omp parallel for private(i) schedule(static)
250 for (i = 0; i < n; ++i) {
251 if (mis_marker[i]) {
252 out->rows_in_F[counter[i]]=i;
253 out->mask_F[i]=counter[i];
254 } else {
255 out->mask_F[i]=-1;
256 }
257 }
258
259 /* Compute row-sum for getting rs(A_FF)^-1*/
260 #pragma omp parallel for private(i,iPtr,j,S) schedule(static)
261 for (i = 0; i < out->n_F; ++i) {
262 S=0;
263 /*printf("[%d ]: [%d] -> ",i, out->rows_in_F[i]);*/
264 for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) {
265 j=A_p->pattern->index[iPtr];
266 /*if (j==out->rows_in_F[i]) printf("diagonal %e",A_p->val[iPtr]);*/
267 if (mis_marker[j])
268 S+=A_p->val[iPtr];
269 }
270 /*printf("-> %e \n",S);*/
271 out->inv_A_FF[i]=1./S;
272 }
273 }
274 }
275
276 /*check whether coarsening process actually makes sense to continue.
277 if coarse matrix at least smaller by 30% then continue, otherwise we stop.*/
278 if ((out->n_F*100/n)<30) {
279 level=1;
280 }
281
282 if ( Paso_noError()) {
283 /* if there are no nodes in the coarse level there is no more work to do */
284 out->n_C=n-out->n_F;
285
286 /*if (out->n_F>500) */
287 out->rows_in_C=MEMALLOC(out->n_C,index_t);
288 out->mask_C=MEMALLOC(n,index_t);
289 if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) {
290 /* creates an index for C from mask */
291 #pragma omp parallel for private(i) schedule(static)
292 for (i = 0; i < n; ++i) counter[i]=! mis_marker[i];
293 Paso_Util_cumsum(n,counter);
294 #pragma omp parallel for private(i) schedule(static)
295 for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1;
296 #pragma omp parallel for private(i) schedule(static)
297 for (i = 0; i < n; ++i) {
298 if (! mis_marker[i]) {
299 out->rows_in_C[counter[i]]=i;
300 out->mask_C[i]=counter[i];
301 } else {
302 out->mask_C[i]=-1;
303 }
304 }
305 }
306 }
307 if ( Paso_noError()) {
308 /* get A_CF block: */
309 out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F);
310 /* get A_FC block: */
311 out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C);
312 /* get A_CC block: */
313 schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C);
314
315 }
316 if ( Paso_noError()) {
317 /*find the pattern of the schur complement with fill in*/
318 temp1=Paso_Pattern_multiply(PATTERN_FORMAT_DEFAULT,out->A_CF->pattern,out->A_FC->pattern);
319 temp2=Paso_Pattern_binop(PATTERN_FORMAT_DEFAULT, schur->pattern, temp1);
320 schur_withFillIn=Paso_SparseMatrix_alloc(A_p->type,temp2,1,1, TRUE);
321 Paso_Pattern_free(temp1);
322 Paso_Pattern_free(temp2);
323 }
324 if ( Paso_noError()) {
325 /* copy values over*/
326 #pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static)
327 for (i = 0; i < schur_withFillIn->numRows; ++i) {
328 for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) {
329 j=schur_withFillIn->pattern->index[iPtr];
330 iPtr_s=schur->pattern->ptr[i];
331 index=&(schur->pattern->index[iPtr_s]);
332 where_p=(index_t*)bsearch(&j,
333 index,
334 schur->pattern->ptr[i + 1]-schur->pattern->ptr[i],
335 sizeof(index_t),
336 Paso_comparIndex);
337 if (where_p!=NULL) {
338 schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)];
339 }
340 }
341 }
342 Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC);
343 out->AMLI_of_Schur=Paso_Solver_getAMLI(schur_withFillIn,level-1,options);
344 }
345 /* allocate work arrays for AMLI application */
346 if (Paso_noError()) {
347 out->x_F=MEMALLOC(n_block*out->n_F,double);
348 out->b_F=MEMALLOC(n_block*out->n_F,double);
349 out->x_C=MEMALLOC(n_block*out->n_C,double);
350 out->b_C=MEMALLOC(n_block*out->n_C,double);
351
352 if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) {
353 #pragma omp parallel for private(i) schedule(static)
354 for (i = 0; i < out->n_F; ++i) {
355 out->x_F[i]=0.;
356 out->b_F[i]=0.;
357 }
358 #pragma omp parallel for private(i) schedule(static)
359 for (i = 0; i < out->n_C; ++i) {
360 out->x_C[i]=0.;
361 out->b_C[i]=0.;
362 }
363 }
364 }
365 Paso_SparseMatrix_free(schur);
366 Paso_SparseMatrix_free(schur_withFillIn);
367 }
368 }
369 }
370 TMPMEMFREE(mis_marker);
371 TMPMEMFREE(counter);
372
373 if (Paso_noError()) {
374 if (verbose && level>0 && !out->coarsest_level) {
375 printf("AMLI: level: %d: %d unknowns eliminated. %d left.\n",level, out->n_F,out->n_C);
376 }
377 return out;
378 } else {
379 Paso_Solver_AMLI_free(out);
380 return NULL;
381 }
382 }
383
384 /**************************************************************/
385
386 /* apply AMLI precondition b-> x
387
388 in fact it solves
389
390 [ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FC ] [ x_F ] = [b_F]
391 [ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C]
392
393 in the form
394
395 b->[b_F,b_C]
396 x_F=invA_FF*b_F
397 b_C=b_C-A_CF*x_F
398 x_C=AMLI(b_C)
399 b_F=b_F-A_FC*x_C
400 x_F=invA_FF*b_F
401 x<-[x_F,x_C]
402
403 should be called within a parallel region
404 barrier synconization should be performed to make sure that the input vector available
405
406 */
407
408 void Paso_Solver_solveAMLI(Paso_Solver_AMLI * amli, double * x, double * b) {
409 dim_t i;
410 double time0=0;
411 double *r=NULL, *x0=NULL,*x_F_temp=NULL;
412 bool_t verbose=0;
413
414 dim_t post_sweeps=amli->post_sweeps;
415 dim_t pre_sweeps=amli->pre_sweeps;
416
417 #ifdef UMFPACK
418 Paso_UMFPACK_Handler * ptr=NULL;
419 #endif
420
421 r=MEMALLOC(amli->n,double);
422 x0=MEMALLOC(amli->n,double);
423 x_F_temp=MEMALLOC(amli->n_F,double);
424
425 if (amli->coarsest_level) {
426
427 time0=Paso_timer();
428 /*If all unknown are eliminated then Jacobi is the best preconditioner*/
429 if (amli->n_F==0 || amli->n_F==amli->n) {
430 Paso_Solver_solveLocalJacobi(amli->A, amli->GS,x,b);
431 }
432 else {
433 #ifdef MKL
434 Paso_MKL1(amli->AOffset1,x,b,verbose);
435 #else
436 #ifdef UMFPACK
437 ptr=(Paso_UMFPACK_Handler *)(amli->solver);
438 Paso_UMFPACK1(&ptr,amli->A,x,b,verbose);
439 amli->solver=(void*) ptr;
440 #else
441 Paso_Solver_solveLocalJacobi(amli->A,amli->GS,x,b);
442 #endif
443 #endif
444 }
445 time0=Paso_timer()-time0;
446 if (verbose) fprintf(stderr,"timing: DIRECT SOLVER: %e\n",time0);
447
448 } else {
449 /* presmoothing */
450 time0=Paso_timer();
451 Paso_Solver_solveLocalJacobi(amli->A,amli->GS,x,b);
452
453 /***************/
454 #pragma omp parallel for private(i) schedule(static)
455 for (i=0;i<amli->n;++i) r[i]=b[i];
456
457 while(pre_sweeps>1) {
458 #pragma omp parallel for private(i) schedule(static)
459 for (i=0;i<amli->n;++i) r[i]+=b[i];
460
461 /* Compute the residual b=b-Ax*/
462 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A,x,1.,r);
463 /* Go round again*/
464 Paso_Solver_solveLocalJacobi(amli->A,amli->GS,x,r);
465 pre_sweeps-=1;
466 }
467 /***************/
468
469 time0=Paso_timer()-time0;
470 if (verbose) fprintf(stderr,"timing: Presmooting: %e\n",time0);
471 /* end of presmoothing */
472
473
474 time0=Paso_timer();
475 #pragma omp parallel for private(i) schedule(static)
476 for (i=0;i<amli->n;++i) r[i]=b[i];
477
478 /*r=b-Ax*/
479 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A,x,1.,r);
480
481 /* r->[b_F,b_C] */
482 #pragma omp parallel for private(i) schedule(static)
483 for (i=0;i<amli->n_F;++i) amli->b_F[i]=r[amli->rows_in_F[i]];
484
485 #pragma omp parallel for private(i) schedule(static)
486 for (i=0;i<amli->n_C;++i) amli->b_C[i]=r[amli->rows_in_C[i]];
487
488 /* x_F=invA_FF*b_F */
489 Paso_Solver_applyBlockDiagonalMatrix(1,amli->n_F,amli->inv_A_FF,amli->A_FF_pivot,amli->x_F,amli->b_F);
490
491 /* b_C=b_C-A_CF*x_F */
492 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A_CF,amli->x_F,1.,amli->b_C);
493
494 time0=Paso_timer()-time0;
495 if (verbose) fprintf(stderr,"timing: Before next level: %e\n",time0);
496
497 /* x_C=AMLI(b_C) */
498 Paso_Solver_solveAMLI(amli->AMLI_of_Schur,amli->x_C,amli->b_C);
499
500 time0=Paso_timer();
501
502 /* b_F=-A_FC*x_C */
503 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A_FC,amli->x_C,0.,amli->b_F);
504 /* x_F_temp=invA_FF*b_F */
505 Paso_Solver_applyBlockDiagonalMatrix(1,amli->n_F,amli->inv_A_FF,amli->A_FF_pivot,x_F_temp,amli->b_F);
506
507 #pragma omp parallel for private(i) schedule(static)
508 for (i=0;i<amli->n_F;++i) {
509 amli->x_F[i]+=x_F_temp[i];
510 }
511
512 /* x<-[x_F,x_C] */
513 #pragma omp parallel for private(i) schedule(static)
514 for (i=0;i<amli->n;++i) {
515 if (amli->mask_C[i]>-1) {
516 x[i]+=amli->x_C[amli->mask_C[i]];
517 } else {
518 x[i]+=amli->x_F[amli->mask_F[i]];
519 }
520 }
521
522 time0=Paso_timer()-time0;
523 if (verbose) fprintf(stderr,"timing: After next level: %e\n",time0);
524
525 /*postsmoothing*/
526 time0=Paso_timer();
527 #pragma omp parallel for private(i) schedule(static)
528 for (i=0;i<amli->n;++i) r[i]=b[i];
529
530 /*r=b-Ax */
531 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A,x,1.,r);
532 Paso_Solver_solveLocalJacobi(amli->A,amli->GS,x0,r);
533
534 #pragma omp parallel for private(i) schedule(static)
535 for (i=0;i<amli->n;++i) x[i]+=x0[i];
536
537
538 /***************/
539 while(post_sweeps>1) {
540
541 #pragma omp parallel for private(i) schedule(static)
542 for (i=0;i<amli->n;++i) r[i]=b[i];
543
544 Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amli->A,x,1.,r);
545 Paso_Solver_solveLocalJacobi(amli->A,amli->GS,x0,r);
546 #pragma omp parallel for private(i) schedule(static)
547 for (i=0;i<amli->n;++i) {
548 x[i]+=x0[i];
549 }
550 post_sweeps-=1;
551 }
552 /**************/
553
554 time0=Paso_timer()-time0;
555 if (verbose) fprintf(stderr,"timing: Postsmoothing: %e\n",time0);
556
557 /*end of postsmoothing*/
558
559 }
560 MEMFREE(r);
561 MEMFREE(x0);
562 MEMFREE(x_F_temp);
563 return;
564 }

  ViewVC Help
Powered by ViewVC 1.1.26