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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3158 - (show annotations)
Mon Sep 6 06:09:11 2010 UTC (8 years, 4 months ago) by gross
File MIME type: text/plain
File size: 14710 byte(s)
GS and Jacobi are now used through the same interface.
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 (! Paso_checkPtr(out)) {
62 out->localSmoother=Paso_Preconditioner_LocalSmoother_alloc(A_p->mainBlock,jacobi,verbose);
63 out->is_local=is_local;
64 }
65 if (Paso_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=Paso_timer();
80 /* allocations: */
81 Paso_Preconditioner_LocalSmoother* out=MEMALLOC(1,Paso_Preconditioner_LocalSmoother);
82 if (! Paso_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 ( ! ( Paso_checkPtr(out->diag) || Paso_checkPtr(out->pivot) ) ) {
90 Paso_SparseMatrix_invMain(A_p, out->diag, out->pivot );
91 }
92
93 }
94 time0=Paso_timer()-time0;
95
96 if (Paso_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, const dim_t sweeps)
125 {
126 const dim_t n= (A_p->mainBlock->numRows) * (A_p->mainBlock->row_block_size);
127
128 double *b_new = smoother->localSmoother->buffer;
129 dim_t nsweeps=sweeps;
130 if (smoother->is_local) {
131 Paso_Preconditioner_LocalSmoother_solve(A_p->mainBlock,smoother->localSmoother,x,b,sweeps);
132 } else {
133 Paso_Copy(n, x, b);
134
135 Paso_Preconditioner_LocalSmoother_Sweep(A_p->mainBlock,smoother->localSmoother,x);
136
137 while (nsweeps > 1 ) {
138 Paso_Copy(n, b_new, b);
139
140 Paso_SystemMatrix_MatrixVector((-1.), A_p, x, 1., b_new); /* b_new = b - A*x */
141
142 Paso_Preconditioner_LocalSmoother_Sweep(A_p->mainBlock,smoother->localSmoother,b_new);
143
144 Paso_AXPY(n, x, 1., b_new);
145 nsweeps--;
146 }
147
148 }
149 }
150 void Paso_Preconditioner_LocalSmoother_solve(Paso_SparseMatrix* A_p, Paso_Preconditioner_LocalSmoother * smoother, double * x, const double * b, const dim_t sweeps)
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 Paso_Copy(n, x, b);
157 Paso_Preconditioner_LocalSmoother_Sweep(A_p,smoother,x);
158
159 while (nsweeps > 1 ) {
160
161 Paso_Copy(n, b_new, b);
162
163 Paso_SparseMatrix_MatrixVector_CSC_OFFSET0((-1.), A_p, x, 1., b_new); /* b_new = b - A*x */
164
165 Paso_Preconditioner_LocalSmoother_Sweep(A_p,smoother,b_new);
166
167 Paso_AXPY(n, x, 1., b_new)
168
169 nsweeps--;
170 }
171 }
172
173 void Paso_Preconditioner_LocalSmoother_Sweep(Paso_SparseMatrix* A, Paso_Preconditioner_LocalSmoother * smoother, double * x)
174 {
175 #ifdef _OPENMP
176 const dim_t nt=omp_get_max_threads();
177 #else
178 const dim_t nt = 1;
179 #endif
180 if (smoother->Jacobi) {
181 Paso_BlockOps_allMV(A->row_block_size,A->numRows,smoother->diag,smoother->pivot,x);
182 } else {
183 if (nt < 2) {
184 Paso_Preconditioner_LocalSmoother_Sweep_sequential(A,smoother,x);
185 } else {
186 Paso_Preconditioner_LocalSmoother_Sweep_colored(A,smoother,x);
187 }
188 }
189 }
190
191 /* inplace Gauss-Seidel sweep in seqential mode: */
192
193 void Paso_Preconditioner_LocalSmoother_Sweep_sequential(Paso_SparseMatrix* A_p, Paso_Preconditioner_LocalSmoother * smoother, double * x)
194 {
195 const dim_t n=A_p->numRows;
196 const dim_t n_block=A_p->row_block_size;
197 const double *diag = smoother->diag;
198 /* const index_t* pivot = smoother->pivot;
199 const dim_t block_size=A_p->block_size; use for block size >3*/
200
201 register dim_t i,k;
202 register index_t iptr_ik, mm;
203
204 const index_t* ptr_main = Paso_SparseMatrix_borrowMainDiagonalPointer(A_p);
205 /* forward substitution */
206
207 if (n_block==1) {
208 Paso_BlockOps_MV_1(&x[0], &diag[0], &x[0]);
209 for (i = 1; i < n; ++i) {
210 mm=ptr_main[i];
211 /* x_i=x_i-a_ik*x_k (with k<i) */
212 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<mm; ++iptr_ik) {
213 k=A_p->pattern->index[iptr_ik];
214 Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
215 }
216 Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
217 }
218 } else if (n_block==2) {
219 Paso_BlockOps_MV_2(&x[0], &diag[0], &x[0]);
220 for (i = 1; i < n; ++i) {
221 mm=ptr_main[i];
222 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<mm; ++iptr_ik) {
223 k=A_p->pattern->index[iptr_ik];
224 Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
225 }
226 Paso_BlockOps_MV_2(&x[2*i], &diag[4*i], &x[2*i]);
227 }
228 } else if (n_block==3) {
229 Paso_BlockOps_MV_3(&x[0], &diag[0], &x[0]);
230 for (i = 1; i < n; ++i) {
231 mm=ptr_main[i];
232 /* x_i=x_i-a_ik*x_k */
233 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<mm; ++iptr_ik) {
234 k=A_p->pattern->index[iptr_ik];
235 Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
236 }
237 Paso_BlockOps_MV_3(&x[3*i], &diag[9*i], &x[3*i]);
238 }
239 } /* add block size >3 */
240
241 /* backward substitution */
242
243 if (n_block==1) {
244 for (i = n-2; i > -1; --i) {
245 mm=ptr_main[i];
246 Paso_BlockOps_MV_1(&x[i], &A_p->val[mm], &x[i]);
247 for (iptr_ik=mm+1; iptr_ik < A_p->pattern->ptr[i+1]; ++iptr_ik) {
248 k=A_p->pattern->index[iptr_ik];
249 Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
250 }
251 Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
252 }
253
254 } else if (n_block==2) {
255 for (i = n-2; i > -1; --i) {
256 mm=ptr_main[i];
257 Paso_BlockOps_MV_2(&x[2*i], &A_p->val[4*mm], &x[2*i]);
258 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_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
261 }
262 Paso_BlockOps_MV_2(&x[2*i], &diag[i*4], &x[2*i]);
263 }
264 } else if (n_block==3) {
265 for (i = n-2; i > -1; --i) {
266 mm=ptr_main[i];
267 Paso_BlockOps_MV_3(&x[3*i], &A_p->val[9*mm], &x[3*i]);
268
269 for (iptr_ik=mm+1; iptr_ik < A_p->pattern->ptr[i+1]; ++iptr_ik) {
270 k=A_p->pattern->index[iptr_ik];
271 Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
272 }
273 Paso_BlockOps_MV_3(&x[3*i], &diag[i*9], &x[3*i]);
274 }
275
276 } /* add block size >3 */
277
278 return;
279 }
280
281 void Paso_Preconditioner_LocalSmoother_Sweep_colored(Paso_SparseMatrix* A_p, Paso_Preconditioner_LocalSmoother * smoother, double * x)
282 {
283 const dim_t n=A_p->numRows;
284 const dim_t n_block=A_p->row_block_size;
285 const double *diag = smoother->diag;
286 index_t* pivot = smoother->pivot;
287 const dim_t block_size=A_p->block_size;
288
289 register dim_t i,k;
290 register index_t color,iptr_ik, mm;
291
292 const index_t* coloring = Paso_Pattern_borrowColoringPointer(A_p->pattern);
293 const dim_t num_colors = Paso_Pattern_getNumColors(A_p->pattern);
294 const index_t* ptr_main = Paso_SparseMatrix_borrowMainDiagonalPointer(A_p);
295
296 /* forward substitution */
297
298
299 /* color = 0 */
300 if (n_block==1) {
301 #pragma omp parallel for schedule(static) private(i)
302 for (i = 0; i < n; ++i) {
303 if (coloring[i]== 0 ) Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
304 }
305 } else if (n_block==2) {
306 #pragma omp parallel for schedule(static) private(i)
307 for (i = 0; i < n; ++i) {
308 if (coloring[i]== 0 ) Paso_BlockOps_MV_2(&x[2*i], &diag[i*4], &x[2*i]);
309 }
310 } else if (n_block==3) {
311 #pragma omp parallel for schedule(static) private(i)
312 for (i = 0; i < n; ++i) {
313 if (coloring[i]==0) Paso_BlockOps_MV_3(&x[3*i], &diag[i*9], &x[3*i]);
314 }
315 } else {
316 #pragma omp parallel for schedule(static) private(i)
317 for (i = 0; i < n; ++i) {
318 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]);
319 }
320 }
321
322 for (color=1;color<num_colors;++color) {
323
324 if (n_block==1) {
325 #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
326 for (i = 0; i < n; ++i) {
327 if (coloring[i]==color) {
328 /* x_i=x_i-a_ik*x_k */
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 if (coloring[k]<color) Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
332 }
333 Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
334 }
335 }
336 } else if (n_block==2) {
337 #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
338 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 if (coloring[k]<color) Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
343 }
344 Paso_BlockOps_MV_2(&x[2*i], &diag[4*i], &x[2*i]);
345 }
346 }
347 } else if (n_block==3) {
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_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
354 }
355 Paso_BlockOps_MV_3(&x[3*i], &diag[9*i], &x[3*i]);
356 }
357 }
358 } else {
359 #pragma omp parallel for schedule(static) private(i,iptr_ik,k)
360 for (i = 0; i < n; ++i) {
361 if (coloring[i] == color) {
362 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
363 k=A_p->pattern->index[iptr_ik];
364 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]);
365 }
366 Paso_BlockOps_Solve_N(n_block, &x[n_block*i], &diag[block_size*i], &pivot[n_block*i], &x[n_block*i]);
367 }
368 }
369 }
370 } /* end of coloring loop */
371
372 /* backward substitution */
373 for (color=(num_colors)-2 ;color>-1;--color) { /* Note: color=(num_colors)-1 is not required */
374 if (n_block==1) {
375 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
376 for (i = 0; i < n; ++i) {
377 if (coloring[i]==color) {
378 mm=ptr_main[i];
379 Paso_BlockOps_MV_1(&x[i], &A_p->val[mm], &x[i]);
380 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
381 k=A_p->pattern->index[iptr_ik];
382 if (coloring[k]>color) Paso_BlockOps_SMV_1(&x[i], &A_p->val[iptr_ik], &x[k]);
383 }
384 Paso_BlockOps_MV_1(&x[i], &diag[i], &x[i]);
385 }
386 }
387 } else if (n_block==2) {
388 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
389 for (i = 0; i < n; ++i) {
390 if (coloring[i]==color) {
391 mm=ptr_main[i];
392 Paso_BlockOps_MV_2(&x[2*i], &A_p->val[4*mm], &x[2*i]);
393 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
394 k=A_p->pattern->index[iptr_ik];
395 if (coloring[k]>color) Paso_BlockOps_SMV_2(&x[2*i], &A_p->val[4*iptr_ik], &x[2*k]);
396 }
397 Paso_BlockOps_MV_2(&x[2*i], &diag[4*i], &x[2*i]);
398 }
399 }
400 } else if (n_block==3) {
401 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
402 for (i = 0; i < n; ++i) {
403 if (coloring[i]==color) {
404 mm=ptr_main[i];
405 Paso_BlockOps_MV_3(&x[3*i], &A_p->val[9*mm], &x[3*i]);
406 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
407 k=A_p->pattern->index[iptr_ik];
408 if (coloring[k]>color) Paso_BlockOps_SMV_3(&x[3*i], &A_p->val[9*iptr_ik], &x[3*k]);
409 }
410 Paso_BlockOps_MV_3(&x[3*i], &diag[9*i], &x[3*i]);
411 }
412 }
413 } else {
414 #pragma omp parallel for schedule(static) private(mm, i,iptr_ik,k)
415 for (i = 0; i < n; ++i) {
416 if (coloring[i]==color) {
417 mm=ptr_main[i];
418 Paso_BlockOps_MV_N( n_block, &x[n_block*i], &A_p->val[block_size*mm], &x[n_block*i] );
419 for (iptr_ik=A_p->pattern->ptr[i];iptr_ik<A_p->pattern->ptr[i+1]; ++iptr_ik) {
420 k=A_p->pattern->index[iptr_ik];
421 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]);
422 }
423 Paso_BlockOps_Solve_N(n_block, &x[n_block*i], &diag[block_size*i], &pivot[n_block*i], &x[n_block*i]);
424 }
425 }
426 }
427 }
428 return;
429 }

  ViewVC Help
Powered by ViewVC 1.1.26