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>
26 #if defined PTHREADS || defined VISUAL_THREAD
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"
52 pthread_mutex_t *amutex;
53 pthread_mutex_t emutex;
57 * the moldyn functions
60 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
62 printf("[moldyn] init\n");
64 /* only needed if compiled without -msse2 (float-store prob!) */
67 memset(moldyn,0,sizeof(t_moldyn));
72 rand_init(&(moldyn->random),NULL,1);
73 moldyn->random.status|=RAND_STAT_VERBOSE;
76 pthread_mutex_init(&emutex,NULL);
82 int moldyn_shutdown(t_moldyn *moldyn) {
88 printf("[moldyn] shutdown\n");
91 for(i=0;i<moldyn->count;i++)
92 pthread_mutex_destroy(&(amutex[i]));
94 pthread_mutex_destroy(&emutex);
97 moldyn_log_shutdown(moldyn);
98 link_cell_shutdown(moldyn);
99 rand_close(&(moldyn->random));
105 int set_int_alg(t_moldyn *moldyn,u8 algo) {
107 printf("[moldyn] integration algorithm: ");
110 case MOLDYN_INTEGRATE_VERLET:
111 moldyn->integrate=velocity_verlet;
112 printf("velocity verlet\n");
115 printf("unknown integration algorithm: %02x\n",algo);
123 int set_cutoff(t_moldyn *moldyn,double cutoff) {
125 moldyn->cutoff=cutoff;
126 moldyn->cutoff_square=cutoff*cutoff;
128 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
133 int set_temperature(t_moldyn *moldyn,double t_ref) {
137 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
142 int set_pressure(t_moldyn *moldyn,double p_ref) {
146 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
151 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
153 moldyn->pt_scale&=(~(P_SCALE_MASK));
154 moldyn->pt_scale|=ptype;
157 printf("[moldyn] p scaling:\n");
159 printf(" p: %s",ptype?"yes":"no ");
161 printf(" | type: %02x | factor: %f",ptype,ptc);
167 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
169 moldyn->pt_scale&=(~(T_SCALE_MASK));
170 moldyn->pt_scale|=ttype;
173 printf("[moldyn] t scaling:\n");
175 printf(" t: %s",ttype?"yes":"no ");
177 printf(" | type: %02x | factor: %f",ttype,ttc);
183 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
185 moldyn->pt_scale=(ptype|ttype);
189 printf("[moldyn] p/t scaling:\n");
191 printf(" p: %s",ptype?"yes":"no ");
193 printf(" | type: %02x | factor: %f",ptype,ptc);
196 printf(" t: %s",ttype?"yes":"no ");
198 printf(" | type: %02x | factor: %f",ttype,ttc);
204 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
210 moldyn->volume=x*y*z;
218 printf("[moldyn] dimensions in A and A^3 respectively:\n");
219 printf(" x: %f\n",moldyn->dim.x);
220 printf(" y: %f\n",moldyn->dim.y);
221 printf(" z: %f\n",moldyn->dim.z);
222 printf(" volume: %f\n",moldyn->volume);
223 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
228 int set_nn_dist(t_moldyn *moldyn,double dist) {
235 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
237 printf("[moldyn] periodic boundary conditions:\n");
240 moldyn->status|=MOLDYN_STAT_PBX;
243 moldyn->status|=MOLDYN_STAT_PBY;
246 moldyn->status|=MOLDYN_STAT_PBZ;
248 printf(" x: %s\n",x?"yes":"no");
249 printf(" y: %s\n",y?"yes":"no");
250 printf(" z: %s\n",z?"yes":"no");
255 int set_potential(t_moldyn *moldyn,u8 type) {
258 case MOLDYN_POTENTIAL_TM:
259 moldyn->func1b=tersoff_mult_1bp;
260 moldyn->func3b_j1=tersoff_mult_3bp_j1;
261 moldyn->func3b_k1=tersoff_mult_3bp_k1;
262 moldyn->func3b_j2=tersoff_mult_3bp_j2;
263 moldyn->func3b_k2=tersoff_mult_3bp_k2;
264 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
266 case MOLDYN_POTENTIAL_AM:
267 moldyn->func3b_j1=albe_mult_3bp_j1;
268 moldyn->func3b_k1=albe_mult_3bp_k1;
269 moldyn->func3b_j2=albe_mult_3bp_j2;
270 moldyn->func3b_k2=albe_mult_3bp_k2;
271 moldyn->check_2b_bond=albe_mult_check_2b_bond;
273 case MOLDYN_POTENTIAL_HO:
274 moldyn->func2b=harmonic_oscillator;
275 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
277 case MOLDYN_POTENTIAL_LJ:
278 moldyn->func2b=lennard_jones;
279 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
282 printf("[moldyn] set potential: unknown type %02x\n",
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->vwrite=timer;
375 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
376 moldyn->vfd=open(filename,
377 O_WRONLY|O_CREAT|O_EXCL,
380 perror("[moldyn] volume log file\n");
383 dprintf(moldyn->vfd,"# volume log file\n");
384 printf("volume (%d)\n",timer);
387 moldyn->swrite=timer;
388 printf("save file (%d)\n",timer);
391 moldyn->awrite=timer;
392 ret=visual_init(moldyn,moldyn->vlsdir);
394 printf("[moldyn] visual init failure\n");
397 printf("visual file (%d)\n",timer);
400 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
401 moldyn->rfd=open(filename,
402 O_WRONLY|O_CREAT|O_EXCL,
405 perror("[moldyn] report fd open");
408 printf("report -> ");
410 snprintf(filename,127,"%s/e_plot.scr",
412 moldyn->epfd=open(filename,
413 O_WRONLY|O_CREAT|O_EXCL,
416 perror("[moldyn] energy plot fd open");
419 dprintf(moldyn->epfd,e_plot_script);
424 snprintf(filename,127,"%s/pressure_plot.scr",
426 moldyn->ppfd=open(filename,
427 O_WRONLY|O_CREAT|O_EXCL,
430 perror("[moldyn] p plot fd open");
433 dprintf(moldyn->ppfd,pressure_plot_script);
438 snprintf(filename,127,"%s/temperature_plot.scr",
440 moldyn->tpfd=open(filename,
441 O_WRONLY|O_CREAT|O_EXCL,
444 perror("[moldyn] t plot fd open");
447 dprintf(moldyn->tpfd,temperature_plot_script);
449 printf("temperature ");
451 dprintf(moldyn->rfd,report_start,
452 moldyn->rauthor,moldyn->rtitle);
456 printf("unknown log type: %02x\n",type);
463 int moldyn_log_shutdown(t_moldyn *moldyn) {
467 printf("[moldyn] log shutdown\n");
471 dprintf(moldyn->rfd,report_energy);
472 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
477 if(moldyn->mfd) close(moldyn->mfd);
481 dprintf(moldyn->rfd,report_pressure);
482 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
489 dprintf(moldyn->rfd,report_temperature);
490 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
495 dprintf(moldyn->rfd,report_end);
497 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
500 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
503 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
512 * creating lattice functions
515 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
516 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin) {
525 pthread_mutex_t *mutex;
531 /* how many atoms do we expect */
534 printf("[moldyn] WARNING: create 'none' lattice called");
536 if(type==CUBIC) new*=1;
537 if(type==FCC) new*=4;
538 if(type==DIAMOND) new*=8;
540 /* allocate space for atoms */
541 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
543 perror("[moldyn] realloc (create lattice)");
547 atom=&(moldyn->atom[count]);
550 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
552 perror("[moldyn] mutex realloc (add atom)");
556 mutex=&(amutex[count]);
559 /* no atoms on the boundaries (only reason: it looks better!) */
573 set_nn_dist(moldyn,lc);
574 ret=cubic_init(a,b,c,lc,atom,&orig);
575 strcpy(name,"cubic");
579 v3_scale(&orig,&orig,0.5);
580 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
581 ret=fcc_init(a,b,c,lc,atom,&orig);
586 v3_scale(&orig,&orig,0.25);
587 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
588 ret=diamond_init(a,b,c,lc,atom,&orig);
589 strcpy(name,"diamond");
592 printf("unknown lattice type (%02x)\n",type);
598 printf("[moldyn] creating lattice failed\n");
599 printf(" amount of atoms\n");
600 printf(" - expected: %d\n",new);
601 printf(" - created: %d\n",ret);
606 printf("[moldyn] created %s lattice with %d atoms\n",name,new);
608 for(ret=0;ret<new;ret++) {
609 atom[ret].element=element;
612 atom[ret].brand=brand;
613 atom[ret].tag=count+ret;
614 check_per_bound(moldyn,&(atom[ret].r));
615 atom[ret].r_0=atom[ret].r;
617 pthread_mutex_init(&(mutex[ret]),NULL);
621 /* update total system mass */
622 total_mass_calc(moldyn);
627 int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
628 t_3dvec *r,t_3dvec *v) {
635 count=(moldyn->count)++; // asshole style!
637 ptr=realloc(atom,(count+1)*sizeof(t_atom));
639 perror("[moldyn] realloc (add atom)");
645 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
647 perror("[moldyn] list realloc (add atom)");
650 moldyn->lc.subcell->list=ptr;
654 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
656 perror("[moldyn] mutex realloc (add atom)");
660 pthread_mutex_init(&(amutex[count]),NULL);
665 /* initialize new atom */
666 memset(&(atom[count]),0,sizeof(t_atom));
669 atom[count].element=element;
670 atom[count].mass=mass;
671 atom[count].brand=brand;
672 atom[count].tag=count;
673 atom[count].attr=attr;
674 check_per_bound(moldyn,&(atom[count].r));
675 atom[count].r_0=atom[count].r;
677 /* update total system mass */
678 total_mass_calc(moldyn);
683 int del_atom(t_moldyn *moldyn,int tag) {
690 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
692 perror("[moldyn]malloc (del atom)");
696 for(cnt=0;cnt<tag;cnt++)
699 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
701 new[cnt-1].tag=cnt-1;
713 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
732 v3_copy(&(atom[count].r),&r);
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 /* fcc lattice init */
751 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
764 /* construct the basis */
765 memset(basis,0,3*sizeof(t_3dvec));
773 /* fill up the room */
781 v3_copy(&(atom[count].r),&r);
784 /* the three face centered atoms */
786 v3_add(&n,&r,&basis[l]);
787 v3_copy(&(atom[count].r),&n);
796 /* coordinate transformation */
797 for(i=0;i<count;i++) {
798 atom[i].r.x-=(a*lc)/2.0;
799 atom[i].r.y-=(b*lc)/2.0;
800 atom[i].r.z-=(c*lc)/2.0;
806 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin) {
811 count=fcc_init(a,b,c,lc,atom,origin);
817 if(origin) v3_add(&o,&o,origin);
819 count+=fcc_init(a,b,c,lc,&atom[count],&o);
824 int destroy_atoms(t_moldyn *moldyn) {
826 if(moldyn->atom) free(moldyn->atom);
831 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
834 * - gaussian distribution of velocities
835 * - zero total momentum
836 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
841 t_3dvec p_total,delta;
846 random=&(moldyn->random);
848 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
850 /* gaussian distribution of velocities */
852 for(i=0;i<moldyn->count;i++) {
853 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
855 v=sigma*rand_get_gauss(random);
857 p_total.x+=atom[i].mass*v;
859 v=sigma*rand_get_gauss(random);
861 p_total.y+=atom[i].mass*v;
863 v=sigma*rand_get_gauss(random);
865 p_total.z+=atom[i].mass*v;
868 /* zero total momentum */
869 v3_scale(&p_total,&p_total,1.0/moldyn->count);
870 for(i=0;i<moldyn->count;i++) {
871 v3_scale(&delta,&p_total,1.0/atom[i].mass);
872 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
875 /* velocity scaling */
876 scale_velocity(moldyn,equi_init);
881 double total_mass_calc(t_moldyn *moldyn) {
887 for(i=0;i<moldyn->count;i++)
888 moldyn->mass+=moldyn->atom[i].mass;
893 double temperature_calc(t_moldyn *moldyn) {
895 /* assume up to date kinetic energy, which is 3/2 N k_B T */
897 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
902 double get_temperature(t_moldyn *moldyn) {
907 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
917 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
920 /* get kinetic energy / temperature & count involved atoms */
923 for(i=0;i<moldyn->count;i++) {
924 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
925 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
930 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
931 else return 0; /* no atoms involved in scaling! */
933 /* (temporary) hack for e,t = 0 */
936 if(moldyn->t_ref!=0.0) {
937 thermal_init(moldyn,equi_init);
941 return 0; /* no scaling needed */
945 /* get scaling factor */
946 scale=moldyn->t_ref/moldyn->t;
950 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
951 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
954 /* velocity scaling */
955 for(i=0;i<moldyn->count;i++) {
956 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
957 v3_scale(&(atom[i].v),&(atom[i].v),scale);
963 double ideal_gas_law_pressure(t_moldyn *moldyn) {
967 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
972 double virial_sum(t_moldyn *moldyn) {
977 /* virial (sum over atom virials) */
985 for(i=0;i<moldyn->count;i++) {
986 virial=&(moldyn->atom[i].virial);
987 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
988 moldyn->vir.xx+=virial->xx;
989 moldyn->vir.yy+=virial->yy;
990 moldyn->vir.zz+=virial->zz;
991 moldyn->vir.xy+=virial->xy;
992 moldyn->vir.xz+=virial->xz;
993 moldyn->vir.yz+=virial->yz;
996 /* global virial (absolute coordinates) */
997 virial=&(moldyn->gvir);
998 moldyn->gv=virial->xx+virial->yy+virial->zz;
1000 return moldyn->virial;
1003 double pressure_calc(t_moldyn *moldyn) {
1007 * with W = 1/3 sum_i f_i r_i (- skipped!)
1008 * virial = sum_i f_i r_i
1010 * => P = (2 Ekin + virial) / (3V)
1013 /* assume up to date virial & up to date kinetic energy */
1015 /* pressure (atom virials) */
1016 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1017 moldyn->p/=(3.0*moldyn->volume);
1019 /* pressure (absolute coordinates) */
1020 moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1021 moldyn->gp/=(3.0*moldyn->volume);
1026 int average_reset(t_moldyn *moldyn) {
1028 printf("[moldyn] average reset\n");
1030 /* update skip value */
1031 moldyn->avg_skip=moldyn->total_steps;
1033 /* kinetic energy */
1037 /* potential energy */
1045 moldyn->virial_sum=0.0;
1056 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1060 if(moldyn->total_steps<moldyn->avg_skip)
1063 denom=moldyn->total_steps+1-moldyn->avg_skip;
1065 /* assume up to date energies, temperature, pressure etc */
1067 /* kinetic energy */
1068 moldyn->k_sum+=moldyn->ekin;
1069 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1070 moldyn->k_avg=moldyn->k_sum/denom;
1071 moldyn->k2_avg=moldyn->k2_sum/denom;
1072 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1074 /* potential energy */
1075 moldyn->v_sum+=moldyn->energy;
1076 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1077 moldyn->v_avg=moldyn->v_sum/denom;
1078 moldyn->v2_avg=moldyn->v2_sum/denom;
1079 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1082 moldyn->t_sum+=moldyn->t;
1083 moldyn->t_avg=moldyn->t_sum/denom;
1086 moldyn->virial_sum+=moldyn->virial;
1087 moldyn->virial_avg=moldyn->virial_sum/denom;
1088 moldyn->gv_sum+=moldyn->gv;
1089 moldyn->gv_avg=moldyn->gv_sum/denom;
1092 moldyn->p_sum+=moldyn->p;
1093 moldyn->p_avg=moldyn->p_sum/denom;
1094 moldyn->gp_sum+=moldyn->gp;
1095 moldyn->gp_avg=moldyn->gp_sum/denom;
1096 moldyn->tp_sum+=moldyn->tp;
1097 moldyn->tp_avg=moldyn->tp_sum/denom;
1102 int get_heat_capacity(t_moldyn *moldyn) {
1106 /* averages needed for heat capacity calc */
1107 if(moldyn->total_steps<moldyn->avg_skip)
1110 /* (temperature average)^2 */
1111 temp2=moldyn->t_avg*moldyn->t_avg;
1112 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1115 /* ideal gas contribution */
1116 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1117 printf(" ideal gas contribution: %f\n",
1118 ighc/moldyn->mass*KILOGRAM/JOULE);
1120 /* specific heat for nvt ensemble */
1121 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1122 moldyn->c_v_nvt/=moldyn->mass;
1124 /* specific heat for nve ensemble */
1125 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1126 moldyn->c_v_nve/=moldyn->mass;
1128 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1129 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1130 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)));
1135 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1151 /* store atomic configuration + dimension */
1152 store=malloc(moldyn->count*sizeof(t_atom));
1154 printf("[moldyn] allocating store mem failed\n");
1157 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1162 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1163 su=pow(2.0-h,ONE_THIRD)-1.0;
1164 dv=(1.0-h)*moldyn->volume;
1166 /* scale up dimension and atom positions */
1167 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1168 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1169 link_cell_shutdown(moldyn);
1170 link_cell_init(moldyn,QUIET);
1171 potential_force_calc(moldyn);
1174 /* restore atomic configuration + dim */
1175 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1178 /* scale down dimension and atom positions */
1179 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1180 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1181 link_cell_shutdown(moldyn);
1182 link_cell_init(moldyn,QUIET);
1183 potential_force_calc(moldyn);
1186 /* calculate pressure */
1187 moldyn->tp=-(y1-y0)/(2.0*dv);
1189 /* restore atomic configuration */
1190 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1192 link_cell_shutdown(moldyn);
1193 link_cell_init(moldyn,QUIET);
1194 //potential_force_calc(moldyn);
1196 /* free store buffer */
1203 double get_pressure(t_moldyn *moldyn) {
1209 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1221 if(x) dim->x*=scale;
1222 if(y) dim->y*=scale;
1223 if(z) dim->z*=scale;
1228 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1239 for(i=0;i<moldyn->count;i++) {
1240 r=&(moldyn->atom[i].r);
1249 int scale_volume(t_moldyn *moldyn) {
1255 vdim=&(moldyn->vis.dim);
1259 /* scaling factor */
1260 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1261 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1262 scale=pow(scale,ONE_THIRD);
1265 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1268 /* scale the atoms and dimensions */
1269 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1270 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1272 /* visualize dimensions */
1279 /* recalculate scaled volume */
1280 moldyn->volume=dim->x*dim->y*dim->z;
1282 /* adjust/reinit linkcell */
1283 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1284 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1285 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1286 link_cell_shutdown(moldyn);
1287 link_cell_init(moldyn,QUIET);
1298 double e_kin_calc(t_moldyn *moldyn) {
1306 for(i=0;i<moldyn->count;i++) {
1307 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1308 moldyn->ekin+=atom[i].ekin;
1311 return moldyn->ekin;
1314 double get_total_energy(t_moldyn *moldyn) {
1316 return(moldyn->ekin+moldyn->energy);
1319 t_3dvec get_total_p(t_moldyn *moldyn) {
1328 for(i=0;i<moldyn->count;i++) {
1329 v3_scale(&p,&(atom[i].v),atom[i].mass);
1330 v3_add(&p_total,&p_total,&p);
1336 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1340 /* nn_dist is the nearest neighbour distance */
1342 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1351 /* linked list / cell method */
1353 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1356 #ifndef LOWMEM_LISTS
1362 /* partitioning the md cell */
1363 lc->nx=moldyn->dim.x/moldyn->cutoff;
1364 lc->x=moldyn->dim.x/lc->nx;
1365 lc->ny=moldyn->dim.y/moldyn->cutoff;
1366 lc->y=moldyn->dim.y/lc->ny;
1367 lc->nz=moldyn->dim.z/moldyn->cutoff;
1368 lc->z=moldyn->dim.z/lc->nz;
1369 lc->cells=lc->nx*lc->ny*lc->nz;
1372 lc->subcell=malloc(lc->cells*sizeof(int*));
1374 lc->subcell=malloc(sizeof(t_lowmem_list));
1376 lc->subcell=malloc(lc->cells*sizeof(t_list));
1379 if(lc->subcell==NULL) {
1380 perror("[moldyn] cell init (malloc)");
1385 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1390 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1393 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1396 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1399 printf(" x: %d x %f A\n",lc->nx,lc->x);
1400 printf(" y: %d x %f A\n",lc->ny,lc->y);
1401 printf(" z: %d x %f A\n",lc->nz,lc->z);
1406 for(i=0;i<lc->cells;i++) {
1407 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1408 if(lc->subcell[i]==NULL) {
1409 perror("[moldyn] list init (malloc)");
1414 printf(" ---> %d malloc %p (%p)\n",
1415 i,lc->subcell[0],lc->subcell);
1419 lc->subcell->head=malloc(lc->cells*sizeof(int));
1420 if(lc->subcell->head==NULL) {
1421 perror("[moldyn] head init (malloc)");
1424 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1425 if(lc->subcell->list==NULL) {
1426 perror("[moldyn] list init (malloc)");
1430 for(i=0;i<lc->cells;i++)
1431 list_init_f(&(lc->subcell[i]));
1434 /* update the list */
1435 link_cell_update(moldyn);
1440 int link_cell_update(t_moldyn *moldyn) {
1458 for(i=0;i<lc->cells;i++)
1460 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1462 lc->subcell->head[i]=-1;
1464 list_destroy_f(&(lc->subcell[i]));
1467 for(count=0;count<moldyn->count;count++) {
1468 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1469 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1470 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1474 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1477 if(p>=MAX_ATOMS_PER_LIST) {
1478 printf("[moldyn] FATAL: amount of atoms too high!\n");
1482 lc->subcell[i+j*nx+k*nxy][p]=count;
1485 lc->subcell->list[count]=lc->subcell->head[p];
1486 lc->subcell->head[p]=count;
1488 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1492 printf(" ---> %d %d malloc %p (%p)\n",
1493 i,count,lc->subcell[i].current,lc->subcell);
1501 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1527 if(i>=nx||j>=ny||k>=nz)
1528 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1531 #ifndef LOWMEM_LISTS
1532 cell[0]=lc->subcell[i+j*nx+k*a];
1534 cell[0]=lc->subcell->head[i+j*nx+k*a];
1536 for(ci=-1;ci<=1;ci++) {
1539 if((x<0)||(x>=nx)) {
1543 for(cj=-1;cj<=1;cj++) {
1546 if((y<0)||(y>=ny)) {
1550 for(ck=-1;ck<=1;ck++) {
1553 if((z<0)||(z>=nz)) {
1557 if(!(ci|cj|ck)) continue;
1559 #ifndef LOWMEM_LISTS
1560 cell[--count2]=lc->subcell[x+y*nx+z*a];
1562 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1567 #ifndef LOWMEM_LISTS
1568 cell[count1++]=lc->subcell[x+y*nx+z*a];
1570 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1582 int link_cell_shutdown(t_moldyn *moldyn) {
1584 #ifndef LOWMEM_LISTS
1592 free(lc->subcell->head);
1593 free(lc->subcell->list);
1597 for(i=0;i<lc->cells;i++) {
1599 free(lc->subcell[i]);
1601 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1602 list_destroy_f(&(lc->subcell[i]));
1612 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1616 t_moldyn_schedule *schedule;
1618 schedule=&(moldyn->schedule);
1619 count=++(schedule->total_sched);
1621 ptr=realloc(schedule->runs,count*sizeof(int));
1623 perror("[moldyn] realloc (runs)");
1627 schedule->runs[count-1]=runs;
1629 ptr=realloc(schedule->tau,count*sizeof(double));
1631 perror("[moldyn] realloc (tau)");
1635 schedule->tau[count-1]=tau;
1637 printf("[moldyn] schedule added:\n");
1638 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1644 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1646 moldyn->schedule.hook=hook;
1647 moldyn->schedule.hook_params=hook_params;
1654 * 'integration of newtons equation' - algorithms
1658 /* start the integration */
1660 int moldyn_integrate(t_moldyn *moldyn) {
1663 unsigned int e,m,s,v,p,t,a;
1665 t_moldyn_schedule *sched;
1670 double energy_scale;
1671 struct timeval t1,t2;
1674 #ifdef VISUAL_THREAD
1676 pthread_t io_thread;
1685 sched=&(moldyn->schedule);
1688 /* initialize linked cell method */
1689 link_cell_init(moldyn,VERBOSE);
1691 /* logging & visualization */
1700 /* sqaure of some variables */
1701 moldyn->tau_square=moldyn->tau*moldyn->tau;
1703 /* get current time */
1704 gettimeofday(&t1,NULL);
1706 /* calculate initial forces */
1707 potential_force_calc(moldyn);
1712 /* some stupid checks before we actually start calculating bullshit */
1713 if(moldyn->cutoff>0.5*moldyn->dim.x)
1714 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
1715 if(moldyn->cutoff>0.5*moldyn->dim.y)
1716 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
1717 if(moldyn->cutoff>0.5*moldyn->dim.z)
1718 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
1720 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
1721 if(ds>0.05*moldyn->nnd)
1722 printf("[moldyn] WARNING: forces too high / tau too small!\n");
1725 /* zero absolute time */
1726 // should have right values!
1728 //moldyn->total_steps=0;
1730 /* debugging, ignore */
1733 /* zero & init moldyn copy */
1734 #ifdef VISUAL_THREAD
1735 memset(&md_copy,0,sizeof(t_moldyn));
1736 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1737 if(atom_copy==NULL) {
1738 perror("[moldyn] malloc atom copy (init)");
1744 printf("##################\n");
1745 printf("# USING PTHREADS #\n");
1746 printf("##################\n");
1748 /* tell the world */
1749 printf("[moldyn] integration start, go get a coffee ...\n");
1751 /* executing the schedule */
1753 while(sched->count<sched->total_sched) {
1755 /* setting amount of runs and finite time step size */
1756 moldyn->tau=sched->tau[sched->count];
1757 moldyn->tau_square=moldyn->tau*moldyn->tau;
1758 moldyn->time_steps=sched->runs[sched->count];
1760 /* energy scaling factor (might change!) */
1761 energy_scale=moldyn->count*EV;
1763 /* integration according to schedule */
1765 for(i=0;i<moldyn->time_steps;i++) {
1767 /* integration step */
1768 moldyn->integrate(moldyn);
1770 /* calculate kinetic energy, temperature and pressure */
1772 temperature_calc(moldyn);
1774 pressure_calc(moldyn);
1776 thermodynamic_pressure_calc(moldyn);
1777 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
1781 /* calculate fluctuations + averages */
1782 average_and_fluctuation_calc(moldyn);
1785 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
1786 scale_velocity(moldyn,FALSE);
1787 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
1788 scale_volume(moldyn);
1790 /* check for log & visualization */
1792 if(!(moldyn->total_steps%e))
1793 dprintf(moldyn->efd,
1795 moldyn->time,moldyn->ekin/energy_scale,
1796 moldyn->energy/energy_scale,
1797 get_total_energy(moldyn)/energy_scale);
1800 if(!(moldyn->total_steps%m)) {
1801 momentum=get_total_p(moldyn);
1802 dprintf(moldyn->mfd,
1803 "%f %f %f %f %f\n",moldyn->time,
1804 momentum.x,momentum.y,momentum.z,
1805 v3_norm(&momentum));
1809 if(!(moldyn->total_steps%p)) {
1810 dprintf(moldyn->pfd,
1811 "%f %f %f %f %f %f %f\n",moldyn->time,
1812 moldyn->p/BAR,moldyn->p_avg/BAR,
1813 moldyn->gp/BAR,moldyn->gp_avg/BAR,
1814 moldyn->tp/BAR,moldyn->tp_avg/BAR);
1818 if(!(moldyn->total_steps%t)) {
1819 dprintf(moldyn->tfd,
1821 moldyn->time,moldyn->t,moldyn->t_avg);
1825 if(!(moldyn->total_steps%v)) {
1826 dprintf(moldyn->vfd,
1827 "%f %f\n",moldyn->time,moldyn->volume);
1831 if(!(moldyn->total_steps%s)) {
1832 snprintf(dir,128,"%s/s-%08.f.save",
1833 moldyn->vlsdir,moldyn->time);
1834 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
1836 if(fd<0) perror("[moldyn] save fd open");
1838 write(fd,moldyn,sizeof(t_moldyn));
1839 write(fd,moldyn->atom,
1840 moldyn->count*sizeof(t_atom));
1846 if(!(moldyn->total_steps%a)) {
1847 #ifdef VISUAL_THREAD
1848 /* check whether thread has not terminated yet */
1850 ret=pthread_join(io_thread,NULL);
1853 /* prepare and start thread */
1854 if(moldyn->count!=md_copy.count) {
1858 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
1860 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1861 if(atom_copy==NULL) {
1862 perror("[moldyn] malloc atom copy (change)");
1866 md_copy.atom=atom_copy;
1867 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
1869 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
1871 perror("[moldyn] create visual atoms thread\n");
1875 visual_atoms(moldyn);
1880 /* display progress */
1882 /* get current time */
1883 gettimeofday(&t2,NULL);
1885 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
1886 sched->count,i,moldyn->total_steps,
1887 moldyn->t,moldyn->t_avg,
1888 moldyn->p/BAR,moldyn->p_avg/BAR,
1889 //moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
1891 (int)(t2.tv_sec-t1.tv_sec));
1895 /* copy over time */
1899 /* increase absolute time */
1900 moldyn->time+=moldyn->tau;
1901 moldyn->total_steps+=1;
1905 /* check for hooks */
1907 printf("\n ## schedule hook %d start ##\n",
1909 sched->hook(moldyn,sched->hook_params);
1910 printf(" ## schedule hook end ##\n");
1913 /* increase the schedule counter */
1921 /* velocity verlet */
1923 int velocity_verlet(t_moldyn *moldyn) {
1926 double tau,tau_square,h;
1931 count=moldyn->count;
1933 tau_square=moldyn->tau_square;
1935 for(i=0;i<count;i++) {
1936 /* check whether fixed atom */
1937 if(atom[i].attr&ATOM_ATTR_FP)
1941 v3_scale(&delta,&(atom[i].v),tau);
1942 v3_add(&(atom[i].r),&(atom[i].r),&delta);
1943 v3_scale(&delta,&(atom[i].f),h*tau_square);
1944 v3_add(&(atom[i].r),&(atom[i].r),&delta);
1945 check_per_bound(moldyn,&(atom[i].r));
1947 /* velocities [actually v(t+tau/2)] */
1948 v3_scale(&delta,&(atom[i].f),h*tau);
1949 v3_add(&(atom[i].v),&(atom[i].v),&delta);
1952 /* criticial check */
1953 moldyn_bc_check(moldyn);
1955 /* neighbour list update */
1956 link_cell_update(moldyn);
1958 /* forces depending on chosen potential */
1960 potential_force_calc(moldyn);
1962 albe_potential_force_calc(moldyn);
1965 for(i=0;i<count;i++) {
1966 /* check whether fixed atom */
1967 if(atom[i].attr&ATOM_ATTR_FP)
1969 /* again velocities [actually v(t+tau)] */
1970 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
1971 v3_add(&(atom[i].v),&(atom[i].v),&delta);
1980 * potentials & corresponding forces & virial routine
1984 /* generic potential and force calculation */
1986 int potential_force_calc(t_moldyn *moldyn) {
1989 t_atom *itom,*jtom,*ktom;
1993 int *neighbour_i[27];
1997 int neighbour_i[27];
2000 t_list neighbour_i[27];
2001 t_list neighbour_i2[27];
2007 count=moldyn->count;
2017 /* reset global virial */
2018 memset(&(moldyn->gvir),0,sizeof(t_virial));
2020 /* reset force, site energy and virial of every atom */
2022 i=omp_get_thread_num();
2023 #pragma omp parallel for private(virial)
2025 for(i=0;i<count;i++) {
2028 v3_zero(&(itom[i].f));
2031 virial=(&(itom[i].virial));
2039 /* reset site energy */
2044 /* get energy, force and virial of every atom */
2046 /* first (and only) loop over atoms i */
2047 for(i=0;i<count;i++) {
2049 /* single particle potential/force */
2050 if(itom[i].attr&ATOM_ATTR_1BP)
2052 moldyn->func1b(moldyn,&(itom[i]));
2054 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2057 /* 2 body pair potential/force */
2059 link_cell_neighbour_index(moldyn,
2060 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2061 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2062 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2067 /* first loop over atoms j */
2068 if(moldyn->func2b) {
2075 while(neighbour_i[j][p]!=-1) {
2077 jtom=&(atom[neighbour_i[j][p]]);
2085 p=lc->subcell->list[p];
2087 this=&(neighbour_i[j]);
2090 if(this->start==NULL)
2094 jtom=this->current->data;
2097 if(jtom==&(itom[i]))
2100 if((jtom->attr&ATOM_ATTR_2BP)&
2101 (itom[i].attr&ATOM_ATTR_2BP)) {
2102 moldyn->func2b(moldyn,
2112 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2118 /* 3 body potential/force */
2120 if(!(itom[i].attr&ATOM_ATTR_3BP))
2123 /* copy the neighbour lists */
2125 /* no copy needed for static lists */
2127 /* no copy needed for lowmem lists */
2129 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2132 /* second loop over atoms j */
2139 while(neighbour_i[j][p]!=-1) {
2141 jtom=&(atom[neighbour_i[j][p]]);
2149 p=lc->subcell->list[p];
2151 this=&(neighbour_i[j]);
2154 if(this->start==NULL)
2159 jtom=this->current->data;
2162 if(jtom==&(itom[i]))
2165 if(!(jtom->attr&ATOM_ATTR_3BP))
2171 if(moldyn->func3b_j1)
2172 moldyn->func3b_j1(moldyn,
2177 /* in first j loop, 3bp run can be skipped */
2178 if(!(moldyn->run3bp))
2181 /* first loop over atoms k */
2182 if(moldyn->func3b_k1) {
2190 while(neighbour_i[k][q]!=-1) {
2192 ktom=&(atom[neighbour_i[k][q]]);
2200 q=lc->subcell->list[q];
2202 that=&(neighbour_i2[k]);
2205 if(that->start==NULL)
2209 ktom=that->current->data;
2212 if(!(ktom->attr&ATOM_ATTR_3BP))
2218 if(ktom==&(itom[i]))
2221 moldyn->func3b_k1(moldyn,
2232 } while(list_next_f(that)!=\
2240 if(moldyn->func3b_j2)
2241 moldyn->func3b_j2(moldyn,
2246 /* second loop over atoms k */
2247 if(moldyn->func3b_k2) {
2255 while(neighbour_i[k][q]!=-1) {
2257 ktom=&(atom[neighbour_i[k][q]]);
2265 q=lc->subcell->list[q];
2267 that=&(neighbour_i2[k]);
2270 if(that->start==NULL)
2274 ktom=that->current->data;
2277 if(!(ktom->attr&ATOM_ATTR_3BP))
2283 if(ktom==&(itom[i]))
2286 moldyn->func3b_k2(moldyn,
2297 } while(list_next_f(that)!=\
2305 /* 2bp post function */
2306 if(moldyn->func3b_j3) {
2307 moldyn->func3b_j3(moldyn,
2316 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2331 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2332 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2334 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2335 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2336 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2340 /* some postprocessing */
2342 #pragma omp parallel for
2344 for(i=0;i<count;i++) {
2345 /* calculate global virial */
2346 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2347 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2348 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2349 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2350 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2351 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2353 /* check forces regarding the given timestep */
2354 if(v3_norm(&(itom[i].f))>\
2355 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2356 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2364 * virial calculation
2367 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2368 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2370 a->virial.xx+=f->x*d->x;
2371 a->virial.yy+=f->y*d->y;
2372 a->virial.zz+=f->z*d->z;
2373 a->virial.xy+=f->x*d->y;
2374 a->virial.xz+=f->x*d->z;
2375 a->virial.yz+=f->y*d->z;
2381 * periodic boundary checking
2384 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2385 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2396 if(moldyn->status&MOLDYN_STAT_PBX) {
2397 if(a->x>=x) a->x-=dim->x;
2398 else if(-a->x>x) a->x+=dim->x;
2400 if(moldyn->status&MOLDYN_STAT_PBY) {
2401 if(a->y>=y) a->y-=dim->y;
2402 else if(-a->y>y) a->y+=dim->y;
2404 if(moldyn->status&MOLDYN_STAT_PBZ) {
2405 if(a->z>=z) a->z-=dim->z;
2406 else if(-a->z>z) a->z+=dim->z;
2413 * debugging / critical check functions
2416 int moldyn_bc_check(t_moldyn *moldyn) {
2429 for(i=0;i<moldyn->count;i++) {
2430 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2431 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2432 i,atom[i].r.x,dim->x/2);
2433 printf("diagnostic:\n");
2434 printf("-----------\natom.r.x:\n");
2436 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2439 ((byte)&(1<<k))?1:0,
2442 printf("---------------\nx=dim.x/2:\n");
2444 memcpy(&byte,(u8 *)(&x)+j,1);
2447 ((byte)&(1<<k))?1:0,
2450 if(atom[i].r.x==x) printf("the same!\n");
2451 else printf("different!\n");
2453 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2454 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2455 i,atom[i].r.y,dim->y/2);
2456 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2457 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2458 i,atom[i].r.z,dim->z/2);
2468 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2475 fd=open(file,O_RDONLY);
2477 perror("[moldyn] load save file open");
2481 fsize=lseek(fd,0,SEEK_END);
2482 lseek(fd,0,SEEK_SET);
2484 size=sizeof(t_moldyn);
2487 cnt=read(fd,moldyn,size);
2489 perror("[moldyn] load save file read (moldyn)");
2495 size=moldyn->count*sizeof(t_atom);
2497 /* correcting possible atom data offset */
2499 if(fsize!=sizeof(t_moldyn)+size) {
2500 corr=fsize-sizeof(t_moldyn)-size;
2501 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2502 printf(" moifying offset:\n");
2503 printf(" - current pos: %d\n",sizeof(t_moldyn));
2504 printf(" - atom size: %d\n",size);
2505 printf(" - file size: %d\n",fsize);
2506 printf(" => correction: %d\n",corr);
2507 lseek(fd,corr,SEEK_CUR);
2510 moldyn->atom=(t_atom *)malloc(size);
2511 if(moldyn->atom==NULL) {
2512 perror("[moldyn] load save file malloc (atoms)");
2517 cnt=read(fd,moldyn->atom,size);
2519 perror("[moldyn] load save file read (atoms)");
2530 int moldyn_free_save_file(t_moldyn *moldyn) {
2537 int moldyn_load(t_moldyn *moldyn) {
2545 * function to find/callback all combinations of 2 body bonds
2548 int process_2b_bonds(t_moldyn *moldyn,void *data,
2549 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2550 void *data,u8 bc)) {
2560 t_list neighbour[27];
2570 for(i=0;i<moldyn->count;i++) {
2571 /* neighbour indexing */
2572 link_cell_neighbour_index(moldyn,
2573 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2574 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2575 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
2580 bc=(j<lc->dnlc)?0:1;
2585 while(neighbour[j][p]!=-1) {
2587 jtom=&(moldyn->atom[neighbour[j][p]]);
2595 p=lc->subcell->list[p];
2597 this=&(neighbour[j]);
2600 if(this->start==NULL)
2605 jtom=this->current->data;
2609 process(moldyn,&(itom[i]),jtom,data,bc);
2616 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2626 * function to find neighboured atoms
2629 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
2630 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
2631 void *data,u8 bc)) {
2641 t_list neighbour[27];
2650 /* neighbour indexing */
2651 link_cell_neighbour_index(moldyn,
2652 (atom->r.x+moldyn->dim.x/2)/lc->x,
2653 (atom->r.y+moldyn->dim.y/2)/lc->x,
2654 (atom->r.z+moldyn->dim.z/2)/lc->x,
2659 bc=(j<lc->dnlc)?0:1;
2664 while(neighbour[j][p]!=-1) {
2666 natom=&(moldyn->atom[neighbour[j][p]]);
2673 natom=&(moldyn->atom[p]);
2674 p=lc->subcell->list[p];
2676 this=&(neighbour[j]);
2679 if(this->start==NULL)
2684 natom=this->current->data;
2688 process(moldyn,atom,natom,data,bc);
2695 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2704 * post processing functions
2707 int get_line(int fd,char *line,int max) {
2714 if(count==max) return count;
2715 ret=read(fd,line+count,1);
2716 if(ret<=0) return ret;
2717 if(line[count]=='\n') {
2718 memset(line+count,0,max-count-1);
2726 int pair_correlation_init(t_moldyn *moldyn,double dr) {
2732 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
2748 for(i=0;i<moldyn->count;i++) {
2750 v3_sub(&dist,&(atom[i].r),&(atom[i].r_0));
2751 check_per_bound(moldyn,&dist);
2752 d2=v3_absolute_square(&dist);
2766 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
2767 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
2768 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
2773 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
2778 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
2779 t_atom *jtom,void *data,u8 bc) {
2786 /* only count pairs once,
2787 * skip same atoms */
2788 if(itom->tag>=jtom->tag)
2792 * pair correlation calc
2799 v3_sub(&dist,&(jtom->r),&(itom->r));
2800 if(bc) check_per_bound(moldyn,&dist);
2801 d=v3_absolute_square(&dist);
2803 /* ignore if greater cutoff */
2804 if(d>moldyn->cutoff_square)
2807 /* fill the slots */
2811 /* should never happen but it does 8) -
2812 * related to -ffloat-store problem! */
2814 printf("[moldyn] WARNING: pcc (%d/%d)",
2820 if(itom->brand!=jtom->brand) {
2825 /* type a - type a bonds */
2827 pcc->stat[s+pcc->o1]+=1;
2829 /* type b - type b bonds */
2830 pcc->stat[s+pcc->o2]+=1;
2836 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
2843 pcc.o1=moldyn->cutoff/dr;
2846 if(pcc.o1*dr<=moldyn->cutoff)
2847 printf("[moldyn] WARNING: pcc (low #slots)\n");
2849 printf("[moldyn] pair correlation calc info:\n");
2850 printf(" time: %f\n",moldyn->time);
2851 printf(" count: %d\n",moldyn->count);
2852 printf(" cutoff: %f\n",moldyn->cutoff);
2853 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
2856 pcc.stat=(double *)ptr;
2859 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
2860 if(pcc.stat==NULL) {
2861 perror("[moldyn] pair correlation malloc");
2866 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
2869 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
2872 for(i=1;i<pcc.o1;i++) {
2873 // normalization: 4 pi r^2 dr
2874 // here: not double counting pairs -> 2 pi r r dr
2875 // ... and actually it's a constant times r^2
2878 pcc.stat[pcc.o1+i]/=norm;
2879 pcc.stat[pcc.o2+i]/=norm;
2884 /* todo: store/print pair correlation function */
2891 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2898 if(itom->tag>=jtom->tag)
2902 v3_sub(&dist,&(jtom->r),&(itom->r));
2903 if(bc) check_per_bound(moldyn,&dist);
2904 d=v3_absolute_square(&dist);
2906 /* ignore if greater or equal cutoff */
2907 if(d>moldyn->cutoff_square)
2910 /* check for potential bond */
2911 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
2914 /* now count this bonding ... */
2917 /* increase total bond counter
2918 * ... double counting!
2923 ba->acnt[jtom->tag]+=1;
2925 ba->bcnt[jtom->tag]+=1;
2928 ba->acnt[itom->tag]+=1;
2930 ba->bcnt[itom->tag]+=1;
2935 int bond_analyze(t_moldyn *moldyn,double *quality) {
2937 // by now: # bonds of type 'a-4b' and 'b-4a' / # bonds total
2945 ba.acnt=malloc(moldyn->count*sizeof(int));
2947 perror("[moldyn] bond analyze malloc (a)");
2950 memset(ba.acnt,0,moldyn->count*sizeof(int));
2952 ba.bcnt=malloc(moldyn->count*sizeof(int));
2954 perror("[moldyn] bond analyze malloc (b)");
2957 memset(ba.bcnt,0,moldyn->count*sizeof(int));
2966 process_2b_bonds(moldyn,&ba,bond_analyze_process);
2968 for(i=0;i<moldyn->count;i++) {
2969 if(atom[i].brand==0) {
2970 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
2974 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
2982 printf("[moldyn] bond analyze: c_cnt=%d | set=%d\n",ccnt,cset);
2983 printf("[moldyn] bond analyze: q_cnt=%d | tot=%d\n",qcnt,ba.tcnt);
2986 quality[0]=1.0*ccnt/cset;
2987 quality[1]=1.0*qcnt/ba.tcnt;
2990 printf("[moldyn] bond analyze: c_bnd_q=%f\n",1.0*qcnt/ba.tcnt);
2991 printf("[moldyn] bond analyze: tot_q=%f\n",1.0*qcnt/ba.tcnt);
2998 * visualization code
3001 int visual_init(t_moldyn *moldyn,char *filebase) {
3003 strncpy(moldyn->vis.fb,filebase,128);
3008 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3015 if(itom->tag>=jtom->tag)
3018 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3021 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3022 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3023 itom->r.x,itom->r.y,itom->r.z,
3024 jtom->r.x,jtom->r.y,jtom->r.z);
3029 #ifdef VISUAL_THREAD
3030 void *visual_atoms(void *ptr) {
3032 int visual_atoms(t_moldyn *moldyn) {
3042 #ifdef VISUAL_THREAD
3056 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3057 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3059 perror("open visual save file fd");
3060 #ifndef VISUAL_THREAD
3065 /* write the actual data file */
3068 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3069 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3071 // atomic configuration
3072 for(i=0;i<moldyn->count;i++)
3073 // atom type, positions, color and kinetic energy
3074 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3078 pse_col[atom[i].element],
3081 // bonds between atoms
3082 #ifndef VISUAL_THREAD
3083 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3088 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3089 -dim.x/2,-dim.y/2,-dim.z/2,
3090 dim.x/2,-dim.y/2,-dim.z/2);
3091 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3092 -dim.x/2,-dim.y/2,-dim.z/2,
3093 -dim.x/2,dim.y/2,-dim.z/2);
3094 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3095 dim.x/2,dim.y/2,-dim.z/2,
3096 dim.x/2,-dim.y/2,-dim.z/2);
3097 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3098 -dim.x/2,dim.y/2,-dim.z/2,
3099 dim.x/2,dim.y/2,-dim.z/2);
3101 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3102 -dim.x/2,-dim.y/2,dim.z/2,
3103 dim.x/2,-dim.y/2,dim.z/2);
3104 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3105 -dim.x/2,-dim.y/2,dim.z/2,
3106 -dim.x/2,dim.y/2,dim.z/2);
3107 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3108 dim.x/2,dim.y/2,dim.z/2,
3109 dim.x/2,-dim.y/2,dim.z/2);
3110 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3111 -dim.x/2,dim.y/2,dim.z/2,
3112 dim.x/2,dim.y/2,dim.z/2);
3114 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3115 -dim.x/2,-dim.y/2,dim.z/2,
3116 -dim.x/2,-dim.y/2,-dim.z/2);
3117 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3118 -dim.x/2,dim.y/2,dim.z/2,
3119 -dim.x/2,dim.y/2,-dim.z/2);
3120 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3121 dim.x/2,-dim.y/2,dim.z/2,
3122 dim.x/2,-dim.y/2,-dim.z/2);
3123 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3124 dim.x/2,dim.y/2,dim.z/2,
3125 dim.x/2,dim.y/2,-dim.z/2);
3130 #ifdef VISUAL_THREAD
3141 * fpu cntrol functions
3144 // set rounding to double (eliminates -ffloat-store!)
3145 int fpu_set_rtd(void) {
3151 ctrl&=~_FPU_EXTENDED;