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

  ViewVC Help
Powered by ViewVC 1.1.26