source: CIVL/examples/omp/helmholtz.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: 15.7 KB
RevLine 
[5f600f6]1# include <stdlib.h>
2# include <stdio.h>
3# include <math.h>
4# include <omp.h>
5
6int main ( int argc, char *argv[] );
7void driver ( int m, int n, int it_max, double alpha, double omega, double tol );
8void error_check ( int m, int n, double alpha, double u[], double f[] );
9void jacobi ( int m, int n, double alpha, double omega, double u[], double f[],
10 double tol, int it_max );
11double *rhs_set ( int m, int n, double alpha );
12double u_exact ( double x, double y );
13double uxx_exact ( double x, double y );
14double uyy_exact ( double x, double y );
15
16/******************************************************************************/
17
18int main ( int argc, char *argv[] )
19
20/******************************************************************************/
21/*
22 Purpose:
23
24 MAIN is the main program for HELMHOLTZ.
25
26 Discussion:
27
28 HELMHOLTZ solves a discretized Helmholtz equation.
29
30 The two dimensional region given is:
31
32 -1 <= X <= +1
33 -1 <= Y <= +1
34
35 The region is discretized by a set of M by N nodes:
36
37 P(I,J) = ( X(I), Y(J) )
38
39 where, for 0 <= I <= M-1, 0 <= J <= N - 1, (C/C++ convention)
40
41 X(I) = ( 2 * I - M + 1 ) / ( M - 1 )
42 Y(J) = ( 2 * J - N + 1 ) / ( N - 1 )
43
44 The Helmholtz equation for the scalar function U(X,Y) is
45
46 - Uxx(X,Y) -Uyy(X,Y) + ALPHA * U(X,Y) = F(X,Y)
47
48 where ALPHA is a positive constant. We suppose that Dirichlet
49 boundary conditions are specified, that is, that the value of
50 U(X,Y) is given for all points along the boundary.
51
52 We suppose that the right hand side function F(X,Y) is specified in
53 such a way that the exact solution is
54
55 U(X,Y) = ( 1 - X^2 ) * ( 1 - Y^2 )
56
57 Using standard finite difference techniques, the second derivatives
58 of U can be approximated by linear combinations of the values
59 of U at neighboring points. Using this fact, the discretized
60 differential equation becomes a set of linear equations of the form:
61
62 A * U = F
63
64 These linear equations are then solved using a form of the Jacobi
65 iterative method with a relaxation factor.
66
67 Directives are used in this code to achieve parallelism.
68 All do loops are parallized with default 'static' scheduling.
69
70 Licensing:
71
72 This code is distributed under the GNU LGPL license.
73
74 Modified:
75
76 19 April 2009
77
78 Author:
79
80 Original FORTRAN77 version by Joseph Robicheaux, Sanjiv Shah.
81 C version by John Burkardt
82*/
83{
84 double alpha = 0.25;
85 int it_max = 100;
86 int m = 500;
87 int n = 500;
88 double omega = 1.1;
89 double tol = 1.0E-08;
90 double wtime;
91
92 printf ( "\n" );
93 printf ( "HELMHOLTZ\n" );
94 printf ( " C/OpenMP version\n" );
95 printf ( "\n" );
96 printf ( " A program which solves the 2D Helmholtz equation.\n" );
97
98 printf ( "\n" );
99 printf ( " This program is being run in parallel.\n" );
100
101 printf ( "\n" );
102 printf ( " Number of processors available = %d\n", omp_get_num_procs ( ) );
103 printf ( " Number of threads = %d\n", omp_get_max_threads ( ) );
104
105 printf ( "\n" );
106 printf ( " The region is [-1,1] x [-1,1].\n" );
107 printf ( " The number of nodes in the X direction is M = %d\n", m );
108 printf ( " The number of nodes in the Y direction is N = %d\n", n );
109 printf ( " Number of variables in linear system M * N = %d\n", m * n );
110 printf ( " The scalar coefficient in the Helmholtz equation is ALPHA = %f\n",
111 alpha );
112 printf ( " The relaxation value is OMEGA = %f\n", omega );
113 printf ( " The error tolerance is TOL = %f\n", tol );
114 printf ( " The maximum number of Jacobi iterations is IT_MAX = %d\n",
115 it_max );
116/*
117 Call the driver routine.
118*/
119 wtime = omp_get_wtime ( );
120
121 driver ( m, n, it_max, alpha, omega, tol );
122
123 wtime = omp_get_wtime ( ) - wtime;
124
125 printf ( "\n" );
126 printf ( " Elapsed wall clock time = %f\n", wtime );
127/*
128 Terminate.
129*/
130 printf ( "\n" );
131 printf ( "HELMHOLTZ\n" );
132 printf ( " Normal end of execution.\n" );
133
134 return 0;
135}
136/******************************************************************************/
137
138void driver ( int m, int n, int it_max, double alpha, double omega, double tol )
139
140/******************************************************************************/
141/*
142 Purpose:
143
144 DRIVER allocates arrays and solves the problem.
145
146 Licensing:
147
148 This code is distributed under the GNU LGPL license.
149
150 Modified:
151
152 21 November 2007
153
154 Author:
155
156 Original FORTRAN77 version by Joseph Robicheaux, Sanjiv Shah.
157 C version by John Burkardt
158
159 Parameters:
160
161 Input, int M, N, the number of grid points in the
162 X and Y directions.
163
164 Input, int IT_MAX, the maximum number of Jacobi
165 iterations allowed.
166
167 Input, double ALPHA, the scalar coefficient in the
168 Helmholtz equation.
169
170 Input, double OMEGA, the relaxation parameter, which
171 should be strictly between 0 and 2. For a pure Jacobi method,
172 use OMEGA = 1.
173
174 Input, double TOL, an error tolerance for the linear
175 equation solver.
176*/
177{
178 double *f;
179 int i;
180 int j;
181 double *u;
182/*
183 Initialize the data.
184*/
185 f = rhs_set ( m, n, alpha );
186
187 u = ( double * ) malloc ( m * n * sizeof ( double ) );
188
189# pragma omp parallel \
190 shared ( m, n, u ) \
191 private ( i, j )
192
193# pragma omp for
194
195 for ( j = 0; j < n; j++ )
196 {
197 for ( i = 0; i < m; i++ )
198 {
199 u[i+j*m] = 0.0;
200 }
201 }
202/*
203 Solve the Helmholtz equation.
204*/
205 jacobi ( m, n, alpha, omega, u, f, tol, it_max );
206/*
207 Determine the error.
208*/
209 error_check ( m, n, alpha, u, f );
210
211 free ( f );
212 free ( u );
213
214 return;
215}
216/******************************************************************************/
217
218void error_check ( int m, int n, double alpha, double u[], double f[] )
219
220/******************************************************************************/
221/*
222 Purpose:
223
224 ERROR_CHECK determines the error in the numerical solution.
225
226 Licensing:
227
228 This code is distributed under the GNU LGPL license.
229
230 Modified:
231
232 21 November 2007
233
234 Author:
235
236 Original FORTRAN77 version by Joseph Robicheaux, Sanjiv Shah.
237 C version by John Burkardt
238
239 Parameters:
240
241 Input, int M, N, the number of grid points in the
242 X and Y directions.
243
244 Input, double ALPHA, the scalar coefficient in the
245 Helmholtz equation. ALPHA should be positive.
246
247 Input, double U[M*N], the solution of the Helmholtz equation
248 at the grid points.
249
250 Input, double F[M*N], values of the right hand side function
251 for the Helmholtz equation at the grid points.
252*/
253{
254 double error_norm;
255 int i;
256 int j;
257 double u_norm;
258 double u_true;
259 double u_true_norm;
260 double x;
261 double y;
262
263 u_norm = 0.0;
264
265# pragma omp parallel \
266 shared ( m, n, u ) \
267 private ( i, j )
268
269# pragma omp for reduction ( + : u_norm )
270
271 for ( j = 0; j < n; j++ )
272 {
273 for ( i = 0; i < m; i++ )
274 {
275 u_norm = u_norm + u[i+j*m] * u[i+j*m];
276 }
277 }
278
279 u_norm = sqrt ( u_norm );
280
281 u_true_norm = 0.0;
282 error_norm = 0.0;
283
284# pragma omp parallel \
285 shared ( m, n, u ) \
286 private ( i, j, u_true, x, y )
287
288# pragma omp for reduction ( + : error_norm, u_true_norm)
289
290 for ( j = 0; j < n; j++ )
291 {
292 for ( i = 0; i < m; i++ )
293 {
294 x = ( double ) ( 2 * i - m + 1 ) / ( double ) ( m - 1 );
295 y = ( double ) ( 2 * j - n + 1 ) / ( double ) ( n - 1 );
296 u_true = u_exact ( x, y );
297 error_norm = error_norm + ( u[i+j*m] - u_true ) * ( u[i+j*m] - u_true );
298 u_true_norm = u_true_norm + u_true * u_true;
299 }
300 }
301
302 error_norm = sqrt ( error_norm );
303 u_true_norm = sqrt ( u_true_norm );
304
305 printf ( "\n" );
306 printf ( " Computed U l2 norm : %f\n", u_norm );
307 printf ( " Computed U_EXACT l2 norm : %f\n", u_true_norm );
308 printf ( " Error l2 norm: %f\n", error_norm );
309
310 return;
311}
312/******************************************************************************/
313
314void jacobi ( int m, int n, double alpha, double omega, double u[], double f[],
315 double tol, int it_max )
316
317/******************************************************************************/
318/*
319 Purpose:
320
321 JACOBI applies the Jacobi iterative method to solve the linear system.
322
323 Licensing:
324
325 This code is distributed under the GNU LGPL license.
326
327 Modified:
328
329 21 November 2007
330
331 Author:
332
333 Original FORTRAN77 version by Joseph Robicheaux, Sanjiv Shah.
334 C version by John Burkardt
335
336 Parameters:
337
338 Input, int M, N, the number of grid points in the
339 X and Y directions.
340
341 Input, double ALPHA, the scalar coefficient in the
342 Helmholtz equation. ALPHA should be positive.
343
344 Input, double OMEGA, the relaxation parameter, which
345 should be strictly between 0 and 2. For a pure Jacobi method,
346 use OMEGA = 1.
347
348 Input/output, double U(M,N), the solution of the Helmholtz
349 equation at the grid points.
350
351 Input, double F(M,N), values of the right hand side function
352 for the Helmholtz equation at the grid points.
353
354 Input, double TOL, an error tolerance for the linear
355 equation solver.
356
357 Input, int IT_MAX, the maximum number of Jacobi
358 iterations allowed.
359*/
360{
361 double ax;
362 double ay;
363 double b;
364 double dx;
365 double dy;
366 double error;
367 double error_norm;
368 int i;
369 int it;
370 int j;
371 double *u_old;
372/*
373 Initialize the coefficients.
374*/
375 dx = 2.0 / ( double ) ( m - 1 );
376 dy = 2.0 / ( double ) ( n - 1 );
377
378 ax = - 1.0 / dx / dx;
379 ay = - 1.0 / dy / dy;
380 b = + 2.0 / dx / dx + 2.0 / dy / dy + alpha;
381
382 u_old = ( double * ) malloc ( m * n * sizeof ( double ) );
383
384 for ( it = 1; it <= it_max; it++ )
385 {
386 error_norm = 0.0;
387/*
388 Copy new solution into old.
389*/
390# pragma omp parallel \
391 shared ( m, n, u, u_old ) \
392 private ( i, j )
393
394# pragma omp for
395 for ( j = 0; j < n; j++ )
396 {
397 for ( i = 0; i < m; i++ )
398 {
399 u_old[i+m*j] = u[i+m*j];
400 }
401 }
402/*
403 Compute stencil, residual, and update.
404*/
405# pragma omp parallel \
406 shared ( ax, ay, b, f, m, n, omega, u, u_old ) \
407 private ( error, i, j )
408
409# pragma omp for reduction ( + : error_norm )
410
411 for ( j = 0; j < n; j++ )
412 {
413 for ( i = 0; i < m; i++ )
414 {
415/*
416 Evaluate the residual.
417*/
418 if ( i == 0 || i == m - 1 || j == 0 || j == n - 1 )
419 {
420 error = u_old[i+j*m] - f[i+j*m];
421 }
422 else
423 {
424 error = ( ax * ( u_old[i-1+j*m] + u_old[i+1+j*m] )
425 + ay * ( u_old[i+(j-1)*m] + u_old[i+(j+1)*m] )
426 + b * u_old[i+j*m] - f[i+j*m] ) / b;
427 }
428/*
429 Update the solution.
430*/
431 u[i+j*m] = u_old[i+j*m] - omega * error;
432/*
433 Accumulate the residual error.
434*/
435 error_norm = error_norm + error * error;
436 }
437 }
438/*
439 Error check.
440*/
441 error_norm = sqrt ( error_norm ) / ( double ) ( m * n );
442
443 printf ( " %d Residual RMS %e\n", it, error_norm );
444
445 if ( error_norm <= tol )
446 {
447 break;
448 }
449
450 }
451
452 printf ( "\n" );
453 printf ( " Total number of iterations %d\n", it );
454
455 free ( u_old );
456
457 return;
458}
459/******************************************************************************/
460
461double *rhs_set ( int m, int n, double alpha )
462
463/******************************************************************************/
464/*
465 Purpose:
466
467 RHS_SET sets the right hand side F(X,Y).
468
469 Discussion:
470
471 The routine assumes that the exact solution and its second
472 derivatives are given by the routine EXACT.
473
474 The appropriate Dirichlet boundary conditions are determined
475 by getting the value of U returned by EXACT.
476
477 The appropriate right hand side function is determined by
478 having EXACT return the values of U, UXX and UYY, and setting
479
480 F = -UXX - UYY + ALPHA * U
481
482 Licensing:
483
484 This code is distributed under the GNU LGPL license.
485
486 Modified:
487
488 21 November 2007
489
490 Author:
491
492 Original FORTRAN77 version by Joseph Robicheaux, Sanjiv Shah.
493 C version by John Burkardt
494
495 Parameters:
496
497 Input, int M, N, the number of grid points in the
498 X and Y directions.
499
500 Input, double ALPHA, the scalar coefficient in the
501 Helmholtz equation. ALPHA should be positive.
502
503 Output, double RHS[M*N], values of the right hand side function
504 for the Helmholtz equation at the grid points.
505*/
506{
507 double *f;
508 double f_norm;
509 int i;
510 int j;
511 double x;
512 double y;
513
514 f = ( double * ) malloc ( m * n * sizeof ( double ) );
515
516# pragma omp parallel \
517 shared ( f, m, n ) \
518 private ( i, j )
519
520# pragma omp for
521
522 for ( j = 0; j < n; j++ )
523 {
524 for ( i = 0; i < m; i++ )
525 {
526 f[i+j*m] = 0.0;
527 }
528 }
529/*
530 Set the boundary conditions.
531*/
532
533# pragma omp parallel \
534 shared ( alpha, f, m, n ) \
535 private ( i, j, x, y )
536 {
537
538# pragma omp for
539 for ( i = 0; i < m; i++ )
540 {
541 j = 0;
542 y = ( double ) ( 2 * j - n + 1 ) / ( double ) ( n - 1 );
543 x = ( double ) ( 2 * i - m + 1 ) / ( double ) ( m - 1 );
544 f[i+j*m] = u_exact ( x, y );
545 }
546
547# pragma omp for
548 for ( i = 0; i < m; i++ )
549 {
550 j = n - 1;
551 y = ( double ) ( 2 * j - n + 1 ) / ( double ) ( n - 1 );
552 x = ( double ) ( 2 * i - m + 1 ) / ( double ) ( m - 1 );
553 f[i+j*m] = u_exact ( x, y );
554 }
555
556# pragma omp for
557 for ( j = 0; j < n; j++ )
558 {
559 i = 0;
560 x = ( double ) ( 2 * i - m + 1 ) / ( double ) ( m - 1 );
561 y = ( double ) ( 2 * j - n + 1 ) / ( double ) ( n - 1 );
562 f[i+j*m] = u_exact ( x, y );
563 }
564
565# pragma omp for
566
567 for ( j = 0; j < n; j++ )
568 {
569 i = m - 1;
570 x = ( double ) ( 2 * i - m + 1 ) / ( double ) ( m - 1 );
571 y = ( double ) ( 2 * j - n + 1 ) / ( double ) ( n - 1 );
572 f[i+j*m] = u_exact ( x, y );
573 }
574/*
575 Set the right hand side F.
576*/
577# pragma omp for
578
579 for ( j = 1; j < n - 1; j++ )
580 {
581 for ( i = 1; i < m - 1; i++ )
582 {
583 x = ( double ) ( 2 * i - m + 1 ) / ( double ) ( m - 1 );
584 y = ( double ) ( 2 * j - n + 1 ) / ( double ) ( n - 1 );
585 f[i+j*m] = - uxx_exact ( x, y ) - uyy_exact ( x, y ) + alpha * u_exact ( x, y );
586 }
587 }
588 }
589
590 f_norm = 0.0;
591
592# pragma omp parallel \
593 shared ( f, m, n ) \
594 private ( i, j )
595
596# pragma omp for reduction ( + : f_norm )
597
598 for ( j = 0; j < n; j++ )
599 {
600 for ( i = 0; i < m; i++ )
601 {
602 f_norm = f_norm + f[i+j*m] * f[i+j*m];
603 }
604 }
605 f_norm = sqrt ( f_norm );
606
607 printf ( "\n" );
608 printf ( " Right hand side l2 norm = %f\n", f_norm );
609
610 return f;
611}
612/******************************************************************************/
613
614double u_exact ( double x, double y )
615
616/******************************************************************************/
617/*
618 Purpose:
619
620 U_EXACT returns the exact value of U(X,Y).
621
622 Licensing:
623
624 This code is distributed under the GNU LGPL license.
625
626 Modified:
627
628 21 November 2007
629
630 Author:
631
632 John Burkardt
633
634 Parameters:
635
636 Input, double X, Y, the point at which the values are needed.
637
638 Output, double U_EXACT, the value of the exact solution.
639*/
640{
641 double value;
642
643 value = ( 1.0 - x * x ) * ( 1.0 - y * y );
644
645 return value;
646}
647/******************************************************************************/
648
649double uxx_exact ( double x, double y )
650
651/******************************************************************************/
652/*
653 Purpose:
654
655 UXX_EXACT returns the exact second X derivative of the solution.
656
657 Licensing:
658
659 This code is distributed under the GNU LGPL license.
660
661 Modified:
662
663 21 November 2007
664
665 Author:
666
667 John Burkardt
668
669 Parameters:
670
671 Input, double X, Y, the point at which the values are needed.
672
673 Output, double UXX_EXACT, the exact second X derivative.
674*/
675{
676 double value;
677
678 value = -2.0 * ( 1.0 + y ) * ( 1.0 - y );
679
680 return value;
681}
682/******************************************************************************/
683
684double uyy_exact ( double x, double y )
685
686/******************************************************************************/
687/*
688 Purpose:
689
690 UYY_EXACT returns the exact second Y derivative of the solution.
691
692 Licensing:
693
694 This code is distributed under the GNU LGPL license.
695
696 Modified:
697
698 21 November 2007
699
700 Author:
701
702 John Burkardt
703
704 Parameters:
705
706 Input, double X, Y, the point at which the values are needed.
707
708 Output, double UYY_EXACT, the exact second Y derivative.
709*/
710{
711 double value;
712
713 value = -2.0 * ( 1.0 + x ) * ( 1.0 - x );
714
715 return value;
716}
Note: See TracBrowser for help on using the repository browser.