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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3259 - (show annotations)
Mon Oct 11 01:48:14 2010 UTC (9 years, 1 month ago) by jfenwick
File MIME type: text/plain
File size: 14858 byte(s)
Merging dudley and scons updates from branches

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

  ViewVC Help
Powered by ViewVC 1.1.26