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 |
|
29 |
/**************************************************************/ |
30 |
|
31 |
/* free all memory used by AMG */ |
32 |
|
33 |
void Paso_Solver_AMG_free(Paso_Solver_AMG * in) { |
34 |
if (in!=NULL) { |
35 |
Paso_Solver_AMG_free(in->AMG_of_Schur); |
36 |
MEMFREE(in->inv_A_FF); |
37 |
MEMFREE(in->A_FF_pivot); |
38 |
Paso_SparseMatrix_free(in->A_FC); |
39 |
Paso_SparseMatrix_free(in->A_CF); |
40 |
MEMFREE(in->rows_in_F); |
41 |
MEMFREE(in->rows_in_C); |
42 |
MEMFREE(in->mask_F); |
43 |
MEMFREE(in->mask_C); |
44 |
MEMFREE(in->x_F); |
45 |
MEMFREE(in->b_F); |
46 |
MEMFREE(in->x_C); |
47 |
MEMFREE(in->b_C); |
48 |
MEMFREE(in); |
49 |
} |
50 |
} |
51 |
|
52 |
/**************************************************************/ |
53 |
|
54 |
/* constructs the block-block factorization of |
55 |
|
56 |
[ A_FF A_FC ] |
57 |
A_p= |
58 |
[ A_CF A_FF ] |
59 |
|
60 |
to |
61 |
|
62 |
[ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] |
63 |
[ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] |
64 |
|
65 |
|
66 |
where S=A_FF-ACF*invA_FF*A_FC within the shape of S |
67 |
|
68 |
then AMG is applied to S again until S becomes empty |
69 |
|
70 |
*/ |
71 |
Paso_Solver_AMG* Paso_Solver_getAMG(Paso_SparseMatrix *A_p,bool_t verbose) { |
72 |
Paso_Solver_AMG* out=NULL; |
73 |
dim_t n=A_p->numRows; |
74 |
dim_t n_block=A_p->row_block_size; |
75 |
index_t* mis_marker=NULL; |
76 |
index_t* counter=NULL; |
77 |
index_t iPtr,*index, *where_p; |
78 |
dim_t i,k; |
79 |
Paso_SparseMatrix * schur=NULL; |
80 |
double A11,A12,A13,A21,A22,A23,A31,A32,A33,D,time0,time1,time2; |
81 |
|
82 |
|
83 |
/* identify independend set of rows/columns */ |
84 |
mis_marker=TMPMEMALLOC(n,index_t); |
85 |
counter=TMPMEMALLOC(n,index_t); |
86 |
out=MEMALLOC(1,Paso_Solver_AMG); |
87 |
out->AMG_of_Schur=NULL; |
88 |
out->inv_A_FF=NULL; |
89 |
out->A_FF_pivot=NULL; |
90 |
out->A_FC=NULL; |
91 |
out->A_CF=NULL; |
92 |
out->rows_in_F=NULL; |
93 |
out->rows_in_C=NULL; |
94 |
out->mask_F=NULL; |
95 |
out->mask_C=NULL; |
96 |
out->x_F=NULL; |
97 |
out->b_F=NULL; |
98 |
out->x_C=NULL; |
99 |
out->b_C=NULL; |
100 |
|
101 |
if ( !(Paso_checkPtr(mis_marker) || Paso_checkPtr(out) || Paso_checkPtr(counter) ) ) { |
102 |
/* identify independend set of rows/columns */ |
103 |
time0=Paso_timer(); |
104 |
#pragma omp parallel for private(i) schedule(static) |
105 |
for (i=0;i<n;++i) mis_marker[i]=-1; |
106 |
Paso_Pattern_mis(A_p->pattern,mis_marker); |
107 |
time2=Paso_timer()-time0; |
108 |
if (Paso_noError()) { |
109 |
#pragma omp parallel for private(i) schedule(static) |
110 |
for (i = 0; i < n; ++i) counter[i]=mis_marker[i]; |
111 |
out->n=n; |
112 |
out->n_block=n_block; |
113 |
out->n_F=Paso_Util_cumsum(n,counter); |
114 |
out->mask_F=MEMALLOC(n,index_t); |
115 |
out->rows_in_F=MEMALLOC(out->n_F,index_t); |
116 |
out->inv_A_FF=MEMALLOC(n_block*n_block*out->n_F,double); |
117 |
out->A_FF_pivot=NULL; /* later use for block size>3 */ |
118 |
if (! (Paso_checkPtr(out->mask_F) || Paso_checkPtr(out->inv_A_FF) || Paso_checkPtr(out->rows_in_F) ) ) { |
119 |
#pragma omp parallel |
120 |
{ |
121 |
/* creates an index for F from mask */ |
122 |
#pragma omp for private(i) schedule(static) |
123 |
for (i = 0; i < out->n_F; ++i) out->rows_in_F[i]=-1; |
124 |
#pragma omp for private(i) schedule(static) |
125 |
for (i = 0; i < n; ++i) { |
126 |
if (mis_marker[i]) { |
127 |
out->rows_in_F[counter[i]]=i; |
128 |
out->mask_F[i]=counter[i]; |
129 |
} else { |
130 |
out->mask_F[i]=-1; |
131 |
} |
132 |
} |
133 |
#pragma omp for private(i, where_p,iPtr,A11,A12,A13,A21,A22,A23,A31,A32,A33,D,index) schedule(static) |
134 |
for (i = 0; i < out->n_F; i++) { |
135 |
/* find main diagonal */ |
136 |
iPtr=A_p->pattern->ptr[out->rows_in_F[i]]; |
137 |
index=&(A_p->pattern->index[iPtr]); |
138 |
where_p=(index_t*)bsearch(&out->rows_in_F[i], |
139 |
index, |
140 |
A_p->pattern->ptr[out->rows_in_F[i] + 1]-A_p->pattern->ptr[out->rows_in_F[i]], |
141 |
sizeof(index_t), |
142 |
Paso_comparIndex); |
143 |
if (where_p==NULL) { |
144 |
Paso_setError(VALUE_ERROR, "Paso_Solver_getAMG: main diagonal element missing."); |
145 |
} else { |
146 |
iPtr+=(index_t)(where_p-index); |
147 |
/* get inverse of A_FF block: */ |
148 |
if (n_block==1) { |
149 |
if (ABS(A_p->val[iPtr])>0.) { |
150 |
out->inv_A_FF[i]=1./A_p->val[iPtr]; |
151 |
} else { |
152 |
Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getAMG: Break-down in AMG decomposition: non-regular main diagonal block."); |
153 |
} |
154 |
} else if (n_block==2) { |
155 |
A11=A_p->val[iPtr*4]; |
156 |
A21=A_p->val[iPtr*4+1]; |
157 |
A12=A_p->val[iPtr*4+2]; |
158 |
A22=A_p->val[iPtr*4+3]; |
159 |
D = A11*A22-A12*A21; |
160 |
if (ABS(D) > 0 ){ |
161 |
D=1./D; |
162 |
out->inv_A_FF[i*4]= A22*D; |
163 |
out->inv_A_FF[i*4+1]=-A21*D; |
164 |
out->inv_A_FF[i*4+2]=-A12*D; |
165 |
out->inv_A_FF[i*4+3]= A11*D; |
166 |
} else { |
167 |
Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getAMG:Break-down in AMG decomposition: non-regular main diagonal block."); |
168 |
} |
169 |
} else if (n_block==3) { |
170 |
A11=A_p->val[iPtr*9 ]; |
171 |
A21=A_p->val[iPtr*9+1]; |
172 |
A31=A_p->val[iPtr*9+2]; |
173 |
A12=A_p->val[iPtr*9+3]; |
174 |
A22=A_p->val[iPtr*9+4]; |
175 |
A32=A_p->val[iPtr*9+5]; |
176 |
A13=A_p->val[iPtr*9+6]; |
177 |
A23=A_p->val[iPtr*9+7]; |
178 |
A33=A_p->val[iPtr*9+8]; |
179 |
D = A11*(A22*A33-A23*A32)+ A12*(A31*A23-A21*A33)+A13*(A21*A32-A31*A22); |
180 |
if (ABS(D) > 0 ){ |
181 |
D=1./D; |
182 |
out->inv_A_FF[i*9 ]=(A22*A33-A23*A32)*D; |
183 |
out->inv_A_FF[i*9+1]=(A31*A23-A21*A33)*D; |
184 |
out->inv_A_FF[i*9+2]=(A21*A32-A31*A22)*D; |
185 |
out->inv_A_FF[i*9+3]=(A13*A32-A12*A33)*D; |
186 |
out->inv_A_FF[i*9+4]=(A11*A33-A31*A13)*D; |
187 |
out->inv_A_FF[i*9+5]=(A12*A31-A11*A32)*D; |
188 |
out->inv_A_FF[i*9+6]=(A12*A23-A13*A22)*D; |
189 |
out->inv_A_FF[i*9+7]=(A13*A21-A11*A23)*D; |
190 |
out->inv_A_FF[i*9+8]=(A11*A22-A12*A21)*D; |
191 |
} else { |
192 |
Paso_setError(ZERO_DIVISION_ERROR, "Paso_Solver_getAMG:Break-down in AMG decomposition: non-regular main diagonal block."); |
193 |
} |
194 |
} |
195 |
} |
196 |
} |
197 |
} /* end parallel region */ |
198 |
|
199 |
if( Paso_noError()) { |
200 |
/* if there are no nodes in the coarse level there is no more work to do */ |
201 |
out->n_C=n-out->n_F; |
202 |
if (out->n_C>0) { |
203 |
out->rows_in_C=MEMALLOC(out->n_C,index_t); |
204 |
out->mask_C=MEMALLOC(n,index_t); |
205 |
if (! (Paso_checkPtr(out->mask_C) || Paso_checkPtr(out->rows_in_C) ) ) { |
206 |
/* creates an index for C from mask */ |
207 |
#pragma omp parallel for private(i) schedule(static) |
208 |
for (i = 0; i < n; ++i) counter[i]=! mis_marker[i]; |
209 |
Paso_Util_cumsum(n,counter); |
210 |
#pragma omp parallel |
211 |
{ |
212 |
#pragma omp for private(i) schedule(static) |
213 |
for (i = 0; i < out->n_C; ++i) out->rows_in_C[i]=-1; |
214 |
#pragma omp for private(i) schedule(static) |
215 |
for (i = 0; i < n; ++i) { |
216 |
if (! mis_marker[i]) { |
217 |
out->rows_in_C[counter[i]]=i; |
218 |
out->mask_C[i]=counter[i]; |
219 |
} else { |
220 |
out->mask_C[i]=-1; |
221 |
} |
222 |
} |
223 |
} /* end parallel region */ |
224 |
/* get A_CF block: */ |
225 |
out->A_CF=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_F,out->rows_in_C,out->mask_F); |
226 |
if (Paso_noError()) { |
227 |
/* get A_FC block: */ |
228 |
out->A_FC=Paso_SparseMatrix_getSubmatrix(A_p,out->n_F,out->n_C,out->rows_in_F,out->mask_C); |
229 |
/* get A_FF block: */ |
230 |
if (Paso_noError()) { |
231 |
schur=Paso_SparseMatrix_getSubmatrix(A_p,out->n_C,out->n_C,out->rows_in_C,out->mask_C); |
232 |
time0=Paso_timer()-time0; |
233 |
if (Paso_noError()) { |
234 |
time1=Paso_timer(); |
235 |
/* update A_CC block to get Schur complement and then apply AMG to it */ |
236 |
Paso_Solver_updateIncompleteSchurComplement(schur,out->A_CF,out->inv_A_FF,out->A_FF_pivot,out->A_FC); |
237 |
time1=Paso_timer()-time1; |
238 |
out->AMG_of_Schur=Paso_Solver_getAMG(schur,verbose); |
239 |
Paso_SparseMatrix_free(schur); |
240 |
} |
241 |
/* allocate work arrays for AMG application */ |
242 |
if (Paso_noError()) { |
243 |
out->x_F=MEMALLOC(n_block*out->n_F,double); |
244 |
out->b_F=MEMALLOC(n_block*out->n_F,double); |
245 |
out->x_C=MEMALLOC(n_block*out->n_C,double); |
246 |
out->b_C=MEMALLOC(n_block*out->n_C,double); |
247 |
if (! (Paso_checkPtr(out->x_F) || Paso_checkPtr(out->b_F) || Paso_checkPtr(out->x_C) || Paso_checkPtr(out->b_C) ) ) { |
248 |
#pragma omp parallel |
249 |
{ |
250 |
#pragma omp for private(i,k) schedule(static) |
251 |
for (i = 0; i < out->n_F; ++i) { |
252 |
for (k=0; k<n_block;++k) { |
253 |
out->x_F[i*n_block+k]=0.; |
254 |
out->b_F[i*n_block+k]=0.; |
255 |
} |
256 |
} |
257 |
#pragma omp for private(i,k) schedule(static) |
258 |
for (i = 0; i < out->n_C; ++i) { |
259 |
for (k=0; k<n_block;++k) { |
260 |
out->x_C[i*n_block+k]=0.; |
261 |
out->b_C[i*n_block+k]=0.; |
262 |
} |
263 |
} |
264 |
} /* end parallel region */ |
265 |
} |
266 |
} |
267 |
} |
268 |
} |
269 |
} |
270 |
} |
271 |
} |
272 |
} |
273 |
} |
274 |
} |
275 |
TMPMEMFREE(mis_marker); |
276 |
TMPMEMFREE(counter); |
277 |
if (Paso_noError()) { |
278 |
if (verbose) { |
279 |
printf("AMG: %d unknowns eliminated. %d left.\n",out->n_F,n-out->n_F); |
280 |
if (out->n_C>0) { |
281 |
printf("timing: AMG: MIS/reordering/elemination : %e/%e/%e\n",time2,time0,time1); |
282 |
} else { |
283 |
printf("timing: AMG: MIS: %e\n",time2); |
284 |
} |
285 |
} |
286 |
return out; |
287 |
} else { |
288 |
Paso_Solver_AMG_free(out); |
289 |
return NULL; |
290 |
} |
291 |
} |
292 |
|
293 |
/**************************************************************/ |
294 |
|
295 |
/* apply AMG precondition b-> x |
296 |
|
297 |
in fact it solves |
298 |
|
299 |
[ I 0 ] [ A_FF 0 ] [ I invA_FF*A_FF ] [ x_F ] = [b_F] |
300 |
[ A_CF*invA_FF I ] [ 0 S ] [ 0 I ] [ x_C ] = [b_C] |
301 |
|
302 |
in the form |
303 |
|
304 |
b->[b_F,b_C] |
305 |
x_F=invA_FF*b_F |
306 |
b_C=b_C-A_CF*x_F |
307 |
x_C=AMG(b_C) |
308 |
b_F=b_F-A_FC*x_C |
309 |
x_F=invA_FF*b_F |
310 |
x<-[x_F,x_C] |
311 |
|
312 |
should be called within a parallel region |
313 |
barrier synconization should be performed to make sure that the input vector available |
314 |
|
315 |
*/ |
316 |
|
317 |
void Paso_Solver_solveAMG(Paso_Solver_AMG * amg, double * x, double * b) { |
318 |
dim_t i,k; |
319 |
dim_t n_block=amg->n_block; |
320 |
|
321 |
if (amg->n_C==0) { |
322 |
/* x=invA_FF*b */ |
323 |
Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,x,b); |
324 |
} else { |
325 |
/* b->[b_F,b_C] */ |
326 |
if (n_block==1) { |
327 |
#pragma omp parallel for private(i) schedule(static) |
328 |
for (i=0;i<amg->n_F;++i) amg->b_F[i]=b[amg->rows_in_F[i]]; |
329 |
#pragma omp parallel for private(i) schedule(static) |
330 |
for (i=0;i<amg->n_C;++i) amg->b_C[i]=b[amg->rows_in_C[i]]; |
331 |
} else { |
332 |
#pragma omp parallel for private(i,k) schedule(static) |
333 |
for (i=0;i<amg->n_F;++i) |
334 |
for (k=0;k<n_block;k++) amg->b_F[amg->n_block*i+k]=b[n_block*amg->rows_in_F[i]+k]; |
335 |
#pragma omp parallel for private(i,k) schedule(static) |
336 |
for (i=0;i<amg->n_C;++i) |
337 |
for (k=0;k<n_block;k++) amg->b_C[amg->n_block*i+k]=b[n_block*amg->rows_in_C[i]+k]; |
338 |
} |
339 |
/* x_F=invA_FF*b_F */ |
340 |
Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F); |
341 |
/* b_C=b_C-A_CF*x_F */ |
342 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_CF,amg->x_F,1.,amg->b_C); |
343 |
/* x_C=AMG(b_C) */ |
344 |
Paso_Solver_solveAMG(amg->AMG_of_Schur,amg->x_C,amg->b_C); |
345 |
/* b_F=b_F-A_FC*x_C */ |
346 |
Paso_SparseMatrix_MatrixVector_CSR_OFFSET0(-1.,amg->A_FC,amg->x_C,1.,amg->b_F); |
347 |
/* x_F=invA_FF*b_F */ |
348 |
Paso_Solver_applyBlockDiagonalMatrix(n_block,amg->n_F,amg->inv_A_FF,amg->A_FF_pivot,amg->x_F,amg->b_F); |
349 |
/* x<-[x_F,x_C] */ |
350 |
if (n_block==1) { |
351 |
#pragma omp parallel for private(i) schedule(static) |
352 |
for (i=0;i<amg->n;++i) { |
353 |
if (amg->mask_C[i]>-1) { |
354 |
x[i]=amg->x_C[amg->mask_C[i]]; |
355 |
} else { |
356 |
x[i]=amg->x_F[amg->mask_F[i]]; |
357 |
} |
358 |
} |
359 |
} else { |
360 |
#pragma omp parallel for private(i,k) schedule(static) |
361 |
for (i=0;i<amg->n;++i) { |
362 |
if (amg->mask_C[i]>-1) { |
363 |
for (k=0;k<n_block;k++) x[n_block*i+k]=amg->x_C[n_block*amg->mask_C[i]+k]; |
364 |
} else { |
365 |
for (k=0;k<n_block;k++) x[n_block*i+k]=amg->x_F[n_block*amg->mask_F[i]+k]; |
366 |
} |
367 |
} |
368 |
} |
369 |
/* all done */ |
370 |
} |
371 |
return; |
372 |
} |
373 |
|
374 |
/* |
375 |
* $Log$ |
376 |
* |
377 |
*/ |