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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1570 - (hide annotations)
Sat May 24 21:31:04 2008 UTC (11 years, 4 months ago) by gross
File MIME type: text/plain
File size: 12931 byte(s)
modifications to PCG to support dynamic scheduling
1 ksteube 1312
2 jgs 150 /* $Id$ */
3    
4 ksteube 1312 /*******************************************************
5     *
6     * Copyright 2003-2007 by ACceSS MNRF
7     * Copyright 2007 by University of Queensland
8     *
9     * http://esscc.uq.edu.au
10     * Primary Business: Queensland, Australia
11     * Licensed under the Open Software License version 3.0
12     * http://www.opensource.org/licenses/osl-3.0.php
13     *
14     *******************************************************/
15 gross 495
16 jgs 150 /* PCG iterations */
17    
18 gross 700 #include "SystemMatrix.h"
19     #include "Paso.h"
20 jgs 150 #include "Solver.h"
21 woo409 757
22 jgs 150 #ifdef _OPENMP
23     #include <omp.h>
24     #endif
25    
26 ksteube 1312 #ifdef PASO_MPI
27     #include <mpi.h>
28     #endif
29    
30 jgs 150 /*
31     *
32     * Purpose
33     * =======
34     *
35     * PCG solves the linear system A*x = b using the
36     * preconditioned conjugate gradient method plus a smoother
37     * A has to be symmetric.
38     *
39     * Convergence test: norm( b - A*x )< TOL.
40     * For other measures, see the above reference.
41     *
42     * Arguments
43     * =========
44     *
45     * r (input) DOUBLE PRECISION array, dimension N.
46     * On entry, residual of inital guess x
47     *
48     * x (input/output) DOUBLE PRECISION array, dimension N.
49     * On input, the initial guess.
50     *
51     * ITER (input/output) INT
52     * On input, the maximum iterations to be performed.
53     * On output, actual number of iterations performed.
54     *
55     * INFO (output) INT
56     *
57     * = SOLVER_NO_ERROR: Successful exit. Iterated approximate solution returned.
58     * = SOLVEr_MAXITER_REACHED
59     * = SOLVER_INPUT_ERROR Illegal parameter:
60     * = SOLVEr_BREAKDOWN: If parameters rHO or OMEGA become smaller
61     * = SOLVER_MEMORY_ERROR : If parameters rHO or OMEGA become smaller
62     *
63     * ==============================================================
64     */
65    
66 gross 1570 /* #define PASO_DYNAMIC_SCHEDULING_MVM */
67    
68     #if defined PASO_DYNAMIC_SCHEDULING_MVM && defined __OPENMP
69     #define USE_DYNAMIC_SCHEDULING
70     #endif
71    
72 jgs 150 err_t Paso_Solver_PCG(
73     Paso_SystemMatrix * A,
74     double * r,
75     double * x,
76     dim_t *iter,
77 gross 584 double * tolerance,
78     Paso_Performance* pp) {
79 jgs 150
80     /* Local variables */
81 gross 1570 dim_t num_iter=0,maxit,num_iter_global, chunk_size=-1, len,rest, n_chunks, np, ipp;
82     register double ss,ss1;
83     dim_t i0, istart, iend;
84 jgs 150 bool_t breakFlag=FALSE, maxIterFlag=FALSE, convergeFlag=FALSE;
85     err_t status = SOLVER_NO_ERROR;
86 ksteube 1312 dim_t n = Paso_SystemMatrix_getTotalNumRows(A);
87 jgs 150 double *resid = tolerance, *rs=NULL, *p=NULL, *v=NULL, *x2=NULL ;
88     double tau_old,tau,beta,delta,gamma_1,gamma_2,alpha,sum_1,sum_2,sum_3,sum_4,sum_5,tol;
89 ksteube 1312 double norm_of_residual,norm_of_residual_global, loc_sum[2], sum[2];
90 gross 495 register double r_tmp,d,rs_tmp,x2_tmp,x_tmp;
91 gross 1570 char* chksz_chr;
92     np=omp_get_max_threads();
93 jgs 150
94 gross 1570 #ifdef USE_DYNAMIC_SCHEDULING
95     chksz_chr=getenv("PASO_CHUNK_SIZE_PCG");
96     if (chksz_chr!=NULL) sscanf(chksz_chr, "%d",&chunk_size);
97     chunk_size=MIN(MAX(1,chunk_size),n/np);
98     n_chunks=n/chunk_size;
99     if (n_chunks*chunk_size<n) n_chunks+=1;
100     #else
101     len=n/np;
102     rest=n-len*np;
103     #endif
104 jgs 150 /* */
105     /*-----------------------------------------------------------------*/
106     /* */
107     /* Start of Calculation : */
108     /* --------------------- */
109     /* */
110     /* */
111     rs=TMPMEMALLOC(n,double);
112     p=TMPMEMALLOC(n,double);
113     v=TMPMEMALLOC(n,double);
114     x2=TMPMEMALLOC(n,double);
115    
116     /* Test the input parameters. */
117    
118     if (n < 0) {
119     status = SOLVER_INPUT_ERROR;
120     } else if (rs==NULL || p==NULL || v==NULL || x2==NULL) {
121     status = SOLVER_MEMORY_ERROR;
122     } else {
123     maxit = *iter;
124     tol = *resid;
125 gross 1556 Performance_startMonitor(pp,PERFORMANCE_SOLVER);
126     /* initialize data */
127 gross 1570 #pragma omp parallel (i0, istart, iend, ipp)
128 jgs 150 {
129 gross 1570 #ifdef USE_DYNAMIC_SCHEDULING
130     #pragma omp for schedule(dynamic, 1)
131     for (ipp=0; ipp < n_chunks; ++ipp) {
132     istart=chunk_size*ipp;
133     iend=MIN(istart+chunk_size,n);
134     #else
135     #pragma omp for schedule(static)
136     for (ipp=0; ipp <np; ++ipp) {
137     istart=len*ipp+MIN(ipp,rest);
138     iend=len*(ipp+1)+MIN(ipp+1,rest);
139     #endif
140     #pragma ivdep
141     for (i0=istart;i0<iend;i0++) {
142     rs[i0]=r[i0];
143     x2[i0]=x[i0];
144     p[i0]=0;
145     v[i0]=0;
146 gross 1556 }
147 gross 1570 #ifdef USE_DYNAMIC_SCHEDULING
148     }
149     #else
150     }
151     #endif
152 gross 1556 }
153     num_iter=0;
154     /* start of iteration */
155     while (!(convergeFlag || maxIterFlag || breakFlag)) {
156 jgs 150 ++(num_iter);
157     /* v=prec(r) */
158 gross 584 Performance_stopMonitor(pp,PERFORMANCE_SOLVER);
159     Performance_startMonitor(pp,PERFORMANCE_PRECONDITIONER);
160 jgs 150 Paso_Solver_solvePreconditioner(A,v,r);
161 gross 584 Performance_stopMonitor(pp,PERFORMANCE_PRECONDITIONER);
162     Performance_startMonitor(pp,PERFORMANCE_SOLVER);
163 jgs 150 /* tau=v*r */
164 gross 1570 sum_1 = 0;
165     #pragma omp parallel (i0, istart, iend, ipp, ss)
166     {
167     #ifdef USE_DYNAMIC_SCHEDULING
168     #pragma omp for schedule(dynamic, 1)
169     for (ipp=0; ipp < n_chunks; ++ipp) {
170     istart=chunk_size*ipp;
171     iend=MIN(istart+chunk_size,n);
172     #else
173     #pragma omp for schedule(static)
174     for (ipp=0; ipp <np; ++ipp) {
175     istart=len*ipp+MIN(ipp,rest);
176     iend=len*(ipp+1)+MIN(ipp+1,rest);
177     #endif
178     ss=0;
179     #pragma ivdep
180     for (i0=istart;i0<iend;i0++) ss+=v[i0]*r[i0];
181     #pragma critical
182     sum_1+=ss;
183     #ifdef USE_DYNAMIC_SCHEDULING
184     }
185     #else
186     }
187     #endif
188     }
189 ksteube 1312 #ifdef PASO_MPI
190     /* In case we have many MPI processes, each of which may have several OMP threads:
191     OMP master participates in an MPI reduction to get global sum_1 */
192 gross 1570 loc_saum[0] = sum_1;
193 gross 1556 MPI_Allreduce(loc_sum, &sum_1, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
194 ksteube 1312 #endif
195 jgs 150 tau_old=tau;
196     tau=sum_1;
197     /* p=v+beta*p */
198 gross 1570 #pragma omp parallel (i0, istart, iend, ipp,beta)
199     {
200     #ifdef USE_DYNAMIC_SCHEDULING
201     #pragma omp for schedule(dynamic, 1)
202     for (ipp=0; ipp < n_chunks; ++ipp) {
203     istart=chunk_size*ipp;
204     iend=MIN(istart+chunk_size,n);
205     #else
206     #pragma omp for schedule(static)
207     for (ipp=0; ipp <np; ++ipp) {
208     istart=len*ipp+MIN(ipp,rest);
209     iend=len*(ipp+1)+MIN(ipp+1,rest);
210     #endif
211     if (num_iter==1) {
212     #pragma ivdep
213     for (i0=istart;i0<iend;i0++) p[i0]=v[i0];
214     } else {
215     beta=tau/tau_old;
216     #pragma ivdep
217     for (i0=istart;i0<iend;i0++) p[i0]=v[i0]+beta*p[i0];
218     }
219     #ifdef USE_DYNAMIC_SCHEDULING
220     }
221     #else
222     }
223     #endif
224 jgs 150 }
225     /* v=A*p */
226 gross 584 Performance_stopMonitor(pp,PERFORMANCE_SOLVER);
227     Performance_startMonitor(pp,PERFORMANCE_MVM);
228 gross 415 Paso_SystemMatrix_MatrixVector_CSR_OFFSET0(ONE, A, p,ZERO,v);
229 gross 584 Performance_stopMonitor(pp,PERFORMANCE_MVM);
230     Performance_startMonitor(pp,PERFORMANCE_SOLVER);
231 gross 1570
232 jgs 150 /* delta=p*v */
233 gross 1570 sum_2 = 0;
234     #pragma omp parallel (i0, istart, iend, ipp,ss)
235     {
236     #ifdef USE_DYNAMIC_SCHEDULING
237     #pragma omp for schedule(dynamic, 1)
238     for (ipp=0; ipp < n_chunks; ++ipp) {
239     istart=chunk_size*ipp;
240     iend=MIN(istart+chunk_size,n);
241     #else
242     #pragma omp for schedule(static)
243     for (ipp=0; ipp <np; ++ipp) {
244     istart=len*ipp+MIN(ipp,rest);
245     iend=len*(ipp+1)+MIN(ipp+1,rest);
246     #endif
247     ss=0;
248     #pragma ivdep
249     for (i0=istart;i0<iend;i0++) ss+=v[i0]*p[i0];
250     #pragma critical
251     sum_2+=ss;
252     #ifdef USE_DYNAMIC_SCHEDULING
253     }
254     #else
255     }
256     #endif
257     }
258 ksteube 1312 #ifdef PASO_MPI
259 gross 1556 loc_sum[0] = sum_2;
260     MPI_Allreduce(loc_sum, &sum_2, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
261 ksteube 1312 #endif
262 jgs 150 delta=sum_2;
263    
264    
265     if (! (breakFlag = (ABS(delta) <= TOLERANCE_FOR_SCALARS))) {
266     /* smoother */
267 gross 1570 sum_3 = 0;
268     sum_4 = 0;
269     #pragma omp parallel (i0, istart, iend, ipp,d, ss, ss1, alpha)
270 gross 1556 {
271 gross 1570 #ifdef USE_DYNAMIC_SCHEDULING
272     #pragma omp for schedule(dynamic, 1)
273     for (ipp=0; ipp < n_chunks; ++ipp) {
274     istart=chunk_size*ipp;
275     iend=MIN(istart+chunk_size,n);
276     #else
277     #pragma omp for schedule(static)
278     for (ipp=0; ipp <np; ++ipp) {
279     istart=len*ipp+MIN(ipp,rest);
280     iend=len*(ipp+1)+MIN(ipp+1,rest);
281     #endif
282     ss=0;
283     ss1=0;
284     alpha=tau/delta;
285     #pragma ivdep
286     for (i0=istart;i0<iend;i0++) {
287     r[i0]-=alpha*v[i0];
288 gross 1556 d=r[i0]-rs[i0];
289 gross 1570 ss+=d*d;
290     ss1+=d*rs[i0];
291 gross 1556 }
292 gross 1570 #pragma critical
293     {
294     sum_3+=ss;
295     sum_4+=ss1;
296     }
297     #ifdef USE_DYNAMIC_SCHEDULING
298     }
299     #else
300     }
301     #endif
302 ksteube 1312 }
303     #ifdef PASO_MPI
304 gross 1556 loc_sum[0] = sum_3;
305     loc_sum[1] = sum_4;
306     MPI_Allreduce(loc_sum, sum, 2, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
307     sum_3=sum[0];
308     sum_4=sum[1];
309 ksteube 1312 #endif
310 gross 1570 sum_5 = 0;
311     #pragma omp parallel (i0, istart, iend, ipp, ss, gamma_1,gamma_2)
312 gross 1556 {
313 gross 1570 gamma_1= ( (ABS(sum_3)<= ZERO) ? 0 : -sum_4/sum_3) ;
314     gamma_2= ONE-gamma_1;
315     #ifdef USE_DYNAMIC_SCHEDULING
316     #pragma omp for schedule(dynamic, 1)
317     for (ipp=0; ipp < n_chunks; ++ipp) {
318     istart=chunk_size*ipp;
319     iend=MIN(istart+chunk_size,n);
320     #else
321     #pragma omp for schedule(static)
322     for (ipp=0; ipp <np; ++ipp) {
323     istart=len*ipp+MIN(ipp,rest);
324     iend=len*(ipp+1)+MIN(ipp+1,rest);
325     #endif
326     ss=0;
327     #pragma ivdep
328     for (i0=istart;i0<iend;i0++) {
329 gross 1556 rs[i0]=gamma_2*rs[i0]+gamma_1*r[i0];
330     x2[i0]+=alpha*p[i0];
331     x[i0]=gamma_2*x[i0]+gamma_1*x2[i0];
332 gross 1570 ss+=rs[i0]*rs[i0];
333     }
334     #pragma omp critical
335     sum_5+=ss;
336     #ifdef USE_DYNAMIC_SCHEDULING
337 gross 1556 }
338 gross 1570 #else
339     }
340     #endif
341 gross 495 }
342 ksteube 1312 #ifdef PASO_MPI
343 gross 1556 loc_sum[0] = sum_5;
344     MPI_Allreduce(loc_sum, &sum_5, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
345 ksteube 1312 #endif
346 gross 495 norm_of_residual=sqrt(sum_5);
347     convergeFlag = norm_of_residual <= tol;
348     maxIterFlag = num_iter == maxit;
349     breakFlag = (ABS(tau) <= TOLERANCE_FOR_SCALARS);
350     }
351 gross 1556 }
352     /* end of iteration */
353     num_iter_global=num_iter;
354     norm_of_residual_global=norm_of_residual;
355     if (maxIterFlag) {
356     status = SOLVER_MAXITER_REACHED;
357     } else if (breakFlag) {
358     status = SOLVER_BREAKDOWN;
359     }
360     Performance_stopMonitor(pp,PERFORMANCE_SOLVER);
361 jgs 150 TMPMEMFREE(rs);
362     TMPMEMFREE(x2);
363     TMPMEMFREE(v);
364     TMPMEMFREE(p);
365     *iter=num_iter_global;
366     *resid=norm_of_residual_global;
367     }
368     /* End of PCG */
369     return status;
370     }

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26