source: CIVL/examples/omp/md_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: 13.8 KB
RevLine 
[20ac35f]1#ifdef _CIVL
2#include <civlc.cvh>
3#endif
4
[71264c4]5# include <stdlib.h>
6# include <stdio.h>
7# include <time.h>
8# include <math.h>
9# include <omp.h>
10
[20ac35f]11#include <assert.h>
12
13#ifdef _CIVL
14$input int ND=1; // Originally 3
15$input int NP=10; // Originally 1000
16$input int NSTEPS=10; // Originally 400
17$input double DT=0.1; // Originally 0.0001
18#else
19#define ND 3
20#define NP 1000
21#define NSTEPS 400
22#define DT 0.0001
23#endif
[7c6a826]24
[71264c4]25int main ( int argc, char *argv[] );
26void compute ( int np, int nd, double pos[], double vel[],
27 double mass, double f[], double *pot, double *kin );
28double dist ( int nd, double r1[], double r2[], double dr[] );
29void initialize ( int np, int nd, double box[], int *seed, double pos[],
30 double vel[], double acc[] );
31double r8_uniform_01 ( int *seed );
32void timestamp ( void );
33void update ( int np, int nd, double pos[], double vel[], double f[],
34 double acc[], double mass, double dt );
35
36/******************************************************************************/
37
38int main ( int argc, char *argv[] )
39
40/******************************************************************************/
41/*
42 Purpose:
43
44 MAIN is the main program for MD_OPENMP.
45
46 Discussion:
47
48 MD implements a simple molecular dynamics simulation.
49
50 The program uses Open MP directives to allow parallel computation.
51
52 The velocity Verlet time integration scheme is used.
53
54 The particles interact with a central pair potential.
55
56 Licensing:
57
58 This code is distributed under the GNU LGPL license.
59
60 Modified:
61
62 30 July 2009
63
64 Author:
65
66 Original FORTRAN77 version by Bill Magro.
67 C version by John Burkardt.
68
69 Parameters:
70
71 None
72*/
73{
74 double *acc;
75 double *box;
[7c6a826]76 double dt = DT;
[71264c4]77 double e0;
78 double *force;
79 int i;
80 int id;
81 double kinetic;
82 double mass = 1.0;
[7c6a826]83 int nd = ND;
84 int np = NP;
[71264c4]85 double *pos;
86 double potential;
87 int proc_num;
88 int seed = 123456789;
89 int step;
[7c6a826]90 int step_num = NSTEPS;
[71264c4]91 int step_print;
92 int step_print_index;
93 int step_print_num;
94 double *vel;
95 double wtime;
96
97 timestamp ( );
98
99 proc_num = omp_get_num_procs ( );
100
101 acc = ( double * ) malloc ( nd * np * sizeof ( double ) );
102 box = ( double * ) malloc ( nd * sizeof ( double ) );
103 force = ( double * ) malloc ( nd * np * sizeof ( double ) );
104 pos = ( double * ) malloc ( nd * np * sizeof ( double ) );
105 vel = ( double * ) malloc ( nd * np * sizeof ( double ) );
106
107 printf ( "\n" );
108 printf ( "MD_OPENMP\n" );
109 printf ( " C/OpenMP version\n" );
110 printf ( "\n" );
111 printf ( " A molecular dynamics program.\n" );
112
113 printf ( "\n" );
114 printf ( " NP, the number of particles in the simulation is %d\n", np );
115 printf ( " STEP_NUM, the number of time steps, is %d\n", step_num );
116 printf ( " DT, the size of each time step, is %f\n", dt );
117
118 printf ( "\n" );
119 printf ( " Number of processors available = %d\n", omp_get_num_procs ( ) );
120 printf ( " Number of threads = %d\n", omp_get_max_threads ( ) );
121/*
122 Set the dimensions of the box.
123*/
124 for ( i = 0; i < nd; i++ )
125 {
126 box[i] = 10.0;
127 }
128
129 printf ( "\n" );
130 printf ( " Initializing positions, velocities, and accelerations.\n" );
131/*
132 Set initial positions, velocities, and accelerations.
133*/
134 initialize ( np, nd, box, &seed, pos, vel, acc );
135/*
136 Compute the forces and energies.
137*/
138 printf ( "\n" );
139 printf ( " Computing initial forces and energies.\n" );
140
141 compute ( np, nd, pos, vel, mass, force, &potential, &kinetic );
142
143 e0 = potential + kinetic;
144/*
145 This is the main time stepping loop:
146 Compute forces and energies,
147 Update positions, velocities, accelerations.
148*/
149 printf ( "\n" );
150 printf ( " At each step, we report the potential and kinetic energies.\n" );
151 printf ( " The sum of these energies should be a constant.\n" );
152 printf ( " As an accuracy check, we also print the relative error\n" );
153 printf ( " in the total energy.\n" );
154 printf ( "\n" );
155 printf ( " Step Potential Kinetic (P+K-E0)/E0\n" );
156 printf ( " Energy P Energy K Relative Energy Error\n" );
157 printf ( "\n" );
158
159 step_print = 0;
160 step_print_index = 0;
161 step_print_num = 10;
162
163 step = 0;
164 printf ( " %8d %14f %14f %14e\n",
165 step, potential, kinetic, ( potential + kinetic - e0 ) / e0 );
166 step_print_index = step_print_index + 1;
167 step_print = ( step_print_index * step_num ) / step_print_num;
168
169 wtime = omp_get_wtime ( );
170
171 for ( step = 1; step <= step_num; step++ )
172 {
173 compute ( np, nd, pos, vel, mass, force, &potential, &kinetic );
174
175 if ( step == step_print )
176 {
177 printf ( " %8d %14f %14f %14e\n", step, potential, kinetic,
178 ( potential + kinetic - e0 ) / e0 );
179 step_print_index = step_print_index + 1;
180 step_print = ( step_print_index * step_num ) / step_print_num;
181 }
182 update ( np, nd, pos, vel, force, acc, mass, dt );
183 }
184 wtime = omp_get_wtime ( ) - wtime;
185
186 printf ( "\n" );
187 printf ( " Elapsed time for main computation:\n" );
188 printf ( " %f seconds.\n", wtime );
189
190 free ( acc );
191 free ( box );
192 free ( force );
193 free ( pos );
194 free ( vel );
195/*
196 Terminate.
197*/
198 printf ( "\n" );
199 printf ( "MD_OPENMP\n" );
200 printf ( " Normal end of execution.\n" );
201
202 printf ( "\n" );
203 timestamp ( );
204
205 return 0;
206}
207/******************************************************************************/
208
209void compute ( int np, int nd, double pos[], double vel[],
210 double mass, double f[], double *pot, double *kin )
211
212/******************************************************************************/
213/*
214 Purpose:
215
216 COMPUTE computes the forces and energies.
217
218 Discussion:
219
220 The computation of forces and energies is fully parallel.
221
222 The potential function V(X) is a harmonic well which smoothly
223 saturates to a maximum value at PI/2:
224
225 v(x) = ( sin ( min ( x, PI2 ) ) )**2
226
227 The derivative of the potential is:
228
229 dv(x) = 2.0 * sin ( min ( x, PI2 ) ) * cos ( min ( x, PI2 ) )
230 = sin ( 2.0 * min ( x, PI2 ) )
231
232 Licensing:
233
234 This code is distributed under the GNU LGPL license.
235
236 Modified:
237
238 21 November 2007
239
240 Author:
241
242 Original FORTRAN77 version by Bill Magro.
243 C version by John Burkardt.
244
245 Parameters:
246
247 Input, int NP, the number of particles.
248
249 Input, int ND, the number of spatial dimensions.
250
251 Input, double POS[ND*NP], the position of each particle.
252
253 Input, double VEL[ND*NP], the velocity of each particle.
254
255 Input, double MASS, the mass of each particle.
256
257 Output, double F[ND*NP], the forces.
258
259 Output, double *POT, the total potential energy.
260
261 Output, double *KIN, the total kinetic energy.
262*/
263{
264 double d;
265 double d2;
266 int i;
267 int j;
268 int k;
269 double ke;
270 double pe;
271 double PI2 = 3.141592653589793 / 2.0;
272 double rij[3];
273
274 pe = 0.0;
275 ke = 0.0;
276
277# pragma omp parallel \
[20ac35f]278 shared ( f, np, pos, vel ) \
279 private ( i, j, k, rij, d, d2, nd )
[71264c4]280
281
282# pragma omp for reduction ( + : pe, ke )
283 for ( k = 0; k < np; k++ )
284 {
[20ac35f]285 nd = ND;
[71264c4]286/*
287 Compute the potential energy and forces.
288*/
289 for ( i = 0; i < nd; i++ )
290 {
291 f[i+k*nd] = 0.0;
292 }
293
294 for ( j = 0; j < np; j++ )
295 {
296 if ( k != j )
297 {
298 d = dist ( nd, pos+k*nd, pos+j*nd, rij );
299/*
300 Attribute half of the potential energy to particle J.
301*/
302 if ( d < PI2 )
303 {
304 d2 = d;
305 }
306 else
307 {
308 d2 = PI2;
309 }
310
311 pe = pe + 0.5 * pow ( sin ( d2 ), 2 );
[20ac35f]312 #ifdef _CIVL
313 $assume(pe != 0);
314 #endif
[71264c4]315
316 for ( i = 0; i < nd; i++ )
317 {
318 f[i+k*nd] = f[i+k*nd] - rij[i] * sin ( 2.0 * d2 ) / d;
319 }
320 }
321 }
322/*
323 Compute the kinetic energy.
324*/
325 for ( i = 0; i < nd; i++ )
326 {
327 ke = ke + vel[i+k*nd] * vel[i+k*nd];
[20ac35f]328 #ifdef _CIVL
329 $assume(ke != 0);
330 #endif
[71264c4]331 }
332 }
333
334 ke = ke * 0.5 * mass;
335
336 *pot = pe;
337 *kin = ke;
338
339 return;
340}
341/******************************************************************************/
342
343double dist ( int nd, double r1[], double r2[], double dr[] )
344
345/******************************************************************************/
346/*
347 Purpose:
348
349 DIST computes the displacement (and its norm) between two particles.
350
351 Licensing:
352
353 This code is distributed under the GNU LGPL license.
354
355 Modified:
356
357 21 November 2007
358
359 Author:
360
361 Original FORTRAN77 version by Bill Magro.
362 C version by John Burkardt.
363
364 Parameters:
365
366 Input, int ND, the number of spatial dimensions.
367
368 Input, double R1[ND], R2[ND], the positions of the particles.
369
370 Output, double DR[ND], the displacement vector.
371
372 Output, double D, the Euclidean norm of the displacement.
373*/
374{
375 double d;
376 int i;
377
378 d = 0.0;
[20ac35f]379 for (int i = 0; i < nd; i++ )
[71264c4]380 {
381 dr[i] = r1[i] - r2[i];
382 d = d + dr[i] * dr[i];
383 }
384 d = sqrt ( d );
385
386 return d;
387}
388/******************************************************************************/
389
390void initialize ( int np, int nd, double box[], int *seed, double pos[],
391 double vel[], double acc[] )
392
393/******************************************************************************/
394/*
395 Purpose:
396
397 INITIALIZE initializes the positions, velocities, and accelerations.
398
399 Licensing:
400
401 This code is distributed under the GNU LGPL license.
402
403 Modified:
404
405 21 November 2007
406
407 Author:
408
409 Original FORTRAN77 version by Bill Magro.
410 C version by John Burkardt.
411
412 Parameters:
413
414 Input, int NP, the number of particles.
415
416 Input, int ND, the number of spatial dimensions.
417
418 Input, double BOX[ND], specifies the maximum position
419 of particles in each dimension.
420
421 Input, int *SEED, a seed for the random number generator.
422
423 Output, double POS[ND*NP], the position of each particle.
424
425 Output, double VEL[ND*NP], the velocity of each particle.
426
427 Output, double ACC[ND*NP], the acceleration of each particle.
428*/
429{
430 int i;
431 int j;
432/*
433 Give the particles random positions within the box.
434*/
435 for ( i = 0; i < nd; i++ )
436 {
437 for ( j = 0; j < np; j++ )
438 {
439 pos[i+j*nd] = box[i] * r8_uniform_01 ( seed );
440 }
441 }
442
443 for ( j = 0; j < np; j++ )
444 {
445 for ( i = 0; i < nd; i++ )
446 {
447 vel[i+j*nd] = 0.0;
448 }
449 }
450 for ( j = 0; j < np; j++ )
451 {
452 for ( i = 0; i < nd; i++ )
453 {
454 acc[i+j*nd] = 0.0;
455 }
456 }
457 return;
458}
459/******************************************************************************/
460
461double r8_uniform_01 ( int *seed )
462
463/******************************************************************************/
464/*
465 Purpose:
466
467 R8_UNIFORM_01 is a unit pseudorandom R8.
468
469 Discussion:
470
471 This routine implements the recursion
472
473 seed = 16807 * seed mod ( 2**31 - 1 )
474 unif = seed / ( 2**31 - 1 )
475
476 The integer arithmetic never requires more than 32 bits,
477 including a sign bit.
478
479 Licensing:
480
481 This code is distributed under the GNU LGPL license.
482
483 Modified:
484
485 11 August 2004
486
487 Author:
488
489 John Burkardt
490
491 Reference:
492
493 Paul Bratley, Bennett Fox, Linus Schrage,
494 A Guide to Simulation,
495 Springer Verlag, pages 201-202, 1983.
496
497 Bennett Fox,
498 Algorithm 647:
499 Implementation and Relative Efficiency of Quasirandom
500 Sequence Generators,
501 ACM Transactions on Mathematical Software,
502 Volume 12, Number 4, pages 362-376, 1986.
503
504 Parameters:
505
506 Input/output, int *SEED, a seed for the random number generator.
507
508 Output, double R8_UNIFORM_01, a new pseudorandom variate, strictly between
509 0 and 1.
510*/
511{
512 int k;
513 double r;
514
515 k = *seed / 127773;
516
517 *seed = 16807 * ( *seed - k * 127773 ) - k * 2836;
518
519 if ( *seed < 0 )
520 {
521 *seed = *seed + 2147483647;
522 }
523
524 r = ( double ) ( *seed ) * 4.656612875E-10;
525
526 return r;
527}
528/******************************************************************************/
529
530void timestamp ( void )
531
532/******************************************************************************/
533/*
534 Purpose:
535
536 TIMESTAMP prints the current YMDHMS date as a time stamp.
537
538 Example:
539
540 31 May 2001 09:45:54 AM
541
542 Licensing:
543
544 This code is distributed under the GNU LGPL license.
545
546 Modified:
547
548 24 September 2003
549
550 Author:
551
552 John Burkardt
553
554 Parameters:
555
556 None
557*/
558{
559# define TIME_SIZE 40
560
561 static char time_buffer[TIME_SIZE];
562 const struct tm *tm;
563 size_t len;
564 time_t now;
565
566 now = time ( NULL );
567 tm = localtime ( &now );
568
569 len = strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm );
570
571 printf ( "%s\n", time_buffer );
572
573 return;
574# undef TIME_SIZE
575}
576/******************************************************************************/
577
578void update ( int np, int nd, double pos[], double vel[], double f[],
579 double acc[], double mass, double dt )
580
581/******************************************************************************/
582/*
583 Purpose:
584
585 UPDATE updates positions, velocities and accelerations.
586
587 Discussion:
588
589 The time integration is fully parallel.
590
591 A velocity Verlet algorithm is used for the updating.
592
593 x(t+dt) = x(t) + v(t) * dt + 0.5 * a(t) * dt * dt
594 v(t+dt) = v(t) + 0.5 * ( a(t) + a(t+dt) ) * dt
595 a(t+dt) = f(t) / m
596
597 Licensing:
598
599 This code is distributed under the GNU LGPL license.
600
601 Modified:
602
603 17 April 2009
604
605 Author:
606
607 Original FORTRAN77 version by Bill Magro.
608 C version by John Burkardt.
609
610 Parameters:
611
612 Input, int NP, the number of particles.
613
614 Input, int ND, the number of spatial dimensions.
615
616 Input/output, double POS[ND*NP], the position of each particle.
617
618 Input/output, double VEL[ND*NP], the velocity of each particle.
619
620 Input, double F[ND*NP], the force on each particle.
621
622 Input/output, double ACC[ND*NP], the acceleration of each particle.
623
624 Input, double MASS, the mass of each particle.
625
626 Input, double DT, the time step.
627*/
628{
629 int i;
630 int j;
631 double rmass;
632
633 rmass = 1.0 / mass;
634
635# pragma omp parallel \
636 shared ( acc, dt, f, nd, np, pos, rmass, vel ) \
637 private ( i, j )
638
639# pragma omp for
640 for ( j = 0; j < np; j++ )
641 {
642 for ( i = 0; i < nd; i++ )
643 {
644 pos[i+j*nd] = pos[i+j*nd] + vel[i+j*nd] * dt + 0.5 * acc[i+j*nd] * dt * dt;
645 vel[i+j*nd] = vel[i+j*nd] + 0.5 * dt * ( f[i+j*nd] * rmass + acc[i+j*nd] );
646 acc[i+j*nd] = f[i+j*nd] * rmass;
647 }
648 }
649
650 return;
651}
Note: See TracBrowser for help on using the repository browser.