source: CIVL/examples/omp/simple/md_openmp.c.s@ 7d77e64

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