2 * moldyn.c - molecular dynamics library main file
4 * author: Frank Zirkelbach <frank.zirkelbach@physik.uni-augsburg.de>
12 #include <sys/types.h>
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"
26 int moldyn_usage(char **argv) {
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("-C <cutoff radius> [m] (%.15f)\n",MOLDYN_CUTOFF);
39 printf("-R <runs> (%d)\n",MOLDYN_RUNS);
40 printf(" -- integration algo --\n");
41 printf(" -I <number> (%d)\n",MOLDYN_INTEGRATE_DEFAULT);
42 printf(" 0: velocity verlet\n");
43 printf(" -- potential --\n");
44 printf(" -P <number> <param1 param2 ...>\n");
45 printf(" 0: harmonic oscillator\n");
46 printf(" param1: spring constant\n");
47 printf(" param2: equilibrium distance\n");
48 printf(" 1: lennard jones\n");
49 printf(" param1: epsilon\n");
50 printf(" param2: sigma\n");
56 int moldyn_parse_argv(t_moldyn *moldyn,int argc,char **argv) {
63 memset(moldyn,0,sizeof(t_moldyn));
66 moldyn->t=MOLDYN_TEMP;
67 moldyn->tau=MOLDYN_TAU;
68 moldyn->time_steps=MOLDYN_RUNS;
69 moldyn->integrate=velocity_verlet;
70 moldyn->potential_force_function=lennard_jones;
77 moldyn->ewrite=atoi(argv[++i]);
78 strncpy(moldyn->efb,argv[++i],64);
81 moldyn->mwrite=atoi(argv[++i]);
82 strncpy(moldyn->mfb,argv[++i],64);
85 moldyn->dwrite=atoi(argv[++i]);
86 strncpy(moldyn->dfb,argv[++i],64);
89 moldyn->swrite=atoi(argv[++i]);
90 strncpy(moldyn->sfb,argv[++i],64);
93 moldyn->vwrite=atoi(argv[++i]);
94 strncpy(moldyn->vfb,argv[++i],64);
97 moldyn->t=atof(argv[++i]);
100 moldyn->tau=atof(argv[++i]);
103 moldyn->cutoff=atof(argv[++i]);
106 moldyn->time_steps=atoi(argv[++i]);
109 /* integration algorithm */
110 switch(atoi(argv[++i])) {
111 case MOLDYN_INTEGRATE_VERLET:
112 moldyn->integrate=velocity_verlet;
115 printf("unknown integration algo %s\n",argv[i]);
121 /* potential + params */
122 switch(atoi(argv[++i])) {
123 case MOLDYN_POTENTIAL_HO:
124 hop.spring_constant=atof(argv[++i]);
125 hop.equilibrium_distance=atof(argv[++i]);
126 moldyn->pot_params=malloc(sizeof(t_ho_params));
127 memcpy(moldyn->pot_params,&hop,sizeof(t_ho_params));
128 moldyn->potential_force_function=harmonic_oscillator;
130 case MOLDYN_POTENTIAL_LJ:
134 ljp.sigma6=s*s*s*s*s*s;
135 ljp.sigma12=ljp.sigma6*ljp.sigma6;
136 moldyn->pot_params=malloc(sizeof(t_lj_params));
137 memcpy(moldyn->pot_params,&ljp,sizeof(t_lj_params));
138 moldyn->potential_force_function=lennard_jones;
141 printf("unknown potential %s\n",argv[i]);
147 printf("unknown option %s\n",argv[i]);
160 int moldyn_log_init(t_moldyn *moldyn) {
168 moldyn->efd=open(moldyn->efb,O_WRONLY|O_CREAT|O_TRUNC);
170 perror("[moldyn] efd open");
173 dprintf(moldyn->efd,"# moldyn total energy logfile\n");
174 moldyn->lvstat|=MOLDYN_LVSTAT_TOTAL_E;
178 moldyn->mfd=open(moldyn->mfb,O_WRONLY|O_CREAT|O_TRUNC);
180 perror("[moldyn] mfd open");
183 dprintf(moldyn->mfd,"# moldyn total momentum logfile\n");
184 moldyn->lvstat|=MOLDYN_LVSTAT_TOTAL_M;
188 moldyn->lvstat|=MOLDYN_LVSTAT_SAVE;
191 moldyn->dfd=open(moldyn->dfb,O_WRONLY|O_CREAT|O_TRUNC);
193 perror("[moldyn] dfd open");
196 write(moldyn->dfd,moldyn,sizeof(t_moldyn));
197 moldyn->lvstat|=MOLDYN_LVSTAT_DUMP;
200 if((moldyn->vwrite)&&(vis)) {
202 visual_init(vis,moldyn->vfb);
203 moldyn->lvstat|=MOLDYN_LVSTAT_VISUAL;
206 moldyn->lvstat|=MOLDYN_LVSTAT_INITIALIZED;
211 int moldyn_log_shutdown(t_moldyn *moldyn) {
213 if(moldyn->efd) close(moldyn->efd);
214 if(moldyn->mfd) close(moldyn->efd);
215 if(moldyn->dfd) close(moldyn->efd);
216 if(moldyn->visual) visual_tini(moldyn->visual);
221 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
225 ret=moldyn_parse_argv(moldyn,argc,argv);
226 if(ret<0) return ret;
228 ret=moldyn_log_init(moldyn);
229 if(ret<0) return ret;
231 rand_init(&(moldyn->random),NULL,1);
232 moldyn->random.status|=RAND_STAT_VERBOSE;
239 int moldyn_shutdown(t_moldyn *moldyn) {
241 moldyn_log_shutdown(moldyn);
242 rand_close(&(moldyn->random));
248 int create_lattice(unsigned char type,int element,double mass,double lc,
249 int a,int b,int c,t_atom **atom) {
257 if(type==FCC) count*=4;
258 if(type==DIAMOND) count*=8;
260 *atom=malloc(count*sizeof(t_atom));
262 perror("malloc (atoms)");
270 ret=fcc_init(a,b,c,lc,*atom,&origin);
273 ret=diamond_init(a,b,c,lc,*atom,&origin);
276 printf("unknown lattice type (%02x)\n",type);
282 printf("ok, there is something wrong ...\n");
283 printf("calculated -> %d atoms\n",count);
284 printf("created -> %d atoms\n",ret);
289 (*atom)[count-1].element=element;
290 (*atom)[count-1].mass=mass;
297 int destroy_lattice(t_atom *atom) {
304 int thermal_init(t_moldyn *moldyn) {
307 * - gaussian distribution of velocities
308 * - zero total momentum
309 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
314 t_3dvec p_total,delta;
319 random=&(moldyn->random);
321 /* gaussian distribution of velocities */
323 for(i=0;i<moldyn->count;i++) {
324 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t/atom[i].mass);
326 v=sigma*rand_get_gauss(random);
328 p_total.x+=atom[i].mass*v;
330 v=sigma*rand_get_gauss(random);
332 p_total.y+=atom[i].mass*v;
334 v=sigma*rand_get_gauss(random);
336 p_total.z+=atom[i].mass*v;
339 /* zero total momentum */
340 v3_scale(&p_total,&p_total,1.0/moldyn->count);
341 for(i=0;i<moldyn->count;i++) {
342 v3_scale(&delta,&p_total,1.0/atom[i].mass);
343 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
346 /* velocity scaling */
347 scale_velocity(moldyn);
352 int scale_velocity(t_moldyn *moldyn) {
361 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
364 for(i=0;i<moldyn->count;i++)
365 e+=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
366 c=sqrt((2.0*e)/(3.0*moldyn->count*K_BOLTZMANN*moldyn->t));
367 for(i=0;i<moldyn->count;i++)
368 v3_scale(&(atom[i].v),&(atom[i].v),(1.0/c));
373 double get_e_kin(t_atom *atom,int count) {
380 for(i=0;i<count;i++) {
381 e+=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
387 double get_e_pot(t_moldyn *moldyn) {
389 return moldyn->energy;
392 double get_total_energy(t_moldyn *moldyn) {
396 e=get_e_kin(moldyn->atom,moldyn->count);
397 e+=get_e_pot(moldyn);
402 t_3dvec get_total_p(t_atom *atom, int count) {
408 for(i=0;i<count;i++) {
409 v3_scale(&p,&(atom[i].v),atom[i].mass);
410 v3_add(&p_total,&p_total,&p);
416 double estimate_time_step(t_moldyn *moldyn,double nn_dist,double t) {
420 tau=0.05*nn_dist/(sqrt(3.0*K_BOLTZMANN*t/moldyn->atom[0].mass));
423 printf("[moldyn] warning: time step (%f > %.15f)\n",
433 /* linked list / cell method */
435 int link_cell_init(t_moldyn *moldyn) {
443 lc->listfd=open("/dev/null",O_WRONLY);
445 /* partitioning the md cell */
446 lc->nx=moldyn->dim.x/moldyn->cutoff;
447 lc->x=moldyn->dim.x/lc->nx;
448 lc->ny=moldyn->dim.y/moldyn->cutoff;
449 lc->y=moldyn->dim.y/lc->ny;
450 lc->nz=moldyn->dim.z/moldyn->cutoff;
451 lc->z=moldyn->dim.z/lc->nz;
453 lc->cells=lc->nx*lc->ny*lc->nz;
454 lc->subcell=malloc(lc->cells*sizeof(t_list));
456 printf("initializing linked cells (%d)\n",lc->cells);
458 for(i=0;i<lc->cells;i++)
459 //list_init(&(lc->subcell[i]),1);
460 list_init(&(lc->subcell[i]));
462 link_cell_update(moldyn);
467 int link_cell_update(t_moldyn *moldyn) {
481 for(i=0;i<lc->cells;i++)
482 list_destroy(&(moldyn->lc.subcell[i]));
484 for(count=0;count<moldyn->count;count++) {
485 i=(atom[count].r.x+(moldyn->dim.x/2))/lc->x;
486 j=(atom[count].r.y+(moldyn->dim.y/2))/lc->y;
487 k=(atom[count].r.z+(moldyn->dim.z/2))/lc->z;
488 list_add_immediate_ptr(&(moldyn->lc.subcell[i+j*nx+k*nx*ny]),
495 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,t_list *cell) {
503 unsigned char bx,by,bz;
514 cell[0]=lc->subcell[i+j*nx+k*a];
515 for(ci=-1;ci<=1;ci++) {
522 for(cj=-1;cj<=1;cj++) {
529 for(ck=-1;ck<=1;ck++) {
536 if(!(ci|cj|ck)) continue;
538 cell[--count2]=lc->subcell[x+y*nx+z*a];
541 cell[count1++]=lc->subcell[x+y*nx+z*a];
550 int link_cell_shutdown(t_moldyn *moldyn) {
557 for(i=0;i<lc->nx*lc->ny*lc->nz;i++)
558 list_shutdown(&(moldyn->lc.subcell[i]));
560 if(lc->listfd) close(lc->listfd);
567 * 'integration of newtons equation' - algorithms
571 /* start the integration */
573 int moldyn_integrate(t_moldyn *moldyn) {
576 unsigned int e,m,s,d,v;
582 /* initialize linked cell method */
583 link_cell_init(moldyn);
585 /* logging & visualization */
592 if(!(moldyn->lvstat&MOLDYN_LVSTAT_INITIALIZED)) {
593 printf("[moldyn] warning, lv system not initialized\n");
597 /* sqaure of some variables */
598 moldyn->tau_square=moldyn->tau*moldyn->tau;
599 moldyn->cutoff_square=moldyn->cutoff*moldyn->cutoff;
601 /* calculate initial forces */
602 moldyn->potential_force_function(moldyn);
604 for(i=0;i<moldyn->time_steps;i++) {
606 /* integration step */
607 moldyn->integrate(moldyn);
609 /* check for log & visualization */
613 "%.15f %.45f\n",i*moldyn->tau,
614 get_total_energy(moldyn));
618 p=get_total_p(moldyn->atom,moldyn->count);
620 "%.15f %.45f\n",i*moldyn->tau,
626 snprintf(fb,128,"%s-%f-%.15f.save",moldyn->sfb,
627 moldyn->t,i*moldyn->tau);
628 fd=open(fb,O_WRONLY|O_TRUNC|O_CREAT);
629 if(fd<0) perror("[moldyn] save fd open");
631 write(fd,moldyn,sizeof(t_moldyn));
632 write(fd,moldyn->atom,
633 moldyn->count*sizeof(t_atom));
639 write(moldyn->dfd,moldyn->atom,
640 moldyn->count*sizeof(t_atom));
645 visual_atoms(moldyn->visual,i*moldyn->tau,
646 moldyn->atom,moldyn->count);
647 printf("\rsteps: %d",i);
656 /* velocity verlet */
658 int velocity_verlet(t_moldyn *moldyn) {
661 double tau,tau_square;
668 tau_square=moldyn->tau_square;
670 for(i=0;i<count;i++) {
672 v3_scale(&delta,&(atom[i].v),tau);
673 v3_add(&(atom[i].r),&(atom[i].r),&delta);
674 v3_scale(&delta,&(atom[i].f),0.5*tau_square/atom[i].mass);
675 v3_add(&(atom[i].r),&(atom[i].r),&delta);
676 v3_per_bound(&(atom[i].r),&(moldyn->dim));
679 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
680 v3_add(&(atom[i].v),&(atom[i].v),&delta);
683 /* neighbour list update */
684 printf("list update ...\n");
685 link_cell_update(moldyn);
688 /* forces depending on chosen potential */
689 printf("calc potential/force ...\n");
690 moldyn->potential_force_function(moldyn);
693 for(i=0;i<count;i++) {
694 /* again velocities */
695 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
696 v3_add(&(atom[i].v),&(atom[i].v),&delta);
705 * potentials & corresponding forces
709 /* harmonic oscillator potential and force */
711 int harmonic_oscillator(t_moldyn *moldyn) {
716 t_list *this,neighbour[27];
719 t_3dvec force,distance;
724 params=moldyn->pot_params;
727 sc=params->spring_constant;
728 equi_dist=params->equilibrium_distance;
731 /* reset energy counter */
734 for(i=0;i<count;i++) {
736 v3_zero(&(atom[i].f));
738 /* determine cell + neighbours */
739 ni=(atom[i].r.x+(moldyn->dim.x/2))/lc->x;
740 nj=(atom[i].r.y+(moldyn->dim.y/2))/lc->y;
741 nk=(atom[i].r.z+(moldyn->dim.z/2))/lc->z;
742 c=link_cell_neighbour_index(moldyn,ni,nj,nk,neighbour);
745 * processing cell of atom i
746 * => no need to check for empty list (1 element at minimum)
748 this=&(neighbour[0]);
751 btom=this->current->data;
754 v3_sub(&distance,&(atom[i].r),&(btom->r));
755 d=v3_norm(&distance);
756 if(d<=moldyn->cutoff) {
757 u+=(0.5*sc*(d-equi_dist)*(d-equi_dist));
758 v3_scale(&force,&distance,
759 -sc*(1.0-(equi_dist/d)));
760 v3_add(&(atom[i].f),&(atom[i].f),&force);
762 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
765 * direct neighbour cells
766 * => no boundary condition check necessary
769 this=&(neighbour[j]);
770 list_reset(this); /* there might not be a single atom */
771 if(this->start!=NULL) {
774 btom=this->current->data;
775 v3_sub(&distance,&(atom[i].r),&(btom->r));
776 d=v3_norm(&distance);
777 if(d<=moldyn->cutoff) {
778 u+=(0.5*sc*(d-equi_dist)*(d-equi_dist));
779 v3_scale(&force,&distance,
780 -sc*(1.0-(equi_dist/d)));
781 v3_add(&(atom[i].f),&(atom[i].f),
784 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
790 * indirect neighbour cells
791 * => check boundary conditions
794 this=&(neighbour[j]);
795 list_reset(this); /* check boundary conditions */
796 if(this->start!=NULL) {
799 btom=this->current->data;
800 v3_sub(&distance,&(atom[i].r),&(btom->r));
801 v3_per_bound(&distance,&(moldyn->dim));
802 d=v3_norm(&distance);
803 if(d<=moldyn->cutoff) {
804 u+=(0.5*sc*(d-equi_dist)*(d-equi_dist));
805 v3_scale(&force,&distance,
806 -sc*(1.0-(equi_dist/d)));
807 v3_add(&(atom[i].f),&(atom[i].f),
810 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
816 moldyn->energy=0.5*u;
821 /* lennard jones potential & force for one sort of atoms */
823 int lennard_jones(t_moldyn *moldyn) {
828 t_list *this,neighbour[27];
831 t_3dvec force,distance;
833 double eps,sig6,sig12;
836 params=moldyn->pot_params;
840 eps=params->epsilon4;
842 sig12=params->sigma12;
844 /* reset energy counter */
847 for(i=0;i<count;i++) {
849 v3_zero(&(atom[i].f));
851 /* determine cell + neighbours */
852 ni=(atom[i].r.x+(moldyn->dim.x/2))/lc->x;
853 nj=(atom[i].r.y+(moldyn->dim.y/2))/lc->y;
854 nk=(atom[i].r.z+(moldyn->dim.z/2))/lc->z;
855 c=link_cell_neighbour_index(moldyn,ni,nj,nk,neighbour);
857 /* processing cell of atom i */
858 this=&(neighbour[0]);
859 list_reset(this); /* list has 1 element at minimum */
861 btom=this->current->data;
864 v3_sub(&distance,&(atom[i].r),&(btom->r));
865 d=v3_absolute_square(&distance); /* 1/r^2 */
866 if(d<=moldyn->cutoff_square) {
870 h1=h2*h2; /* 1/r^12 */
871 u+=eps*(sig12*h1-sig6*h2);
878 v3_scale(&force,&distance,d);
879 v3_add(&(atom[i].f),&(atom[i].f),&force);
881 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
883 /* neighbours not doing boundary condition overflow */
885 this=&(neighbour[j]);
886 list_reset(this); /* there might not be a single atom */
887 if(this->start!=NULL) {
890 btom=this->current->data;
891 v3_sub(&distance,&(atom[i].r),&(btom->r));
892 d=v3_absolute_square(&distance); /* r^2 */
893 if(d<=moldyn->cutoff_square) {
897 h1=h2*h2; /* 1/r^12 */
898 u+=eps*(sig12*h1-sig6*h2);
905 v3_scale(&force,&distance,d);
906 v3_add(&(atom[i].f),&(atom[i].f),
909 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
914 /* neighbours due to boundary conditions */
916 this=&(neighbour[j]);
917 list_reset(this); /* check boundary conditions */
918 if(this->start!=NULL) {
921 btom=this->current->data;
922 v3_sub(&distance,&(atom[i].r),&(btom->r));
923 v3_per_bound(&distance,&(moldyn->dim));
924 d=v3_absolute_square(&distance); /* r^2 */
925 if(d<=moldyn->cutoff_square) {
929 h1=h2*h2; /* 1/r^12 */
930 u+=eps*(sig12*h1-sig6*h2);
937 v3_scale(&force,&distance,d);
938 v3_add(&(atom[i].f),&(atom[i].f),
941 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
947 moldyn->energy=0.5*u;
952 /* tersoff potential & force for 2 sorts of atoms */
954 int tersoff(t_moldyn *moldyn) {
956 t_tersoff_params *params;
957 t_atom *atom,*btom,*ktom;
959 t_list *this,*thisk,neighbour[27],neighbourk[27];
967 params=moldyn->pot_params;
972 /* reset energy counter */
975 for(i=0;i<count;i++) {
977 v3_zero(&(atom[i].f));
979 /* determin cell neighbours */
980 ni=(atom[i].r.x+(moldyn->dim.x/2))/lc->x;
981 nj=(atom[i].r.y+(moldyn->dim.y/2))/lc->y;
982 nk=(atom[i].r.z+(moldyn->dim.z/2))/lc->z;
983 c=link_cell_neighbour_index(moldyn,ni,nj,nk,neighbour);
986 * processing cell of atom i
987 * => no need to check for empty list (1 element at minimum)
989 this=&(neighbour[0]);
992 btom=this->current->data;
998 v3_sub(&dist_ij,btom,&(atom[i]));
999 d_ij=v3_norm(&dist_ij);
1015 arg1=PI*(d_ij-R)/s_r;
1016 f_c=0.5+0.5*cos(arg1);
1017 df_c=-0.5*sin(arg1)*(PI/(s_r*d_ij));
1018 f_r=A*exp(-lambda*d_ij);
1019 f_a=-B*exp(-mu*d_ij);
1026 /* end 2 body stuff */
1028 /* determine cell neighbours of btom */
1029 ki=(btom->r.x+(moldyn->dim.x/2))/lc->x;
1030 kj=(btom->r.y+(moldyn->dim.y/2))/lc->y;
1031 kk=(btom->r.z+(moldyn->dim.z/2))/lc->z;
1032 ck=link_cell_neighbour_index(moldyn,ki,kj,kk,
1036 thisk=&(neighbourk[0]);
1039 ktom=thisk->current->data;
1042 if(ktom==&(atom[i]))
1045 /* 3 body stuff (1) */
1054 /* end 3 body stuff (1) */
1057 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1059 /* direct neighbours of btom cell */
1061 thisk=&(neighbourk[k]);
1063 if(thisk->start!=NULL) {
1066 ktom=thisk->current->data;
1067 if(ktom==&(atom[i]))
1070 /* 3 body stuff (2) */
1072 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1077 /* indirect neighbours of btom cell */
1078 for(k=ck;k<27;k++) {
1079 thisk=&(neighbourk[k]);
1081 if(thisk->start!=NULL) {
1084 ktom=thisk->current->data;
1085 if(ktom==&(atom[i]))
1090 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1096 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
1099 * direct neighbour cells of atom i
1102 this=&(neighbour[j]);
1104 if(this->start!=NULL) {
1107 btom=this->current->data;
1112 /* determine cell neighbours of btom */
1113 ki=(btom->r.x+(moldyn->dim.x/2))/lc->x;
1114 kj=(btom->r.y+(moldyn->dim.y/2))/lc->y;
1115 kk=(btom->r.z+(moldyn->dim.z/2))/lc->z;
1116 ck=link_cell_neighbour_index(moldyn,ki,kj,kk,
1120 thisk=&(neighbourk[0]);
1123 ktom=thisk->current->data;
1126 if(ktom==&(atom[i]))
1129 /* 3 body stuff (1) */
1131 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1133 /* direct neighbours of btom cell */
1135 thisk=&(neighbourk[k]);
1137 if(thisk->start!=NULL) {
1140 ktom=thisk->current->data;
1141 if(ktom==&(atom[i]))
1144 /* 3 body stuff (2) */
1146 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1151 /* indirect neighbours of btom cell */
1152 for(k=ck;k<27;k++) {
1153 thisk=&(neighbourk[k]);
1155 if(thisk->start!=NULL) {
1158 ktom=thisk->current->data;
1159 if(ktom==&(atom[i]))
1162 /* 3 body stuff (3) */
1164 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1170 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
1176 * indirect neighbour cells of atom i
1179 this=&(neighbour[j]);
1181 if(this->start!=NULL) {
1184 btom=this->current->data;
1189 /* determine cell neighbours of btom */
1190 ki=(btom->r.x+(moldyn->dim.x/2))/lc->x;
1191 kj=(btom->r.y+(moldyn->dim.y/2))/lc->y;
1192 kk=(btom->r.z+(moldyn->dim.z/2))/lc->z;
1193 ck=link_cell_neighbour_index(moldyn,ki,kj,kk,
1197 thisk=&(neighbourk[0]);
1200 ktom=thisk->current->data;
1203 if(ktom==&(atom[i]))
1206 /* 3 body stuff (1) */
1208 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1210 /* direct neighbours of btom cell */
1212 thisk=&(neighbourk[k]);
1214 if(thisk->start!=NULL) {
1217 ktom=thisk->current->data;
1218 if(ktom==&(atom[i]))
1221 /* 3 body stuff (2) */
1223 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1228 /* indirect neighbours of btom cell */
1229 for(k=ck;k<27;k++) {
1230 thisk=&(neighbourk[k]);
1232 if(thisk->start!=NULL) {
1235 ktom=thisk->current->data;
1236 if(ktom==&(atom[i]))
1239 /* 3 body stuff (3) */
1241 } while(list_next(thisk)!=L_NO_NEXT_ELEMENT);
1247 } while(list_next(this)!=L_NO_NEXT_ELEMENT);
1254 moldyn->energy=0.5*u;