basic linked list / cell method support implemented
[physik/posic.git] / moldyn.c
1 /*
2  * moldyn.c - molecular dynamics library main file
3  *
4  * author: Frank Zirkelbach <frank.zirkelbach@physik.uni-augsburg.de>
5  *
6  */
7
8 #define _GNU_SOURCE
9 #include <stdio.h>
10 #include <stdlib.h>
11 #include <string.h>
12 #include <sys/types.h>
13 #include <sys/stat.h>
14 #include <fcntl.h>
15 #include <unistd.h>
16 #include <math.h>
17
18 #include "moldyn.h"
19
20 #include "math/math.h"
21 #include "init/init.h"
22 #include "random/random.h"
23 #include "visual/visual.h"
24 #include "list/list.h"
25
26 int moldyn_usage(char **argv) {
27
28         printf("\n%s usage:\n\n",argv[0]);
29         printf("--- general options ---\n");
30         printf("-E <steps> <file> (log total energy)\n");
31         printf("-M <steps> <file> (log total momentum)\n");
32         printf("-D <steps> <file> (dump total information)\n");
33         printf("-S <steps> <filebase> (single save file)\n");
34         printf("-V <steps> <filebase> (rasmol file)\n");
35         printf("--- physics options ---\n");
36         printf("-T <temperature> [K] (%f)\n",MOLDYN_TEMP);
37         printf("-t <timestep tau> [s] (%.15f)\n",MOLDYN_TAU);
38         printf("-R <runs> (%d)\n",MOLDYN_RUNS);
39         printf("\n");
40
41         return 0;
42 }
43
44 int moldyn_parse_argv(t_moldyn *moldyn,int argc,char **argv) {
45
46         int i;
47
48         memset(moldyn,0,sizeof(t_moldyn));
49
50         /* default values */
51         moldyn->t=MOLDYN_TEMP;
52         moldyn->tau=MOLDYN_TAU;
53         moldyn->time_steps=MOLDYN_RUNS;
54
55         /* parse argv */
56         for(i=1;i<argc;i++) {
57                 if(argv[i][0]=='-') {
58                         switch(argv[i][1]){
59                                 case 'E':
60                                         moldyn->ewrite=atoi(argv[++i]);
61                                         strncpy(moldyn->efb,argv[++i],64);
62                                         break;
63                                 case 'M':
64                                         moldyn->mwrite=atoi(argv[++i]);
65                                         strncpy(moldyn->mfb,argv[++i],64);
66                                         break;
67                                 case 'D':
68                                         moldyn->dwrite=atoi(argv[++i]);
69                                         strncpy(moldyn->dfb,argv[++i],64);
70                                         break;
71                                 case 'S':
72                                         moldyn->swrite=atoi(argv[++i]);
73                                         strncpy(moldyn->sfb,argv[++i],64);
74                                         break;
75                                 case 'V':
76                                         moldyn->vwrite=atoi(argv[++i]);
77                                         strncpy(moldyn->vfb,argv[++i],64);
78                                         break;
79                                 case 'T':
80                                         moldyn->t=atof(argv[++i]);
81                                         break;
82                                 case 't':
83                                         moldyn->tau=atof(argv[++i]);
84                                         break;
85                                 case 'R':
86                                         moldyn->time_steps=atoi(argv[++i]);
87                                         break;
88                                 default:
89                                         printf("unknown option %s\n",argv[i]);
90                                         moldyn_usage(argv);
91                                         return -1;
92                         }
93                 } else {
94                         moldyn_usage(argv);
95                         return -1;
96                 }
97         }
98
99         return 0;
100 }
101
102 int moldyn_log_init(t_moldyn *moldyn,void *v) {
103
104         moldyn->lvstat=0;
105         t_visual *vis;
106
107         vis=v;
108
109         if(moldyn->ewrite) {
110                 moldyn->efd=open(moldyn->efb,O_WRONLY|O_CREAT|O_TRUNC);
111                 if(moldyn->efd<0) {
112                         perror("[moldyn] efd open");
113                         return moldyn->efd;
114                 }
115                 dprintf(moldyn->efd,"# moldyn total energy logfile\n");
116                 moldyn->lvstat|=MOLDYN_LVSTAT_TOTAL_E;
117         }
118
119         if(moldyn->mwrite) {
120                 moldyn->mfd=open(moldyn->mfb,O_WRONLY|O_CREAT|O_TRUNC);
121                 if(moldyn->mfd<0) {
122                         perror("[moldyn] mfd open");
123                         return moldyn->mfd;
124                 }
125                 dprintf(moldyn->mfd,"# moldyn total momentum logfile\n");
126                 moldyn->lvstat|=MOLDYN_LVSTAT_TOTAL_M;
127         }
128
129         if(moldyn->swrite)
130                 moldyn->lvstat|=MOLDYN_LVSTAT_SAVE;
131
132         if(moldyn->dwrite) {
133                 moldyn->dfd=open(moldyn->dfb,O_WRONLY|O_CREAT|O_TRUNC);
134                 if(moldyn->dfd<0) {
135                         perror("[moldyn] dfd open");
136                         return moldyn->dfd;
137                 }
138                 write(moldyn->dfd,moldyn,sizeof(t_moldyn));
139                 moldyn->lvstat|=MOLDYN_LVSTAT_DUMP;
140         }
141
142         if((moldyn->vwrite)&&(vis)) {
143                 moldyn->visual=vis;
144                 visual_init(vis,moldyn->vfb);
145                 moldyn->lvstat|=MOLDYN_LVSTAT_VISUAL;
146         }
147
148         moldyn->lvstat|=MOLDYN_LVSTAT_INITIALIZED;
149
150         return 0;
151 }
152
153 int moldyn_shutdown(t_moldyn *moldyn) {
154
155         if(moldyn->efd) close(moldyn->efd);
156         if(moldyn->mfd) close(moldyn->efd);
157         if(moldyn->dfd) close(moldyn->efd);
158         if(moldyn->visual) visual_tini(moldyn->visual);
159
160         return 0;
161 }
162
163 int create_lattice(unsigned char type,int element,double mass,double lc,
164                    int a,int b,int c,t_atom **atom) {
165
166         int count;
167         int ret;
168         t_3dvec origin;
169
170         count=a*b*c;
171
172         if(type==FCC) count*=4;
173         if(type==DIAMOND) count*=8;
174
175         *atom=malloc(count*sizeof(t_atom));
176         if(*atom==NULL) {
177                 perror("malloc (atoms)");
178                 return -1;
179         }
180
181         v3_zero(&origin);
182
183         switch(type) {
184                 case FCC:
185                         ret=fcc_init(a,b,c,lc,*atom,&origin);
186                         break;
187                 case DIAMOND:
188                         ret=diamond_init(a,b,c,lc,*atom,&origin);
189                         break;
190                 default:
191                         printf("unknown lattice type (%02x)\n",type);
192                         return -1;
193         }
194
195         /* debug */
196         if(ret!=count) {
197                 printf("ok, there is something wrong ...\n");
198                 printf("calculated -> %d atoms\n",count);
199                 printf("created -> %d atoms\n",ret);
200                 return -1;
201         }
202
203         while(count) {
204                 (*atom)[count-1].element=element;
205                 (*atom)[count-1].mass=mass;
206                 count-=1;
207         }
208
209         return ret;
210 }
211
212 int destroy_lattice(t_atom *atom) {
213
214         if(atom) free(atom);
215
216         return 0;
217 }
218
219 int thermal_init(t_moldyn *moldyn,t_random *random,int count) {
220
221         /*
222          * - gaussian distribution of velocities
223          * - zero total momentum
224          * - velocity scaling (E = 3/2 N k T), E: kinetic energy
225          */
226
227         int i;
228         double v,sigma;
229         t_3dvec p_total,delta;
230         t_atom *atom;
231
232         atom=moldyn->atom;
233
234         /* gaussian distribution of velocities */
235         v3_zero(&p_total);
236         for(i=0;i<count;i++) {
237                 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t/atom[i].mass);
238                 /* x direction */
239                 v=sigma*rand_get_gauss(random);
240                 atom[i].v.x=v;
241                 p_total.x+=atom[i].mass*v;
242                 /* y direction */
243                 v=sigma*rand_get_gauss(random);
244                 atom[i].v.y=v;
245                 p_total.y+=atom[i].mass*v;
246                 /* z direction */
247                 v=sigma*rand_get_gauss(random);
248                 atom[i].v.z=v;
249                 p_total.z+=atom[i].mass*v;
250         }
251
252         /* zero total momentum */
253         v3_scale(&p_total,&p_total,1.0/count);
254         for(i=0;i<count;i++) {
255                 v3_scale(&delta,&p_total,1.0/atom[i].mass);
256                 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
257         }
258
259         /* velocity scaling */
260         scale_velocity(moldyn,count);
261
262         return 0;
263 }
264
265 int scale_velocity(t_moldyn *moldyn,int count) {
266
267         int i;
268         double e,c;
269         t_atom *atom;
270
271         atom=moldyn->atom;
272
273         /*
274          * - velocity scaling (E = 3/2 N k T), E: kinetic energy
275          */
276         e=0.0;
277         for(i=0;i<count;i++)
278                 e+=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
279         c=sqrt((2.0*e)/(3.0*count*K_BOLTZMANN*moldyn->t));
280         for(i=0;i<count;i++)
281                 v3_scale(&(atom[i].v),&(atom[i].v),(1.0/c));
282
283         return 0;
284 }
285
286 double get_e_kin(t_atom *atom,int count) {
287
288         int i;
289         double e;
290
291         e=0.0;
292
293         for(i=0;i<count;i++) {
294                 e+=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
295         }
296
297         return e;
298 }
299
300 double get_e_pot(t_moldyn *moldyn) {
301
302         return(moldyn->potential(moldyn));
303 }
304
305 double get_total_energy(t_moldyn *moldyn) {
306
307         double e;
308
309         e=get_e_kin(moldyn->atom,moldyn->count);
310         e+=get_e_pot(moldyn);
311
312         return e;
313 }
314
315 t_3dvec get_total_p(t_atom *atom, int count) {
316
317         t_3dvec p,p_total;
318         int i;
319
320         v3_zero(&p_total);
321         for(i=0;i<count;i++) {
322                 v3_scale(&p,&(atom[i].v),atom[i].mass);
323                 v3_add(&p_total,&p_total,&p);
324         }
325
326         return p_total;
327 }
328
329 double estimate_time_step(t_moldyn *moldyn,double nn_dist,double t) {
330
331         double tau;
332
333         tau=0.05*nn_dist/(sqrt(3.0*K_BOLTZMANN*t/moldyn->atom[0].mass));
334         tau*=1.0E-9;
335         if(tau<moldyn->tau)
336                 printf("[moldyn] warning: time step  (%f > %.15f)\n",
337                        moldyn->tau,tau);
338
339         return tau;     
340 }
341
342 /*
343  * numerical tricks
344  */
345
346 /* verlet list */
347
348 int verlet_list_init(t_moldyn *moldyn) {
349
350         int i,fd;
351
352         fd=open("/dev/null",O_WRONLY);
353
354         for(i=0;i<moldyn->count;i++)
355                 list_init(&(moldyn->atom[i].verlet),fd);
356
357         moldyn->r_verlet=1.1*moldyn->cutoff;
358         /* +moldyn->tau*\
359                 sqrt(3.0*K_BOLTZMANN*moldyn->t/moldyn->atom[0].mass); */
360
361         printf("debug: r verlet = %.15f\n",moldyn->r_verlet);
362         printf("       r cutoff = %.15f\n",moldyn->cutoff);
363         printf("       dim      = %.15f\n",moldyn->dim.x);
364
365         /* make sure to update the list in the beginning */
366         moldyn->dr_max1=moldyn->r_verlet;
367         moldyn->dr_max2=moldyn->r_verlet;
368
369         return 0;
370 }
371
372 int link_cell_init(t_moldyn *moldyn) {
373
374         t_linkcell *lc;
375
376         lc=&(moldyn->lc);
377
378         /* partitioning the md cell */
379         lc->nx=moldyn->dim.x/moldyn->cutoff;
380         lc->x=moldyn->dim.x/lc->nx;
381         lc->ny=moldyn->dim.y/moldyn->cutoff;
382         lc->y=moldyn->dim.y/lc->ny;
383         lc->nz=moldyn->dim.z/moldyn->cutoff;
384         lc->z=moldyn->dim.z/lc->nz;
385
386         lc->subcell=malloc(lc->nx*lc->ny*lc->nz*sizeof(t_list));
387
388         link_cell_update(moldyn);
389         
390         return 0;
391 }
392
393 int verlet_list_update(t_moldyn *moldyn) {
394
395         int i,j;
396         t_3dvec d;
397         t_atom *atom;
398
399         atom=moldyn->atom;
400
401         puts("debug: list update start");
402
403         for(i=0;i<moldyn->count;i++) {
404                 list_destroy(&(atom[i].verlet));
405                 for(j=0;j<moldyn->count;j++) {
406                         if(i!=j) {
407                                 v3_sub(&d,&(atom[i].r),&(atom[j].r));
408                                 v3_per_bound(&d,&(moldyn->dim));
409                                 if(v3_norm(&d)<=moldyn->r_verlet)
410                                         list_add_immediate_ptr(&(atom[i].verlet),&(atom[j]));
411                         }
412                 }
413         }
414
415         moldyn->dr_max1=0.0;
416         moldyn->dr_max2=0.0;
417
418         puts("debug: list update end");
419
420         return 0;
421 }
422
423 int link_cell_update(t_moldyn *moldyn) {
424
425         int count,i,j,k;
426         int nx,ny,nz;
427         t_atom *atom;
428
429         atom=moldyn->atom;
430         nx=moldyn->lc.nx; ny=moldyn->lc.ny; nz=moldyn->lc.nz;
431
432         for(i=0;i<nx*ny*nz;i++)
433                 list_destroy(&(moldyn->lc.subcell[i]));
434         
435         for(count=0;count<moldyn->count;count++) {
436                 for(i=0;i<nx;i++) {
437                         if((atom[count].r.x>=i*moldyn->dim.x) && \
438                            (atom[count].r.x<(i+1)*moldyn->dim.x))
439                                 break;
440                 }
441                 for(j=0;j<ny;j++) {
442                         if((atom[count].r.y>=j*moldyn->dim.y) && \
443                            (atom[count].r.y<(j+1)*moldyn->dim.y))
444                                 break;
445                 }
446                 for(k=0;k<nz;k++) {
447                         if((atom[count].r.z>=k*moldyn->dim.z) && \
448                            (atom[count].r.z<(k+1)*moldyn->dim.z))
449                                 break;
450                 }
451                 list_add_immediate_ptr(&(moldyn->lc.subcell[i+j*nx+k*nx*ny]),
452                                        &(atom[count]));
453         }
454
455         return 0;
456 }
457
458 int verlet_list_shutdown(t_moldyn *moldyn) {
459
460         int i;
461
462         for(i=0;i<moldyn->count;i++)
463                 list_shutdown(&(moldyn->atom[i].verlet));
464
465         return 0;
466 }
467
468 int link_cell_shutdown(t_moldyn *moldyn) {
469
470         int i;
471         t_linkcell *lc;
472
473         lc=&(moldyn->lc);
474
475         for(i=0;i<lc->nx*lc->ny*lc->nz;i++)
476                 list_shutdown(&(moldyn->lc.subcell[i]));
477
478         return 0;
479 }
480
481 /*
482  *
483  * 'integration of newtons equation' - algorithms
484  *
485  */
486
487 /* start the integration */
488
489 int moldyn_integrate(t_moldyn *moldyn) {
490
491         int i;
492         unsigned int e,m,s,d,v;
493         t_3dvec p;
494         double rlc;
495
496         int fd;
497         char fb[128];
498
499         /* logging & visualization */
500         e=moldyn->ewrite;
501         m=moldyn->mwrite;
502         s=moldyn->swrite;
503         d=moldyn->dwrite;
504         v=moldyn->vwrite;
505
506         /* verlet list */       
507         rlc=moldyn->r_verlet-moldyn->cutoff;
508
509         if(!(moldyn->lvstat&MOLDYN_LVSTAT_INITIALIZED)) {
510                 printf("[moldyn] warning, lv system not initialized\n");
511                 return -1;
512         }
513
514         /* create the verlet list */
515         verlet_list_update(moldyn);
516
517         /* calculate initial forces */
518         moldyn->force(moldyn);
519
520         for(i=0;i<moldyn->time_steps;i++) {
521                 /* check for velet list update */
522                 printf(".");
523                 if(moldyn->dr_max1+moldyn->dr_max2>rlc) {
524                         printf("\n");
525                         verlet_list_update(moldyn);
526                 }
527
528                 /* integration step */
529                 moldyn->integrate(moldyn);
530
531                 /* check for log & visualiziation */
532                 if(e) {
533                         if(!(i%e))
534                                 dprintf(moldyn->efd,
535                                         "%.15f %.45f\n",i*moldyn->tau,
536                                         get_total_energy(moldyn));
537                 }
538                 if(m) {
539                         if(!(i%m)) {
540                                 p=get_total_p(moldyn->atom,moldyn->count);
541                                 dprintf(moldyn->mfd,
542                                         "%.15f %.45f\n",i*moldyn->tau,
543                                         v3_norm(&p));
544                         }
545                 }
546                 if(s) {
547                         if(!(i%s)) {
548                                 snprintf(fb,128,"%s-%f-%.15f.save",moldyn->sfb,
549                                          moldyn->t,i*moldyn->tau);
550                                 fd=open(fb,O_WRONLY|O_TRUNC|O_CREAT);
551                                 if(fd<0) perror("[moldyn] save fd open");
552                                 else {
553                                         write(fd,moldyn,sizeof(t_moldyn));
554                                         write(fd,moldyn->atom,
555                                               moldyn->count*sizeof(t_atom));
556                                 }
557                         }       
558                 }
559                 if(d) {
560                         if(!(i%d))
561                                 write(moldyn->dfd,moldyn->atom,
562                                       moldyn->count*sizeof(t_atom));
563
564                 }
565                 if(v) {
566                         if(!(i%v))
567                                 visual_atoms(moldyn->visual,i*moldyn->tau,
568                                              moldyn->atom,moldyn->count);
569                 }
570         }
571
572         return 0;
573 }
574
575 /* velocity verlet */
576
577 int velocity_verlet(t_moldyn *moldyn) {
578
579         int i,count;
580         double tau,tau_square,dr;
581         t_3dvec delta;
582         t_atom *atom;
583
584         atom=moldyn->atom;
585         count=moldyn->count;
586         tau=moldyn->tau;
587
588         tau_square=tau*tau;
589
590         for(i=0;i<count;i++) {
591                 /* new positions */
592                 v3_scale(&delta,&(atom[i].v),tau);
593                 v3_add(&(atom[i].r),&(atom[i].r),&delta);
594                 v3_add(&(atom[i].dr),&(atom[i].dr),&delta);
595                 v3_scale(&delta,&(atom[i].f),0.5*tau_square/atom[i].mass);
596                 v3_add(&(atom[i].r),&(atom[i].r),&delta);
597                 v3_add(&(atom[i].dr),&(atom[i].dr),&delta);
598                 v3_per_bound(&(atom[i].r),&(moldyn->dim));
599
600                 /* set maximum dr (possible list update) */
601                 dr=v3_norm(&(atom[i].dr));
602                 if(dr>moldyn->dr_max1) {
603                         moldyn->dr_max2=moldyn->dr_max1;
604                         moldyn->dr_max1=dr;
605                 }
606                 else if(dr>moldyn->dr_max2) moldyn->dr_max2=dr;
607
608                 /* velocities */
609                 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
610                 v3_add(&(atom[i].v),&(atom[i].v),&delta);
611         }
612
613         /* forces depending on chosen potential */
614         moldyn->force(moldyn);
615
616         for(i=0;i<count;i++) {
617                 /* again velocities */
618                 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
619                 v3_add(&(atom[i].v),&(atom[i].v),&delta);
620         }
621
622         return 0;
623 }
624
625
626 /*
627  *
628  * potentials & corresponding forces
629  * 
630  */
631
632 /* harmonic oscillator potential and force */
633
634 double potential_harmonic_oscillator(t_moldyn *moldyn) {
635
636         t_ho_params *params;
637         t_atom *atom;
638         int i,j;
639         int count;
640         t_3dvec distance;
641         double d,u;
642         double sc,equi_dist;
643
644         params=moldyn->pot_params;
645         atom=moldyn->atom;
646         sc=params->spring_constant;
647         equi_dist=params->equilibrium_distance;
648         count=moldyn->count;
649
650         u=0.0;
651         for(i=0;i<count;i++) {
652                 for(j=0;j<i;j++) {
653                         v3_sub(&distance,&(atom[i].r),&(atom[j].r));
654                         d=v3_norm(&distance);
655                         u+=(0.5*sc*(d-equi_dist)*(d-equi_dist));
656                 }
657         }
658
659         return u;
660 }
661
662 int force_harmonic_oscillator(t_moldyn *moldyn) {
663
664         t_ho_params *params;
665         int i,j,count;
666         t_atom *atom;
667         t_3dvec distance;
668         t_3dvec force;
669         double d;
670         double sc,equi_dist;
671
672         atom=moldyn->atom;      
673         count=moldyn->count;
674         params=moldyn->pot_params;
675         sc=params->spring_constant;
676         equi_dist=params->equilibrium_distance;
677
678         for(i=0;i<count;i++) v3_zero(&(atom[i].f));
679
680         for(i=0;i<count;i++) {
681                 for(j=0;j<i;j++) {
682                         v3_sub(&distance,&(atom[i].r),&(atom[j].r));
683                         v3_per_bound(&distance,&(moldyn->dim));
684                         d=v3_norm(&distance);
685                         if(d<=moldyn->cutoff) {
686                                 v3_scale(&force,&distance,
687                                          -sc*(1.0-(equi_dist/d)));
688                                 v3_add(&(atom[i].f),&(atom[i].f),&force);
689                                 v3_sub(&(atom[j].f),&(atom[j].f),&force);
690                         }
691                 }
692         }
693
694         return 0;
695 }
696
697
698 /* lennard jones potential & force for one sort of atoms */
699  
700 double potential_lennard_jones(t_moldyn *moldyn) {
701
702         t_lj_params *params;
703         t_atom *atom,*btom;
704         int i;
705         int count;
706         t_3dvec distance;
707         double d,help;
708         double u;
709         double eps,sig6,sig12;
710
711         params=moldyn->pot_params;
712         atom=moldyn->atom;
713         count=moldyn->count;
714         eps=params->epsilon4;
715         sig6=params->sigma6;
716         sig12=params->sigma12;
717
718         u=0.0;
719         for(i=0;i<count;i++) {
720                 list_reset(&(atom[i].verlet));
721                 if(atom[i].verlet.current==NULL) continue;
722                 while(1) {
723                         btom=atom[i].verlet.current->data;
724                         v3_sub(&distance,&(atom[i].r),&(btom->r));
725                         v3_per_bound(&distance,&(moldyn->dim));
726                         d=1.0/v3_absolute_square(&distance);    /* 1/r^2 */
727                         help=d*d;                               /* 1/r^4 */
728                         help*=d;                                /* 1/r^6 */
729                         d=help*help;                            /* 1/r^12 */
730                         u+=eps*(sig12*d-sig6*help);
731                         if(list_next(&(atom[i].verlet))==L_NO_NEXT_ELEMENT)
732                                 break;
733                 }
734         }
735         
736         return u;
737 }
738
739 int force_lennard_jones(t_moldyn *moldyn) {
740
741         t_lj_params *params;
742         int i,count;
743         t_atom *atom,*btom;
744         t_3dvec distance;
745         t_3dvec force;
746         double d,h1,h2;
747         double eps,sig6,sig12;
748
749         atom=moldyn->atom;      
750         count=moldyn->count;
751         params=moldyn->pot_params;
752         eps=params->epsilon4;
753         sig6=6*params->sigma6;
754         sig12=12*params->sigma12;
755
756         for(i=0;i<count;i++) v3_zero(&(atom[i].f));
757
758         for(i=0;i<count;i++) {
759                 list_reset(&(atom[i].verlet));
760                 if(atom[i].verlet.current==NULL) continue;
761                 while(1) {
762                         btom=atom[i].verlet.current->data;
763                         v3_sub(&distance,&(atom[i].r),&(btom->r));
764                         v3_per_bound(&distance,&(moldyn->dim));
765                         d=v3_absolute_square(&distance);
766                         if(d<=moldyn->cutoff_square) {
767                                 h1=1.0/d;                       /* 1/r^2 */
768                                 d=h1*h1;                        /* 1/r^4 */
769                                 h2=d*d;                         /* 1/r^8 */
770                                 h1*=d;                          /* 1/r^6 */
771                                 h1*=h2;                         /* 1/r^14 */
772                                 h1*=sig12;
773                                 h2*=sig6;
774                                 /* actually there would be a '-',       *
775                                  * but f=-d/dr potential                */
776                                 d=h1+h2;
777                                 d*=eps;
778                                 v3_scale(&force,&distance,d);
779                                 v3_add(&(atom[i].f),&(atom[i].f),&force);
780                                 //v3_sub(&(atom[j].f),&(atom[j].f),&force);
781                         }
782                         if(list_next(&(atom[i].verlet))==L_NO_NEXT_ELEMENT)
783                                 break;
784                 }
785         }
786
787         return 0;
788 }
789