/[escript]/branches/domexper/paso/src/AMLI.c
ViewVC logotype

Contents of /branches/domexper/paso/src/AMLI.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3234 - (show annotations)
Mon Oct 4 01:46:30 2010 UTC (9 years, 2 months ago) by jfenwick
File MIME type: text/plain
File size: 19099 byte(s)
Some subdirs need to have changes pulled over but all of the unit tests 
except for modellib appear to work

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

  ViewVC Help
Powered by ViewVC 1.1.26