2 * moldyn.c - molecular dynamics library main file
4 * author: Frank Zirkelbach <frank.zirkelbach@physik.uni-augsburg.de>
12 #include <sys/types.h>
21 #include "report/report.h"
24 * global variables, pse and atom colors (only needed here)
27 static char *pse_name[]={
49 static char *pse_col[]={
72 * the moldyn functions
75 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
77 printf("[moldyn] init\n");
79 memset(moldyn,0,sizeof(t_moldyn));
84 rand_init(&(moldyn->random),NULL,1);
85 moldyn->random.status|=RAND_STAT_VERBOSE;
90 int moldyn_shutdown(t_moldyn *moldyn) {
92 printf("[moldyn] shutdown\n");
94 moldyn_log_shutdown(moldyn);
95 link_cell_shutdown(moldyn);
96 rand_close(&(moldyn->random));
102 int set_int_alg(t_moldyn *moldyn,u8 algo) {
104 printf("[moldyn] integration algorithm: ");
107 case MOLDYN_INTEGRATE_VERLET:
108 moldyn->integrate=velocity_verlet;
109 printf("velocity verlet\n");
112 printf("unknown integration algorithm: %02x\n",algo);
120 int set_cutoff(t_moldyn *moldyn,double cutoff) {
122 moldyn->cutoff=cutoff;
124 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
129 int set_bondlen(t_moldyn *moldyn,double b0,double b1,double bm) {
131 moldyn->bondlen[0]=b0*b0;
132 moldyn->bondlen[1]=b1*b1;
134 moldyn->bondlen[2]=b0*b1;
136 moldyn->bondlen[2]=bm*bm;
141 int set_temperature(t_moldyn *moldyn,double t_ref) {
145 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
150 int set_pressure(t_moldyn *moldyn,double p_ref) {
154 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
159 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
161 moldyn->pt_scale=(ptype|ttype);
165 printf("[moldyn] p/t scaling:\n");
167 printf(" p: %s",ptype?"yes":"no ");
169 printf(" | type: %02x | factor: %f",ptype,ptc);
172 printf(" t: %s",ttype?"yes":"no ");
174 printf(" | type: %02x | factor: %f",ttype,ttc);
180 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
186 moldyn->volume=x*y*z;
194 moldyn->dv=0.000001*moldyn->volume;
196 printf("[moldyn] dimensions in A and A^3 respectively:\n");
197 printf(" x: %f\n",moldyn->dim.x);
198 printf(" y: %f\n",moldyn->dim.y);
199 printf(" z: %f\n",moldyn->dim.z);
200 printf(" volume: %f\n",moldyn->volume);
201 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
202 printf(" delta volume (pressure calc): %f\n",moldyn->dv);
207 int set_nn_dist(t_moldyn *moldyn,double dist) {
214 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
216 printf("[moldyn] periodic boundary conditions:\n");
219 moldyn->status|=MOLDYN_STAT_PBX;
222 moldyn->status|=MOLDYN_STAT_PBY;
225 moldyn->status|=MOLDYN_STAT_PBZ;
227 printf(" x: %s\n",x?"yes":"no");
228 printf(" y: %s\n",y?"yes":"no");
229 printf(" z: %s\n",z?"yes":"no");
234 int set_potential1b(t_moldyn *moldyn,pf_func1b func) {
241 int set_potential2b(t_moldyn *moldyn,pf_func2b func) {
248 int set_potential3b_j1(t_moldyn *moldyn,pf_func2b func) {
250 moldyn->func3b_j1=func;
255 int set_potential3b_j2(t_moldyn *moldyn,pf_func2b func) {
257 moldyn->func3b_j2=func;
262 int set_potential3b_j3(t_moldyn *moldyn,pf_func2b func) {
264 moldyn->func3b_j3=func;
269 int set_potential3b_k1(t_moldyn *moldyn,pf_func3b func) {
271 moldyn->func3b_k1=func;
276 int set_potential3b_k2(t_moldyn *moldyn,pf_func3b func) {
278 moldyn->func3b_k2=func;
283 int set_potential_params(t_moldyn *moldyn,void *params) {
285 moldyn->pot_params=params;
290 int set_avg_skip(t_moldyn *moldyn,int skip) {
292 printf("[moldyn] skip %d steps before starting average calc\n",skip);
293 moldyn->avg_skip=skip;
298 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
300 strncpy(moldyn->vlsdir,dir,127);
305 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
307 strncpy(moldyn->rauthor,author,63);
308 strncpy(moldyn->rtitle,title,63);
313 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
318 printf("[moldyn] set log: ");
321 case LOG_TOTAL_ENERGY:
322 moldyn->ewrite=timer;
323 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
324 moldyn->efd=open(filename,
325 O_WRONLY|O_CREAT|O_EXCL,
328 perror("[moldyn] energy log fd open");
331 dprintf(moldyn->efd,"# total energy log file\n");
332 printf("total energy (%d)\n",timer);
334 case LOG_TOTAL_MOMENTUM:
335 moldyn->mwrite=timer;
336 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
337 moldyn->mfd=open(filename,
338 O_WRONLY|O_CREAT|O_EXCL,
341 perror("[moldyn] momentum log fd open");
344 dprintf(moldyn->efd,"# total momentum log file\n");
345 printf("total momentum (%d)\n",timer);
348 moldyn->pwrite=timer;
349 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
350 moldyn->pfd=open(filename,
351 O_WRONLY|O_CREAT|O_EXCL,
354 perror("[moldyn] pressure log file\n");
357 dprintf(moldyn->pfd,"# pressure log file\n");
358 printf("pressure (%d)\n",timer);
360 case LOG_TEMPERATURE:
361 moldyn->twrite=timer;
362 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
363 moldyn->tfd=open(filename,
364 O_WRONLY|O_CREAT|O_EXCL,
367 perror("[moldyn] temperature log file\n");
370 dprintf(moldyn->tfd,"# temperature log file\n");
371 printf("temperature (%d)\n",timer);
374 moldyn->swrite=timer;
375 printf("save file (%d)\n",timer);
378 moldyn->vwrite=timer;
379 ret=visual_init(moldyn,moldyn->vlsdir);
381 printf("[moldyn] visual init failure\n");
384 printf("visual file (%d)\n",timer);
387 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
388 moldyn->rfd=open(filename,
389 O_WRONLY|O_CREAT|O_EXCL,
392 perror("[moldyn] report fd open");
395 printf("report -> ");
397 snprintf(filename,127,"%s/e_plot.scr",
399 moldyn->epfd=open(filename,
400 O_WRONLY|O_CREAT|O_EXCL,
403 perror("[moldyn] energy plot fd open");
406 dprintf(moldyn->epfd,e_plot_script);
411 snprintf(filename,127,"%s/pressure_plot.scr",
413 moldyn->ppfd=open(filename,
414 O_WRONLY|O_CREAT|O_EXCL,
417 perror("[moldyn] p plot fd open");
420 dprintf(moldyn->ppfd,pressure_plot_script);
425 snprintf(filename,127,"%s/temperature_plot.scr",
427 moldyn->tpfd=open(filename,
428 O_WRONLY|O_CREAT|O_EXCL,
431 perror("[moldyn] t plot fd open");
434 dprintf(moldyn->tpfd,temperature_plot_script);
436 printf("temperature ");
438 dprintf(moldyn->rfd,report_start,
439 moldyn->rauthor,moldyn->rtitle);
443 printf("unknown log type: %02x\n",type);
450 int moldyn_log_shutdown(t_moldyn *moldyn) {
454 printf("[moldyn] log shutdown\n");
458 dprintf(moldyn->rfd,report_energy);
459 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
464 if(moldyn->mfd) close(moldyn->mfd);
468 dprintf(moldyn->rfd,report_pressure);
469 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
476 dprintf(moldyn->rfd,report_temperature);
477 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
482 dprintf(moldyn->rfd,report_end);
484 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
487 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
490 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
499 * creating lattice functions
502 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
503 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin) {
514 /* how many atoms do we expect */
515 if(type==CUBIC) new*=1;
516 if(type==FCC) new*=4;
517 if(type==DIAMOND) new*=8;
519 /* allocate space for atoms */
520 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
522 perror("[moldyn] realloc (create lattice)");
526 atom=&(moldyn->atom[count]);
528 /* no atoms on the boundaries (only reason: it looks better!) */
542 set_nn_dist(moldyn,lc);
543 ret=cubic_init(a,b,c,lc,atom,&orig);
547 v3_scale(&orig,&orig,0.5);
548 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
549 ret=fcc_init(a,b,c,lc,atom,&orig);
553 v3_scale(&orig,&orig,0.25);
554 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
555 ret=diamond_init(a,b,c,lc,atom,&orig);
558 printf("unknown lattice type (%02x)\n",type);
564 printf("[moldyn] creating lattice failed\n");
565 printf(" amount of atoms\n");
566 printf(" - expected: %d\n",new);
567 printf(" - created: %d\n",ret);
572 printf("[moldyn] created lattice with %d atoms\n",new);
574 for(ret=0;ret<new;ret++) {
575 atom[ret].element=element;
578 atom[ret].brand=brand;
579 atom[ret].tag=count+ret;
580 check_per_bound(moldyn,&(atom[ret].r));
581 atom[ret].r_0=atom[ret].r;
584 /* update total system mass */
585 total_mass_calc(moldyn);
590 int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
591 t_3dvec *r,t_3dvec *v) {
598 count=(moldyn->count)++; // asshole style!
600 ptr=realloc(atom,(count+1)*sizeof(t_atom));
602 perror("[moldyn] realloc (add atom)");
609 /* initialize new atom */
610 memset(&(atom[count]),0,sizeof(t_atom));
613 atom[count].element=element;
614 atom[count].mass=mass;
615 atom[count].brand=brand;
616 atom[count].tag=count;
617 atom[count].attr=attr;
618 check_per_bound(moldyn,&(atom[count].r));
619 atom[count].r_0=atom[count].r;
621 /* update total system mass */
622 total_mass_calc(moldyn);
627 int del_atom(t_moldyn *moldyn,int tag) {
634 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
636 perror("[moldyn]malloc (del atom)");
640 for(cnt=0;cnt<tag;cnt++)
643 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
645 new[cnt-1].tag=cnt-1;
657 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
676 v3_copy(&(atom[count].r),&r);
685 for(i=0;i<count;i++) {
686 atom[i].r.x-=(a*lc)/2.0;
687 atom[i].r.y-=(b*lc)/2.0;
688 atom[i].r.z-=(c*lc)/2.0;
694 /* fcc lattice init */
695 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
708 /* construct the basis */
709 memset(basis,0,3*sizeof(t_3dvec));
717 /* fill up the room */
725 v3_copy(&(atom[count].r),&r);
728 /* the three face centered atoms */
730 v3_add(&n,&r,&basis[l]);
731 v3_copy(&(atom[count].r),&n);
740 /* coordinate transformation */
741 for(i=0;i<count;i++) {
742 atom[i].r.x-=(a*lc)/2.0;
743 atom[i].r.y-=(b*lc)/2.0;
744 atom[i].r.z-=(c*lc)/2.0;
750 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
755 count=fcc_init(a,b,c,lc,atom,origin);
761 if(origin) v3_add(&o,&o,origin);
763 count+=fcc_init(a,b,c,lc,&atom[count],&o);
768 int destroy_atoms(t_moldyn *moldyn) {
770 if(moldyn->atom) free(moldyn->atom);
775 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
778 * - gaussian distribution of velocities
779 * - zero total momentum
780 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
785 t_3dvec p_total,delta;
790 random=&(moldyn->random);
792 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
794 /* gaussian distribution of velocities */
796 for(i=0;i<moldyn->count;i++) {
797 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
799 v=sigma*rand_get_gauss(random);
801 p_total.x+=atom[i].mass*v;
803 v=sigma*rand_get_gauss(random);
805 p_total.y+=atom[i].mass*v;
807 v=sigma*rand_get_gauss(random);
809 p_total.z+=atom[i].mass*v;
812 /* zero total momentum */
813 v3_scale(&p_total,&p_total,1.0/moldyn->count);
814 for(i=0;i<moldyn->count;i++) {
815 v3_scale(&delta,&p_total,1.0/atom[i].mass);
816 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
819 /* velocity scaling */
820 scale_velocity(moldyn,equi_init);
825 double total_mass_calc(t_moldyn *moldyn) {
831 for(i=0;i<moldyn->count;i++)
832 moldyn->mass+=moldyn->atom[i].mass;
837 double temperature_calc(t_moldyn *moldyn) {
839 /* assume up to date kinetic energy, which is 3/2 N k_B T */
841 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
846 double get_temperature(t_moldyn *moldyn) {
851 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
861 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
864 /* get kinetic energy / temperature & count involved atoms */
867 for(i=0;i<moldyn->count;i++) {
868 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
869 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
874 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
875 else return 0; /* no atoms involved in scaling! */
877 /* (temporary) hack for e,t = 0 */
880 if(moldyn->t_ref!=0.0) {
881 thermal_init(moldyn,equi_init);
885 return 0; /* no scaling needed */
889 /* get scaling factor */
890 scale=moldyn->t_ref/moldyn->t;
894 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
895 scale=1.0+(scale-1.0)/moldyn->t_tc;
898 /* velocity scaling */
899 for(i=0;i<moldyn->count;i++) {
900 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
901 v3_scale(&(atom[i].v),&(atom[i].v),scale);
907 double ideal_gas_law_pressure(t_moldyn *moldyn) {
911 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
916 double virial_sum(t_moldyn *moldyn) {
922 /* virial (sum over atom virials) */
924 for(i=0;i<moldyn->count;i++) {
925 virial=&(moldyn->atom[i].virial);
926 v+=(virial->xx+virial->yy+virial->zz);
930 /* global virial (absolute coordinates) */
931 virial=&(moldyn->gvir);
932 moldyn->gv=virial->xx+virial->yy+virial->zz;
934 return moldyn->virial;
937 double pressure_calc(t_moldyn *moldyn) {
941 * with W = 1/3 sum_i f_i r_i (- skipped!)
942 * virial = sum_i f_i r_i
944 * => P = (2 Ekin + virial) / (3V)
947 /* assume up to date virial & up to date kinetic energy */
949 /* pressure (atom virials) */
950 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
951 moldyn->p/=(3.0*moldyn->volume);
953 /* pressure (absolute coordinates) */
954 moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
955 moldyn->gp/=(3.0*moldyn->volume);
960 int average_reset(t_moldyn *moldyn) {
962 printf("[moldyn] average reset\n");
964 /* update skip value */
965 moldyn->avg_skip=moldyn->total_steps;
971 /* potential energy */
979 moldyn->virial_sum=0.0;
989 int average_and_fluctuation_calc(t_moldyn *moldyn) {
993 if(moldyn->total_steps<moldyn->avg_skip)
996 denom=moldyn->total_steps+1-moldyn->avg_skip;
998 /* assume up to date energies, temperature, pressure etc */
1000 /* kinetic energy */
1001 moldyn->k_sum+=moldyn->ekin;
1002 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1003 moldyn->k_avg=moldyn->k_sum/denom;
1004 moldyn->k2_avg=moldyn->k2_sum/denom;
1005 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1007 /* potential energy */
1008 moldyn->v_sum+=moldyn->energy;
1009 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1010 moldyn->v_avg=moldyn->v_sum/denom;
1011 moldyn->v2_avg=moldyn->v2_sum/denom;
1012 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1015 moldyn->t_sum+=moldyn->t;
1016 moldyn->t_avg=moldyn->t_sum/denom;
1019 moldyn->virial_sum+=moldyn->virial;
1020 moldyn->virial_avg=moldyn->virial_sum/denom;
1021 moldyn->gv_sum+=moldyn->gv;
1022 moldyn->gv_avg=moldyn->gv_sum/denom;
1025 moldyn->p_sum+=moldyn->p;
1026 moldyn->p_avg=moldyn->p_sum/denom;
1027 moldyn->gp_sum+=moldyn->gp;
1028 moldyn->gp_avg=moldyn->gp_sum/denom;
1033 int get_heat_capacity(t_moldyn *moldyn) {
1037 /* averages needed for heat capacity calc */
1038 if(moldyn->total_steps<moldyn->avg_skip)
1041 /* (temperature average)^2 */
1042 temp2=moldyn->t_avg*moldyn->t_avg;
1043 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1046 /* ideal gas contribution */
1047 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1048 printf(" ideal gas contribution: %f\n",
1049 ighc/moldyn->mass*KILOGRAM/JOULE);
1051 /* specific heat for nvt ensemble */
1052 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1053 moldyn->c_v_nvt/=moldyn->mass;
1055 /* specific heat for nve ensemble */
1056 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1057 moldyn->c_v_nve/=moldyn->mass;
1059 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1060 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1061 printf(" --> <dV2> sim: %f experimental: %f\n",moldyn->dv2_avg,1.5*moldyn->count*K_B2*moldyn->t_avg*moldyn->t_avg*(1.0-1.5*moldyn->count*K_BOLTZMANN/(700*moldyn->mass*JOULE/KILOGRAM)));
1066 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1070 double u_up,u_down,dv;
1082 dv=8*scale*scale*scale*moldyn->volume;
1084 store=malloc(moldyn->count*sizeof(t_atom));
1086 printf("[moldyn] allocating store mem failed\n");
1090 /* save unscaled potential energy + atom/dim configuration */
1091 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1094 /* scale up dimension and atom positions */
1095 scale_dim(moldyn,SCALE_UP,scale,TRUE,TRUE,TRUE);
1096 scale_atoms(moldyn,SCALE_UP,scale,TRUE,TRUE,TRUE);
1097 link_cell_shutdown(moldyn);
1098 link_cell_init(moldyn,QUIET);
1099 potential_force_calc(moldyn);
1100 u_up=moldyn->energy;
1102 /* restore atomic configuration + dim */
1103 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1106 /* scale down dimension and atom positions */
1107 scale_dim(moldyn,SCALE_DOWN,scale,TRUE,TRUE,TRUE);
1108 scale_atoms(moldyn,SCALE_DOWN,scale,TRUE,TRUE,TRUE);
1109 link_cell_shutdown(moldyn);
1110 link_cell_init(moldyn,QUIET);
1111 potential_force_calc(moldyn);
1112 u_down=moldyn->energy;
1114 /* calculate pressure */
1115 p=-(u_up-u_down)/dv;
1116 printf("-------> %.10f %.10f %f\n",u_up/EV/moldyn->count,u_down/EV/moldyn->count,p/BAR);
1118 /* restore atomic configuration + dim */
1119 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1122 /* restore energy */
1123 potential_force_calc(moldyn);
1125 link_cell_shutdown(moldyn);
1126 link_cell_init(moldyn,QUIET);
1131 double get_pressure(t_moldyn *moldyn) {
1137 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1149 if(x) dim->x*=scale;
1150 if(y) dim->y*=scale;
1151 if(z) dim->z*=scale;
1156 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1167 for(i=0;i<moldyn->count;i++) {
1168 r=&(moldyn->atom[i].r);
1177 int scale_volume(t_moldyn *moldyn) {
1183 vdim=&(moldyn->vis.dim);
1187 /* scaling factor */
1188 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1189 scale=1.0-(moldyn->p_ref-moldyn->p)/moldyn->p_tc;
1190 scale=pow(scale,ONE_THIRD);
1193 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1195 moldyn->debug=scale;
1197 /* scale the atoms and dimensions */
1198 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1199 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1201 /* visualize dimensions */
1208 /* recalculate scaled volume */
1209 moldyn->volume=dim->x*dim->y*dim->z;
1211 /* adjust/reinit linkcell */
1212 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1213 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1214 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1215 link_cell_shutdown(moldyn);
1216 link_cell_init(moldyn,QUIET);
1227 double e_kin_calc(t_moldyn *moldyn) {
1235 for(i=0;i<moldyn->count;i++) {
1236 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1237 moldyn->ekin+=atom[i].ekin;
1240 return moldyn->ekin;
1243 double get_total_energy(t_moldyn *moldyn) {
1245 return(moldyn->ekin+moldyn->energy);
1248 t_3dvec get_total_p(t_moldyn *moldyn) {
1257 for(i=0;i<moldyn->count;i++) {
1258 v3_scale(&p,&(atom[i].v),atom[i].mass);
1259 v3_add(&p_total,&p_total,&p);
1265 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1269 /* nn_dist is the nearest neighbour distance */
1271 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1280 /* linked list / cell method */
1282 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1289 /* partitioning the md cell */
1290 lc->nx=moldyn->dim.x/moldyn->cutoff;
1291 lc->x=moldyn->dim.x/lc->nx;
1292 lc->ny=moldyn->dim.y/moldyn->cutoff;
1293 lc->y=moldyn->dim.y/lc->ny;
1294 lc->nz=moldyn->dim.z/moldyn->cutoff;
1295 lc->z=moldyn->dim.z/lc->nz;
1296 lc->cells=lc->nx*lc->ny*lc->nz;
1299 lc->subcell=malloc(lc->cells*sizeof(int*));
1301 lc->subcell=malloc(lc->cells*sizeof(t_list));
1304 if(lc->subcell==NULL) {
1305 perror("[moldyn] cell init (malloc)");
1310 printf("[moldyn] FATAL: less then 27 subcells!\n");
1314 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1317 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1320 printf(" x: %d x %f A\n",lc->nx,lc->x);
1321 printf(" y: %d x %f A\n",lc->ny,lc->y);
1322 printf(" z: %d x %f A\n",lc->nz,lc->z);
1327 for(i=0;i<lc->cells;i++) {
1328 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1329 if(lc->subcell[i]==NULL) {
1330 perror("[moldyn] list init (malloc)");
1335 printf(" ---> %d malloc %p (%p)\n",
1336 i,lc->subcell[0],lc->subcell);
1340 for(i=0;i<lc->cells;i++)
1341 list_init_f(&(lc->subcell[i]));
1344 /* update the list */
1345 link_cell_update(moldyn);
1350 int link_cell_update(t_moldyn *moldyn) {
1366 for(i=0;i<lc->cells;i++)
1368 memset(lc->subcell[i],0,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1370 list_destroy_f(&(lc->subcell[i]));
1373 for(count=0;count<moldyn->count;count++) {
1374 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1375 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1376 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1380 while(lc->subcell[i+j*nx+k*nx*ny][p]!=0)
1383 if(p>=MAX_ATOMS_PER_LIST) {
1384 printf("[moldyn] FATAL: amount of atoms too high!\n");
1388 lc->subcell[i+j*nx+k*nx*ny][p]=count;
1390 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nx*ny]),
1394 printf(" ---> %d %d malloc %p (%p)\n",
1395 i,count,lc->subcell[i].current,lc->subcell);
1403 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1427 if(i>=nx||j>=ny||k>=nz)
1428 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1431 cell[0]=lc->subcell[i+j*nx+k*a];
1432 for(ci=-1;ci<=1;ci++) {
1435 if((x<0)||(x>=nx)) {
1439 for(cj=-1;cj<=1;cj++) {
1442 if((y<0)||(y>=ny)) {
1446 for(ck=-1;ck<=1;ck++) {
1449 if((z<0)||(z>=nz)) {
1453 if(!(ci|cj|ck)) continue;
1455 cell[--count2]=lc->subcell[x+y*nx+z*a];
1458 cell[count1++]=lc->subcell[x+y*nx+z*a];
1469 int link_cell_shutdown(t_moldyn *moldyn) {
1476 for(i=0;i<lc->cells;i++) {
1478 free(lc->subcell[i]);
1480 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1481 list_destroy_f(&(lc->subcell[i]));
1490 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1494 t_moldyn_schedule *schedule;
1496 schedule=&(moldyn->schedule);
1497 count=++(schedule->total_sched);
1499 ptr=realloc(schedule->runs,count*sizeof(int));
1501 perror("[moldyn] realloc (runs)");
1505 schedule->runs[count-1]=runs;
1507 ptr=realloc(schedule->tau,count*sizeof(double));
1509 perror("[moldyn] realloc (tau)");
1513 schedule->tau[count-1]=tau;
1515 printf("[moldyn] schedule added:\n");
1516 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1522 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1524 moldyn->schedule.hook=hook;
1525 moldyn->schedule.hook_params=hook_params;
1532 * 'integration of newtons equation' - algorithms
1536 /* start the integration */
1538 int moldyn_integrate(t_moldyn *moldyn) {
1541 unsigned int e,m,s,v,p,t;
1543 t_moldyn_schedule *sched;
1548 double energy_scale;
1549 struct timeval t1,t2;
1552 sched=&(moldyn->schedule);
1555 /* initialize linked cell method */
1556 link_cell_init(moldyn,VERBOSE);
1558 /* logging & visualization */
1566 /* sqaure of some variables */
1567 moldyn->tau_square=moldyn->tau*moldyn->tau;
1568 moldyn->cutoff_square=moldyn->cutoff*moldyn->cutoff;
1570 /* get current time */
1571 gettimeofday(&t1,NULL);
1573 /* calculate initial forces */
1574 potential_force_calc(moldyn);
1579 /* some stupid checks before we actually start calculating bullshit */
1580 if(moldyn->cutoff>0.5*moldyn->dim.x)
1581 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
1582 if(moldyn->cutoff>0.5*moldyn->dim.y)
1583 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
1584 if(moldyn->cutoff>0.5*moldyn->dim.z)
1585 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
1586 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
1587 if(ds>0.05*moldyn->nnd)
1588 printf("[moldyn] WARNING: forces too high / tau too small!\n");
1590 /* zero absolute time */
1592 moldyn->total_steps=0;
1594 /* debugging, ignore */
1597 /* tell the world */
1598 printf("[moldyn] integration start, go get a coffee ...\n");
1600 /* executing the schedule */
1602 while(sched->count<sched->total_sched) {
1604 /* setting amount of runs and finite time step size */
1605 moldyn->tau=sched->tau[sched->count];
1606 moldyn->tau_square=moldyn->tau*moldyn->tau;
1607 moldyn->time_steps=sched->runs[sched->count];
1609 /* energy scaling factor (might change!) */
1610 energy_scale=moldyn->count*EV;
1612 /* integration according to schedule */
1614 for(i=0;i<moldyn->time_steps;i++) {
1616 /* integration step */
1617 moldyn->integrate(moldyn);
1619 /* calculate kinetic energy, temperature and pressure */
1621 temperature_calc(moldyn);
1623 pressure_calc(moldyn);
1624 average_and_fluctuation_calc(moldyn);
1627 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
1628 scale_velocity(moldyn,FALSE);
1629 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
1630 scale_volume(moldyn);
1632 /* check for log & visualization */
1634 if(!(moldyn->total_steps%e))
1635 dprintf(moldyn->efd,
1637 moldyn->time,moldyn->ekin/energy_scale,
1638 moldyn->energy/energy_scale,
1639 get_total_energy(moldyn)/energy_scale);
1642 if(!(moldyn->total_steps%m)) {
1643 momentum=get_total_p(moldyn);
1644 dprintf(moldyn->mfd,
1645 "%f %f %f %f %f\n",moldyn->time,
1646 momentum.x,momentum.y,momentum.z,
1647 v3_norm(&momentum));
1651 if(!(moldyn->total_steps%p)) {
1652 dprintf(moldyn->pfd,
1653 "%f %f %f %f %f\n",moldyn->time,
1654 moldyn->p/BAR,moldyn->p_avg/BAR,
1655 moldyn->gp/BAR,moldyn->gp_avg/BAR);
1659 if(!(moldyn->total_steps%t)) {
1660 dprintf(moldyn->tfd,
1662 moldyn->time,moldyn->t,moldyn->t_avg);
1666 if(!(moldyn->total_steps%s)) {
1667 snprintf(dir,128,"%s/s-%07.f.save",
1668 moldyn->vlsdir,moldyn->time);
1669 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
1671 if(fd<0) perror("[moldyn] save fd open");
1673 write(fd,moldyn,sizeof(t_moldyn));
1674 write(fd,moldyn->atom,
1675 moldyn->count*sizeof(t_atom));
1681 if(!(moldyn->total_steps%v)) {
1682 visual_atoms(moldyn);
1686 /* display progress */
1687 //if(!(moldyn->total_steps%10)) {
1688 /* get current time */
1689 gettimeofday(&t2,NULL);
1691 printf("\rsched:%d, steps:%d/%d, T:%3.1f/%3.1f P:%4.1f/%4.1f V:%6.1f (%d)",
1692 sched->count,i,moldyn->total_steps,
1693 moldyn->t,moldyn->t_avg,
1694 moldyn->p_avg/BAR,moldyn->gp_avg/BAR,
1696 (int)(t2.tv_sec-t1.tv_sec));
1700 /* copy over time */
1704 /* increase absolute time */
1705 moldyn->time+=moldyn->tau;
1706 moldyn->total_steps+=1;
1710 /* check for hooks */
1712 printf("\n ## schedule hook %d start ##\n",
1714 sched->hook(moldyn,sched->hook_params);
1715 printf(" ## schedule hook end ##\n");
1718 /* increase the schedule counter */
1726 /* velocity verlet */
1728 int velocity_verlet(t_moldyn *moldyn) {
1731 double tau,tau_square,h;
1736 count=moldyn->count;
1738 tau_square=moldyn->tau_square;
1740 for(i=0;i<count;i++) {
1741 /* check whether fixed atom */
1742 if(atom[i].attr&ATOM_ATTR_FP)
1746 v3_scale(&delta,&(atom[i].v),tau);
1747 v3_add(&(atom[i].r),&(atom[i].r),&delta);
1748 v3_scale(&delta,&(atom[i].f),h*tau_square);
1749 v3_add(&(atom[i].r),&(atom[i].r),&delta);
1750 check_per_bound(moldyn,&(atom[i].r));
1752 /* velocities [actually v(t+tau/2)] */
1753 v3_scale(&delta,&(atom[i].f),h*tau);
1754 v3_add(&(atom[i].v),&(atom[i].v),&delta);
1757 /* criticial check */
1758 moldyn_bc_check(moldyn);
1760 /* neighbour list update */
1761 link_cell_update(moldyn);
1763 /* forces depending on chosen potential */
1764 potential_force_calc(moldyn);
1766 for(i=0;i<count;i++) {
1767 /* check whether fixed atom */
1768 if(atom[i].attr&ATOM_ATTR_FP)
1770 /* again velocities [actually v(t+tau)] */
1771 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
1772 v3_add(&(atom[i].v),&(atom[i].v),&delta);
1781 * potentials & corresponding forces & virial routine
1785 /* generic potential and force calculation */
1787 int potential_force_calc(t_moldyn *moldyn) {
1790 t_atom *itom,*jtom,*ktom;
1794 int *neighbour_i[27];
1798 t_list neighbour_i[27];
1799 t_list neighbour_i2[27];
1805 count=moldyn->count;
1815 /* reset global virial */
1816 memset(&(moldyn->gvir),0,sizeof(t_virial));
1818 /* reset force, site energy and virial of every atom */
1819 for(i=0;i<count;i++) {
1822 v3_zero(&(itom[i].f));
1825 virial=(&(itom[i].virial));
1833 /* reset site energy */
1838 /* get energy, force and virial of every atom */
1840 /* first (and only) loop over atoms i */
1841 for(i=0;i<count;i++) {
1843 /* single particle potential/force */
1844 if(itom[i].attr&ATOM_ATTR_1BP)
1846 moldyn->func1b(moldyn,&(itom[i]));
1848 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
1851 /* 2 body pair potential/force */
1853 link_cell_neighbour_index(moldyn,
1854 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
1855 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
1856 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
1861 /* first loop over atoms j */
1862 if(moldyn->func2b) {
1869 while(neighbour_i[j][p]!=0) {
1871 jtom=&(atom[neighbour_i[j][p]]);
1874 if(jtom==&(itom[i]))
1877 if((jtom->attr&ATOM_ATTR_2BP)&
1878 (itom[i].attr&ATOM_ATTR_2BP)) {
1879 moldyn->func2b(moldyn,
1886 this=&(neighbour_i[j]);
1889 if(this->start==NULL)
1893 jtom=this->current->data;
1895 if(jtom==&(itom[i]))
1898 if((jtom->attr&ATOM_ATTR_2BP)&
1899 (itom[i].attr&ATOM_ATTR_2BP)) {
1900 moldyn->func2b(moldyn,
1905 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
1911 /* 3 body potential/force */
1913 if(!(itom[i].attr&ATOM_ATTR_3BP))
1916 /* copy the neighbour lists */
1918 /* no copy needed for static lists */
1920 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
1923 /* second loop over atoms j */
1930 while(neighbour_i[j][p]!=0) {
1932 jtom=&(atom[neighbour_i[j][p]]);
1935 this=&(neighbour_i[j]);
1938 if(this->start==NULL)
1943 jtom=this->current->data;
1946 if(jtom==&(itom[i]))
1949 if(!(jtom->attr&ATOM_ATTR_3BP))
1955 if(moldyn->func3b_j1)
1956 moldyn->func3b_j1(moldyn,
1961 /* in first j loop, 3bp run can be skipped */
1962 if(!(moldyn->run3bp))
1965 /* first loop over atoms k */
1966 if(moldyn->func3b_k1) {
1974 while(neighbour_i[j][q]!=0) {
1976 ktom=&(atom[neighbour_i[k][q]]);
1979 that=&(neighbour_i2[k]);
1982 if(that->start==NULL)
1986 ktom=that->current->data;
1989 if(!(ktom->attr&ATOM_ATTR_3BP))
1995 if(ktom==&(itom[i]))
1998 moldyn->func3b_k1(moldyn,
2006 } while(list_next_f(that)!=\
2014 if(moldyn->func3b_j2)
2015 moldyn->func3b_j2(moldyn,
2020 /* second loop over atoms k */
2021 if(moldyn->func3b_k2) {
2029 while(neighbour_i[j][q]!=0) {
2031 ktom=&(atom[neighbour_i[k][q]]);
2034 that=&(neighbour_i2[k]);
2037 if(that->start==NULL)
2041 ktom=that->current->data;
2044 if(!(ktom->attr&ATOM_ATTR_3BP))
2050 if(ktom==&(itom[i]))
2053 moldyn->func3b_k2(moldyn,
2062 } while(list_next_f(that)!=\
2070 /* 2bp post function */
2071 if(moldyn->func3b_j3) {
2072 moldyn->func3b_j3(moldyn,
2079 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2094 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2095 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2097 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2098 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2099 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2103 /* some postprocessing */
2104 for(i=0;i<count;i++) {
2105 /* calculate global virial */
2106 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2107 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2108 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2109 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2110 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2111 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2113 /* check forces regarding the given timestep */
2114 if(v3_norm(&(itom[i].f))>\
2115 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2116 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2124 * virial calculation
2127 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2128 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2130 a->virial.xx+=f->x*d->x;
2131 a->virial.yy+=f->y*d->y;
2132 a->virial.zz+=f->z*d->z;
2133 a->virial.xy+=f->x*d->y;
2134 a->virial.xz+=f->x*d->z;
2135 a->virial.yz+=f->y*d->z;
2141 * periodic boundary checking
2144 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2145 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2156 if(moldyn->status&MOLDYN_STAT_PBX) {
2157 if(a->x>=x) a->x-=dim->x;
2158 else if(-a->x>x) a->x+=dim->x;
2160 if(moldyn->status&MOLDYN_STAT_PBY) {
2161 if(a->y>=y) a->y-=dim->y;
2162 else if(-a->y>y) a->y+=dim->y;
2164 if(moldyn->status&MOLDYN_STAT_PBZ) {
2165 if(a->z>=z) a->z-=dim->z;
2166 else if(-a->z>z) a->z+=dim->z;
2173 * debugging / critical check functions
2176 int moldyn_bc_check(t_moldyn *moldyn) {
2189 for(i=0;i<moldyn->count;i++) {
2190 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2191 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2192 i,atom[i].r.x,dim->x/2);
2193 printf("diagnostic:\n");
2194 printf("-----------\natom.r.x:\n");
2196 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2199 ((byte)&(1<<k))?1:0,
2202 printf("---------------\nx=dim.x/2:\n");
2204 memcpy(&byte,(u8 *)(&x)+j,1);
2207 ((byte)&(1<<k))?1:0,
2210 if(atom[i].r.x==x) printf("the same!\n");
2211 else printf("different!\n");
2213 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2214 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2215 i,atom[i].r.y,dim->y/2);
2216 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2217 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2218 i,atom[i].r.z,dim->z/2);
2228 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2235 fd=open(file,O_RDONLY);
2237 perror("[moldyn] load save file open");
2241 fsize=lseek(fd,0,SEEK_END);
2242 lseek(fd,0,SEEK_SET);
2244 size=sizeof(t_moldyn);
2247 cnt=read(fd,moldyn,size);
2249 perror("[moldyn] load save file read (moldyn)");
2255 size=moldyn->count*sizeof(t_atom);
2257 /* correcting possible atom data offset */
2259 if(fsize!=sizeof(t_moldyn)+size) {
2260 corr=fsize-sizeof(t_moldyn)-size;
2261 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2262 printf(" moifying offset:\n");
2263 printf(" - current pos: %d\n",sizeof(t_moldyn));
2264 printf(" - atom size: %d\n",size);
2265 printf(" - file size: %d\n",fsize);
2266 printf(" => correction: %d\n",corr);
2267 lseek(fd,corr,SEEK_CUR);
2270 moldyn->atom=(t_atom *)malloc(size);
2271 if(moldyn->atom==NULL) {
2272 perror("[moldyn] load save file malloc (atoms)");
2277 cnt=read(fd,moldyn->atom,size);
2279 perror("[moldyn] load save file read (atoms)");
2290 int moldyn_load(t_moldyn *moldyn) {
2298 * post processing functions
2301 int get_line(int fd,char *line,int max) {
2308 if(count==max) return count;
2309 ret=read(fd,line+count,1);
2310 if(ret<=0) return ret;
2311 if(line[count]=='\n') {
2319 int pair_correlation_init(t_moldyn *moldyn,double dr) {
2325 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
2341 for(i=0;i<moldyn->count;i++) {
2343 v3_sub(&dist,&(atom[i].r),&(atom[i].r_0));
2344 check_per_bound(moldyn,&dist);
2345 d2=v3_absolute_square(&dist);
2361 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
2371 t_list neighbour[27];
2380 unsigned char ibrand;
2384 slots=2.0*moldyn->cutoff/dr;
2387 if(slots*dr<=moldyn->cutoff)
2388 printf("[moldyn] WARNING: pcc (low #slots)\n");
2390 printf("[moldyn] pair correlation calc info:\n");
2391 printf(" time: %f\n",moldyn->time);
2392 printf(" count: %d\n",moldyn->count);
2393 printf(" cutoff: %f\n",moldyn->cutoff);
2394 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
2400 stat=(double *)malloc(3*slots*sizeof(double));
2402 perror("[moldyn] pair correlation malloc");
2407 memset(stat,0,3*slots*sizeof(double));
2409 link_cell_init(moldyn,VERBOSE);
2413 for(i=0;i<moldyn->count;i++) {
2414 /* neighbour indexing */
2415 link_cell_neighbour_index(moldyn,
2416 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2417 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2418 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
2421 /* brand of atom i */
2422 ibrand=itom[i].brand;
2426 bc=(j<lc->dnlc)?0:1;
2431 while(neighbour[j][p]!=0) {
2433 jtom=&(moldyn->atom[neighbour[j][p]]);
2436 this=&(neighbour[j]);
2439 if(this->start==NULL)
2444 jtom=this->current->data;
2446 /* only count pairs once,
2447 * skip same atoms */
2448 if(itom[i].tag>=jtom->tag)
2452 * pair correlation calc
2456 v3_sub(&dist,&(jtom->r),&(itom[i].r));
2457 if(bc) check_per_bound(moldyn,&dist);
2458 d=v3_absolute_square(&dist);
2460 /* ignore if greater or equal cutoff */
2461 if(d>=4.0*moldyn->cutoff_square)
2464 /* fill the slots */
2468 /* should never happen but it does 8) -
2469 * related to -ffloat-store problem! */
2471 printf("[moldyn] WARNING: pcc (%d/%d)",
2477 if(ibrand!=jtom->brand) {
2482 /* type a - type a bonds */
2486 /* type b - type b bonds */
2492 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2498 for(i=1;i<slots;i++) {
2499 // normalization: 4 pi r^2 dr
2500 // here: not double counting pairs -> 2 pi r r dr
2501 // ... and actually it's a constant times r^2
2504 stat[slots+i]/=norm;
2510 /* todo: store/print pair correlation function */
2516 link_cell_shutdown(moldyn);
2521 int analyze_bonds(t_moldyn *moldyn) {
2530 * visualization code
2533 int visual_init(t_moldyn *moldyn,char *filebase) {
2535 strncpy(moldyn->vis.fb,filebase,128);
2540 int visual_atoms(t_moldyn *moldyn) {
2554 t_list neighbour[27];
2570 sprintf(file,"%s/atomic_conf_%07.f.xyz",v->fb,moldyn->time);
2571 fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
2573 perror("open visual save file fd");
2577 /* write the actual data file */
2580 dprintf(fd,"# [P] %d %07.f <%f,%f,%f>\n",
2581 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
2583 // atomic configuration
2584 for(i=0;i<moldyn->count;i++) {
2585 // atom type, positions, color and kinetic energy
2586 dprintf(fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
2590 pse_col[atom[i].element],
2594 * bond detection should usually be done by potential
2595 * functions. brrrrr! EVIL!
2597 * todo: potentials need to export a 'find_bonds' function!
2600 // bonds between atoms
2601 if(!(atom[i].attr&ATOM_ATTR_VB))
2603 link_cell_neighbour_index(moldyn,
2604 (atom[i].r.x+moldyn->dim.x/2)/lc->x,
2605 (atom[i].r.y+moldyn->dim.y/2)/lc->y,
2606 (atom[i].r.z+moldyn->dim.z/2)/lc->z,
2612 while(neighbour[j][p]!=0) {
2613 btom=&(atom[neighbour[j][p]]);
2616 list_reset_f(&neighbour[j]);
2617 if(neighbour[j].start==NULL)
2620 btom=neighbour[j].current->data;
2622 if(btom==&atom[i]) // skip identical atoms
2624 //if(btom<&atom[i]) // skip half of them
2626 v3_sub(&dist,&(atom[i].r),&(btom->r));
2627 if(bc) check_per_bound(moldyn,&dist);
2628 d2=v3_absolute_square(&dist);
2629 brand=atom[i].brand;
2630 if(brand==btom->brand) {
2631 if(d2>moldyn->bondlen[brand])
2635 if(d2>moldyn->bondlen[2])
2638 dprintf(fd,"# [B] %f %f %f %f %f %f\n",
2639 atom[i].r.x,atom[i].r.y,atom[i].r.z,
2640 btom->r.x,btom->r.y,btom->r.z);
2644 } while(list_next_f(&neighbour[j])!=L_NO_NEXT_ELEMENT);
2651 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2652 -dim.x/2,-dim.y/2,-dim.z/2,
2653 dim.x/2,-dim.y/2,-dim.z/2);
2654 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2655 -dim.x/2,-dim.y/2,-dim.z/2,
2656 -dim.x/2,dim.y/2,-dim.z/2);
2657 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2658 dim.x/2,dim.y/2,-dim.z/2,
2659 dim.x/2,-dim.y/2,-dim.z/2);
2660 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2661 -dim.x/2,dim.y/2,-dim.z/2,
2662 dim.x/2,dim.y/2,-dim.z/2);
2664 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2665 -dim.x/2,-dim.y/2,dim.z/2,
2666 dim.x/2,-dim.y/2,dim.z/2);
2667 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2668 -dim.x/2,-dim.y/2,dim.z/2,
2669 -dim.x/2,dim.y/2,dim.z/2);
2670 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2671 dim.x/2,dim.y/2,dim.z/2,
2672 dim.x/2,-dim.y/2,dim.z/2);
2673 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2674 -dim.x/2,dim.y/2,dim.z/2,
2675 dim.x/2,dim.y/2,dim.z/2);
2677 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2678 -dim.x/2,-dim.y/2,dim.z/2,
2679 -dim.x/2,-dim.y/2,-dim.z/2);
2680 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2681 -dim.x/2,dim.y/2,dim.z/2,
2682 -dim.x/2,dim.y/2,-dim.z/2);
2683 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2684 dim.x/2,-dim.y/2,dim.z/2,
2685 dim.x/2,-dim.y/2,-dim.z/2);
2686 dprintf(fd,"# [D] %f %f %f %f %f %f\n",
2687 dim.x/2,dim.y/2,dim.z/2,
2688 dim.x/2,dim.y/2,-dim.z/2);