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 <fpu_control.h>
31 #include "report/report.h"
33 /* potential includes */
34 #include "potentials/harmonic_oscillator.h"
35 #include "potentials/lennard_jones.h"
36 #include "potentials/albe.h"
38 #include "potentials/tersoff_orig.h"
40 #include "potentials/tersoff.h"
51 * the moldyn functions
54 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
56 printf("[moldyn] init\n");
58 /* only needed if compiled without -msse2 (float-store prob!) */
61 memset(moldyn,0,sizeof(t_moldyn));
66 rand_init(&(moldyn->random),NULL,1);
67 moldyn->random.status|=RAND_STAT_VERBOSE;
72 int moldyn_shutdown(t_moldyn *moldyn) {
74 printf("[moldyn] shutdown\n");
76 moldyn_log_shutdown(moldyn);
77 link_cell_shutdown(moldyn);
78 rand_close(&(moldyn->random));
84 int set_int_alg(t_moldyn *moldyn,u8 algo) {
86 printf("[moldyn] integration algorithm: ");
89 case MOLDYN_INTEGRATE_VERLET:
90 moldyn->integrate=velocity_verlet;
91 printf("velocity verlet\n");
94 printf("unknown integration algorithm: %02x\n",algo);
102 int set_cutoff(t_moldyn *moldyn,double cutoff) {
104 moldyn->cutoff=cutoff;
105 moldyn->cutoff_square=cutoff*cutoff;
107 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
112 int set_temperature(t_moldyn *moldyn,double t_ref) {
116 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
121 int set_pressure(t_moldyn *moldyn,double p_ref) {
125 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
130 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
132 moldyn->pt_scale&=(~(P_SCALE_MASK));
133 moldyn->pt_scale|=ptype;
136 printf("[moldyn] p scaling:\n");
138 printf(" p: %s",ptype?"yes":"no ");
140 printf(" | type: %02x | factor: %f",ptype,ptc);
146 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
148 moldyn->pt_scale&=(~(T_SCALE_MASK));
149 moldyn->pt_scale|=ttype;
152 printf("[moldyn] t scaling:\n");
154 printf(" t: %s",ttype?"yes":"no ");
156 printf(" | type: %02x | factor: %f",ttype,ttc);
162 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
164 moldyn->pt_scale=(ptype|ttype);
168 printf("[moldyn] p/t scaling:\n");
170 printf(" p: %s",ptype?"yes":"no ");
172 printf(" | type: %02x | factor: %f",ptype,ptc);
175 printf(" t: %s",ttype?"yes":"no ");
177 printf(" | type: %02x | factor: %f",ttype,ttc);
183 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
189 moldyn->volume=x*y*z;
197 printf("[moldyn] dimensions in A and A^3 respectively:\n");
198 printf(" x: %f\n",moldyn->dim.x);
199 printf(" y: %f\n",moldyn->dim.y);
200 printf(" z: %f\n",moldyn->dim.z);
201 printf(" volume: %f\n",moldyn->volume);
202 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
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_potential(t_moldyn *moldyn,u8 type) {
237 case MOLDYN_POTENTIAL_TM:
238 moldyn->func1b=tersoff_mult_1bp;
239 moldyn->func3b_j1=tersoff_mult_3bp_j1;
240 moldyn->func3b_k1=tersoff_mult_3bp_k1;
241 moldyn->func3b_j2=tersoff_mult_3bp_j2;
242 moldyn->func3b_k2=tersoff_mult_3bp_k2;
243 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
245 case MOLDYN_POTENTIAL_AM:
246 moldyn->func3b_j1=albe_mult_3bp_j1;
247 moldyn->func3b_k1=albe_mult_3bp_k1;
248 moldyn->func3b_j2=albe_mult_3bp_j2;
249 moldyn->func3b_k2=albe_mult_3bp_k2;
250 moldyn->check_2b_bond=albe_mult_check_2b_bond;
252 case MOLDYN_POTENTIAL_HO:
253 moldyn->func2b=harmonic_oscillator;
254 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
256 case MOLDYN_POTENTIAL_LJ:
257 moldyn->func2b=lennard_jones;
258 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
261 printf("[moldyn] set potential: unknown type %02x\n",
269 int set_avg_skip(t_moldyn *moldyn,int skip) {
271 printf("[moldyn] skip %d steps before starting average calc\n",skip);
272 moldyn->avg_skip=skip;
277 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
279 strncpy(moldyn->vlsdir,dir,127);
284 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
286 strncpy(moldyn->rauthor,author,63);
287 strncpy(moldyn->rtitle,title,63);
292 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
297 printf("[moldyn] set log: ");
300 case LOG_TOTAL_ENERGY:
301 moldyn->ewrite=timer;
302 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
303 moldyn->efd=open(filename,
304 O_WRONLY|O_CREAT|O_EXCL,
307 perror("[moldyn] energy log fd open");
310 dprintf(moldyn->efd,"# total energy log file\n");
311 printf("total energy (%d)\n",timer);
313 case LOG_TOTAL_MOMENTUM:
314 moldyn->mwrite=timer;
315 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
316 moldyn->mfd=open(filename,
317 O_WRONLY|O_CREAT|O_EXCL,
320 perror("[moldyn] momentum log fd open");
323 dprintf(moldyn->efd,"# total momentum log file\n");
324 printf("total momentum (%d)\n",timer);
327 moldyn->pwrite=timer;
328 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
329 moldyn->pfd=open(filename,
330 O_WRONLY|O_CREAT|O_EXCL,
333 perror("[moldyn] pressure log file\n");
336 dprintf(moldyn->pfd,"# pressure log file\n");
337 printf("pressure (%d)\n",timer);
339 case LOG_TEMPERATURE:
340 moldyn->twrite=timer;
341 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
342 moldyn->tfd=open(filename,
343 O_WRONLY|O_CREAT|O_EXCL,
346 perror("[moldyn] temperature log file\n");
349 dprintf(moldyn->tfd,"# temperature log file\n");
350 printf("temperature (%d)\n",timer);
353 moldyn->vwrite=timer;
354 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
355 moldyn->vfd=open(filename,
356 O_WRONLY|O_CREAT|O_EXCL,
359 perror("[moldyn] volume log file\n");
362 dprintf(moldyn->vfd,"# volume log file\n");
363 printf("volume (%d)\n",timer);
366 moldyn->swrite=timer;
367 printf("save file (%d)\n",timer);
370 moldyn->awrite=timer;
371 ret=visual_init(moldyn,moldyn->vlsdir);
373 printf("[moldyn] visual init failure\n");
376 printf("visual file (%d)\n",timer);
379 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
380 moldyn->rfd=open(filename,
381 O_WRONLY|O_CREAT|O_EXCL,
384 perror("[moldyn] report fd open");
387 printf("report -> ");
389 snprintf(filename,127,"%s/e_plot.scr",
391 moldyn->epfd=open(filename,
392 O_WRONLY|O_CREAT|O_EXCL,
395 perror("[moldyn] energy plot fd open");
398 dprintf(moldyn->epfd,e_plot_script);
403 snprintf(filename,127,"%s/pressure_plot.scr",
405 moldyn->ppfd=open(filename,
406 O_WRONLY|O_CREAT|O_EXCL,
409 perror("[moldyn] p plot fd open");
412 dprintf(moldyn->ppfd,pressure_plot_script);
417 snprintf(filename,127,"%s/temperature_plot.scr",
419 moldyn->tpfd=open(filename,
420 O_WRONLY|O_CREAT|O_EXCL,
423 perror("[moldyn] t plot fd open");
426 dprintf(moldyn->tpfd,temperature_plot_script);
428 printf("temperature ");
430 dprintf(moldyn->rfd,report_start,
431 moldyn->rauthor,moldyn->rtitle);
435 printf("unknown log type: %02x\n",type);
442 int moldyn_log_shutdown(t_moldyn *moldyn) {
446 printf("[moldyn] log shutdown\n");
450 dprintf(moldyn->rfd,report_energy);
451 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
456 if(moldyn->mfd) close(moldyn->mfd);
460 dprintf(moldyn->rfd,report_pressure);
461 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
468 dprintf(moldyn->rfd,report_temperature);
469 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
474 dprintf(moldyn->rfd,report_end);
476 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
479 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
482 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
491 * creating lattice functions
494 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
495 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin) {
507 /* how many atoms do we expect */
510 printf("[moldyn] WARNING: create 'none' lattice called");
512 if(type==CUBIC) new*=1;
513 if(type==FCC) new*=4;
514 if(type==DIAMOND) new*=8;
516 /* allocate space for atoms */
517 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
519 perror("[moldyn] realloc (create lattice)");
523 atom=&(moldyn->atom[count]);
525 /* no atoms on the boundaries (only reason: it looks better!) */
539 set_nn_dist(moldyn,lc);
540 ret=cubic_init(a,b,c,lc,atom,&orig);
541 strcpy(name,"cubic");
545 v3_scale(&orig,&orig,0.5);
546 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
547 ret=fcc_init(a,b,c,lc,atom,&orig);
552 v3_scale(&orig,&orig,0.25);
553 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
554 ret=diamond_init(a,b,c,lc,atom,&orig);
555 strcpy(name,"diamond");
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 %s lattice with %d atoms\n",name,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)");
608 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
610 perror("[moldyn] list realloc (add atom)");
613 moldyn->lc.subcell->list=ptr;
618 /* initialize new atom */
619 memset(&(atom[count]),0,sizeof(t_atom));
622 atom[count].element=element;
623 atom[count].mass=mass;
624 atom[count].brand=brand;
625 atom[count].tag=count;
626 atom[count].attr=attr;
627 check_per_bound(moldyn,&(atom[count].r));
628 atom[count].r_0=atom[count].r;
630 /* update total system mass */
631 total_mass_calc(moldyn);
636 int del_atom(t_moldyn *moldyn,int tag) {
643 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
645 perror("[moldyn]malloc (del atom)");
649 for(cnt=0;cnt<tag;cnt++)
652 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
654 new[cnt-1].tag=cnt-1;
666 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
685 v3_copy(&(atom[count].r),&r);
694 for(i=0;i<count;i++) {
695 atom[i].r.x-=(a*lc)/2.0;
696 atom[i].r.y-=(b*lc)/2.0;
697 atom[i].r.z-=(c*lc)/2.0;
703 /* fcc lattice init */
704 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
717 /* construct the basis */
718 memset(basis,0,3*sizeof(t_3dvec));
726 /* fill up the room */
734 v3_copy(&(atom[count].r),&r);
737 /* the three face centered atoms */
739 v3_add(&n,&r,&basis[l]);
740 v3_copy(&(atom[count].r),&n);
749 /* coordinate transformation */
750 for(i=0;i<count;i++) {
751 atom[i].r.x-=(a*lc)/2.0;
752 atom[i].r.y-=(b*lc)/2.0;
753 atom[i].r.z-=(c*lc)/2.0;
759 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
764 count=fcc_init(a,b,c,lc,atom,origin);
770 if(origin) v3_add(&o,&o,origin);
772 count+=fcc_init(a,b,c,lc,&atom[count],&o);
777 int destroy_atoms(t_moldyn *moldyn) {
779 if(moldyn->atom) free(moldyn->atom);
784 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
787 * - gaussian distribution of velocities
788 * - zero total momentum
789 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
794 t_3dvec p_total,delta;
799 random=&(moldyn->random);
801 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
803 /* gaussian distribution of velocities */
805 for(i=0;i<moldyn->count;i++) {
806 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
808 v=sigma*rand_get_gauss(random);
810 p_total.x+=atom[i].mass*v;
812 v=sigma*rand_get_gauss(random);
814 p_total.y+=atom[i].mass*v;
816 v=sigma*rand_get_gauss(random);
818 p_total.z+=atom[i].mass*v;
821 /* zero total momentum */
822 v3_scale(&p_total,&p_total,1.0/moldyn->count);
823 for(i=0;i<moldyn->count;i++) {
824 v3_scale(&delta,&p_total,1.0/atom[i].mass);
825 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
828 /* velocity scaling */
829 scale_velocity(moldyn,equi_init);
834 double total_mass_calc(t_moldyn *moldyn) {
840 for(i=0;i<moldyn->count;i++)
841 moldyn->mass+=moldyn->atom[i].mass;
846 double temperature_calc(t_moldyn *moldyn) {
848 /* assume up to date kinetic energy, which is 3/2 N k_B T */
850 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
855 double get_temperature(t_moldyn *moldyn) {
860 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
870 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
873 /* get kinetic energy / temperature & count involved atoms */
876 for(i=0;i<moldyn->count;i++) {
877 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
878 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
883 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
884 else return 0; /* no atoms involved in scaling! */
886 /* (temporary) hack for e,t = 0 */
889 if(moldyn->t_ref!=0.0) {
890 thermal_init(moldyn,equi_init);
894 return 0; /* no scaling needed */
898 /* get scaling factor */
899 scale=moldyn->t_ref/moldyn->t;
903 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
904 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
907 /* velocity scaling */
908 for(i=0;i<moldyn->count;i++) {
909 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
910 v3_scale(&(atom[i].v),&(atom[i].v),scale);
916 double ideal_gas_law_pressure(t_moldyn *moldyn) {
920 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
925 double virial_sum(t_moldyn *moldyn) {
930 /* virial (sum over atom virials) */
938 for(i=0;i<moldyn->count;i++) {
939 virial=&(moldyn->atom[i].virial);
940 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
941 moldyn->vir.xx+=virial->xx;
942 moldyn->vir.yy+=virial->yy;
943 moldyn->vir.zz+=virial->zz;
944 moldyn->vir.xy+=virial->xy;
945 moldyn->vir.xz+=virial->xz;
946 moldyn->vir.yz+=virial->yz;
949 /* global virial (absolute coordinates) */
950 virial=&(moldyn->gvir);
951 moldyn->gv=virial->xx+virial->yy+virial->zz;
953 return moldyn->virial;
956 double pressure_calc(t_moldyn *moldyn) {
960 * with W = 1/3 sum_i f_i r_i (- skipped!)
961 * virial = sum_i f_i r_i
963 * => P = (2 Ekin + virial) / (3V)
966 /* assume up to date virial & up to date kinetic energy */
968 /* pressure (atom virials) */
969 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
970 moldyn->p/=(3.0*moldyn->volume);
972 /* pressure (absolute coordinates) */
973 moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
974 moldyn->gp/=(3.0*moldyn->volume);
979 int average_reset(t_moldyn *moldyn) {
981 printf("[moldyn] average reset\n");
983 /* update skip value */
984 moldyn->avg_skip=moldyn->total_steps;
990 /* potential energy */
998 moldyn->virial_sum=0.0;
1009 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1013 if(moldyn->total_steps<moldyn->avg_skip)
1016 denom=moldyn->total_steps+1-moldyn->avg_skip;
1018 /* assume up to date energies, temperature, pressure etc */
1020 /* kinetic energy */
1021 moldyn->k_sum+=moldyn->ekin;
1022 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1023 moldyn->k_avg=moldyn->k_sum/denom;
1024 moldyn->k2_avg=moldyn->k2_sum/denom;
1025 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1027 /* potential energy */
1028 moldyn->v_sum+=moldyn->energy;
1029 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1030 moldyn->v_avg=moldyn->v_sum/denom;
1031 moldyn->v2_avg=moldyn->v2_sum/denom;
1032 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1035 moldyn->t_sum+=moldyn->t;
1036 moldyn->t_avg=moldyn->t_sum/denom;
1039 moldyn->virial_sum+=moldyn->virial;
1040 moldyn->virial_avg=moldyn->virial_sum/denom;
1041 moldyn->gv_sum+=moldyn->gv;
1042 moldyn->gv_avg=moldyn->gv_sum/denom;
1045 moldyn->p_sum+=moldyn->p;
1046 moldyn->p_avg=moldyn->p_sum/denom;
1047 moldyn->gp_sum+=moldyn->gp;
1048 moldyn->gp_avg=moldyn->gp_sum/denom;
1049 moldyn->tp_sum+=moldyn->tp;
1050 moldyn->tp_avg=moldyn->tp_sum/denom;
1055 int get_heat_capacity(t_moldyn *moldyn) {
1059 /* averages needed for heat capacity calc */
1060 if(moldyn->total_steps<moldyn->avg_skip)
1063 /* (temperature average)^2 */
1064 temp2=moldyn->t_avg*moldyn->t_avg;
1065 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1068 /* ideal gas contribution */
1069 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1070 printf(" ideal gas contribution: %f\n",
1071 ighc/moldyn->mass*KILOGRAM/JOULE);
1073 /* specific heat for nvt ensemble */
1074 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1075 moldyn->c_v_nvt/=moldyn->mass;
1077 /* specific heat for nve ensemble */
1078 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1079 moldyn->c_v_nve/=moldyn->mass;
1081 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1082 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1083 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)));
1088 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1104 /* store atomic configuration + dimension */
1105 store=malloc(moldyn->count*sizeof(t_atom));
1107 printf("[moldyn] allocating store mem failed\n");
1110 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1115 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1116 su=pow(2.0-h,ONE_THIRD)-1.0;
1117 dv=(1.0-h)*moldyn->volume;
1119 /* scale up dimension and atom positions */
1120 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1121 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1122 link_cell_shutdown(moldyn);
1123 link_cell_init(moldyn,QUIET);
1124 potential_force_calc(moldyn);
1127 /* restore atomic configuration + dim */
1128 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1131 /* scale down dimension and atom positions */
1132 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1133 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1134 link_cell_shutdown(moldyn);
1135 link_cell_init(moldyn,QUIET);
1136 potential_force_calc(moldyn);
1139 /* calculate pressure */
1140 moldyn->tp=-(y1-y0)/(2.0*dv);
1142 /* restore atomic configuration */
1143 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1145 link_cell_shutdown(moldyn);
1146 link_cell_init(moldyn,QUIET);
1147 //potential_force_calc(moldyn);
1149 /* free store buffer */
1156 double get_pressure(t_moldyn *moldyn) {
1162 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1174 if(x) dim->x*=scale;
1175 if(y) dim->y*=scale;
1176 if(z) dim->z*=scale;
1181 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1192 for(i=0;i<moldyn->count;i++) {
1193 r=&(moldyn->atom[i].r);
1202 int scale_volume(t_moldyn *moldyn) {
1208 vdim=&(moldyn->vis.dim);
1212 /* scaling factor */
1213 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1214 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1215 scale=pow(scale,ONE_THIRD);
1218 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1221 /* scale the atoms and dimensions */
1222 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1223 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1225 /* visualize dimensions */
1232 /* recalculate scaled volume */
1233 moldyn->volume=dim->x*dim->y*dim->z;
1235 /* adjust/reinit linkcell */
1236 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1237 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1238 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1239 link_cell_shutdown(moldyn);
1240 link_cell_init(moldyn,QUIET);
1251 double e_kin_calc(t_moldyn *moldyn) {
1259 for(i=0;i<moldyn->count;i++) {
1260 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1261 moldyn->ekin+=atom[i].ekin;
1264 return moldyn->ekin;
1267 double get_total_energy(t_moldyn *moldyn) {
1269 return(moldyn->ekin+moldyn->energy);
1272 t_3dvec get_total_p(t_moldyn *moldyn) {
1281 for(i=0;i<moldyn->count;i++) {
1282 v3_scale(&p,&(atom[i].v),atom[i].mass);
1283 v3_add(&p_total,&p_total,&p);
1289 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1293 /* nn_dist is the nearest neighbour distance */
1295 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1304 /* linked list / cell method */
1306 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1309 #ifndef LOWMEM_LISTS
1315 /* partitioning the md cell */
1316 lc->nx=moldyn->dim.x/moldyn->cutoff;
1317 lc->x=moldyn->dim.x/lc->nx;
1318 lc->ny=moldyn->dim.y/moldyn->cutoff;
1319 lc->y=moldyn->dim.y/lc->ny;
1320 lc->nz=moldyn->dim.z/moldyn->cutoff;
1321 lc->z=moldyn->dim.z/lc->nz;
1322 lc->cells=lc->nx*lc->ny*lc->nz;
1325 lc->subcell=malloc(lc->cells*sizeof(int*));
1327 lc->subcell=malloc(sizeof(t_lowmem_list));
1329 lc->subcell=malloc(lc->cells*sizeof(t_list));
1332 if(lc->subcell==NULL) {
1333 perror("[moldyn] cell init (malloc)");
1338 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1343 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1346 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1349 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1352 printf(" x: %d x %f A\n",lc->nx,lc->x);
1353 printf(" y: %d x %f A\n",lc->ny,lc->y);
1354 printf(" z: %d x %f A\n",lc->nz,lc->z);
1359 for(i=0;i<lc->cells;i++) {
1360 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1361 if(lc->subcell[i]==NULL) {
1362 perror("[moldyn] list init (malloc)");
1367 printf(" ---> %d malloc %p (%p)\n",
1368 i,lc->subcell[0],lc->subcell);
1372 lc->subcell->head=malloc(lc->cells*sizeof(int));
1373 if(lc->subcell->head==NULL) {
1374 perror("[moldyn] head init (malloc)");
1377 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1378 if(lc->subcell->list==NULL) {
1379 perror("[moldyn] list init (malloc)");
1383 for(i=0;i<lc->cells;i++)
1384 list_init_f(&(lc->subcell[i]));
1387 /* update the list */
1388 link_cell_update(moldyn);
1393 int link_cell_update(t_moldyn *moldyn) {
1411 for(i=0;i<lc->cells;i++)
1413 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1415 lc->subcell->head[i]=-1;
1417 list_destroy_f(&(lc->subcell[i]));
1420 for(count=0;count<moldyn->count;count++) {
1421 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1422 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1423 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1427 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1430 if(p>=MAX_ATOMS_PER_LIST) {
1431 printf("[moldyn] FATAL: amount of atoms too high!\n");
1435 lc->subcell[i+j*nx+k*nxy][p]=count;
1438 lc->subcell->list[count]=lc->subcell->head[p];
1439 lc->subcell->head[p]=count;
1441 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1445 printf(" ---> %d %d malloc %p (%p)\n",
1446 i,count,lc->subcell[i].current,lc->subcell);
1454 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1480 if(i>=nx||j>=ny||k>=nz)
1481 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1484 #ifndef LOWMEM_LISTS
1485 cell[0]=lc->subcell[i+j*nx+k*a];
1487 cell[0]=lc->subcell->head[i+j*nx+k*a];
1489 for(ci=-1;ci<=1;ci++) {
1492 if((x<0)||(x>=nx)) {
1496 for(cj=-1;cj<=1;cj++) {
1499 if((y<0)||(y>=ny)) {
1503 for(ck=-1;ck<=1;ck++) {
1506 if((z<0)||(z>=nz)) {
1510 if(!(ci|cj|ck)) continue;
1512 #ifndef LOWMEM_LISTS
1513 cell[--count2]=lc->subcell[x+y*nx+z*a];
1515 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1520 #ifndef LOWMEM_LISTS
1521 cell[count1++]=lc->subcell[x+y*nx+z*a];
1523 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1535 int link_cell_shutdown(t_moldyn *moldyn) {
1537 #ifndef LOWMEM_LISTS
1545 free(lc->subcell->head);
1546 free(lc->subcell->list);
1550 for(i=0;i<lc->cells;i++) {
1552 free(lc->subcell[i]);
1554 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1555 list_destroy_f(&(lc->subcell[i]));
1565 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1569 t_moldyn_schedule *schedule;
1571 schedule=&(moldyn->schedule);
1572 count=++(schedule->total_sched);
1574 ptr=realloc(schedule->runs,count*sizeof(int));
1576 perror("[moldyn] realloc (runs)");
1580 schedule->runs[count-1]=runs;
1582 ptr=realloc(schedule->tau,count*sizeof(double));
1584 perror("[moldyn] realloc (tau)");
1588 schedule->tau[count-1]=tau;
1590 printf("[moldyn] schedule added:\n");
1591 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1597 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1599 moldyn->schedule.hook=hook;
1600 moldyn->schedule.hook_params=hook_params;
1607 * 'integration of newtons equation' - algorithms
1611 /* start the integration */
1613 int moldyn_integrate(t_moldyn *moldyn) {
1616 unsigned int e,m,s,v,p,t,a;
1618 t_moldyn_schedule *sched;
1623 double energy_scale;
1624 struct timeval t1,t2;
1629 pthread_t io_thread;
1638 sched=&(moldyn->schedule);
1641 /* initialize linked cell method */
1642 link_cell_init(moldyn,VERBOSE);
1644 /* logging & visualization */
1653 /* sqaure of some variables */
1654 moldyn->tau_square=moldyn->tau*moldyn->tau;
1656 /* get current time */
1657 gettimeofday(&t1,NULL);
1659 /* calculate initial forces */
1660 potential_force_calc(moldyn);
1665 /* some stupid checks before we actually start calculating bullshit */
1666 if(moldyn->cutoff>0.5*moldyn->dim.x)
1667 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
1668 if(moldyn->cutoff>0.5*moldyn->dim.y)
1669 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
1670 if(moldyn->cutoff>0.5*moldyn->dim.z)
1671 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
1673 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
1674 if(ds>0.05*moldyn->nnd)
1675 printf("[moldyn] WARNING: forces too high / tau too small!\n");
1678 /* zero absolute time */
1679 // should have right values!
1681 //moldyn->total_steps=0;
1683 /* debugging, ignore */
1686 /* zero & init moldyn copy */
1688 memset(&md_copy,0,sizeof(t_moldyn));
1689 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1690 if(atom_copy==NULL) {
1691 perror("[moldyn] malloc atom copy (init)");
1696 /* tell the world */
1697 printf("[moldyn] integration start, go get a coffee ...\n");
1699 /* executing the schedule */
1701 while(sched->count<sched->total_sched) {
1703 /* setting amount of runs and finite time step size */
1704 moldyn->tau=sched->tau[sched->count];
1705 moldyn->tau_square=moldyn->tau*moldyn->tau;
1706 moldyn->time_steps=sched->runs[sched->count];
1708 /* energy scaling factor (might change!) */
1709 energy_scale=moldyn->count*EV;
1711 /* integration according to schedule */
1713 for(i=0;i<moldyn->time_steps;i++) {
1715 /* integration step */
1716 moldyn->integrate(moldyn);
1718 /* calculate kinetic energy, temperature and pressure */
1720 temperature_calc(moldyn);
1722 pressure_calc(moldyn);
1724 thermodynamic_pressure_calc(moldyn);
1725 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
1729 /* calculate fluctuations + averages */
1730 average_and_fluctuation_calc(moldyn);
1733 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
1734 scale_velocity(moldyn,FALSE);
1735 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
1736 scale_volume(moldyn);
1738 /* check for log & visualization */
1740 if(!(moldyn->total_steps%e))
1741 dprintf(moldyn->efd,
1743 moldyn->time,moldyn->ekin/energy_scale,
1744 moldyn->energy/energy_scale,
1745 get_total_energy(moldyn)/energy_scale);
1748 if(!(moldyn->total_steps%m)) {
1749 momentum=get_total_p(moldyn);
1750 dprintf(moldyn->mfd,
1751 "%f %f %f %f %f\n",moldyn->time,
1752 momentum.x,momentum.y,momentum.z,
1753 v3_norm(&momentum));
1757 if(!(moldyn->total_steps%p)) {
1758 dprintf(moldyn->pfd,
1759 "%f %f %f %f %f %f %f\n",moldyn->time,
1760 moldyn->p/BAR,moldyn->p_avg/BAR,
1761 moldyn->gp/BAR,moldyn->gp_avg/BAR,
1762 moldyn->tp/BAR,moldyn->tp_avg/BAR);
1766 if(!(moldyn->total_steps%t)) {
1767 dprintf(moldyn->tfd,
1769 moldyn->time,moldyn->t,moldyn->t_avg);
1773 if(!(moldyn->total_steps%v)) {
1774 dprintf(moldyn->vfd,
1775 "%f %f\n",moldyn->time,moldyn->volume);
1779 if(!(moldyn->total_steps%s)) {
1780 snprintf(dir,128,"%s/s-%08.f.save",
1781 moldyn->vlsdir,moldyn->time);
1782 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
1784 if(fd<0) perror("[moldyn] save fd open");
1786 write(fd,moldyn,sizeof(t_moldyn));
1787 write(fd,moldyn->atom,
1788 moldyn->count*sizeof(t_atom));
1794 if(!(moldyn->total_steps%a)) {
1796 /* check whether thread has not terminated yet */
1798 ret=pthread_join(io_thread,NULL);
1801 /* prepare and start thread */
1802 if(moldyn->count!=md_copy.count) {
1806 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
1808 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1809 if(atom_copy==NULL) {
1810 perror("[moldyn] malloc atom copy (change)");
1814 md_copy.atom=atom_copy;
1815 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
1817 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
1819 perror("[moldyn] create visual atoms thread\n");
1823 visual_atoms(moldyn);
1828 /* display progress */
1830 /* get current time */
1831 gettimeofday(&t2,NULL);
1833 printf("\rsched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)",
1834 sched->count,i,moldyn->total_steps,
1835 moldyn->t,moldyn->t_avg,
1836 moldyn->p/BAR,moldyn->p_avg/BAR,
1837 //moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
1839 (int)(t2.tv_sec-t1.tv_sec));
1843 /* copy over time */
1847 /* increase absolute time */
1848 moldyn->time+=moldyn->tau;
1849 moldyn->total_steps+=1;
1853 /* check for hooks */
1855 printf("\n ## schedule hook %d start ##\n",
1857 sched->hook(moldyn,sched->hook_params);
1858 printf(" ## schedule hook end ##\n");
1861 /* increase the schedule counter */
1869 /* velocity verlet */
1871 int velocity_verlet(t_moldyn *moldyn) {
1874 double tau,tau_square,h;
1879 count=moldyn->count;
1881 tau_square=moldyn->tau_square;
1883 for(i=0;i<count;i++) {
1884 /* check whether fixed atom */
1885 if(atom[i].attr&ATOM_ATTR_FP)
1889 v3_scale(&delta,&(atom[i].v),tau);
1890 v3_add(&(atom[i].r),&(atom[i].r),&delta);
1891 v3_scale(&delta,&(atom[i].f),h*tau_square);
1892 v3_add(&(atom[i].r),&(atom[i].r),&delta);
1893 check_per_bound(moldyn,&(atom[i].r));
1895 /* velocities [actually v(t+tau/2)] */
1896 v3_scale(&delta,&(atom[i].f),h*tau);
1897 v3_add(&(atom[i].v),&(atom[i].v),&delta);
1900 /* criticial check */
1901 moldyn_bc_check(moldyn);
1903 /* neighbour list update */
1904 link_cell_update(moldyn);
1906 /* forces depending on chosen potential */
1908 potential_force_calc(moldyn);
1910 albe_potential_force_calc(moldyn);
1913 for(i=0;i<count;i++) {
1914 /* check whether fixed atom */
1915 if(atom[i].attr&ATOM_ATTR_FP)
1917 /* again velocities [actually v(t+tau)] */
1918 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
1919 v3_add(&(atom[i].v),&(atom[i].v),&delta);
1928 * potentials & corresponding forces & virial routine
1932 /* generic potential and force calculation */
1934 int potential_force_calc(t_moldyn *moldyn) {
1937 t_atom *itom,*jtom,*ktom;
1941 int *neighbour_i[27];
1945 int neighbour_i[27];
1948 t_list neighbour_i[27];
1949 t_list neighbour_i2[27];
1955 count=moldyn->count;
1965 /* reset global virial */
1966 memset(&(moldyn->gvir),0,sizeof(t_virial));
1968 /* reset force, site energy and virial of every atom */
1970 i=omp_get_thread_num();
1971 #pragma omp parallel for private(virial)
1973 for(i=0;i<count;i++) {
1976 v3_zero(&(itom[i].f));
1979 virial=(&(itom[i].virial));
1987 /* reset site energy */
1992 /* get energy, force and virial of every atom */
1994 /* first (and only) loop over atoms i */
1995 for(i=0;i<count;i++) {
1997 /* single particle potential/force */
1998 if(itom[i].attr&ATOM_ATTR_1BP)
2000 moldyn->func1b(moldyn,&(itom[i]));
2002 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2005 /* 2 body pair potential/force */
2007 link_cell_neighbour_index(moldyn,
2008 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2009 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2010 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2015 /* first loop over atoms j */
2016 if(moldyn->func2b) {
2023 while(neighbour_i[j][p]!=-1) {
2025 jtom=&(atom[neighbour_i[j][p]]);
2033 p=lc->subcell->list[p];
2035 this=&(neighbour_i[j]);
2038 if(this->start==NULL)
2042 jtom=this->current->data;
2045 if(jtom==&(itom[i]))
2048 if((jtom->attr&ATOM_ATTR_2BP)&
2049 (itom[i].attr&ATOM_ATTR_2BP)) {
2050 moldyn->func2b(moldyn,
2060 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2066 /* 3 body potential/force */
2068 if(!(itom[i].attr&ATOM_ATTR_3BP))
2071 /* copy the neighbour lists */
2073 /* no copy needed for static lists */
2075 /* no copy needed for lowmem lists */
2077 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2080 /* second loop over atoms j */
2087 while(neighbour_i[j][p]!=-1) {
2089 jtom=&(atom[neighbour_i[j][p]]);
2097 p=lc->subcell->list[p];
2099 this=&(neighbour_i[j]);
2102 if(this->start==NULL)
2107 jtom=this->current->data;
2110 if(jtom==&(itom[i]))
2113 if(!(jtom->attr&ATOM_ATTR_3BP))
2119 if(moldyn->func3b_j1)
2120 moldyn->func3b_j1(moldyn,
2125 /* in first j loop, 3bp run can be skipped */
2126 if(!(moldyn->run3bp))
2129 /* first loop over atoms k */
2130 if(moldyn->func3b_k1) {
2138 while(neighbour_i[k][q]!=-1) {
2140 ktom=&(atom[neighbour_i[k][q]]);
2148 q=lc->subcell->list[q];
2150 that=&(neighbour_i2[k]);
2153 if(that->start==NULL)
2157 ktom=that->current->data;
2160 if(!(ktom->attr&ATOM_ATTR_3BP))
2166 if(ktom==&(itom[i]))
2169 moldyn->func3b_k1(moldyn,
2179 } while(list_next_f(that)!=\
2187 if(moldyn->func3b_j2)
2188 moldyn->func3b_j2(moldyn,
2193 /* second loop over atoms k */
2194 if(moldyn->func3b_k2) {
2202 while(neighbour_i[k][q]!=-1) {
2204 ktom=&(atom[neighbour_i[k][q]]);
2212 q=lc->subcell->list[q];
2214 that=&(neighbour_i2[k]);
2217 if(that->start==NULL)
2221 ktom=that->current->data;
2224 if(!(ktom->attr&ATOM_ATTR_3BP))
2230 if(ktom==&(itom[i]))
2233 moldyn->func3b_k2(moldyn,
2244 } while(list_next_f(that)!=\
2252 /* 2bp post function */
2253 if(moldyn->func3b_j3) {
2254 moldyn->func3b_j3(moldyn,
2263 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2278 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2279 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2281 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2282 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2283 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2287 /* some postprocessing */
2289 #pragma omp parallel for
2291 for(i=0;i<count;i++) {
2292 /* calculate global virial */
2293 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2294 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2295 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2296 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2297 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2298 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2300 /* check forces regarding the given timestep */
2301 if(v3_norm(&(itom[i].f))>\
2302 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2303 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2311 * virial calculation
2314 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2315 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2317 a->virial.xx+=f->x*d->x;
2318 a->virial.yy+=f->y*d->y;
2319 a->virial.zz+=f->z*d->z;
2320 a->virial.xy+=f->x*d->y;
2321 a->virial.xz+=f->x*d->z;
2322 a->virial.yz+=f->y*d->z;
2328 * periodic boundary checking
2331 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2332 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2343 if(moldyn->status&MOLDYN_STAT_PBX) {
2344 if(a->x>=x) a->x-=dim->x;
2345 else if(-a->x>x) a->x+=dim->x;
2347 if(moldyn->status&MOLDYN_STAT_PBY) {
2348 if(a->y>=y) a->y-=dim->y;
2349 else if(-a->y>y) a->y+=dim->y;
2351 if(moldyn->status&MOLDYN_STAT_PBZ) {
2352 if(a->z>=z) a->z-=dim->z;
2353 else if(-a->z>z) a->z+=dim->z;
2360 * debugging / critical check functions
2363 int moldyn_bc_check(t_moldyn *moldyn) {
2376 for(i=0;i<moldyn->count;i++) {
2377 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2378 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2379 i,atom[i].r.x,dim->x/2);
2380 printf("diagnostic:\n");
2381 printf("-----------\natom.r.x:\n");
2383 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2386 ((byte)&(1<<k))?1:0,
2389 printf("---------------\nx=dim.x/2:\n");
2391 memcpy(&byte,(u8 *)(&x)+j,1);
2394 ((byte)&(1<<k))?1:0,
2397 if(atom[i].r.x==x) printf("the same!\n");
2398 else printf("different!\n");
2400 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2401 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2402 i,atom[i].r.y,dim->y/2);
2403 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2404 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2405 i,atom[i].r.z,dim->z/2);
2415 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2422 fd=open(file,O_RDONLY);
2424 perror("[moldyn] load save file open");
2428 fsize=lseek(fd,0,SEEK_END);
2429 lseek(fd,0,SEEK_SET);
2431 size=sizeof(t_moldyn);
2434 cnt=read(fd,moldyn,size);
2436 perror("[moldyn] load save file read (moldyn)");
2442 size=moldyn->count*sizeof(t_atom);
2444 /* correcting possible atom data offset */
2446 if(fsize!=sizeof(t_moldyn)+size) {
2447 corr=fsize-sizeof(t_moldyn)-size;
2448 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2449 printf(" moifying offset:\n");
2450 printf(" - current pos: %d\n",sizeof(t_moldyn));
2451 printf(" - atom size: %d\n",size);
2452 printf(" - file size: %d\n",fsize);
2453 printf(" => correction: %d\n",corr);
2454 lseek(fd,corr,SEEK_CUR);
2457 moldyn->atom=(t_atom *)malloc(size);
2458 if(moldyn->atom==NULL) {
2459 perror("[moldyn] load save file malloc (atoms)");
2464 cnt=read(fd,moldyn->atom,size);
2466 perror("[moldyn] load save file read (atoms)");
2477 int moldyn_free_save_file(t_moldyn *moldyn) {
2484 int moldyn_load(t_moldyn *moldyn) {
2492 * function to find/callback all combinations of 2 body bonds
2495 int process_2b_bonds(t_moldyn *moldyn,void *data,
2496 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2497 void *data,u8 bc)) {
2507 t_list neighbour[27];
2517 for(i=0;i<moldyn->count;i++) {
2518 /* neighbour indexing */
2519 link_cell_neighbour_index(moldyn,
2520 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2521 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2522 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
2527 bc=(j<lc->dnlc)?0:1;
2532 while(neighbour[j][p]!=-1) {
2534 jtom=&(moldyn->atom[neighbour[j][p]]);
2542 p=lc->subcell->list[p];
2544 this=&(neighbour[j]);
2547 if(this->start==NULL)
2552 jtom=this->current->data;
2556 process(moldyn,&(itom[i]),jtom,data,bc);
2563 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2573 * function to find neighboured atoms
2576 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
2577 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
2578 void *data,u8 bc)) {
2588 t_list neighbour[27];
2597 /* neighbour indexing */
2598 link_cell_neighbour_index(moldyn,
2599 (atom->r.x+moldyn->dim.x/2)/lc->x,
2600 (atom->r.y+moldyn->dim.y/2)/lc->x,
2601 (atom->r.z+moldyn->dim.z/2)/lc->x,
2606 bc=(j<lc->dnlc)?0:1;
2611 while(neighbour[j][p]!=-1) {
2613 natom=&(moldyn->atom[neighbour[j][p]]);
2620 natom=&(moldyn->atom[p]);
2621 p=lc->subcell->list[p];
2623 this=&(neighbour[j]);
2626 if(this->start==NULL)
2631 natom=this->current->data;
2635 process(moldyn,atom,natom,data,bc);
2642 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2651 * post processing functions
2654 int get_line(int fd,char *line,int max) {
2661 if(count==max) return count;
2662 ret=read(fd,line+count,1);
2663 if(ret<=0) return ret;
2664 if(line[count]=='\n') {
2665 memset(line+count,0,max-count-1);
2673 int pair_correlation_init(t_moldyn *moldyn,double dr) {
2679 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
2695 for(i=0;i<moldyn->count;i++) {
2697 v3_sub(&dist,&(atom[i].r),&(atom[i].r_0));
2698 check_per_bound(moldyn,&dist);
2699 d2=v3_absolute_square(&dist);
2713 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
2714 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
2715 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
2720 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
2725 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
2726 t_atom *jtom,void *data,u8 bc) {
2733 /* only count pairs once,
2734 * skip same atoms */
2735 if(itom->tag>=jtom->tag)
2739 * pair correlation calc
2746 v3_sub(&dist,&(jtom->r),&(itom->r));
2747 if(bc) check_per_bound(moldyn,&dist);
2748 d=v3_absolute_square(&dist);
2750 /* ignore if greater cutoff */
2751 if(d>moldyn->cutoff_square)
2754 /* fill the slots */
2758 /* should never happen but it does 8) -
2759 * related to -ffloat-store problem! */
2761 printf("[moldyn] WARNING: pcc (%d/%d)",
2767 if(itom->brand!=jtom->brand) {
2772 /* type a - type a bonds */
2774 pcc->stat[s+pcc->o1]+=1;
2776 /* type b - type b bonds */
2777 pcc->stat[s+pcc->o2]+=1;
2783 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
2790 pcc.o1=moldyn->cutoff/dr;
2793 if(pcc.o1*dr<=moldyn->cutoff)
2794 printf("[moldyn] WARNING: pcc (low #slots)\n");
2796 printf("[moldyn] pair correlation calc info:\n");
2797 printf(" time: %f\n",moldyn->time);
2798 printf(" count: %d\n",moldyn->count);
2799 printf(" cutoff: %f\n",moldyn->cutoff);
2800 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
2803 pcc.stat=(double *)ptr;
2806 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
2807 if(pcc.stat==NULL) {
2808 perror("[moldyn] pair correlation malloc");
2813 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
2816 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
2819 for(i=1;i<pcc.o1;i++) {
2820 // normalization: 4 pi r^2 dr
2821 // here: not double counting pairs -> 2 pi r r dr
2822 // ... and actually it's a constant times r^2
2825 pcc.stat[pcc.o1+i]/=norm;
2826 pcc.stat[pcc.o2+i]/=norm;
2831 /* todo: store/print pair correlation function */
2838 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2845 if(itom->tag>=jtom->tag)
2849 v3_sub(&dist,&(jtom->r),&(itom->r));
2850 if(bc) check_per_bound(moldyn,&dist);
2851 d=v3_absolute_square(&dist);
2853 /* ignore if greater or equal cutoff */
2854 if(d>moldyn->cutoff_square)
2857 /* check for potential bond */
2858 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
2861 /* now count this bonding ... */
2864 /* increase total bond counter
2865 * ... double counting!
2870 ba->acnt[jtom->tag]+=1;
2872 ba->bcnt[jtom->tag]+=1;
2875 ba->acnt[itom->tag]+=1;
2877 ba->bcnt[itom->tag]+=1;
2882 int bond_analyze(t_moldyn *moldyn,double *quality) {
2884 // by now: # bonds of type 'a-4b' and 'b-4a' / # bonds total
2892 ba.acnt=malloc(moldyn->count*sizeof(int));
2894 perror("[moldyn] bond analyze malloc (a)");
2897 memset(ba.acnt,0,moldyn->count*sizeof(int));
2899 ba.bcnt=malloc(moldyn->count*sizeof(int));
2901 perror("[moldyn] bond analyze malloc (b)");
2904 memset(ba.bcnt,0,moldyn->count*sizeof(int));
2913 process_2b_bonds(moldyn,&ba,bond_analyze_process);
2915 for(i=0;i<moldyn->count;i++) {
2916 if(atom[i].brand==0) {
2917 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
2921 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
2929 printf("[moldyn] bond analyze: c_cnt=%d | set=%d\n",ccnt,cset);
2930 printf("[moldyn] bond analyze: q_cnt=%d | tot=%d\n",qcnt,ba.tcnt);
2933 quality[0]=1.0*ccnt/cset;
2934 quality[1]=1.0*qcnt/ba.tcnt;
2937 printf("[moldyn] bond analyze: c_bnd_q=%f\n",1.0*qcnt/ba.tcnt);
2938 printf("[moldyn] bond analyze: tot_q=%f\n",1.0*qcnt/ba.tcnt);
2945 * visualization code
2948 int visual_init(t_moldyn *moldyn,char *filebase) {
2950 strncpy(moldyn->vis.fb,filebase,128);
2955 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2962 if(itom->tag>=jtom->tag)
2965 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
2968 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
2969 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
2970 itom->r.x,itom->r.y,itom->r.z,
2971 jtom->r.x,jtom->r.y,jtom->r.z);
2977 void *visual_atoms(void *ptr) {
2979 int visual_atoms(t_moldyn *moldyn) {
3003 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3004 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3006 perror("open visual save file fd");
3012 /* write the actual data file */
3015 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3016 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3018 // atomic configuration
3019 for(i=0;i<moldyn->count;i++)
3020 // atom type, positions, color and kinetic energy
3021 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3025 pse_col[atom[i].element],
3028 // bonds between atoms
3029 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3033 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3034 -dim.x/2,-dim.y/2,-dim.z/2,
3035 dim.x/2,-dim.y/2,-dim.z/2);
3036 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3037 -dim.x/2,-dim.y/2,-dim.z/2,
3038 -dim.x/2,dim.y/2,-dim.z/2);
3039 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3040 dim.x/2,dim.y/2,-dim.z/2,
3041 dim.x/2,-dim.y/2,-dim.z/2);
3042 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3043 -dim.x/2,dim.y/2,-dim.z/2,
3044 dim.x/2,dim.y/2,-dim.z/2);
3046 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3047 -dim.x/2,-dim.y/2,dim.z/2,
3048 dim.x/2,-dim.y/2,dim.z/2);
3049 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3050 -dim.x/2,-dim.y/2,dim.z/2,
3051 -dim.x/2,dim.y/2,dim.z/2);
3052 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3053 dim.x/2,dim.y/2,dim.z/2,
3054 dim.x/2,-dim.y/2,dim.z/2);
3055 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3056 -dim.x/2,dim.y/2,dim.z/2,
3057 dim.x/2,dim.y/2,dim.z/2);
3059 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3060 -dim.x/2,-dim.y/2,dim.z/2,
3061 -dim.x/2,-dim.y/2,-dim.z/2);
3062 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3063 -dim.x/2,dim.y/2,dim.z/2,
3064 -dim.x/2,dim.y/2,-dim.z/2);
3065 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3066 dim.x/2,-dim.y/2,dim.z/2,
3067 dim.x/2,-dim.y/2,-dim.z/2);
3068 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3069 dim.x/2,dim.y/2,dim.z/2,
3070 dim.x/2,dim.y/2,-dim.z/2);
3086 * fpu cntrol functions
3089 // set rounding to double (eliminates -ffloat-store!)
3090 int fpu_set_rtd(void) {
3096 ctrl&=~_FPU_EXTENDED;