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,
517 u8 p_type,t_part_vals *p_vals) {
526 pthread_mutex_t *mutex;
532 /* how many atoms do we expect */
535 printf("[moldyn] WARNING: create 'none' lattice called");
537 if(type==CUBIC) new*=1;
538 if(type==FCC) new*=4;
539 if(type==DIAMOND) new*=8;
541 /* allocate space for atoms */
542 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
544 perror("[moldyn] realloc (create lattice)");
548 atom=&(moldyn->atom[count]);
551 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
553 perror("[moldyn] mutex realloc (add atom)");
557 mutex=&(amutex[count]);
560 /* no atoms on the boundaries (only reason: it looks better!) */
574 set_nn_dist(moldyn,lc);
575 ret=cubic_init(a,b,c,lc,atom,&orig,p_type,p_vals);
576 strcpy(name,"cubic");
580 v3_scale(&orig,&orig,0.5);
581 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
582 ret=fcc_init(a,b,c,lc,atom,&orig,p_type,p_vals);
587 v3_scale(&orig,&orig,0.25);
588 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
589 ret=diamond_init(a,b,c,lc,atom,&orig,p_type,p_vals);
590 strcpy(name,"diamond");
593 printf("unknown lattice type (%02x)\n",type);
599 printf("[moldyn] creating %s lattice (lc=%f) incomplete\n",
601 printf(" (ignore in case of partial lattice creation)\n");
602 printf(" amount of atoms\n");
603 printf(" - expected: %d\n",new);
604 printf(" - created: %d\n",ret);
609 printf("[moldyn] created %s lattice with %d atoms\n",name,ret);
611 for(new=0;new<ret;new++) {
612 atom[new].element=element;
615 atom[new].brand=brand;
616 atom[new].tag=count+new;
617 check_per_bound(moldyn,&(atom[new].r));
618 atom[new].r_0=atom[new].r;
620 pthread_mutex_init(&(mutex[new]),NULL);
625 ptr=realloc(moldyn->atom,moldyn->count*sizeof(t_atom));
627 perror("[moldyn] realloc (create lattice - alloc fix)");
633 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
635 perror("[moldyn] list realloc (create lattice)");
638 moldyn->lc.subcell->list=ptr;
641 /* update total system mass */
642 total_mass_calc(moldyn);
647 int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
648 t_3dvec *r,t_3dvec *v) {
655 count=(moldyn->count)++; // asshole style!
657 ptr=realloc(atom,(count+1)*sizeof(t_atom));
659 perror("[moldyn] realloc (add atom)");
665 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
667 perror("[moldyn] list realloc (add atom)");
670 moldyn->lc.subcell->list=ptr;
674 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
676 perror("[moldyn] mutex realloc (add atom)");
680 pthread_mutex_init(&(amutex[count]),NULL);
685 /* initialize new atom */
686 memset(&(atom[count]),0,sizeof(t_atom));
689 atom[count].element=element;
690 atom[count].mass=mass;
691 atom[count].brand=brand;
692 atom[count].tag=count;
693 atom[count].attr=attr;
694 check_per_bound(moldyn,&(atom[count].r));
695 atom[count].r_0=atom[count].r;
697 /* update total system mass */
698 total_mass_calc(moldyn);
703 int del_atom(t_moldyn *moldyn,int tag) {
710 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
712 perror("[moldyn]malloc (del atom)");
716 for(cnt=0;cnt<tag;cnt++)
719 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
721 new[cnt-1].tag=cnt-1;
733 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
734 u8 p_type,t_part_vals *p_vals) {
751 /* shift partition values */
753 p.x=p_vals->p.x+(a*lc)/2.0;
754 p.y=p_vals->p.y+(b*lc)/2.0;
755 p.z=p_vals->p.z+(c*lc)/2.0;
767 if(v3_absolute_square(&dist)<(p_vals->r*p_vals->r)) {
768 v3_copy(&(atom[count].r),&r);
774 if(v3_absolute_square(&dist)>=(p_vals->r*p_vals->r)) {
775 v3_copy(&(atom[count].r),&r);
781 if((fabs(dist.x)<p_vals->d.x)&&
782 (fabs(dist.y)<p_vals->d.y)&&
783 (fabs(dist.z)<p_vals->d.z)) {
784 v3_copy(&(atom[count].r),&r);
790 if((fabs(dist.x)>=p_vals->d.x)||
791 (fabs(dist.y)>=p_vals->d.y)||
792 (fabs(dist.z)>=p_vals->d.z)) {
793 v3_copy(&(atom[count].r),&r);
798 v3_copy(&(atom[count].r),&r);
809 for(i=0;i<count;i++) {
810 atom[i].r.x-=(a*lc)/2.0;
811 atom[i].r.y-=(b*lc)/2.0;
812 atom[i].r.z-=(c*lc)/2.0;
818 /* fcc lattice init */
819 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
820 u8 p_type,t_part_vals *p_vals) {
837 /* construct the basis */
838 memset(basis,0,3*sizeof(t_3dvec));
846 /* shift partition values */
848 p.x=p_vals->p.x+(a*lc)/2.0;
849 p.y=p_vals->p.y+(b*lc)/2.0;
850 p.z=p_vals->p.z+(c*lc)/2.0;
853 /* fill up the room */
864 if(v3_absolute_square(&dist)<(p_vals->r*p_vals->r)) {
865 v3_copy(&(atom[count].r),&r);
871 if(v3_absolute_square(&dist)>=(p_vals->r*p_vals->r)) {
872 v3_copy(&(atom[count].r),&r);
878 if((fabs(dist.x)<p_vals->d.x)&&
879 (fabs(dist.y)<p_vals->d.y)&&
880 (fabs(dist.z)<p_vals->d.z)) {
881 v3_copy(&(atom[count].r),&r);
887 if((fabs(dist.x)>=p_vals->d.x)||
888 (fabs(dist.y)>=p_vals->d.y)||
889 (fabs(dist.z)>=p_vals->d.z)) {
890 v3_copy(&(atom[count].r),&r);
895 v3_copy(&(atom[count].r),&r);
899 /* the three face centered atoms */
901 v3_add(&n,&r,&basis[l]);
905 if(v3_absolute_square(&dist)<(p_vals->r*p_vals->r)) {
906 v3_copy(&(atom[count].r),&n);
912 if(v3_absolute_square(&dist)>=(p_vals->r*p_vals->r)) {
913 v3_copy(&(atom[count].r),&n);
919 if((fabs(dist.x)<p_vals->d.x)&&
920 (fabs(dist.y)<p_vals->d.y)&&
921 (fabs(dist.z)<p_vals->d.z)) {
922 v3_copy(&(atom[count].r),&n);
928 if((fabs(dist.x)>=p_vals->d.x)||
929 (fabs(dist.y)>=p_vals->d.y)||
930 (fabs(dist.z)>=p_vals->d.z)) {
931 v3_copy(&(atom[count].r),&n);
936 v3_copy(&(atom[count].r),&n);
948 /* coordinate transformation */
949 for(i=0;i<count;i++) {
950 atom[i].r.x-=(a*lc)/2.0;
951 atom[i].r.y-=(b*lc)/2.0;
952 atom[i].r.z-=(c*lc)/2.0;
958 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
959 u8 p_type,t_part_vals *p_vals) {
964 count=fcc_init(a,b,c,lc,atom,origin,p_type,p_vals);
970 if(origin) v3_add(&o,&o,origin);
972 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_type,p_vals);
977 int destroy_atoms(t_moldyn *moldyn) {
979 if(moldyn->atom) free(moldyn->atom);
984 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
987 * - gaussian distribution of velocities
988 * - zero total momentum
989 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
994 t_3dvec p_total,delta;
999 random=&(moldyn->random);
1001 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
1003 /* gaussian distribution of velocities */
1005 for(i=0;i<moldyn->count;i++) {
1006 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
1008 v=sigma*rand_get_gauss(random);
1010 p_total.x+=atom[i].mass*v;
1012 v=sigma*rand_get_gauss(random);
1014 p_total.y+=atom[i].mass*v;
1016 v=sigma*rand_get_gauss(random);
1018 p_total.z+=atom[i].mass*v;
1021 /* zero total momentum */
1022 v3_scale(&p_total,&p_total,1.0/moldyn->count);
1023 for(i=0;i<moldyn->count;i++) {
1024 v3_scale(&delta,&p_total,1.0/atom[i].mass);
1025 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
1028 /* velocity scaling */
1029 scale_velocity(moldyn,equi_init);
1034 double total_mass_calc(t_moldyn *moldyn) {
1040 for(i=0;i<moldyn->count;i++)
1041 moldyn->mass+=moldyn->atom[i].mass;
1043 return moldyn->mass;
1046 double temperature_calc(t_moldyn *moldyn) {
1048 /* assume up to date kinetic energy, which is 3/2 N k_B T */
1051 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
1057 double get_temperature(t_moldyn *moldyn) {
1062 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
1072 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1075 /* get kinetic energy / temperature & count involved atoms */
1078 for(i=0;i<moldyn->count;i++) {
1079 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
1080 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1085 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1086 else return 0; /* no atoms involved in scaling! */
1088 /* (temporary) hack for e,t = 0 */
1091 if(moldyn->t_ref!=0.0) {
1092 thermal_init(moldyn,equi_init);
1096 return 0; /* no scaling needed */
1100 /* get scaling factor */
1101 scale=moldyn->t_ref/moldyn->t;
1105 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1106 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1109 /* velocity scaling */
1110 for(i=0;i<moldyn->count;i++) {
1111 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1112 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1118 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1122 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1127 double virial_sum(t_moldyn *moldyn) {
1132 /* virial (sum over atom virials) */
1140 for(i=0;i<moldyn->count;i++) {
1141 virial=&(moldyn->atom[i].virial);
1142 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1143 moldyn->vir.xx+=virial->xx;
1144 moldyn->vir.yy+=virial->yy;
1145 moldyn->vir.zz+=virial->zz;
1146 moldyn->vir.xy+=virial->xy;
1147 moldyn->vir.xz+=virial->xz;
1148 moldyn->vir.yz+=virial->yz;
1151 /* global virial (absolute coordinates) */
1152 //virial=&(moldyn->gvir);
1153 //moldyn->gv=virial->xx+virial->yy+virial->zz;
1155 return moldyn->virial;
1158 double pressure_calc(t_moldyn *moldyn) {
1162 * with W = 1/3 sum_i f_i r_i (- skipped!)
1163 * virial = sum_i f_i r_i
1165 * => P = (2 Ekin + virial) / (3V)
1168 /* assume up to date virial & up to date kinetic energy */
1170 /* pressure (atom virials) */
1171 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1172 moldyn->p/=(3.0*moldyn->volume);
1174 moldyn->px=2.0*moldyn->ekinx+moldyn->vir.xx;
1175 moldyn->px/=moldyn->volume;
1176 moldyn->py=2.0*moldyn->ekiny+moldyn->vir.yy;
1177 moldyn->py/=moldyn->volume;
1178 moldyn->pz=2.0*moldyn->ekinz+moldyn->vir.zz;
1179 moldyn->pz/=moldyn->volume;
1181 /* pressure (absolute coordinates) */
1182 //moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1183 //moldyn->gp/=(3.0*moldyn->volume);
1188 int average_reset(t_moldyn *moldyn) {
1190 printf("[moldyn] average reset\n");
1192 /* update skip value */
1193 moldyn->avg_skip=moldyn->total_steps;
1195 /* kinetic energy */
1199 /* potential energy */
1207 moldyn->virial_sum=0.0;
1208 //moldyn->gv_sum=0.0;
1212 //moldyn->gp_sum=0.0;
1218 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1222 if(moldyn->total_steps<moldyn->avg_skip)
1225 denom=moldyn->total_steps+1-moldyn->avg_skip;
1227 /* assume up to date energies, temperature, pressure etc */
1229 /* kinetic energy */
1230 moldyn->k_sum+=moldyn->ekin;
1231 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1232 moldyn->k_avg=moldyn->k_sum/denom;
1233 moldyn->k2_avg=moldyn->k2_sum/denom;
1234 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1236 /* potential energy */
1237 moldyn->v_sum+=moldyn->energy;
1238 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1239 moldyn->v_avg=moldyn->v_sum/denom;
1240 moldyn->v2_avg=moldyn->v2_sum/denom;
1241 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1244 moldyn->t_sum+=moldyn->t;
1245 moldyn->t_avg=moldyn->t_sum/denom;
1248 moldyn->virial_sum+=moldyn->virial;
1249 moldyn->virial_avg=moldyn->virial_sum/denom;
1250 //moldyn->gv_sum+=moldyn->gv;
1251 //moldyn->gv_avg=moldyn->gv_sum/denom;
1254 moldyn->p_sum+=moldyn->p;
1255 moldyn->p_avg=moldyn->p_sum/denom;
1256 //moldyn->gp_sum+=moldyn->gp;
1257 //moldyn->gp_avg=moldyn->gp_sum/denom;
1258 moldyn->tp_sum+=moldyn->tp;
1259 moldyn->tp_avg=moldyn->tp_sum/denom;
1264 int get_heat_capacity(t_moldyn *moldyn) {
1268 /* averages needed for heat capacity calc */
1269 if(moldyn->total_steps<moldyn->avg_skip)
1272 /* (temperature average)^2 */
1273 temp2=moldyn->t_avg*moldyn->t_avg;
1274 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1277 /* ideal gas contribution */
1278 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1279 printf(" ideal gas contribution: %f\n",
1280 ighc/moldyn->mass*KILOGRAM/JOULE);
1282 /* specific heat for nvt ensemble */
1283 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1284 moldyn->c_v_nvt/=moldyn->mass;
1286 /* specific heat for nve ensemble */
1287 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1288 moldyn->c_v_nve/=moldyn->mass;
1290 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1291 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1292 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)));
1297 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1313 /* store atomic configuration + dimension */
1314 store=malloc(moldyn->count*sizeof(t_atom));
1316 printf("[moldyn] allocating store mem failed\n");
1319 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1324 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1325 su=pow(2.0-h,ONE_THIRD)-1.0;
1326 dv=(1.0-h)*moldyn->volume;
1328 /* scale up dimension and atom positions */
1329 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1330 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1331 link_cell_shutdown(moldyn);
1332 link_cell_init(moldyn,QUIET);
1333 potential_force_calc(moldyn);
1336 /* restore atomic configuration + dim */
1337 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1340 /* scale down dimension and atom positions */
1341 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1342 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1343 link_cell_shutdown(moldyn);
1344 link_cell_init(moldyn,QUIET);
1345 potential_force_calc(moldyn);
1348 /* calculate pressure */
1349 moldyn->tp=-(y1-y0)/(2.0*dv);
1351 /* restore atomic configuration */
1352 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1354 link_cell_shutdown(moldyn);
1355 link_cell_init(moldyn,QUIET);
1356 //potential_force_calc(moldyn);
1358 /* free store buffer */
1365 double get_pressure(t_moldyn *moldyn) {
1371 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1383 if(x) dim->x*=scale;
1384 if(y) dim->y*=scale;
1385 if(z) dim->z*=scale;
1390 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1401 for(i=0;i<moldyn->count;i++) {
1402 r=&(moldyn->atom[i].r);
1411 int scale_atoms_ind(t_moldyn *moldyn,double x,double y,double z) {
1416 for(i=0;i<moldyn->count;i++) {
1417 r=&(moldyn->atom[i].r);
1426 int scale_dim_ind(t_moldyn *moldyn,double x,double y,double z) {
1439 int scale_volume(t_moldyn *moldyn) {
1446 vdim=&(moldyn->vis.dim);
1450 /* scaling factor */
1452 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1453 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1454 scale=pow(scale,ONE_THIRD);
1457 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1461 sx=1.0-(moldyn->p_ref-moldyn->px)*moldyn->p_tc*moldyn->tau;
1462 sy=1.0-(moldyn->p_ref-moldyn->py)*moldyn->p_tc*moldyn->tau;
1463 sz=1.0-(moldyn->p_ref-moldyn->pz)*moldyn->p_tc*moldyn->tau;
1464 sx=pow(sx,ONE_THIRD);
1465 sy=pow(sy,ONE_THIRD);
1466 sz=pow(sz,ONE_THIRD);
1468 /* scale the atoms and dimensions */
1469 //scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1470 //scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1471 scale_atoms_ind(moldyn,sx,sy,sz);
1472 scale_dim_ind(moldyn,sx,sy,sz);
1474 /* visualize dimensions */
1481 /* recalculate scaled volume */
1482 moldyn->volume=dim->x*dim->y*dim->z;
1484 /* adjust/reinit linkcell */
1485 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1486 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1487 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1488 link_cell_shutdown(moldyn);
1489 link_cell_init(moldyn,QUIET);
1503 double e_kin_calc(t_moldyn *moldyn) {
1514 for(i=0;i<moldyn->count;i++) {
1515 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1516 moldyn->ekin+=atom[i].ekin;
1517 moldyn->ekinx+=0.5*atom[i].mass*atom[i].v.x*atom[i].v.x;
1518 moldyn->ekiny+=0.5*atom[i].mass*atom[i].v.y*atom[i].v.y;
1519 moldyn->ekinz+=0.5*atom[i].mass*atom[i].v.z*atom[i].v.z;
1522 return moldyn->ekin;
1525 double get_total_energy(t_moldyn *moldyn) {
1527 return(moldyn->ekin+moldyn->energy);
1530 t_3dvec get_total_p(t_moldyn *moldyn) {
1539 for(i=0;i<moldyn->count;i++) {
1540 v3_scale(&p,&(atom[i].v),atom[i].mass);
1541 v3_add(&p_total,&p_total,&p);
1547 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1551 /* nn_dist is the nearest neighbour distance */
1553 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1562 /* linked list / cell method */
1564 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1567 #ifndef LOWMEM_LISTS
1573 /* partitioning the md cell */
1574 lc->nx=moldyn->dim.x/moldyn->cutoff;
1575 lc->x=moldyn->dim.x/lc->nx;
1576 lc->ny=moldyn->dim.y/moldyn->cutoff;
1577 lc->y=moldyn->dim.y/lc->ny;
1578 lc->nz=moldyn->dim.z/moldyn->cutoff;
1579 lc->z=moldyn->dim.z/lc->nz;
1580 lc->cells=lc->nx*lc->ny*lc->nz;
1583 lc->subcell=malloc(lc->cells*sizeof(int*));
1585 lc->subcell=malloc(sizeof(t_lowmem_list));
1587 lc->subcell=malloc(lc->cells*sizeof(t_list));
1590 if(lc->subcell==NULL) {
1591 perror("[moldyn] cell init (malloc)");
1596 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1601 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1604 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1607 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1610 printf(" x: %d x %f A\n",lc->nx,lc->x);
1611 printf(" y: %d x %f A\n",lc->ny,lc->y);
1612 printf(" z: %d x %f A\n",lc->nz,lc->z);
1617 for(i=0;i<lc->cells;i++) {
1618 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1619 if(lc->subcell[i]==NULL) {
1620 perror("[moldyn] list init (malloc)");
1625 printf(" ---> %d malloc %p (%p)\n",
1626 i,lc->subcell[0],lc->subcell);
1630 lc->subcell->head=malloc(lc->cells*sizeof(int));
1631 if(lc->subcell->head==NULL) {
1632 perror("[moldyn] head init (malloc)");
1635 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1636 if(lc->subcell->list==NULL) {
1637 perror("[moldyn] list init (malloc)");
1641 for(i=0;i<lc->cells;i++)
1642 list_init_f(&(lc->subcell[i]));
1645 /* update the list */
1646 link_cell_update(moldyn);
1651 int link_cell_update(t_moldyn *moldyn) {
1669 for(i=0;i<lc->cells;i++)
1671 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1673 lc->subcell->head[i]=-1;
1675 list_destroy_f(&(lc->subcell[i]));
1678 for(count=0;count<moldyn->count;count++) {
1679 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1680 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1681 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1685 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1688 if(p>=MAX_ATOMS_PER_LIST) {
1689 printf("[moldyn] FATAL: amount of atoms too high!\n");
1693 lc->subcell[i+j*nx+k*nxy][p]=count;
1696 lc->subcell->list[count]=lc->subcell->head[p];
1697 lc->subcell->head[p]=count;
1699 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1703 printf(" ---> %d %d malloc %p (%p)\n",
1704 i,count,lc->subcell[i].current,lc->subcell);
1712 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1738 if(i>=nx||j>=ny||k>=nz)
1739 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1742 #ifndef LOWMEM_LISTS
1743 cell[0]=lc->subcell[i+j*nx+k*a];
1745 cell[0]=lc->subcell->head[i+j*nx+k*a];
1747 for(ci=-1;ci<=1;ci++) {
1750 if((x<0)||(x>=nx)) {
1754 for(cj=-1;cj<=1;cj++) {
1757 if((y<0)||(y>=ny)) {
1761 for(ck=-1;ck<=1;ck++) {
1764 if((z<0)||(z>=nz)) {
1768 if(!(ci|cj|ck)) continue;
1770 #ifndef LOWMEM_LISTS
1771 cell[--count2]=lc->subcell[x+y*nx+z*a];
1773 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1778 #ifndef LOWMEM_LISTS
1779 cell[count1++]=lc->subcell[x+y*nx+z*a];
1781 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1793 int link_cell_shutdown(t_moldyn *moldyn) {
1795 #ifndef LOWMEM_LISTS
1803 free(lc->subcell->head);
1804 free(lc->subcell->list);
1808 for(i=0;i<lc->cells;i++) {
1810 free(lc->subcell[i]);
1812 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1813 list_destroy_f(&(lc->subcell[i]));
1823 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1827 t_moldyn_schedule *schedule;
1829 schedule=&(moldyn->schedule);
1830 count=++(schedule->total_sched);
1832 ptr=realloc(schedule->runs,count*sizeof(int));
1834 perror("[moldyn] realloc (runs)");
1838 schedule->runs[count-1]=runs;
1840 ptr=realloc(schedule->tau,count*sizeof(double));
1842 perror("[moldyn] realloc (tau)");
1846 schedule->tau[count-1]=tau;
1848 printf("[moldyn] schedule added:\n");
1849 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1855 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1857 moldyn->schedule.hook=hook;
1858 moldyn->schedule.hook_params=hook_params;
1865 * 'integration of newtons equation' - algorithms
1869 /* start the integration */
1871 int moldyn_integrate(t_moldyn *moldyn) {
1874 unsigned int e,m,s,v,p,t,a;
1876 t_moldyn_schedule *sched;
1881 double energy_scale;
1882 struct timeval t1,t2;
1885 #ifdef VISUAL_THREAD
1887 pthread_t io_thread;
1896 sched=&(moldyn->schedule);
1899 /* initialize linked cell method */
1900 link_cell_init(moldyn,VERBOSE);
1902 /* logging & visualization */
1911 /* sqaure of some variables */
1912 moldyn->tau_square=moldyn->tau*moldyn->tau;
1914 /* get current time */
1915 gettimeofday(&t1,NULL);
1917 /* calculate initial forces */
1918 potential_force_calc(moldyn);
1923 /* some stupid checks before we actually start calculating bullshit */
1924 if(moldyn->cutoff>0.5*moldyn->dim.x)
1925 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
1926 if(moldyn->cutoff>0.5*moldyn->dim.y)
1927 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
1928 if(moldyn->cutoff>0.5*moldyn->dim.z)
1929 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
1931 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
1932 if(ds>0.05*moldyn->nnd)
1933 printf("[moldyn] WARNING: forces too high / tau too small!\n");
1936 /* zero absolute time */
1937 // should have right values!
1939 //moldyn->total_steps=0;
1941 /* debugging, ignore */
1944 /* zero & init moldyn copy */
1945 #ifdef VISUAL_THREAD
1946 memset(&md_copy,0,sizeof(t_moldyn));
1947 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1948 if(atom_copy==NULL) {
1949 perror("[moldyn] malloc atom copy (init)");
1955 printf("##################\n");
1956 printf("# USING PTHREADS #\n");
1957 printf("##################\n");
1959 /* tell the world */
1960 printf("[moldyn] integration start, go get a coffee ...\n");
1962 /* executing the schedule */
1964 while(sched->count<sched->total_sched) {
1966 /* setting amount of runs and finite time step size */
1967 moldyn->tau=sched->tau[sched->count];
1968 moldyn->tau_square=moldyn->tau*moldyn->tau;
1969 moldyn->time_steps=sched->runs[sched->count];
1971 /* energy scaling factor (might change!) */
1972 energy_scale=moldyn->count*EV;
1974 /* integration according to schedule */
1976 for(i=0;i<moldyn->time_steps;i++) {
1978 /* integration step */
1979 moldyn->integrate(moldyn);
1981 /* calculate kinetic energy, temperature and pressure */
1983 temperature_calc(moldyn);
1985 pressure_calc(moldyn);
1987 thermodynamic_pressure_calc(moldyn);
1988 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
1992 /* calculate fluctuations + averages */
1993 average_and_fluctuation_calc(moldyn);
1996 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
1997 scale_velocity(moldyn,FALSE);
1998 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
1999 scale_volume(moldyn);
2001 /* check for log & visualization */
2003 if(!(moldyn->total_steps%e))
2004 dprintf(moldyn->efd,
2006 moldyn->time,moldyn->ekin/energy_scale,
2007 moldyn->energy/energy_scale,
2008 get_total_energy(moldyn)/energy_scale);
2011 if(!(moldyn->total_steps%m)) {
2012 momentum=get_total_p(moldyn);
2013 dprintf(moldyn->mfd,
2014 "%f %f %f %f %f\n",moldyn->time,
2015 momentum.x,momentum.y,momentum.z,
2016 v3_norm(&momentum));
2020 if(!(moldyn->total_steps%p)) {
2021 dprintf(moldyn->pfd,
2022 "%f %f %f %f %f %f %f\n",moldyn->time,
2023 moldyn->p/BAR,moldyn->p_avg/BAR,
2024 moldyn->p/BAR,moldyn->p_avg/BAR,
2025 moldyn->tp/BAR,moldyn->tp_avg/BAR);
2029 if(!(moldyn->total_steps%t)) {
2030 dprintf(moldyn->tfd,
2032 moldyn->time,moldyn->t,moldyn->t_avg);
2036 if(!(moldyn->total_steps%v)) {
2037 dprintf(moldyn->vfd,
2038 "%f %f %f %f %f\n",moldyn->time,
2046 if(!(moldyn->total_steps%s)) {
2047 snprintf(dir,128,"%s/s-%08.f.save",
2048 moldyn->vlsdir,moldyn->time);
2049 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
2051 if(fd<0) perror("[moldyn] save fd open");
2053 write(fd,moldyn,sizeof(t_moldyn));
2054 write(fd,moldyn->atom,
2055 moldyn->count*sizeof(t_atom));
2061 if(!(moldyn->total_steps%a)) {
2062 #ifdef VISUAL_THREAD
2063 /* check whether thread has not terminated yet */
2065 ret=pthread_join(io_thread,NULL);
2068 /* prepare and start thread */
2069 if(moldyn->count!=md_copy.count) {
2073 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
2075 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2076 if(atom_copy==NULL) {
2077 perror("[moldyn] malloc atom copy (change)");
2081 md_copy.atom=atom_copy;
2082 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
2084 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
2086 perror("[moldyn] create visual atoms thread\n");
2090 visual_atoms(moldyn);
2095 /* display progress */
2097 /* get current time */
2098 gettimeofday(&t2,NULL);
2100 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
2101 sched->count,i,moldyn->total_steps,
2102 moldyn->t,moldyn->t_avg,
2103 moldyn->p/BAR,moldyn->p_avg/BAR,
2104 //moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
2106 (int)(t2.tv_sec-t1.tv_sec));
2110 /* copy over time */
2114 /* increase absolute time */
2115 moldyn->time+=moldyn->tau;
2116 moldyn->total_steps+=1;
2120 /* check for hooks */
2122 printf("\n ## schedule hook %d start ##\n",
2124 sched->hook(moldyn,sched->hook_params);
2125 printf(" ## schedule hook end ##\n");
2128 /* increase the schedule counter */
2136 /* velocity verlet */
2138 int velocity_verlet(t_moldyn *moldyn) {
2141 double tau,tau_square,h;
2146 count=moldyn->count;
2148 tau_square=moldyn->tau_square;
2150 for(i=0;i<count;i++) {
2151 /* check whether fixed atom */
2152 if(atom[i].attr&ATOM_ATTR_FP)
2156 v3_scale(&delta,&(atom[i].v),tau);
2157 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2158 v3_scale(&delta,&(atom[i].f),h*tau_square);
2159 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2160 check_per_bound(moldyn,&(atom[i].r));
2162 /* velocities [actually v(t+tau/2)] */
2163 v3_scale(&delta,&(atom[i].f),h*tau);
2164 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2167 /* criticial check */
2168 moldyn_bc_check(moldyn);
2170 /* neighbour list update */
2171 link_cell_update(moldyn);
2173 /* forces depending on chosen potential */
2175 potential_force_calc(moldyn);
2177 albe_potential_force_calc(moldyn);
2180 for(i=0;i<count;i++) {
2181 /* check whether fixed atom */
2182 if(atom[i].attr&ATOM_ATTR_FP)
2184 /* again velocities [actually v(t+tau)] */
2185 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2186 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2195 * potentials & corresponding forces & virial routine
2199 /* generic potential and force calculation */
2201 int potential_force_calc(t_moldyn *moldyn) {
2204 t_atom *itom,*jtom,*ktom;
2208 int *neighbour_i[27];
2212 int neighbour_i[27];
2215 t_list neighbour_i[27];
2216 t_list neighbour_i2[27];
2222 count=moldyn->count;
2232 /* reset global virial */
2233 memset(&(moldyn->gvir),0,sizeof(t_virial));
2235 /* reset force, site energy and virial of every atom */
2237 i=omp_get_thread_num();
2238 #pragma omp parallel for private(virial)
2240 for(i=0;i<count;i++) {
2243 v3_zero(&(itom[i].f));
2246 virial=(&(itom[i].virial));
2254 /* reset site energy */
2259 /* get energy, force and virial of every atom */
2261 /* first (and only) loop over atoms i */
2262 for(i=0;i<count;i++) {
2264 /* single particle potential/force */
2265 if(itom[i].attr&ATOM_ATTR_1BP)
2267 moldyn->func1b(moldyn,&(itom[i]));
2269 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2272 /* 2 body pair potential/force */
2274 link_cell_neighbour_index(moldyn,
2275 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2276 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2277 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2282 /* first loop over atoms j */
2283 if(moldyn->func2b) {
2290 while(neighbour_i[j][p]!=-1) {
2292 jtom=&(atom[neighbour_i[j][p]]);
2300 p=lc->subcell->list[p];
2302 this=&(neighbour_i[j]);
2305 if(this->start==NULL)
2309 jtom=this->current->data;
2312 if(jtom==&(itom[i]))
2315 if((jtom->attr&ATOM_ATTR_2BP)&
2316 (itom[i].attr&ATOM_ATTR_2BP)) {
2317 moldyn->func2b(moldyn,
2327 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2333 /* 3 body potential/force */
2335 if(!(itom[i].attr&ATOM_ATTR_3BP))
2338 /* copy the neighbour lists */
2340 /* no copy needed for static lists */
2342 /* no copy needed for lowmem lists */
2344 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2347 /* second loop over atoms j */
2354 while(neighbour_i[j][p]!=-1) {
2356 jtom=&(atom[neighbour_i[j][p]]);
2364 p=lc->subcell->list[p];
2366 this=&(neighbour_i[j]);
2369 if(this->start==NULL)
2374 jtom=this->current->data;
2377 if(jtom==&(itom[i]))
2380 if(!(jtom->attr&ATOM_ATTR_3BP))
2386 if(moldyn->func3b_j1)
2387 moldyn->func3b_j1(moldyn,
2392 /* in first j loop, 3bp run can be skipped */
2393 if(!(moldyn->run3bp))
2396 /* first loop over atoms k */
2397 if(moldyn->func3b_k1) {
2405 while(neighbour_i[k][q]!=-1) {
2407 ktom=&(atom[neighbour_i[k][q]]);
2415 q=lc->subcell->list[q];
2417 that=&(neighbour_i2[k]);
2420 if(that->start==NULL)
2424 ktom=that->current->data;
2427 if(!(ktom->attr&ATOM_ATTR_3BP))
2433 if(ktom==&(itom[i]))
2436 moldyn->func3b_k1(moldyn,
2447 } while(list_next_f(that)!=\
2455 if(moldyn->func3b_j2)
2456 moldyn->func3b_j2(moldyn,
2461 /* second loop over atoms k */
2462 if(moldyn->func3b_k2) {
2470 while(neighbour_i[k][q]!=-1) {
2472 ktom=&(atom[neighbour_i[k][q]]);
2480 q=lc->subcell->list[q];
2482 that=&(neighbour_i2[k]);
2485 if(that->start==NULL)
2489 ktom=that->current->data;
2492 if(!(ktom->attr&ATOM_ATTR_3BP))
2498 if(ktom==&(itom[i]))
2501 moldyn->func3b_k2(moldyn,
2512 } while(list_next_f(that)!=\
2520 /* 2bp post function */
2521 if(moldyn->func3b_j3) {
2522 moldyn->func3b_j3(moldyn,
2531 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2546 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2547 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2549 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2550 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2551 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2555 /* some postprocessing */
2557 #pragma omp parallel for
2559 for(i=0;i<count;i++) {
2560 /* calculate global virial */
2561 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2562 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2563 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2564 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2565 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2566 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2568 /* check forces regarding the given timestep */
2569 if(v3_norm(&(itom[i].f))>\
2570 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2571 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2579 * virial calculation
2582 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2583 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2585 a->virial.xx+=f->x*d->x;
2586 a->virial.yy+=f->y*d->y;
2587 a->virial.zz+=f->z*d->z;
2588 a->virial.xy+=f->x*d->y;
2589 a->virial.xz+=f->x*d->z;
2590 a->virial.yz+=f->y*d->z;
2596 * periodic boundary checking
2599 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2600 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2611 if(moldyn->status&MOLDYN_STAT_PBX) {
2612 if(a->x>=x) a->x-=dim->x;
2613 else if(-a->x>x) a->x+=dim->x;
2615 if(moldyn->status&MOLDYN_STAT_PBY) {
2616 if(a->y>=y) a->y-=dim->y;
2617 else if(-a->y>y) a->y+=dim->y;
2619 if(moldyn->status&MOLDYN_STAT_PBZ) {
2620 if(a->z>=z) a->z-=dim->z;
2621 else if(-a->z>z) a->z+=dim->z;
2628 * debugging / critical check functions
2631 int moldyn_bc_check(t_moldyn *moldyn) {
2644 for(i=0;i<moldyn->count;i++) {
2645 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2646 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2647 i,atom[i].r.x,dim->x/2);
2648 printf("diagnostic:\n");
2649 printf("-----------\natom.r.x:\n");
2651 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2654 ((byte)&(1<<k))?1:0,
2657 printf("---------------\nx=dim.x/2:\n");
2659 memcpy(&byte,(u8 *)(&x)+j,1);
2662 ((byte)&(1<<k))?1:0,
2665 if(atom[i].r.x==x) printf("the same!\n");
2666 else printf("different!\n");
2668 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2669 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2670 i,atom[i].r.y,dim->y/2);
2671 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2672 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2673 i,atom[i].r.z,dim->z/2);
2683 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2690 fd=open(file,O_RDONLY);
2692 perror("[moldyn] load save file open");
2696 fsize=lseek(fd,0,SEEK_END);
2697 lseek(fd,0,SEEK_SET);
2699 size=sizeof(t_moldyn);
2702 cnt=read(fd,moldyn,size);
2704 perror("[moldyn] load save file read (moldyn)");
2710 size=moldyn->count*sizeof(t_atom);
2712 /* correcting possible atom data offset */
2714 if(fsize!=sizeof(t_moldyn)+size) {
2715 corr=fsize-sizeof(t_moldyn)-size;
2716 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2717 printf(" moifying offset:\n");
2718 printf(" - current pos: %d\n",sizeof(t_moldyn));
2719 printf(" - atom size: %d\n",size);
2720 printf(" - file size: %d\n",fsize);
2721 printf(" => correction: %d\n",corr);
2722 lseek(fd,corr,SEEK_CUR);
2725 moldyn->atom=(t_atom *)malloc(size);
2726 if(moldyn->atom==NULL) {
2727 perror("[moldyn] load save file malloc (atoms)");
2732 cnt=read(fd,moldyn->atom,size);
2734 perror("[moldyn] load save file read (atoms)");
2741 amutex=malloc(moldyn->count*sizeof(pthread_mutex_t));
2743 perror("[moldyn] load save file (mutexes)");
2746 for(cnt=0;cnt<moldyn->count;cnt++)
2747 pthread_mutex_init(&(amutex[cnt]),NULL);
2755 int moldyn_free_save_file(t_moldyn *moldyn) {
2762 int moldyn_load(t_moldyn *moldyn) {
2770 * function to find/callback all combinations of 2 body bonds
2773 int process_2b_bonds(t_moldyn *moldyn,void *data,
2774 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2775 void *data,u8 bc)) {
2785 t_list neighbour[27];
2795 for(i=0;i<moldyn->count;i++) {
2796 /* neighbour indexing */
2797 link_cell_neighbour_index(moldyn,
2798 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2799 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2800 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
2805 bc=(j<lc->dnlc)?0:1;
2810 while(neighbour[j][p]!=-1) {
2812 jtom=&(moldyn->atom[neighbour[j][p]]);
2820 p=lc->subcell->list[p];
2822 this=&(neighbour[j]);
2825 if(this->start==NULL)
2830 jtom=this->current->data;
2834 process(moldyn,&(itom[i]),jtom,data,bc);
2841 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2851 * function to find neighboured atoms
2854 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
2855 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
2856 void *data,u8 bc)) {
2866 t_list neighbour[27];
2875 /* neighbour indexing */
2876 link_cell_neighbour_index(moldyn,
2877 (atom->r.x+moldyn->dim.x/2)/lc->x,
2878 (atom->r.y+moldyn->dim.y/2)/lc->x,
2879 (atom->r.z+moldyn->dim.z/2)/lc->x,
2884 bc=(j<lc->dnlc)?0:1;
2889 while(neighbour[j][p]!=-1) {
2891 natom=&(moldyn->atom[neighbour[j][p]]);
2898 natom=&(moldyn->atom[p]);
2899 p=lc->subcell->list[p];
2901 this=&(neighbour[j]);
2904 if(this->start==NULL)
2909 natom=this->current->data;
2913 process(moldyn,atom,natom,data,bc);
2920 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2929 * post processing functions
2932 int get_line(int fd,char *line,int max) {
2939 if(count==max) return count;
2940 ret=read(fd,line+count,1);
2941 if(ret<=0) return ret;
2942 if(line[count]=='\n') {
2943 memset(line+count,0,max-count-1);
2951 int pair_correlation_init(t_moldyn *moldyn,double dr) {
2957 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
2973 for(i=0;i<moldyn->count;i++) {
2975 v3_sub(&dist,&(atom[i].r),&(atom[i].r_0));
2976 check_per_bound(moldyn,&dist);
2977 d2=v3_absolute_square(&dist);
2991 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
2992 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
2993 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
2998 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
3003 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
3004 t_atom *jtom,void *data,u8 bc) {
3011 /* only count pairs once,
3012 * skip same atoms */
3013 if(itom->tag>=jtom->tag)
3017 * pair correlation calc
3024 v3_sub(&dist,&(jtom->r),&(itom->r));
3025 if(bc) check_per_bound(moldyn,&dist);
3026 d=v3_absolute_square(&dist);
3028 /* ignore if greater cutoff */
3029 if(d>moldyn->cutoff_square)
3032 /* fill the slots */
3036 /* should never happen but it does 8) -
3037 * related to -ffloat-store problem! */
3039 printf("[moldyn] WARNING: pcc (%d/%d)",
3045 if(itom->brand!=jtom->brand) {
3050 /* type a - type a bonds */
3052 pcc->stat[s+pcc->o1]+=1;
3054 /* type b - type b bonds */
3055 pcc->stat[s+pcc->o2]+=1;
3061 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
3068 pcc.o1=moldyn->cutoff/dr;
3071 if(pcc.o1*dr<=moldyn->cutoff)
3072 printf("[moldyn] WARNING: pcc (low #slots)\n");
3074 printf("[moldyn] pair correlation calc info:\n");
3075 printf(" time: %f\n",moldyn->time);
3076 printf(" count: %d\n",moldyn->count);
3077 printf(" cutoff: %f\n",moldyn->cutoff);
3078 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
3081 pcc.stat=(double *)ptr;
3084 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
3085 if(pcc.stat==NULL) {
3086 perror("[moldyn] pair correlation malloc");
3091 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
3094 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
3097 for(i=1;i<pcc.o1;i++) {
3098 // normalization: 4 pi r^2 dr
3099 // here: not double counting pairs -> 2 pi r r dr
3100 // ... and actually it's a constant times r^2
3103 pcc.stat[pcc.o1+i]/=norm;
3104 pcc.stat[pcc.o2+i]/=norm;
3109 /* todo: store/print pair correlation function */
3116 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3123 if(itom->tag>=jtom->tag)
3127 v3_sub(&dist,&(jtom->r),&(itom->r));
3128 if(bc) check_per_bound(moldyn,&dist);
3129 d=v3_absolute_square(&dist);
3131 /* ignore if greater or equal cutoff */
3132 if(d>moldyn->cutoff_square)
3135 /* check for potential bond */
3136 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3139 /* now count this bonding ... */
3142 /* increase total bond counter
3143 * ... double counting!
3148 ba->acnt[jtom->tag]+=1;
3150 ba->bcnt[jtom->tag]+=1;
3153 ba->acnt[itom->tag]+=1;
3155 ba->bcnt[itom->tag]+=1;
3160 int bond_analyze(t_moldyn *moldyn,double *quality) {
3162 // by now: # bonds of type 'a-4b' and 'b-4a' / # bonds total
3170 ba.acnt=malloc(moldyn->count*sizeof(int));
3172 perror("[moldyn] bond analyze malloc (a)");
3175 memset(ba.acnt,0,moldyn->count*sizeof(int));
3177 ba.bcnt=malloc(moldyn->count*sizeof(int));
3179 perror("[moldyn] bond analyze malloc (b)");
3182 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3191 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3193 for(i=0;i<moldyn->count;i++) {
3194 if(atom[i].brand==0) {
3195 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3199 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3207 printf("[moldyn] bond analyze: c_cnt=%d | set=%d\n",ccnt,cset);
3208 printf("[moldyn] bond analyze: q_cnt=%d | tot=%d\n",qcnt,ba.tcnt);
3211 quality[0]=1.0*ccnt/cset;
3212 quality[1]=1.0*qcnt/ba.tcnt;
3215 printf("[moldyn] bond analyze: c_bnd_q=%f\n",1.0*qcnt/ba.tcnt);
3216 printf("[moldyn] bond analyze: tot_q=%f\n",1.0*qcnt/ba.tcnt);
3223 * visualization code
3226 int visual_init(t_moldyn *moldyn,char *filebase) {
3228 strncpy(moldyn->vis.fb,filebase,128);
3233 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3240 if(itom->tag>=jtom->tag)
3243 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3246 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3247 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3248 itom->r.x,itom->r.y,itom->r.z,
3249 jtom->r.x,jtom->r.y,jtom->r.z);
3254 #ifdef VISUAL_THREAD
3255 void *visual_atoms(void *ptr) {
3257 int visual_atoms(t_moldyn *moldyn) {
3267 #ifdef VISUAL_THREAD
3281 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3282 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3284 perror("open visual save file fd");
3285 #ifndef VISUAL_THREAD
3290 /* write the actual data file */
3293 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3294 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3296 // atomic configuration
3297 for(i=0;i<moldyn->count;i++)
3298 // atom type, positions, color and kinetic energy
3299 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3303 pse_col[atom[i].element],
3306 // bonds between atoms
3307 #ifndef VISUAL_THREAD
3308 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3313 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3314 -dim.x/2,-dim.y/2,-dim.z/2,
3315 dim.x/2,-dim.y/2,-dim.z/2);
3316 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3317 -dim.x/2,-dim.y/2,-dim.z/2,
3318 -dim.x/2,dim.y/2,-dim.z/2);
3319 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3320 dim.x/2,dim.y/2,-dim.z/2,
3321 dim.x/2,-dim.y/2,-dim.z/2);
3322 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3323 -dim.x/2,dim.y/2,-dim.z/2,
3324 dim.x/2,dim.y/2,-dim.z/2);
3326 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3327 -dim.x/2,-dim.y/2,dim.z/2,
3328 dim.x/2,-dim.y/2,dim.z/2);
3329 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3330 -dim.x/2,-dim.y/2,dim.z/2,
3331 -dim.x/2,dim.y/2,dim.z/2);
3332 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3333 dim.x/2,dim.y/2,dim.z/2,
3334 dim.x/2,-dim.y/2,dim.z/2);
3335 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3336 -dim.x/2,dim.y/2,dim.z/2,
3337 dim.x/2,dim.y/2,dim.z/2);
3339 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3340 -dim.x/2,-dim.y/2,dim.z/2,
3341 -dim.x/2,-dim.y/2,-dim.z/2);
3342 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3343 -dim.x/2,dim.y/2,dim.z/2,
3344 -dim.x/2,dim.y/2,-dim.z/2);
3345 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3346 dim.x/2,-dim.y/2,dim.z/2,
3347 dim.x/2,-dim.y/2,-dim.z/2);
3348 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3349 dim.x/2,dim.y/2,dim.z/2,
3350 dim.x/2,dim.y/2,-dim.z/2);
3355 #ifdef VISUAL_THREAD
3366 * fpu cntrol functions
3369 // set rounding to double (eliminates -ffloat-store!)
3370 int fpu_set_rtd(void) {
3376 ctrl&=~_FPU_EXTENDED;