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