source: CIVL/examples/omp/poisson_openmp.c

main
Last change on this file was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 12.4 KB
Line 
1# include <stdlib.h>
2# include <stdio.h>
3# include <math.h>
4# include <time.h>
5# include <omp.h>
6
7# define NX 10
8# define NY 10
9
10int main ( int argc, char *argv[] );
11double r8mat_rms ( int nx, int ny, double a[NX][NY] );
12void rhs ( int nx, int ny, double f[NX][NY] );
13void sweep ( int nx, int ny, double dx, double dy, double f[NX][NY],
14 int itold, int itnew, double u[NX][NY], double unew[NX][NY] );
15void timestamp ( void );
16double u_exact ( double x, double y );
17double uxxyy_exact ( double x, double y );
18
19/******************************************************************************/
20
21int main ( int argc, char *argv[] )
22
23/******************************************************************************/
24/*
25 Purpose:
26
27 MAIN is the main program for POISSON_OPENMP.
28
29 Discussion:
30
31 POISSON_OPENMP is a program for solving the Poisson problem.
32
33 This program uses OpenMP for parallel execution.
34
35 The Poisson equation
36
37 - DEL^2 U(X,Y) = F(X,Y)
38
39 is solved on the unit square [0,1] x [0,1] using a grid of NX by
40 NX evenly spaced points. The first and last points in each direction
41 are boundary points.
42
43 The boundary conditions and F are set so that the exact solution is
44
45 U(x,y) = sin ( pi * x * y )
46
47 so that
48
49 - DEL^2 U(x,y) = pi^2 * ( x^2 + y^2 ) * sin ( pi * x * y )
50
51 The Jacobi iteration is repeatedly applied until convergence is detected.
52
53 For convenience in writing the discretized equations, we assume that NX = NY.
54
55 Licensing:
56
57 This code is distributed under the GNU LGPL license.
58
59 Modified:
60
61 14 December 2011
62
63 Author:
64
65 John Burkardt
66*/
67{
68 int converged;
69 double diff;
70 double dx;
71 double dy;
72 double error;
73 double f[NX][NY];
74 int i;
75 int id;
76 int itnew;
77 int itold;
78 int j;
79 int jt;
80 int jt_max = 20;
81 int nx = NX;
82 int ny = NY;
83 double tolerance = 0.000001;
84 double u[NX][NY];
85 double u_norm;
86 double udiff[NX][NY];
87 double uexact[NX][NY];
88 double unew[NX][NY];
89 double unew_norm;
90 double wtime;
91 double x;
92 double y;
93
94 dx = 1.0 / ( double ) ( nx - 1 );
95 dy = 1.0 / ( double ) ( ny - 1 );
96/*
97 Print a message.
98*/
99 timestamp ( );
100 printf ( "\n" );
101 printf ( "POISSON_OPENMP:\n" );
102 printf ( " C version\n" );
103 printf ( " A program for solving the Poisson equation.\n" );
104 printf ( "\n" );
105 printf ( " Use OpenMP for parallel execution.\n" );
106 printf ( " The number of processors is %d\n", omp_get_num_procs ( ) );
107# pragma omp parallel
108{
109 id = omp_get_thread_num ( );
110 if ( id == 0 )
111 {
112 printf ( " The maximum number of threads is %d\n", omp_get_num_threads ( ) );
113 }
114}
115 printf ( "\n" );
116 printf ( " -DEL^2 U = F(X,Y)\n" );
117 printf ( "\n" );
118 printf ( " on the rectangle 0 <= X <= 1, 0 <= Y <= 1.\n" );
119 printf ( "\n" );
120 printf ( " F(X,Y) = pi^2 * ( x^2 + y^2 ) * sin ( pi * x * y )\n" );
121 printf ( "\n" );
122 printf ( " The number of interior X grid points is %d\n", nx );
123 printf ( " The number of interior Y grid points is %d\n", ny );
124 printf ( " The X grid spacing is %f\n", dx );
125 printf ( " The Y grid spacing is %f\n", dy );
126/*
127 Set the right hand side array F.
128*/
129 rhs ( nx, ny, f );
130/*
131 Set the initial solution estimate UNEW.
132 We are "allowed" to pick up the boundary conditions exactly.
133*/
134 for ( j = 0; j < ny; j++ )
135 {
136 for ( i = 0; i < nx; i++ )
137 {
138 if ( i == 0 || i == nx - 1 || j == 0 || j == ny - 1 )
139 {
140 unew[i][j] = f[i][j];
141 }
142 else
143 {
144 unew[i][j] = 0.0;
145 }
146 }
147 }
148 unew_norm = r8mat_rms ( nx, ny, unew );
149/*
150 Set up the exact solution UEXACT.
151*/
152 for ( j = 0; j < ny; j++ )
153 {
154 y = ( double ) ( j ) / ( double ) ( ny - 1 );
155 for ( i = 0; i < nx; i++ )
156 {
157 x = ( double ) ( i ) / ( double ) ( nx - 1 );
158 uexact[i][j] = u_exact ( x, y );
159 }
160 }
161 u_norm = r8mat_rms ( nx, ny, uexact );
162 printf ( " RMS of exact solution = %g\n", u_norm );
163/*
164 Do the iteration.
165*/
166 converged = 0;
167
168 printf ( "\n" );
169 printf ( " Step ||Unew|| ||Unew-U|| ||Unew-Exact||\n" );
170 printf ( "\n" );
171
172 for ( j = 0; j < ny; j++ )
173 {
174 for ( i = 0; i < nx; i++ )
175 {
176 udiff[i][j] = unew[i][j] - uexact[i][j];
177 }
178 }
179 error = r8mat_rms ( nx, ny, udiff );
180 printf ( " %4d %14g %14g\n", 0, unew_norm, error );
181
182 wtime = omp_get_wtime ( );
183
184 itnew = 0;
185
186 for ( ; ; )
187 {
188 itold = itnew;
189 itnew = itold + 500;
190/*
191 SWEEP carries out 500 Jacobi steps in parallel before we come
192 back to check for convergence.
193*/
194 sweep ( nx, ny, dx, dy, f, itold, itnew, u, unew );
195/*
196 Check for convergence.
197*/
198 u_norm = unew_norm;
199 unew_norm = r8mat_rms ( nx, ny, unew );
200
201 for ( j = 0; j < ny; j++ )
202 {
203 for ( i = 0; i < nx; i++ )
204 {
205 udiff[i][j] = unew[i][j] - u[i][j];
206 }
207 }
208 diff = r8mat_rms ( nx, ny, udiff );
209
210 for ( j = 0; j < ny; j++ )
211 {
212 for ( i = 0; i < nx; i++ )
213 {
214 udiff[i][j] = unew[i][j] - uexact[i][j];
215 }
216 }
217 error = r8mat_rms ( nx, ny, udiff );
218
219 printf ( " %4d %14g %14g %14g\n", itnew, unew_norm, diff, error );
220
221 if ( diff <= tolerance )
222 {
223 converged = 1;
224 break;
225 }
226
227 }
228
229 if ( converged )
230 {
231 printf ( " The iteration has converged.\n" );
232 }
233 else
234 {
235 printf ( " The iteration has NOT converged.\n" );
236 }
237
238 wtime = omp_get_wtime ( ) - wtime;
239 printf ( "\n" );
240 printf ( " Elapsed seconds = %g\n", wtime );
241/*
242 Terminate.
243*/
244 printf ( "\n" );
245 printf ( "POISSON_OPENMP:\n" );
246 printf ( " Normal end of execution.\n" );
247 printf ( "\n" );
248 timestamp ( );
249
250 return 0;
251}
252/******************************************************************************/
253
254double r8mat_rms ( int nx, int ny, double a[NX][NY] )
255
256/******************************************************************************/
257/*
258 Purpose:
259
260 R8MAT_RMS returns the RMS norm of a vector stored as a matrix.
261
262 Licensing:
263
264 This code is distributed under the GNU LGPL license.
265
266 Modified:
267
268 01 March 2003
269
270 Author:
271
272 John Burkardt
273
274 Parameters:
275
276 Input, int NX, NY, the number of rows and columns in A.
277
278 Input, double A[NX][NY], the vector.
279
280 Output, double R8MAT_RMS, the root mean square of the entries of A.
281*/
282{
283 int i;
284 int j;
285 double v;
286
287 v = 0.0;
288
289 for ( j = 0; j < ny; j++ )
290 {
291 for ( i = 0; i < nx; i++ )
292 {
293 v = v + a[i][j] * a[i][j];
294 }
295 }
296
297 v = sqrt ( v / ( double ) ( nx * ny ) );
298
299 return v;
300}
301/******************************************************************************/
302
303void rhs ( int nx, int ny, double f[NX][NY] )
304
305/******************************************************************************/
306/*
307 Purpose:
308
309 RHS initializes the right hand side "vector".
310
311 Discussion:
312
313 It is convenient for us to set up RHS as a 2D array. However, each
314 entry of RHS is really the right hand side of a linear system of the
315 form
316
317 A * U = F
318
319 In cases where U(I,J) is a boundary value, then the equation is simply
320
321 U(I,J) = F(i,j)
322
323 and F(I,J) holds the boundary data.
324
325 Otherwise, the equation has the form
326
327 (1/DX^2) * ( U(I+1,J)+U(I-1,J)+U(I,J-1)+U(I,J+1)-4*U(I,J) ) = F(I,J)
328
329 where DX is the spacing and F(I,J) is the value at X(I), Y(J) of
330
331 pi^2 * ( x^2 + y^2 ) * sin ( pi * x * y )
332
333 Licensing:
334
335 This code is distributed under the GNU LGPL license.
336
337 Modified:
338
339 28 October 2011
340
341 Author:
342
343 John Burkardt
344
345 Parameters:
346
347 Input, int NX, NY, the X and Y grid dimensions.
348
349 Output, double F[NX][NY], the initialized right hand side data.
350*/
351{
352 double fnorm;
353 int i;
354 int j;
355 double x;
356 double y;
357/*
358 The "boundary" entries of F store the boundary values of the solution.
359 The "interior" entries of F store the right hand sides of the Poisson equation.
360*/
361 for ( j = 0; j < ny; j++ )
362 {
363 y = ( double ) ( j ) / ( double ) ( ny - 1 );
364 for ( i = 0; i < nx; i++ )
365 {
366 x = ( double ) ( i ) / ( double ) ( nx - 1 );
367 if ( i == 0 || i == nx - 1 || j == 0 || j == ny - 1 )
368 {
369 f[i][j] = u_exact ( x, y );
370 }
371 else
372 {
373 f[i][j] = - uxxyy_exact ( x, y );
374 }
375 }
376 }
377
378 fnorm = r8mat_rms ( nx, ny, f );
379
380 printf ( " RMS of F = %g\n", fnorm );
381
382 return;
383}
384/******************************************************************************/
385
386void sweep ( int nx, int ny, double dx, double dy, double f[NX][NY],
387 int itold, int itnew, double u[NX][NY], double unew[NX][NY] )
388
389/******************************************************************************/
390/*
391 Purpose:
392
393 SWEEP carries out one step of the Jacobi iteration.
394
395 Discussion:
396
397 Assuming DX = DY, we can approximate
398
399 - ( d/dx d/dx + d/dy d/dy ) U(X,Y)
400
401 by
402
403 ( U(i-1,j) + U(i+1,j) + U(i,j-1) + U(i,j+1) - 4*U(i,j) ) / dx / dy
404
405 The discretization employed below will not be correct in the general
406 case where DX and DY are not equal. It's only a little more complicated
407 to allow DX and DY to be different, but we're not going to worry about
408 that right now.
409
410 Licensing:
411
412 This code is distributed under the GNU LGPL license.
413
414 Modified:
415
416 14 December 2011
417
418 Author:
419
420 John Burkardt
421
422 Parameters:
423
424 Input, int NX, NY, the X and Y grid dimensions.
425
426 Input, double DX, DY, the spacing between grid points.
427
428 Input, double F[NX][NY], the right hand side data.
429
430 Input, int ITOLD, the iteration index on input.
431
432 Input, int ITNEW, the desired iteration index
433 on output.
434
435 Input, double U[NX][NY], the solution estimate on
436 iteration ITNEW-1.
437
438 Input/output, double UNEW[NX][NY], on input, the solution
439 estimate on iteration ITOLD. On output, the solution estimate on
440 iteration ITNEW.
441*/
442{
443 int i;
444 int it;
445 int j;
446
447# pragma omp parallel shared ( dx, dy, f, itnew, itold, nx, ny, u, unew ) private ( i, it, j )
448
449 for ( it = itold + 1; it <= itnew; it++ )
450 {
451/*
452 Save the current estimate.
453*/
454# pragma omp for
455 for ( j = 0; j < ny; j++ )
456 {
457 for ( i = 0; i < nx; i++ )
458 {
459 u[i][j] = unew[i][j];
460 }
461 }
462/*
463 Compute a new estimate.
464*/
465# pragma omp for
466 for ( j = 0; j < ny; j++ )
467 {
468 for ( i = 0; i < nx; i++ )
469 {
470 if ( i == 0 || j == 0 || i == nx - 1 || j == ny - 1 )
471 {
472 unew[i][j] = f[i][j];
473 }
474 else
475 {
476 unew[i][j] = 0.25 * (
477 u[i-1][j] + u[i][j+1] + u[i][j-1] + u[i+1][j] + f[i][j] * dx * dy );
478 }
479 }
480 }
481
482 }
483 return;
484}
485/******************************************************************************/
486
487void timestamp ( void )
488
489/******************************************************************************/
490/*
491 Purpose:
492
493 TIMESTAMP prints the current YMDHMS date as a time stamp.
494
495 Example:
496
497 31 May 2001 09:45:54 AM
498
499 Licensing:
500
501 This code is distributed under the GNU LGPL license.
502
503 Modified:
504
505 24 September 2003
506
507 Author:
508
509 John Burkardt
510
511 Parameters:
512
513 None
514*/
515{
516# define TIME_SIZE 40
517
518 static char time_buffer[TIME_SIZE];
519 const struct tm *tm;
520 time_t now;
521
522 now = time ( NULL );
523 tm = localtime ( &now );
524
525 strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm );
526
527 printf ( "%s\n", time_buffer );
528
529 return;
530# undef TIME_SIZE
531}
532/******************************************************************************/
533
534double u_exact ( double x, double y )
535
536/******************************************************************************/
537/*
538 Purpose:
539
540 U_EXACT evaluates the exact solution.
541
542 Licensing:
543
544 This code is distributed under the GNU LGPL license.
545
546 Modified:
547
548 25 October 2011
549
550 Author:
551
552 John Burkardt
553
554 Parameters:
555
556 Input, double X, Y, the coordinates of a point.
557
558 Output, double U_EXACT, the value of the exact solution
559 at (X,Y).
560*/
561{
562 double pi = 3.141592653589793;
563 double value;
564
565 value = sin ( pi * x * y );
566
567 return value;
568}
569/******************************************************************************/
570
571double uxxyy_exact ( double x, double y )
572
573/******************************************************************************/
574/*
575 Purpose:
576
577 UXXYY_EXACT evaluates ( d/dx d/dx + d/dy d/dy ) of the exact solution.
578
579 Licensing:
580
581 This code is distributed under the GNU LGPL license.
582
583 Modified:
584
585 25 October 2011
586
587 Author:
588
589 John Burkardt
590
591 Parameters:
592
593 Input, double X, Y, the coordinates of a point.
594
595 Output, double UXXYY_EXACT, the value of
596 ( d/dx d/dx + d/dy d/dy ) of the exact solution at (X,Y).
597*/
598{
599 double pi = 3.141592653589793;
600 double value;
601
602 value = - pi * pi * ( x * x + y * y ) * sin ( pi * x * y );
603
604 return value;
605}
606# undef NX
607# undef NY
Note: See TracBrowser for help on using the repository browser.