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

  ViewVC Help
Powered by ViewVC 1.1.26