source: CIVL/mods/dev.civl.abc/examples/omp/poisson_openmp.c

main
Last change on this file was aad342c, checked in by Stephen Siegel <siegel@…>, 3 years ago

Performing huge refactor to incorporate ABC, GMC, and SARL into CIVL repo and use Java modules.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5664 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 161
8# define NY 161
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 v = sqrt ( v / ( double ) ( nx * ny ) );
297
298 return v;
299}
300/******************************************************************************/
301
302void rhs ( int nx, int ny, double f[NX][NY] )
303
304/******************************************************************************/
305/*
306 Purpose:
307
308 RHS initializes the right hand side "vector".
309
310 Discussion:
311
312 It is convenient for us to set up RHS as a 2D array. However, each
313 entry of RHS is really the right hand side of a linear system of the
314 form
315
316 A * U = F
317
318 In cases where U(I,J) is a boundary value, then the equation is simply
319
320 U(I,J) = F(i,j)
321
322 and F(I,J) holds the boundary data.
323
324 Otherwise, the equation has the form
325
326 (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)
327
328 where DX is the spacing and F(I,J) is the value at X(I), Y(J) of
329
330 pi^2 * ( x^2 + y^2 ) * sin ( pi * x * y )
331
332 Licensing:
333
334 This code is distributed under the GNU LGPL license.
335
336 Modified:
337
338 28 October 2011
339
340 Author:
341
342 John Burkardt
343
344 Parameters:
345
346 Input, int NX, NY, the X and Y grid dimensions.
347
348 Output, double F[NX][NY], the initialized right hand side data.
349*/
350{
351 double fnorm;
352 int i;
353 int j;
354 double x;
355 double y;
356/*
357 The "boundary" entries of F store the boundary values of the solution.
358 The "interior" entries of F store the right hand sides of the Poisson equation.
359*/
360 for ( j = 0; j < ny; j++ )
361 {
362 y = ( double ) ( j ) / ( double ) ( ny - 1 );
363 for ( i = 0; i < nx; i++ )
364 {
365 x = ( double ) ( i ) / ( double ) ( nx - 1 );
366 if ( i == 0 || i == nx - 1 || j == 0 || j == ny - 1 )
367 {
368 f[i][j] = u_exact ( x, y );
369 }
370 else
371 {
372 f[i][j] = - uxxyy_exact ( x, y );
373 }
374 }
375 }
376
377 fnorm = r8mat_rms ( nx, ny, f );
378
379 printf ( " RMS of F = %g\n", fnorm );
380
381 return;
382}
383/******************************************************************************/
384
385void sweep ( int nx, int ny, double dx, double dy, double f[NX][NY],
386 int itold, int itnew, double u[NX][NY], double unew[NX][NY] )
387
388/******************************************************************************/
389/*
390 Purpose:
391
392 SWEEP carries out one step of the Jacobi iteration.
393
394 Discussion:
395
396 Assuming DX = DY, we can approximate
397
398 - ( d/dx d/dx + d/dy d/dy ) U(X,Y)
399
400 by
401
402 ( U(i-1,j) + U(i+1,j) + U(i,j-1) + U(i,j+1) - 4*U(i,j) ) / dx / dy
403
404 The discretization employed below will not be correct in the general
405 case where DX and DY are not equal. It's only a little more complicated
406 to allow DX and DY to be different, but we're not going to worry about
407 that right now.
408
409 Licensing:
410
411 This code is distributed under the GNU LGPL license.
412
413 Modified:
414
415 14 December 2011
416
417 Author:
418
419 John Burkardt
420
421 Parameters:
422
423 Input, int NX, NY, the X and Y grid dimensions.
424
425 Input, double DX, DY, the spacing between grid points.
426
427 Input, double F[NX][NY], the right hand side data.
428
429 Input, int ITOLD, the iteration index on input.
430
431 Input, int ITNEW, the desired iteration index
432 on output.
433
434 Input, double U[NX][NY], the solution estimate on
435 iteration ITNEW-1.
436
437 Input/output, double UNEW[NX][NY], on input, the solution
438 estimate on iteration ITOLD. On output, the solution estimate on
439 iteration ITNEW.
440*/
441{
442 int i;
443 int it;
444 int j;
445
446# pragma omp parallel shared ( dx, dy, f, itnew, itold, nx, ny, u, unew ) private ( i, it, j )
447
448 for ( it = itold + 1; it <= itnew; it++ )
449 {
450/*
451 Save the current estimate.
452*/
453# pragma omp for
454 for ( j = 0; j < ny; j++ )
455 {
456 for ( i = 0; i < nx; i++ )
457 {
458 u[i][j] = unew[i][j];
459 }
460 }
461/*
462 Compute a new estimate.
463*/
464# pragma omp for
465 for ( j = 0; j < ny; j++ )
466 {
467 for ( i = 0; i < nx; i++ )
468 {
469 if ( i == 0 || j == 0 || i == nx - 1 || j == ny - 1 )
470 {
471 unew[i][j] = f[i][j];
472 }
473 else
474 {
475 unew[i][j] = 0.25 * (
476 u[i-1][j] + u[i][j+1] + u[i][j-1] + u[i+1][j] + f[i][j] * dx * dy );
477 }
478 }
479 }
480
481 }
482 return;
483}
484/******************************************************************************/
485
486void timestamp ( void )
487
488/******************************************************************************/
489/*
490 Purpose:
491
492 TIMESTAMP prints the current YMDHMS date as a time stamp.
493
494 Example:
495
496 31 May 2001 09:45:54 AM
497
498 Licensing:
499
500 This code is distributed under the GNU LGPL license.
501
502 Modified:
503
504 24 September 2003
505
506 Author:
507
508 John Burkardt
509
510 Parameters:
511
512 None
513*/
514{
515# define TIME_SIZE 40
516
517 static char time_buffer[TIME_SIZE];
518 const struct tm *tm;
519 time_t now;
520
521 now = time ( NULL );
522 tm = localtime ( &now );
523
524 strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm );
525
526 printf ( "%s\n", time_buffer );
527
528 return;
529# undef TIME_SIZE
530}
531/******************************************************************************/
532
533double u_exact ( double x, double y )
534
535/******************************************************************************/
536/*
537 Purpose:
538
539 U_EXACT evaluates the exact solution.
540
541 Licensing:
542
543 This code is distributed under the GNU LGPL license.
544
545 Modified:
546
547 25 October 2011
548
549 Author:
550
551 John Burkardt
552
553 Parameters:
554
555 Input, double X, Y, the coordinates of a point.
556
557 Output, double U_EXACT, the value of the exact solution
558 at (X,Y).
559*/
560{
561 double pi = 3.141592653589793;
562 double value;
563
564 value = sin ( pi * x * y );
565
566 return value;
567}
568/******************************************************************************/
569
570double uxxyy_exact ( double x, double y )
571
572/******************************************************************************/
573/*
574 Purpose:
575
576 UXXYY_EXACT evaluates ( d/dx d/dx + d/dy d/dy ) of the exact solution.
577
578 Licensing:
579
580 This code is distributed under the GNU LGPL license.
581
582 Modified:
583
584 25 October 2011
585
586 Author:
587
588 John Burkardt
589
590 Parameters:
591
592 Input, double X, Y, the coordinates of a point.
593
594 Output, double UXXYY_EXACT, the value of
595 ( d/dx d/dx + d/dy d/dy ) of the exact solution at (X,Y).
596*/
597{
598 double pi = 3.141592653589793;
599 double value;
600
601 value = - pi * pi * ( x * x + y * y ) * sin ( pi * x * y );
602
603 return value;
604}
605# undef NX
606# undef NY
Note: See TracBrowser for help on using the repository browser.