/[escript]/trunk/paso/src/Smoother.c
ViewVC logotype

Annotation of /trunk/paso/src/Smoother.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3379 - (hide annotations)
Wed Nov 24 04:48:49 2010 UTC (8 years, 3 months ago) by gross
File MIME type: text/plain
File size: 14515 byte(s)
some clarification on lumping
1 jfenwick 1851
2     /*******************************************************
3     *
4 jfenwick 2881 * Copyright (c) 2003-2010 by University of Queensland
5 jfenwick 1851 * 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 gross 3097 /* Paso: Gauss-Seidel */
18 jfenwick 1851
19     /**************************************************************/
20    
21     /* Author: artak@uq.edu.au */
22    
23     /**************************************************************/
24    
25     #include "Paso.h"
26 gross 3120 #include "Preconditioner.h"
27 jfenwick 1851 #include "PasoUtil.h"
28 gross 3097 #include "BlockOps.h"
29 jfenwick 1851
30     #include <stdio.h>
31    
32 gross 3094
33 jfenwick 1851 /**************************************************************/
34    
35 gross 3158 /* free all memory used by Smoother */
36 jfenwick 1851
37 gross 3158 void Paso_Preconditioner_Smoother_free(Paso_Preconditioner_Smoother * in) {
38 jfenwick 1851 if (in!=NULL) {
39 gross 3158 Paso_Preconditioner_LocalSmoother_free(in->localSmoother);
40 jfenwick 1851 MEMFREE(in);
41     }
42     }
43 gross 3158 void Paso_Preconditioner_LocalSmoother_free(Paso_Preconditioner_LocalSmoother * in) {
44 gross 3094 if (in!=NULL) {
45     MEMFREE(in->diag);
46     MEMFREE(in->pivot);
47 gross 3097 MEMFREE(in->buffer);
48 gross 3094 MEMFREE(in);
49     }
50     }
51 jfenwick 1851 /**************************************************************/
52    
53 gross 3094 /* constructs the symmetric Gauss-Seidel preconditioner
54 jfenwick 1851
55     */
56 gross 3158 Paso_Preconditioner_Smoother* Paso_Preconditioner_Smoother_alloc(Paso_SystemMatrix * A_p, const bool_t jacobi, const bool_t is_local, const bool_t verbose)
57 gross 3094 {
58    
59 jfenwick 1851 /* allocations: */
60 gross 3158 Paso_Preconditioner_Smoother* out=MEMALLOC(1,Paso_Preconditioner_Smoother);
61 jfenwick 3259 if (! Esys_checkPtr(out)) {
62 gross 3158 out->localSmoother=Paso_Preconditioner_LocalSmoother_alloc(A_p->mainBlock,jacobi,verbose);
63 gross 3094 out->is_local=is_local;
64 jfenwick 1851 }
65 jfenwick 3259 if (Esys_MPIInfo_noError(A_p->mpi_info)) {
66 jfenwick 1851 return out;
67 gross 3094 } else {
68 gross 3158 Paso_Preconditioner_Smoother_free(out);
69 jfenwick 1851 return NULL;
70     }
71     }
72 jfenwick 3259 Paso_Preconditioner_LocalSmoother* Paso_Preconditioner_LocalSmoother_alloc(Paso_SparseMatrix * A_p, const bool_t jacobi, bool_t verbose)
73 gross 3158 {
74 gross 3094
75     dim_t n=A_p->numRows;
76     dim_t n_block=A_p->row_block_size;
77     dim_t block_size=A_p->block_size;
78    
79 jfenwick 3259 double time0=Esys_timer();
80 gross 3094 /* allocations: */
81 gross 3158 Paso_Preconditioner_LocalSmoother* out=MEMALLOC(1,Paso_Preconditioner_LocalSmoother);
82 jfenwick 3259 if (! Esys_checkPtr(out)) {
83 gross 3094
84     out->diag=MEMALLOC( ((size_t) n) * ((size_t) block_size),double);
85     out->pivot=MEMALLOC( ((size_t) n) * ((size_t) n_block), index_t);
86 gross 3097 out->buffer=MEMALLOC( ((size_t) n) * ((size_t) n_block), double);
87 gross 3158 out->Jacobi=jacobi;
88 gross 3094
89 jfenwick 3259 if ( ! ( Esys_checkPtr(out->diag) || Esys_checkPtr(out->pivot) ) ) {
90 gross 3094 Paso_SparseMatrix_invMain(A_p, out->diag, out->pivot );
91     }
92    
93     }
94 jfenwick 3259 time0=Esys_timer()-time0;
95 gross 3094
96 jfenwick 3259 if (Esys_noError()) {
97 gross 3094 return out;
98     } else {
99 gross 3158 Paso_Preconditioner_LocalSmoother_free(out);
100 gross 3094 return NULL;
101     }
102     }
103 jfenwick 1851
104 gross 3097 /*
105 jfenwick 1851
106 gross 3097 performs a few sweeps of the from
107 jfenwick 1851
108 gross 3097 S (x_{k} - x_{k-1}) = b - A x_{k-1}
109 jfenwick 1851
110 gross 3097 where x_{0}=0 and S provides some approximatioon of A.
111 jfenwick 1851
112 gross 3097 Under MPI S is build on using A_p->mainBlock only.
113 gross 3158 if Smoother is local the defect b - A x_{k-1} is calculated using A_p->mainBlock only.
114 gross 3094
115 jfenwick 1851 */
116    
117 gross 3159 void Paso_Preconditioner_Smoother_solve(Paso_SystemMatrix* A_p, Paso_Preconditioner_Smoother * smoother, double * x, const double * b,
118     const dim_t sweeps, const bool_t x_is_initial)
119 gross 3094 {
120 gross 3097 const dim_t n= (A_p->mainBlock->numRows) * (A_p->mainBlock->row_block_size);
121    
122 gross 3158 double *b_new = smoother->localSmoother->buffer;
123     dim_t nsweeps=sweeps;
124     if (smoother->is_local) {
125 gross 3159 Paso_Preconditioner_LocalSmoother_solve(A_p->mainBlock,smoother->localSmoother,x,b,sweeps,x_is_initial);
126 gross 3097 } else {
127 gross 3159 if (! x_is_initial) {
128     Paso_Copy(n, x, b);
129     Paso_Preconditioner_LocalSmoother_Sweep(A_p->mainBlock,smoother->localSmoother,x);
130     nsweeps--;
131     }
132     while (nsweeps > 0 ) {
133 gross 3158 Paso_Copy(n, b_new, b);
134 gross 3097 Paso_SystemMatrix_MatrixVector((-1.), A_p, x, 1., b_new); /* b_new = b - A*x */
135 gross 3159 Paso_Preconditioner_LocalSmoother_Sweep(A_p->mainBlock,smoother->localSmoother,b_new);
136 gross 3158 Paso_AXPY(n, x, 1., b_new);
137     nsweeps--;
138 gross 3097 }
139    
140 gross 3094 }
141 gross 3097 }
142 gross 3159 void Paso_Preconditioner_LocalSmoother_solve(Paso_SparseMatrix* A_p, Paso_Preconditioner_LocalSmoother * smoother, double * x, const double * b,
143     const dim_t sweeps, const bool_t x_is_initial)
144 gross 3097 {
145     const dim_t n= (A_p->numRows) * (A_p->row_block_size);
146 gross 3158 double *b_new = smoother->buffer;
147     dim_t nsweeps=sweeps;
148 gross 3094
149 gross 3159 if (! x_is_initial) {
150     Paso_Copy(n, x, b);
151     Paso_Preconditioner_LocalSmoother_Sweep(A_p,smoother,x);
152     nsweeps--;
153     }
154 gross 3105
155 gross 3159 while (nsweeps > 0 ) {
156 gross 3158 Paso_Copy(n, b_new, b);
157 gross 3283
158     Paso_SparseMatrix_MatrixVector_CSR_OFFSET0((-1.), A_p, x, 1., b_new); /* b_new = b - A*x */
159 gross 3158 Paso_Preconditioner_LocalSmoother_Sweep(A_p,smoother,b_new);
160 gross 3159 Paso_AXPY(n, x, 1., b_new);
161 gross 3158 nsweeps--;
162 gross 3094 }
163     }
164 jfenwick 1851
165 gross 3158 void Paso_Preconditioner_LocalSmoother_Sweep(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * smoother, double * x)
166 gross 3094 {
167     const dim_t nt=omp_get_max_threads();
168 gross 3158 if (smoother->Jacobi) {
169     Paso_BlockOps_allMV(A->row_block_size,A->numRows,smoother->diag,smoother->pivot,x);
170 gross 3097 } else {
171 gross 3158 if (nt < 2) {
172     Paso_Preconditioner_LocalSmoother_Sweep_sequential(A,smoother,x);
173     } else {
174     Paso_Preconditioner_LocalSmoother_Sweep_colored(A,smoother,x);
175     }
176 gross 3094 }
177 jfenwick 1851 }
178    
179 gross 3097 /* inplace Gauss-Seidel sweep in seqential mode: */
180 gross 3094
181 gross 3158 void Paso_Preconditioner_LocalSmoother_Sweep_sequential(Paso_SparseMatrix* A_p, Paso_Preconditioner_LocalSmoother * smoother, double * x)
182 gross 3097 {
183     const dim_t n=A_p->numRows;
184     const dim_t n_block=A_p->row_block_size;
185 gross 3158 const double *diag = smoother->diag;
186     /* const index_t* pivot = smoother->pivot;
187 gross 3097 const dim_t block_size=A_p->block_size; use for block size >3*/
188 gross 3094
189 gross 3097 register dim_t i,k;
190     register index_t iptr_ik, mm;
191    
192     const index_t* ptr_main = Paso_SparseMatrix_borrowMainDiagonalPointer(A_p);
193     /* forward substitution */
194    
195     if (n_block==1) {
196 gross 3379 x[0]*=diag[0];
197 gross 3098 for (i = 1; i < n; ++i) {
198 gross 3097 mm=ptr_main[i];
199 gross 3105 /* x_i=x_i-a_ik*x_k (with k<i) */
200 gross 3097 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<mm; ++iptr_ik) {
201 gross 3379 k=A_p->pattern->index[iptr_ik];
202     x[i]-=A_p->val[iptr_ik]*x[k];
203 gross 3097 }
204 gross 3379 x[i]*=diag[i];
205 gross 3097 }
206     } else if (n_block==2) {
207 gross 3098 Paso_BlockOps_MV_2(&x[0], &diag[0], &x[0]);
208     for (i = 1; i < n; ++i) {
209 gross 3097 mm=ptr_main[i];
210     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<mm; ++iptr_ik) {
211     k=A_p->pattern->index[iptr_ik];
212     Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
213     }
214 gross 3098 Paso_BlockOps_MV_2(&x[2*i], &diag[4*i], &x[2*i]);
215 gross 3097 }
216     } else if (n_block==3) {
217 gross 3098 Paso_BlockOps_MV_3(&x[0], &diag[0], &x[0]);
218     for (i = 1; i < n; ++i) {
219 gross 3097 mm=ptr_main[i];
220     /* x_i=x_i-a_ik*x_k */
221     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<mm; ++iptr_ik) {
222     k=A_p->pattern->index[iptr_ik];
223     Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
224     }
225 gross 3098 Paso_BlockOps_MV_3(&x[3*i], &diag[9*i], &x[3*i]);
226 gross 3097 }
227     } /* add block size >3 */
228 gross 3105
229 gross 3097 /* backward substitution */
230    
231     if (n_block==1) {
232 gross 3105 for (i = n-2; i > -1; --i) {
233 gross 3098 mm=ptr_main[i];
234     Paso_BlockOps_MV_1(&x[i], &A_p->val[mm], &x[i]);
235     for (iptr_ik=mm+1; iptr_ik < A_p->pattern->ptr[i+1]; ++iptr_ik) {
236     k=A_p->pattern->index[iptr_ik];
237     Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
238     }
239     Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
240 gross 3097 }
241 gross 3098
242 gross 3097 } else if (n_block==2) {
243 gross 3105 for (i = n-2; i > -1; --i) {
244 gross 3098 mm=ptr_main[i];
245     Paso_BlockOps_MV_2(&x[2*i], &A_p->val[4*mm], &x[2*i]);
246     for (iptr_ik=mm+1; iptr_ik < A_p->pattern->ptr[i+1]; ++iptr_ik) {
247     k=A_p->pattern->index[iptr_ik];
248     Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
249     }
250     Paso_BlockOps_MV_2(&x[2*i], &diag[i*4], &x[2*i]);
251 gross 3283
252 gross 3097 }
253     } else if (n_block==3) {
254 gross 3105 for (i = n-2; i > -1; --i) {
255 gross 3098 mm=ptr_main[i];
256     Paso_BlockOps_MV_3(&x[3*i], &A_p->val[9*mm], &x[3*i]);
257 gross 3097
258 gross 3098 for (iptr_ik=mm+1; iptr_ik < A_p->pattern->ptr[i+1]; ++iptr_ik) {
259     k=A_p->pattern->index[iptr_ik];
260     Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
261     }
262     Paso_BlockOps_MV_3(&x[3*i], &diag[i*9], &x[3*i]);
263 gross 3097 }
264    
265     } /* add block size >3 */
266    
267     return;
268     }
269 gross 3094
270 gross 3158 void Paso_Preconditioner_LocalSmoother_Sweep_colored(Paso_SparseMatrix* A_p, Paso_Preconditioner_LocalSmoother * smoother, double * x)
271 gross 3094 {
272     const dim_t n=A_p->numRows;
273     const dim_t n_block=A_p->row_block_size;
274 gross 3158 const double *diag = smoother->diag;
275     index_t* pivot = smoother->pivot;
276 gross 3097 const dim_t block_size=A_p->block_size;
277 gross 3094
278     register dim_t i,k;
279     register index_t color,iptr_ik, mm;
280    
281     const index_t* coloring = Paso_Pattern_borrowColoringPointer(A_p->pattern);
282     const dim_t num_colors = Paso_Pattern_getNumColors(A_p->pattern);
283     const index_t* ptr_main = Paso_SparseMatrix_borrowMainDiagonalPointer(A_p);
284    
285     /* forward substitution */
286 gross 3105
287 gross 3094
288     /* color = 0 */
289     if (n_block==1) {
290 gross 3097 #pragma omp parallel for schedule(static) private(i)
291 gross 3094 for (i = 0; i < n; ++i) {
292 gross 3097 if (coloring[i]== 0 ) Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
293     }
294 gross 3094 } else if (n_block==2) {
295 gross 3097 #pragma omp parallel for schedule(static) private(i)
296 gross 3094 for (i = 0; i < n; ++i) {
297 gross 3097 if (coloring[i]== 0 ) Paso_BlockOps_MV_2(&x[2*i], &diag[i*4], &x[2*i]);
298 gross 3094 }
299 gross 3097 } else if (n_block==3) {
300     #pragma omp parallel for schedule(static) private(i)
301 gross 3094 for (i = 0; i < n; ++i) {
302 gross 3097 if (coloring[i]==0) Paso_BlockOps_MV_3(&x[3*i], &diag[i*9], &x[3*i]);
303 gross 3094 }
304 gross 3097 } else {
305     #pragma omp parallel for schedule(static) private(i)
306     for (i = 0; i < n; ++i) {
307     if (coloring[i]==0) Paso_BlockOps_Solve_N(n_block, &x[n_block*i], &diag[block_size*i], &pivot[n_block*i], &x[n_block*i]);
308     }
309     }
310 gross 3094
311     for (color=1;color<num_colors;++color) {
312 gross 3105
313 gross 3094 if (n_block==1) {
314 gross 3097 #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
315 gross 3094 for (i = 0; i < n; ++i) {
316     if (coloring[i]==color) {
317     /* x_i=x_i-a_ik*x_k */
318     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
319 gross 3105 k=A_p->pattern->index[iptr_ik];
320 gross 3097 if (coloring[k]<color) Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
321 gross 3094 }
322 gross 3097 Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
323 gross 3094 }
324     }
325     } else if (n_block==2) {
326 gross 3097 #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
327 gross 3094 for (i = 0; i < n; ++i) {
328     if (coloring[i]==color) {
329     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
330     k=A_p->pattern->index[iptr_ik];
331 gross 3097 if (coloring[k]<color) Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
332 gross 3094 }
333 gross 3097 Paso_BlockOps_MV_2(&x[2*i], &diag[4*i], &x[2*i]);
334 gross 3094 }
335     }
336     } else if (n_block==3) {
337 gross 3097 #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
338 gross 3094 for (i = 0; i < n; ++i) {
339     if (coloring[i]==color) {
340     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
341     k=A_p->pattern->index[iptr_ik];
342 gross 3097 if (coloring[k]<color) Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
343 gross 3094 }
344 gross 3097 Paso_BlockOps_MV_3(&x[3*i], &diag[9*i], &x[3*i]);
345 gross 3094 }
346     }
347 gross 3097 } else {
348     #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
349     for (i = 0; i < n; ++i) {
350     if (coloring[i] == color) {
351     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
352     k=A_p->pattern->index[iptr_ik];
353     if (coloring[k]<color) Paso_BlockOps_SMV_N(n_block, &x[n_block*i], &A_p->val[block_size*iptr_ik], &x[n_block*k]);
354     }
355     Paso_BlockOps_Solve_N(n_block, &x[n_block*i], &diag[block_size*i], &pivot[n_block*i], &x[n_block*i]);
356     }
357     }
358     }
359 gross 3094 } /* end of coloring loop */
360    
361     /* backward substitution */
362     for (color=(num_colors)-2 ;color>-1;--color) { /* Note: color=(num_colors)-1 is not required */
363     if (n_block==1) {
364 gross 3097 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
365 gross 3094 for (i = 0; i < n; ++i) {
366     if (coloring[i]==color) {
367     mm=ptr_main[i];
368 gross 3097 Paso_BlockOps_MV_1(&x[i], &A_p->val[mm], &x[i]);
369 gross 3094 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
370     k=A_p->pattern->index[iptr_ik];
371 gross 3097 if (coloring[k]>color) Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
372 gross 3094 }
373 gross 3097 Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
374 gross 3094 }
375     }
376     } else if (n_block==2) {
377 gross 3097 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
378 gross 3094 for (i = 0; i < n; ++i) {
379     if (coloring[i]==color) {
380     mm=ptr_main[i];
381 gross 3097 Paso_BlockOps_MV_2(&x[2*i], &A_p->val[4*mm], &x[2*i]);
382 gross 3094 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
383     k=A_p->pattern->index[iptr_ik];
384 gross 3097 if (coloring[k]>color) Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
385 gross 3094 }
386 gross 3097 Paso_BlockOps_MV_2(&x[2*i], &diag[4*i], &x[2*i]);
387 gross 3094 }
388     }
389     } else if (n_block==3) {
390 gross 3097 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
391 gross 3094 for (i = 0; i < n; ++i) {
392     if (coloring[i]==color) {
393     mm=ptr_main[i];
394 gross 3097 Paso_BlockOps_MV_3(&x[3*i], &A_p->val[9*mm], &x[3*i]);
395 gross 3094 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
396     k=A_p->pattern->index[iptr_ik];
397 gross 3097 if (coloring[k]>color) Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
398 gross 3094 }
399 gross 3097 Paso_BlockOps_MV_3(&x[3*i], &diag[9*i], &x[3*i]);
400 gross 3094 }
401     }
402 gross 3097 } else {
403     #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
404     for (i = 0; i < n; ++i) {
405     if (coloring[i]==color) {
406     mm=ptr_main[i];
407     Paso_BlockOps_MV_N( n_block, &x[n_block*i], &A_p->val[block_size*mm], &x[n_block*i] );
408     for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
409     k=A_p->pattern->index[iptr_ik];
410     if (coloring[k]>color) Paso_BlockOps_SMV_N(n_block, &x[n_block*i], &A_p->val[block_size*iptr_ik], &x[n_block*k]);
411     }
412     Paso_BlockOps_Solve_N(n_block, &x[n_block*i], &diag[block_size*i], &pivot[n_block*i], &x[n_block*i]);
413     }
414     }
415     }
416 gross 3094 }
417     return;
418 jfenwick 3259 }

  ViewVC Help
Powered by ViewVC 1.1.26