1 |
|
2 |
/******************************************************* |
3 |
* |
4 |
* Copyright (c) 2003-2008 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 with reordering */ |
18 |
|
19 |
/**************************************************************/ |
20 |
|
21 |
/* Author: artak@access.edu.au */ |
22 |
|
23 |
/**************************************************************/ |
24 |
|
25 |
#include "Paso.h" |
26 |
#include "Solver.h" |
27 |
#include "PasoUtil.h" |
28 |
#include "UMFPACK.h" |
29 |
#include "Pattern_coupling.h" |
30 |
|
31 |
/**************************************************************/ |
32 |
|
33 |
/* free all memory used by AMG */ |
34 |
|
35 |
void Paso_Solver_AMG_free(Paso_Solver_AMG * in) { |
36 |
if (in!=NULL) { |
37 |
Paso_Solver_AMG_free(in->AMG_of_Schur); |
38 |
Paso_Solver_Jacobi_free(in->GS); |
39 |
MEMFREE(in->inv_A_FF); |
40 |
MEMFREE(in->A_FF_pivot); |
41 |
Paso_SparseMatrix_free(in->A_FC); |
42 |
Paso_SparseMatrix_free(in->A_CF); |
43 |
MEMFREE(in->rows_in_F); |
44 |
MEMFREE(in->rows_in_C); |
45 |
MEMFREE(in->mask_F); |
46 |
MEMFREE(in->mask_C); |
47 |
MEMFREE(in->x_F); |
48 |
MEMFREE(in->b_F); |
49 |
MEMFREE(in->x_C); |
50 |
MEMFREE(in->b_C); |
51 |
MEMFREE(in); |
52 |
} |
53 |
} |
54 |
|
55 |
/**************************************************************/ |
56 |
|
57 |
/* constructs the block-block factorization of |
58 |
|
59 |
[ A_FF A_FC ] |
60 |
A_p= |
61 |
[ A_CF A_FF ] |
62 |
|
63 |
to |
64 |
|
65 |
[ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] |
66 |
[ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] |
67 |
|
68 |
|
69 |
where S=A_FF-ACF*invA_FF*A_FC within the shape of S |
70 |
|
71 |
then AMG is applied to S again until S becomes empty |
72 |
|
73 |
*/ |
74 |
Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose,dim_t level, double couplingParam) { |
75 |
Paso_Solver_AMG* out=NULL; |
76 |
dim_t n=A_p->numRows; |
77 |
dim_t n_block=A_p->row_block_size; |
78 |
index_t* mis_marker=NULL; |
79 |
index_t* counter=NULL; |
80 |
index_t iPtr,*index, *where_p, iPtr_s; |
81 |
dim_t i,k,j; |
82 |
Paso_SparseMatrix * schur=NULL; |
83 |
Paso_SparseMatrix * schur_withFillIn=NULL; |
84 |
double time0=0,time1=0,time2=0,S=0; |
85 |
/*Paso_Pattern* test;*/ |
86 |
|
87 |
/* identify independend set of rows/columns */ |
88 |
mis_marker=TMPMEMALLOC(n,index_t); |
89 |
counter=TMPMEMALLOC(n,index_t); |
90 |
out=MEMALLOC(1,Paso_Solver_AMG); |
91 |
out->AMG_of_Schur=NULL; |
92 |
out->inv_A_FF=NULL; |
93 |
out->A_FF_pivot=NULL; |
94 |
out->A_FC=NULL; |
95 |
out->A_CF=NULL; |
96 |
out->rows_in_F=NULL; |
97 |
out->rows_in_C=NULL; |
98 |
out->mask_F=NULL; |
99 |
out->mask_C=NULL; |
100 |
out->x_F=NULL; |
101 |
out->b_F=NULL; |
102 |
out->x_C=NULL; |
103 |
out->b_C=NULL; |
104 |
out->GS=NULL; |
105 |
out->A=Paso_SparseMatrix_getReference(A_p); |
106 |
/*out->GS=Paso_Solver_getGS(A_p,verbose);*/ |
107 |
out->GS=Paso_Solver_getJacobi(A_p); |
108 |
/*out->GS->sweeps=2;*/ |
109 |
out->level=level; |
110 |
|
111 |
if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) { |
112 |
/* identify independend set of rows/columns */ |
113 |
#pragma omp parallel for private(i) schedule(static) |
114 |
for (i=0;i<n;++i) mis_marker[i]=-1; |
115 |
/*Paso_Pattern_RS(A_p,mis_marker,0.25);*/ |
116 |
/*Paso_Pattern_Aggregiation(A_p,mis_marker,0.5);*/ |
117 |
Paso_Pattern_coup(A_p,mis_marker,couplingParam); |
118 |
if (Paso_noError()) { |
119 |
#pragma omp parallel for private(i) schedule(static) |
120 |
for (i = 0; i < n; ++i) counter[i]=mis_marker[i]; |
121 |
out->n=n; |
122 |
out->n_block=n_block; |
123 |
out->n_F=Paso_Util_cumsum(n,counter); |
124 |
out->mask_F=MEMALLOC(n,index_t); |
125 |
out->rows_in_F=MEMALLOC(out->n_F,index_t); |
126 |
out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double); |
127 |
out->A_FF_pivot=NULL; /* later use for block size>3 */ |
128 |
if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) { |
129 |
/* creates an index for F from mask */ |
130 |
#pragma omp parallel for private(i) schedule(static) |
131 |
for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1; |
132 |
#pragma omp parallel for private(i) schedule(static) |
133 |
for (i = 0; i < n; ++i) { |
134 |
if (mis_marker[i]) { |
135 |
out->rows_in_F[counter[i]]=i; |
136 |
out->mask_F[i]=counter[i]; |
137 |
} else { |
138 |
out->mask_F[i]=-1; |
139 |
} |
140 |
} |
141 |
|
142 |
/* Compute row-sum for getting rs(A_FF)^-1*/ |
143 |
#pragma omp parallel for private(i,iPtr,j,S) schedule(static) |
144 |
for (i = 0; i < out->n_F; ++i) { |
145 |
S=0; |
146 |
for (iPtr=A_p->pattern->ptr[out->rows_in_F[i]];iPtr<A_p->pattern->ptr[out->rows_in_F[i] + 1]; ++iPtr) { |
147 |
j=A_p->pattern->index[iPtr]; |
148 |
if (mis_marker[j]) |
149 |
S+=A_p->val[iPtr]; |
150 |
} |
151 |
out->inv_A_FF[i]=1./S; |
152 |
} |
153 |
|
154 |
if( Paso_noError()) { |
155 |
/* if there are no nodes in the coarse level there is no more work to do */ |
156 |
out->n_C=n-out->n_F; |
157 |
if (level<3) { |
158 |
/*if (out->n_F>500) {*/ |
159 |
out->rows_in_C=MEMALLOC(out->n_C,index_t); |
160 |
out->mask_C=MEMALLOC(n,index_t); |
161 |
if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) { |
162 |
/* creates an index for C from mask */ |
163 |
#pragma omp parallel for private(i) schedule(static) |
164 |
for (i = 0; i < n; ++i) counter[i]=! mis_marker[i]; |
165 |
Paso_Util_cumsum(n,counter); |
166 |
#pragma omp parallel for private(i) schedule(static) |
167 |
for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1; |
168 |
#pragma omp parallel for private(i) schedule(static) |
169 |
for (i = 0; i < n; ++i) { |
170 |
if (! mis_marker[i]) { |
171 |
out->rows_in_C[counter[i]]=i; |
172 |
out->mask_C[i]=counter[i]; |
173 |
} else { |
174 |
out->mask_C[i]=-1; |
175 |
} |
176 |
} |
177 |
/* get A_CF block: */ |
178 |
out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F); |
179 |
if (Paso_noError()) { |
180 |
/* get A_FC block: */ |
181 |
out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C); |
182 |
/* get A_CC block: */ |
183 |
if (Paso_noError()) { |
184 |
schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C); |
185 |
/*find the pattern of the schur complement with fill in*/ |
186 |
|
187 |
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); |
188 |
|
189 |
fprintf(stderr,"Sparsity of Schure: %dx%d LEN %d Percentage %f\n",schur_withFillIn->pattern->numOutput,schur_withFillIn->pattern->numInput,schur_withFillIn->len,schur_withFillIn->len/(1.*schur_withFillIn->pattern->numOutput*schur_withFillIn->pattern->numInput)); |
190 |
|
191 |
/* copy values over*/ |
192 |
#pragma omp parallel for private(i,iPtr,j,iPtr_s,index,where_p) schedule(static) |
193 |
for (i = 0; i < schur_withFillIn->numRows; ++i) { |
194 |
for (iPtr=schur_withFillIn->pattern->ptr[i];iPtr<schur_withFillIn->pattern->ptr[i + 1]; ++iPtr) { |
195 |
j=schur_withFillIn->pattern->index[iPtr]; |
196 |
iPtr_s=schur->pattern->ptr[i]; |
197 |
schur_withFillIn->val[iPtr]=0.; |
198 |
index=&(schur->pattern->index[iPtr_s]); |
199 |
where_p=(index_t*)bsearch(&j, |
200 |
index, |
201 |
schur->pattern->ptr[i + 1]-schur->pattern->ptr[i], |
202 |
sizeof(index_t), |
203 |
Paso_comparIndex); |
204 |
if (where_p!=NULL) { |
205 |
schur_withFillIn->val[iPtr]=schur->val[iPtr_s+(index_t)(where_p-index)]; |
206 |
} |
207 |
} |
208 |
} |
209 |
|
210 |
if (Paso_noError()) { |
211 |
Paso_Solver_updateIncompleteSchurComplement(schur_withFillIn,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC); |
212 |
out->AMG_of_Schur=Paso_Solver_getAMG(schur_withFillIn,verbose,level+1,couplingParam); |
213 |
Paso_SparseMatrix_free(schur); |
214 |
} |
215 |
/* allocate work arrays for AMG application */ |
216 |
if (Paso_noError()) { |
217 |
out->x_F=MEMALLOC(n_block*out->n_F,double); |
218 |
out->b_F=MEMALLOC(n_block*out->n_F,double); |
219 |
out->x_C=MEMALLOC(n_block*out->n_C,double); |
220 |
out->b_C=MEMALLOC(n_block*out->n_C,double); |
221 |
|
222 |
if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) { |
223 |
#pragma omp parallel for private(i,k) schedule(static) |
224 |
for (i = 0; i < out->n_F; ++i) { |
225 |
for (k=0; k<n_block;++k) { |
226 |
out->x_F[i*n_block+k]=0.; |
227 |
out->b_F[i*n_block+k]=0.; |
228 |
} |
229 |
} |
230 |
#pragma omp parallel for private(i,k) schedule(static) |
231 |
for (i = 0; i < out->n_C; ++i) { |
232 |
for (k=0; k<n_block;++k) { |
233 |
out->x_C[i*n_block+k]=0.; |
234 |
out->b_C[i*n_block+k]=0.; |
235 |
} |
236 |
} |
237 |
} |
238 |
} |
239 |
} |
240 |
} |
241 |
} |
242 |
} |
243 |
} |
244 |
} |
245 |
} |
246 |
} |
247 |
TMPMEMFREE(mis_marker); |
248 |
TMPMEMFREE(counter); |
249 |
if (Paso_noError()) { |
250 |
if (verbose) { |
251 |
printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F); |
252 |
if (level<3) { |
253 |
/*if (out->n_F<500) {*/ |
254 |
printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1); |
255 |
} else { |
256 |
printf("timing: AMG: MIS: %e\n",time2); |
257 |
} |
258 |
} |
259 |
return out; |
260 |
} else { |
261 |
Paso_Solver_AMG_free(out); |
262 |
return NULL; |
263 |
} |
264 |
} |
265 |
|
266 |
/**************************************************************/ |
267 |
|
268 |
/* apply AMG precondition b-> x |
269 |
|
270 |
in fact it solves |
271 |
|
272 |
[ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] [ x_F ] = [b_F] |
273 |
[ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C] |
274 |
|
275 |
in the form |
276 |
|
277 |
b->[b_F,b_C] |
278 |
x_F=invA_FF*b_F |
279 |
b_C=b_C-A_CF*x_F |
280 |
x_C=AMG(b_C) |
281 |
b_F=b_F-A_FC*x_C |
282 |
x_F=invA_FF*b_F |
283 |
x<-[x_F,x_C] |
284 |
|
285 |
should be called within a parallel region |
286 |
barrier synconization should be performed to make sure that the input vector available |
287 |
|
288 |
*/ |
289 |
|
290 |
void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) { |
291 |
dim_t i; |
292 |
double *r=MEMALLOC(amg->n,double); |
293 |
/*Paso_Solver_GS* GS=NULL;*/ |
294 |
double *x0=MEMALLOC(amg->n,double); |
295 |
double time0=0; |
296 |
|
297 |
if (amg->level==3) { |
298 |
/*if (amg->n_F<=500) {*/ |
299 |
time0=Paso_timer(); |
300 |
|
301 |
Paso_Solver_solveJacobi(amg->GS,x,b); |
302 |
|
303 |
/* #pragma omp parallel for private(i) schedule(static) |
304 |
for (i=0;i<amg->n;++i) r[i]=b[i]; |
305 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r); |
306 |
Paso_Solver_solveGS(amg->GS,x0,r); |
307 |
#pragma omp parallel for private(i) schedule(static) |
308 |
for (i=0;i<amg->n;++i) { |
309 |
x[i]+=x0[i]; |
310 |
} |
311 |
*/ |
312 |
/*Paso_UMFPACK1(amg->A,x,b,0);*/ |
313 |
|
314 |
time0=Paso_timer()-time0; |
315 |
/*fprintf(stderr,"timing: DIRECT SOLVER: %e/\n",time0);*/ |
316 |
} else { |
317 |
|
318 |
/* presmoothing */ |
319 |
Paso_Solver_solveJacobi(amg->GS,x,b); |
320 |
/* |
321 |
#pragma omp parallel for private(i) schedule(static) |
322 |
for (i=0;i<amg->n;++i) r[i]=b[i]; |
323 |
|
324 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r); |
325 |
Paso_Solver_solveGS(amg->GS,x0,r); |
326 |
|
327 |
#pragma omp parallel for private(i) schedule(static) |
328 |
for (i=0;i<amg->n;++i) { |
329 |
x[i]+=x0[i]; |
330 |
} |
331 |
*/ |
332 |
/* end of presmoothing */ |
333 |
|
334 |
#pragma omp parallel for private(i) schedule(static) |
335 |
for (i=0;i<amg->n;++i) r[i]=b[i]; |
336 |
|
337 |
/*r=b-Ax*/ |
338 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r); |
339 |
|
340 |
/* b->[b_F,b_C] */ |
341 |
#pragma omp parallel for private(i) schedule(static) |
342 |
for (i=0;i<amg->n_F;++i) amg->b_F[i]=r[amg->rows_in_F[i]]; |
343 |
|
344 |
#pragma omp parallel for private(i) schedule(static) |
345 |
for (i=0;i<amg->n_C;++i) amg->b_C[i]=r[amg->rows_in_C[i]]; |
346 |
|
347 |
/* x_F=invA_FF*b_F */ |
348 |
Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F); |
349 |
|
350 |
/* b_C=b_C-A_CF*x_F */ |
351 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C); |
352 |
|
353 |
/* x_C=AMG(b_C) */ |
354 |
Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C); |
355 |
|
356 |
/* b_F=b_F-A_FC*x_C */ |
357 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F); |
358 |
/* x_F=invA_FF*b_F */ |
359 |
Paso_Solver_applyBlockDiagonalMatrix(1,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F); |
360 |
/* x<-[x_F,x_C] */ |
361 |
|
362 |
#pragma omp parallel for private(i) schedule(static) |
363 |
for (i=0;i<amg->n;++i) { |
364 |
if (amg->mask_C[i]>-1) { |
365 |
x[i]=amg->x_C[amg->mask_C[i]]; |
366 |
} else { |
367 |
x[i]=amg->x_F[amg->mask_F[i]]; |
368 |
} |
369 |
} |
370 |
|
371 |
/*postsmoothing*/ |
372 |
|
373 |
#pragma omp parallel for private(i) schedule(static) |
374 |
for (i=0;i<amg->n;++i) r[i]=b[i]; |
375 |
|
376 |
/*r=b-Ax*/ |
377 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r); |
378 |
Paso_Solver_solveJacobi(amg->GS,x0,r); |
379 |
|
380 |
#pragma omp parallel for private(i) schedule(static) |
381 |
for (i=0;i<amg->n;++i) x[i]+=x0[i]; |
382 |
|
383 |
|
384 |
/*#pragma omp parallel for private(i) schedule(static) |
385 |
for (i=0;i<amg->n;++i) r[i]=b[i]; |
386 |
|
387 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A,x,1.,r); |
388 |
Paso_Solver_solveGS(amg->GS,x0,r); |
389 |
|
390 |
#pragma omp parallel for private(i) schedule(static) |
391 |
for (i=0;i<amg->n;++i) { |
392 |
x[i]+=x0[i]; |
393 |
} |
394 |
*/ |
395 |
/*end of postsmoothing*/ |
396 |
|
397 |
} |
398 |
MEMFREE(r); |
399 |
MEMFREE(x0); |
400 |
return; |
401 |
} |
402 |
|
403 |
/* |
404 |
* $Log$ |
405 |
* |
406 |
*/ |