source: CIVL/examples/omp/md_openmp.c@ d48b845

1.23 2.0 main test-branch
Last change on this file since d48b845 was 7c6a826, checked in by Matthew B. Dwyer <matthewbdwyer@…>, 11 years ago

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

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