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"
54 pthread_mutex_t *amutex;
55 pthread_mutex_t emutex;
59 * the moldyn functions
62 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
64 printf("[moldyn] init\n");
66 /* only needed if compiled without -msse2 (float-store prob!) */
69 memset(moldyn,0,sizeof(t_moldyn));
74 rand_init(&(moldyn->random),NULL,1);
75 moldyn->random.status|=RAND_STAT_VERBOSE;
78 pthread_mutex_init(&emutex,NULL);
81 #ifdef CONSTRAINT_110_5832
82 printf("\n\n\nWARNING! WARNING! WARNING!\n\n\n");
83 printf("\n\n\n!! -- constraints enabled -- !!\n\n\n");
85 #ifdef CONSTRAINT_11X_5832
86 printf("\n\n\nWARNING! WARNING! WARNING!\n\n\n");
87 printf("\n\n\n!! -- constraints enabled -- !!\n\n\n");
92 int moldyn_shutdown(t_moldyn *moldyn) {
98 printf("[moldyn] shutdown\n");
101 for(i=0;i<moldyn->count;i++)
102 pthread_mutex_destroy(&(amutex[i]));
104 pthread_mutex_destroy(&emutex);
107 moldyn_log_shutdown(moldyn);
108 link_cell_shutdown(moldyn);
109 rand_close(&(moldyn->random));
115 int set_int_alg(t_moldyn *moldyn,u8 algo) {
117 printf("[moldyn] integration algorithm: ");
120 case MOLDYN_INTEGRATE_VERLET:
121 moldyn->integrate=velocity_verlet;
122 printf("velocity verlet\n");
125 printf("unknown integration algorithm: %02x\n",algo);
133 int set_cutoff(t_moldyn *moldyn,double cutoff) {
135 moldyn->cutoff=cutoff;
136 moldyn->cutoff_square=cutoff*cutoff;
138 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
143 int set_temperature(t_moldyn *moldyn,double t_ref) {
147 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
152 int set_pressure(t_moldyn *moldyn,double p_ref) {
156 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
161 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
163 moldyn->pt_scale&=(~(P_SCALE_MASK));
164 moldyn->pt_scale|=ptype;
167 printf("[moldyn] p scaling:\n");
169 printf(" p: %s",ptype?"yes":"no ");
171 printf(" | type: %02x | factor: %f",ptype,ptc);
177 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
179 moldyn->pt_scale&=(~(T_SCALE_MASK));
180 moldyn->pt_scale|=ttype;
183 printf("[moldyn] t scaling:\n");
185 printf(" t: %s",ttype?"yes":"no ");
187 printf(" | type: %02x | factor: %f",ttype,ttc);
193 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
195 moldyn->pt_scale=(ptype|ttype);
199 printf("[moldyn] p/t scaling:\n");
201 printf(" p: %s",ptype?"yes":"no ");
203 printf(" | type: %02x | factor: %f",ptype,ptc);
206 printf(" t: %s",ttype?"yes":"no ");
208 printf(" | type: %02x | factor: %f",ttype,ttc);
214 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
220 moldyn->volume=x*y*z;
228 printf("[moldyn] dimensions in A and A^3 respectively:\n");
229 printf(" x: %f\n",moldyn->dim.x);
230 printf(" y: %f\n",moldyn->dim.y);
231 printf(" z: %f\n",moldyn->dim.z);
232 printf(" volume: %f\n",moldyn->volume);
233 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
238 int set_nn_dist(t_moldyn *moldyn,double dist) {
245 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
247 printf("[moldyn] periodic boundary conditions:\n");
250 moldyn->status|=MOLDYN_STAT_PBX;
253 moldyn->status|=MOLDYN_STAT_PBY;
256 moldyn->status|=MOLDYN_STAT_PBZ;
258 printf(" x: %s\n",x?"yes":"no");
259 printf(" y: %s\n",y?"yes":"no");
260 printf(" z: %s\n",z?"yes":"no");
265 int set_potential(t_moldyn *moldyn,u8 type) {
268 case MOLDYN_POTENTIAL_TM:
269 //moldyn->func1b=tersoff_mult_1bp;
270 moldyn->func3b_j1=tersoff_mult_3bp_j1;
271 moldyn->func3b_k1=tersoff_mult_3bp_k1;
272 moldyn->func3b_j2=tersoff_mult_3bp_j2;
273 moldyn->func3b_k2=tersoff_mult_3bp_k2;
274 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
276 case MOLDYN_POTENTIAL_AM:
277 moldyn->func3b_j1=albe_mult_3bp_j1;
278 moldyn->func3b_k1=albe_mult_3bp_k1;
279 moldyn->func3b_j2=albe_mult_3bp_j2;
280 moldyn->func3b_k2=albe_mult_3bp_k2;
281 moldyn->check_2b_bond=albe_mult_check_2b_bond;
283 case MOLDYN_POTENTIAL_HO:
284 moldyn->func2b=harmonic_oscillator;
285 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
287 case MOLDYN_POTENTIAL_LJ:
288 moldyn->func2b=lennard_jones;
289 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
292 printf("[moldyn] set potential: unknown type %02x\n",
300 int set_avg_skip(t_moldyn *moldyn,int skip) {
302 printf("[moldyn] skip %d steps before starting average calc\n",skip);
303 moldyn->avg_skip=skip;
308 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
310 strncpy(moldyn->vlsdir,dir,127);
315 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
317 strncpy(moldyn->rauthor,author,63);
318 strncpy(moldyn->rtitle,title,63);
323 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
328 printf("[moldyn] set log: ");
331 case LOG_TOTAL_ENERGY:
332 moldyn->ewrite=timer;
333 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
334 moldyn->efd=open(filename,
335 O_WRONLY|O_CREAT|O_EXCL,
338 perror("[moldyn] energy log fd open");
341 dprintf(moldyn->efd,"# total energy log file\n");
342 printf("total energy (%d)\n",timer);
344 case LOG_TOTAL_MOMENTUM:
345 moldyn->mwrite=timer;
346 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
347 moldyn->mfd=open(filename,
348 O_WRONLY|O_CREAT|O_EXCL,
351 perror("[moldyn] momentum log fd open");
354 dprintf(moldyn->efd,"# total momentum log file\n");
355 printf("total momentum (%d)\n",timer);
358 moldyn->pwrite=timer;
359 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
360 moldyn->pfd=open(filename,
361 O_WRONLY|O_CREAT|O_EXCL,
364 perror("[moldyn] pressure log file\n");
367 dprintf(moldyn->pfd,"# pressure log file\n");
368 printf("pressure (%d)\n",timer);
370 case LOG_TEMPERATURE:
371 moldyn->twrite=timer;
372 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
373 moldyn->tfd=open(filename,
374 O_WRONLY|O_CREAT|O_EXCL,
377 perror("[moldyn] temperature log file\n");
380 dprintf(moldyn->tfd,"# temperature log file\n");
381 printf("temperature (%d)\n",timer);
384 moldyn->vwrite=timer;
385 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
386 moldyn->vfd=open(filename,
387 O_WRONLY|O_CREAT|O_EXCL,
390 perror("[moldyn] volume log file\n");
393 dprintf(moldyn->vfd,"# volume log file\n");
394 printf("volume (%d)\n",timer);
397 moldyn->swrite=timer;
398 printf("save file (%d)\n",timer);
401 moldyn->awrite=timer;
402 ret=visual_init(moldyn,moldyn->vlsdir);
404 printf("[moldyn] visual init failure\n");
407 printf("visual file (%d)\n",timer);
410 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
411 moldyn->rfd=open(filename,
412 O_WRONLY|O_CREAT|O_EXCL,
415 perror("[moldyn] report fd open");
418 printf("report -> ");
420 snprintf(filename,127,"%s/e_plot.scr",
422 moldyn->epfd=open(filename,
423 O_WRONLY|O_CREAT|O_EXCL,
426 perror("[moldyn] energy plot fd open");
429 dprintf(moldyn->epfd,e_plot_script);
434 snprintf(filename,127,"%s/pressure_plot.scr",
436 moldyn->ppfd=open(filename,
437 O_WRONLY|O_CREAT|O_EXCL,
440 perror("[moldyn] p plot fd open");
443 dprintf(moldyn->ppfd,pressure_plot_script);
448 snprintf(filename,127,"%s/temperature_plot.scr",
450 moldyn->tpfd=open(filename,
451 O_WRONLY|O_CREAT|O_EXCL,
454 perror("[moldyn] t plot fd open");
457 dprintf(moldyn->tpfd,temperature_plot_script);
459 printf("temperature ");
461 dprintf(moldyn->rfd,report_start,
462 moldyn->rauthor,moldyn->rtitle);
466 printf("unknown log type: %02x\n",type);
473 int moldyn_log_shutdown(t_moldyn *moldyn) {
477 printf("[moldyn] log shutdown\n");
481 dprintf(moldyn->rfd,report_energy);
482 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
487 if(moldyn->mfd) close(moldyn->mfd);
491 dprintf(moldyn->rfd,report_pressure);
492 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
499 dprintf(moldyn->rfd,report_temperature);
500 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
505 dprintf(moldyn->rfd,report_end);
507 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
510 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
513 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
522 * creating lattice functions
525 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,
526 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin,
527 t_part_params *p_params,t_defect_params *d_params,
528 t_offset_params *o_params) {
537 pthread_mutex_t *mutex;
543 /* how many atoms do we expect */
546 printf("[moldyn] WARNING: create 'none' lattice called");
548 if(type==CUBIC) new*=1;
549 if(type==FCC) new*=4;
550 if(type==DIAMOND) new*=8;
554 switch(d_params->stype) {
555 case DEFECT_STYPE_DB_X:
556 case DEFECT_STYPE_DB_Y:
557 case DEFECT_STYPE_DB_Z:
558 case DEFECT_STYPE_DB_R:
562 printf("[moldyn] WARNING: cl unknown defect\n");
567 /* allocate space for atoms */
568 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
570 perror("[moldyn] realloc (create lattice)");
574 atom=&(moldyn->atom[count]);
577 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
579 perror("[moldyn] mutex realloc (add atom)");
583 mutex=&(amutex[count]);
586 /* no atoms on the boundaries (only reason: it looks better!) */
601 v3_add(&orig,&orig,&(o_params->o));
602 set_nn_dist(moldyn,lc);
603 ret=cubic_init(a,b,c,lc,atom,&orig,p_params,d_params);
604 strcpy(name,"cubic");
608 v3_scale(&orig,&orig,0.5);
610 v3_add(&orig,&orig,&(o_params->o));
611 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
612 ret=fcc_init(a,b,c,lc,atom,&orig,p_params,d_params);
617 v3_scale(&orig,&orig,0.25);
619 v3_add(&orig,&orig,&(o_params->o));
620 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
621 ret=diamond_init(a,b,c,lc,atom,&orig,p_params,d_params);
622 strcpy(name,"diamond");
625 printf("unknown lattice type (%02x)\n",type);
631 printf("[moldyn] creating %s lattice (lc=%f) incomplete\n",
633 printf(" (ignore for partial lattice creation)\n");
634 printf(" amount of atoms\n");
635 printf(" - expected: %d\n",new);
636 printf(" - created: %d\n",ret);
641 printf("[moldyn] created %s lattice with %d atoms\n",name,ret);
643 for(new=0;new<ret;new++) {
644 atom[new].element=element;
645 atom[new].mass=pse_mass[element];
647 atom[new].brand=brand;
648 atom[new].tag=count+new;
649 check_per_bound(moldyn,&(atom[new].r));
650 atom[new].r_0=atom[new].r;
652 pthread_mutex_init(&(mutex[new]),NULL);
656 atom[new].element=d_params->element;
657 atom[new].mass=pse_mass[d_params->element];
658 atom[new].attr=d_params->attr;
659 atom[new].brand=d_params->brand;
660 atom[new].tag=count+new;
661 check_per_bound(moldyn,&(atom[new].r));
662 atom[new].r_0=atom[new].r;
664 pthread_mutex_init(&(mutex[new]),NULL);
670 ptr=realloc(moldyn->atom,moldyn->count*sizeof(t_atom));
672 perror("[moldyn] realloc (create lattice - alloc fix)");
677 // WHAT ABOUT AMUTEX !!!!
680 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
682 perror("[moldyn] list realloc (create lattice)");
685 moldyn->lc.subcell->list=ptr;
688 /* update total system mass */
689 total_mass_calc(moldyn);
694 int add_atom(t_moldyn *moldyn,int element,u8 brand,u8 attr,
695 t_3dvec *r,t_3dvec *v) {
702 count=(moldyn->count)++; // asshole style!
704 ptr=realloc(atom,(count+1)*sizeof(t_atom));
706 perror("[moldyn] realloc (add atom)");
712 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
714 perror("[moldyn] list realloc (add atom)");
717 moldyn->lc.subcell->list=ptr;
721 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
723 perror("[moldyn] mutex realloc (add atom)");
727 pthread_mutex_init(&(amutex[count]),NULL);
732 /* initialize new atom */
733 memset(&(atom[count]),0,sizeof(t_atom));
736 atom[count].element=element;
737 atom[count].mass=pse_mass[element];
738 atom[count].brand=brand;
739 atom[count].tag=count;
740 atom[count].attr=attr;
741 check_per_bound(moldyn,&(atom[count].r));
742 atom[count].r_0=atom[count].r;
744 /* update total system mass */
745 total_mass_calc(moldyn);
750 int del_atom(t_moldyn *moldyn,int tag) {
754 #if defined LOWMEM_LISTS || defined PTHREADS
760 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
762 perror("[moldyn]malloc (del atom)");
766 for(cnt=0;cnt<tag;cnt++)
769 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
771 new[cnt-1].tag=cnt-1;
780 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
782 perror("[moldyn] list realloc (del atom)");
785 moldyn->lc.subcell->list=ptr;
787 link_cell_update(moldyn);
791 ptr=realloc(amutex,moldyn->count*sizeof(pthread_mutex_t));
793 perror("[moldyn] mutex realloc (add atom)");
797 pthread_mutex_destroy(&(amutex[moldyn->count+1]));
804 #define set_atom_positions(pos) \
805 if(d_params->type) {\
806 d_o.x=0; d_o.y=0; d_o.z=0;\
807 d_d.x=0; d_d.y=0; d_d.z=0;\
808 switch(d_params->stype) {\
809 case DEFECT_STYPE_DB_X:\
813 case DEFECT_STYPE_DB_Y:\
817 case DEFECT_STYPE_DB_Z:\
821 case DEFECT_STYPE_DB_R:\
824 printf("[moldyn] WARNING: unknown defect\n");\
827 v3_add(&dr,&pos,&d_o);\
828 v3_copy(&(atom[count].r),&dr);\
830 v3_add(&dr,&pos,&d_d);\
831 v3_copy(&(atom[count].r),&dr);\
835 v3_copy(&(atom[count].r),&pos);\
840 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
841 t_part_params *p_params,t_defect_params *d_params) {
861 /* shift partition values */
863 p.x=p_params->p.x+(a*lc)/2.0;
864 p.y=p_params->p.y+(b*lc)/2.0;
865 p.z=p_params->p.z+(c*lc)/2.0;
874 switch(p_params->type) {
877 if(v3_absolute_square(&dist)<
878 (p_params->r*p_params->r)) {
879 set_atom_positions(r);
884 if(v3_absolute_square(&dist)>=
885 (p_params->r*p_params->r)) {
886 set_atom_positions(r);
891 if((fabs(dist.x)<p_params->d.x)&&
892 (fabs(dist.y)<p_params->d.y)&&
893 (fabs(dist.z)<p_params->d.z)) {
894 set_atom_positions(r);
899 if((fabs(dist.x)>=p_params->d.x)||
900 (fabs(dist.y)>=p_params->d.y)||
901 (fabs(dist.z)>=p_params->d.z)) {
902 set_atom_positions(r);
906 set_atom_positions(r);
916 for(i=0;i<count;i++) {
917 atom[i].r.x-=(a*lc)/2.0;
918 atom[i].r.y-=(b*lc)/2.0;
919 atom[i].r.z-=(c*lc)/2.0;
925 /* fcc lattice init */
926 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
927 t_part_params *p_params,t_defect_params *d_params) {
945 /* construct the basis */
946 memset(basis,0,3*sizeof(t_3dvec));
954 /* shift partition values */
956 p.x=p_params->p.x+(a*lc)/2.0;
957 p.y=p_params->p.y+(b*lc)/2.0;
958 p.z=p_params->p.z+(c*lc)/2.0;
961 /* fill up the room */
969 switch(p_params->type) {
972 if(v3_absolute_square(&dist)<
973 (p_params->r*p_params->r)) {
974 set_atom_positions(r);
979 if(v3_absolute_square(&dist)>=
980 (p_params->r*p_params->r)) {
981 set_atom_positions(r);
986 if((fabs(dist.x)<p_params->d.x)&&
987 (fabs(dist.y)<p_params->d.y)&&
988 (fabs(dist.z)<p_params->d.z)) {
989 set_atom_positions(r);
994 if((fabs(dist.x)>=p_params->d.x)||
995 (fabs(dist.y)>=p_params->d.y)||
996 (fabs(dist.z)>=p_params->d.z)) {
997 set_atom_positions(r);
1001 set_atom_positions(r);
1004 /* the three face centered atoms */
1006 v3_add(&n,&r,&basis[l]);
1007 switch(p_params->type) {
1009 v3_sub(&dist,&n,&p);
1010 if(v3_absolute_square(&dist)<
1011 (p_params->r*p_params->r)) {
1012 set_atom_positions(n);
1015 case PART_OUTSIDE_R:
1016 v3_sub(&dist,&n,&p);
1017 if(v3_absolute_square(&dist)>=
1018 (p_params->r*p_params->r)) {
1019 set_atom_positions(n);
1023 v3_sub(&dist,&n,&p);
1024 if((fabs(dist.x)<p_params->d.x)&&
1025 (fabs(dist.y)<p_params->d.y)&&
1026 (fabs(dist.z)<p_params->d.z)) {
1027 set_atom_positions(n);
1030 case PART_OUTSIDE_D:
1031 v3_sub(&dist,&n,&p);
1032 if((fabs(dist.x)>=p_params->d.x)||
1033 (fabs(dist.y)>=p_params->d.y)||
1034 (fabs(dist.z)>=p_params->d.z)) {
1035 set_atom_positions(n);
1039 set_atom_positions(n);
1050 /* coordinate transformation */
1051 for(i=0;i<count;i++) {
1052 atom[i].r.x-=(a*lc)/2.0;
1053 atom[i].r.y-=(b*lc)/2.0;
1054 atom[i].r.z-=(c*lc)/2.0;
1060 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
1061 t_part_params *p_params,t_defect_params *d_params) {
1066 count=fcc_init(a,b,c,lc,atom,origin,p_params,d_params);
1072 if(origin) v3_add(&o,&o,origin);
1074 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_params,d_params);
1079 int destroy_atoms(t_moldyn *moldyn) {
1081 if(moldyn->atom) free(moldyn->atom);
1086 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
1089 * - gaussian distribution of velocities
1090 * - zero total momentum
1091 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1096 t_3dvec p_total,delta;
1101 random=&(moldyn->random);
1103 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
1105 /* gaussian distribution of velocities */
1107 for(i=0;i<moldyn->count;i++) {
1108 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
1110 v=sigma*rand_get_gauss(random);
1112 p_total.x+=atom[i].mass*v;
1114 v=sigma*rand_get_gauss(random);
1116 p_total.y+=atom[i].mass*v;
1118 v=sigma*rand_get_gauss(random);
1120 p_total.z+=atom[i].mass*v;
1123 /* zero total momentum */
1124 v3_scale(&p_total,&p_total,1.0/moldyn->count);
1125 for(i=0;i<moldyn->count;i++) {
1126 v3_scale(&delta,&p_total,1.0/atom[i].mass);
1127 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
1130 /* velocity scaling */
1131 scale_velocity(moldyn,equi_init);
1136 double total_mass_calc(t_moldyn *moldyn) {
1142 for(i=0;i<moldyn->count;i++)
1143 moldyn->mass+=moldyn->atom[i].mass;
1145 return moldyn->mass;
1148 double temperature_calc(t_moldyn *moldyn) {
1150 /* assume up to date kinetic energy, which is 3/2 N k_B T */
1153 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
1159 double get_temperature(t_moldyn *moldyn) {
1164 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
1174 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1177 /* get kinetic energy / temperature & count involved atoms */
1180 for(i=0;i<moldyn->count;i++) {
1181 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
1182 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1187 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1188 else return 0; /* no atoms involved in scaling! */
1190 /* (temporary) hack for e,t = 0 */
1193 if(moldyn->t_ref!=0.0) {
1194 thermal_init(moldyn,equi_init);
1198 return 0; /* no scaling needed */
1202 /* get scaling factor */
1203 scale=moldyn->t_ref/moldyn->t;
1207 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1208 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1211 /* velocity scaling */
1212 for(i=0;i<moldyn->count;i++) {
1213 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1214 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1220 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1224 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1229 double virial_sum(t_moldyn *moldyn) {
1234 /* virial (sum over atom virials) */
1242 for(i=0;i<moldyn->count;i++) {
1243 virial=&(moldyn->atom[i].virial);
1244 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1245 moldyn->vir.xx+=virial->xx;
1246 moldyn->vir.yy+=virial->yy;
1247 moldyn->vir.zz+=virial->zz;
1248 moldyn->vir.xy+=virial->xy;
1249 moldyn->vir.xz+=virial->xz;
1250 moldyn->vir.yz+=virial->yz;
1253 /* global virial (absolute coordinates) */
1254 //virial=&(moldyn->gvir);
1255 //moldyn->gv=virial->xx+virial->yy+virial->zz;
1257 return moldyn->virial;
1260 double pressure_calc(t_moldyn *moldyn) {
1264 * with W = 1/3 sum_i f_i r_i (- skipped!)
1265 * virial = sum_i f_i r_i
1267 * => P = (2 Ekin + virial) / (3V)
1270 /* assume up to date virial & up to date kinetic energy */
1272 /* pressure (atom virials) */
1273 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1274 moldyn->p/=(3.0*moldyn->volume);
1276 //moldyn->px=2.0*moldyn->ekinx+moldyn->vir.xx;
1277 //moldyn->px/=moldyn->volume;
1278 //moldyn->py=2.0*moldyn->ekiny+moldyn->vir.yy;
1279 //moldyn->py/=moldyn->volume;
1280 //moldyn->pz=2.0*moldyn->ekinz+moldyn->vir.zz;
1281 //moldyn->pz/=moldyn->volume;
1283 /* pressure (absolute coordinates) */
1284 //moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1285 //moldyn->gp/=(3.0*moldyn->volume);
1290 int average_reset(t_moldyn *moldyn) {
1292 printf("[moldyn] average reset\n");
1294 /* update skip value */
1295 moldyn->avg_skip=moldyn->total_steps;
1297 /* kinetic energy */
1301 /* potential energy */
1309 moldyn->virial_sum=0.0;
1310 //moldyn->gv_sum=0.0;
1314 //moldyn->gp_sum=0.0;
1320 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1324 if(moldyn->total_steps<moldyn->avg_skip)
1327 denom=moldyn->total_steps+1-moldyn->avg_skip;
1329 /* assume up to date energies, temperature, pressure etc */
1331 /* kinetic energy */
1332 moldyn->k_sum+=moldyn->ekin;
1333 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1334 moldyn->k_avg=moldyn->k_sum/denom;
1335 moldyn->k2_avg=moldyn->k2_sum/denom;
1336 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1338 /* potential energy */
1339 moldyn->v_sum+=moldyn->energy;
1340 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1341 moldyn->v_avg=moldyn->v_sum/denom;
1342 moldyn->v2_avg=moldyn->v2_sum/denom;
1343 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1346 moldyn->t_sum+=moldyn->t;
1347 moldyn->t_avg=moldyn->t_sum/denom;
1350 moldyn->virial_sum+=moldyn->virial;
1351 moldyn->virial_avg=moldyn->virial_sum/denom;
1352 //moldyn->gv_sum+=moldyn->gv;
1353 //moldyn->gv_avg=moldyn->gv_sum/denom;
1356 moldyn->p_sum+=moldyn->p;
1357 moldyn->p_avg=moldyn->p_sum/denom;
1358 //moldyn->gp_sum+=moldyn->gp;
1359 //moldyn->gp_avg=moldyn->gp_sum/denom;
1360 moldyn->tp_sum+=moldyn->tp;
1361 moldyn->tp_avg=moldyn->tp_sum/denom;
1366 int get_heat_capacity(t_moldyn *moldyn) {
1370 /* averages needed for heat capacity calc */
1371 if(moldyn->total_steps<moldyn->avg_skip)
1374 /* (temperature average)^2 */
1375 temp2=moldyn->t_avg*moldyn->t_avg;
1376 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1379 /* ideal gas contribution */
1380 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1381 printf(" ideal gas contribution: %f\n",
1382 ighc/moldyn->mass*KILOGRAM/JOULE);
1384 /* specific heat for nvt ensemble */
1385 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1386 moldyn->c_v_nvt/=moldyn->mass;
1388 /* specific heat for nve ensemble */
1389 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1390 moldyn->c_v_nve/=moldyn->mass;
1392 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1393 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1394 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)));
1399 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1415 /* store atomic configuration + dimension */
1416 store=malloc(moldyn->count*sizeof(t_atom));
1418 printf("[moldyn] allocating store mem failed\n");
1421 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1426 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1427 su=pow(2.0-h,ONE_THIRD)-1.0;
1428 dv=(1.0-h)*moldyn->volume;
1430 /* scale up dimension and atom positions */
1431 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1432 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1433 link_cell_shutdown(moldyn);
1434 link_cell_init(moldyn,QUIET);
1435 potential_force_calc(moldyn);
1438 /* restore atomic configuration + dim */
1439 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1442 /* scale down dimension and atom positions */
1443 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1444 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1445 link_cell_shutdown(moldyn);
1446 link_cell_init(moldyn,QUIET);
1447 potential_force_calc(moldyn);
1450 /* calculate pressure */
1451 moldyn->tp=-(y1-y0)/(2.0*dv);
1453 /* restore atomic configuration */
1454 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1456 link_cell_shutdown(moldyn);
1457 link_cell_init(moldyn,QUIET);
1458 //potential_force_calc(moldyn);
1460 /* free store buffer */
1467 double get_pressure(t_moldyn *moldyn) {
1473 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1485 if(x) dim->x*=scale;
1486 if(y) dim->y*=scale;
1487 if(z) dim->z*=scale;
1492 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1503 for(i=0;i<moldyn->count;i++) {
1504 r=&(moldyn->atom[i].r);
1513 int scale_atoms_ind(t_moldyn *moldyn,double x,double y,double z) {
1518 for(i=0;i<moldyn->count;i++) {
1519 r=&(moldyn->atom[i].r);
1528 int scale_dim_ind(t_moldyn *moldyn,double x,double y,double z) {
1541 int scale_volume(t_moldyn *moldyn) {
1548 vdim=&(moldyn->vis.dim);
1552 /* scaling factor */
1553 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1554 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1555 scale=pow(scale,ONE_THIRD);
1558 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1563 sx=1.0-(moldyn->p_ref-moldyn->px)*moldyn->p_tc*moldyn->tau;
1564 sy=1.0-(moldyn->p_ref-moldyn->py)*moldyn->p_tc*moldyn->tau;
1565 sz=1.0-(moldyn->p_ref-moldyn->pz)*moldyn->p_tc*moldyn->tau;
1566 sx=pow(sx,ONE_THIRD);
1567 sy=pow(sy,ONE_THIRD);
1568 sz=pow(sz,ONE_THIRD);
1571 /* scale the atoms and dimensions */
1572 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1573 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1574 //scale_atoms_ind(moldyn,sx,sy,sz);
1575 //scale_dim_ind(moldyn,sx,sy,sz);
1577 /* visualize dimensions */
1584 /* recalculate scaled volume */
1585 moldyn->volume=dim->x*dim->y*dim->z;
1587 /* adjust/reinit linkcell */
1588 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1589 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1590 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1591 link_cell_shutdown(moldyn);
1592 link_cell_init(moldyn,QUIET);
1606 double e_kin_calc(t_moldyn *moldyn) {
1613 //moldyn->ekinx=0.0;
1614 //moldyn->ekiny=0.0;
1615 //moldyn->ekinz=0.0;
1617 for(i=0;i<moldyn->count;i++) {
1618 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1619 moldyn->ekin+=atom[i].ekin;
1620 //moldyn->ekinx+=0.5*atom[i].mass*atom[i].v.x*atom[i].v.x;
1621 //moldyn->ekiny+=0.5*atom[i].mass*atom[i].v.y*atom[i].v.y;
1622 //moldyn->ekinz+=0.5*atom[i].mass*atom[i].v.z*atom[i].v.z;
1625 return moldyn->ekin;
1628 double get_total_energy(t_moldyn *moldyn) {
1630 return(moldyn->ekin+moldyn->energy);
1633 t_3dvec get_total_p(t_moldyn *moldyn) {
1642 for(i=0;i<moldyn->count;i++) {
1643 v3_scale(&p,&(atom[i].v),atom[i].mass);
1644 v3_add(&p_total,&p_total,&p);
1650 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1654 /* nn_dist is the nearest neighbour distance */
1656 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1665 /* linked list / cell method */
1667 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1670 #ifndef LOWMEM_LISTS
1676 /* partitioning the md cell */
1677 lc->nx=moldyn->dim.x/moldyn->cutoff;
1678 lc->x=moldyn->dim.x/lc->nx;
1679 lc->ny=moldyn->dim.y/moldyn->cutoff;
1680 lc->y=moldyn->dim.y/lc->ny;
1681 lc->nz=moldyn->dim.z/moldyn->cutoff;
1682 lc->z=moldyn->dim.z/lc->nz;
1683 lc->cells=lc->nx*lc->ny*lc->nz;
1686 lc->subcell=malloc(lc->cells*sizeof(int*));
1688 lc->subcell=malloc(sizeof(t_lowmem_list));
1690 lc->subcell=malloc(lc->cells*sizeof(t_list));
1693 if(lc->subcell==NULL) {
1694 perror("[moldyn] cell init (malloc)");
1699 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1704 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1707 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1710 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1713 printf(" x: %d x %f A\n",lc->nx,lc->x);
1714 printf(" y: %d x %f A\n",lc->ny,lc->y);
1715 printf(" z: %d x %f A\n",lc->nz,lc->z);
1720 for(i=0;i<lc->cells;i++) {
1721 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1722 if(lc->subcell[i]==NULL) {
1723 perror("[moldyn] list init (malloc)");
1728 printf(" ---> %d malloc %p (%p)\n",
1729 i,lc->subcell[0],lc->subcell);
1733 lc->subcell->head=malloc(lc->cells*sizeof(int));
1734 if(lc->subcell->head==NULL) {
1735 perror("[moldyn] head init (malloc)");
1738 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1739 if(lc->subcell->list==NULL) {
1740 perror("[moldyn] list init (malloc)");
1744 for(i=0;i<lc->cells;i++)
1745 list_init_f(&(lc->subcell[i]));
1748 /* update the list */
1749 link_cell_update(moldyn);
1754 int link_cell_update(t_moldyn *moldyn) {
1772 for(i=0;i<lc->cells;i++)
1774 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1776 lc->subcell->head[i]=-1;
1778 list_destroy_f(&(lc->subcell[i]));
1781 for(count=0;count<moldyn->count;count++) {
1782 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1783 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1784 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1788 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1791 if(p>=MAX_ATOMS_PER_LIST) {
1792 printf("[moldyn] FATAL: amount of atoms too high!\n");
1796 lc->subcell[i+j*nx+k*nxy][p]=count;
1799 lc->subcell->list[count]=lc->subcell->head[p];
1800 lc->subcell->head[p]=count;
1802 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1806 printf(" ---> %d %d malloc %p (%p)\n",
1807 i,count,lc->subcell[i].current,lc->subcell);
1815 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1841 if(i>=nx||j>=ny||k>=nz)
1842 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1845 #ifndef LOWMEM_LISTS
1846 cell[0]=lc->subcell[i+j*nx+k*a];
1848 cell[0]=lc->subcell->head[i+j*nx+k*a];
1850 for(ci=-1;ci<=1;ci++) {
1853 if((x<0)||(x>=nx)) {
1857 for(cj=-1;cj<=1;cj++) {
1860 if((y<0)||(y>=ny)) {
1864 for(ck=-1;ck<=1;ck++) {
1867 if((z<0)||(z>=nz)) {
1871 if(!(ci|cj|ck)) continue;
1873 #ifndef LOWMEM_LISTS
1874 cell[--count2]=lc->subcell[x+y*nx+z*a];
1876 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1881 #ifndef LOWMEM_LISTS
1882 cell[count1++]=lc->subcell[x+y*nx+z*a];
1884 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1896 int link_cell_shutdown(t_moldyn *moldyn) {
1898 #ifndef LOWMEM_LISTS
1906 free(lc->subcell->head);
1907 free(lc->subcell->list);
1911 for(i=0;i<lc->cells;i++) {
1913 free(lc->subcell[i]);
1915 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1916 list_destroy_f(&(lc->subcell[i]));
1926 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1930 t_moldyn_schedule *schedule;
1932 schedule=&(moldyn->schedule);
1933 count=++(schedule->total_sched);
1935 ptr=realloc(schedule->runs,count*sizeof(int));
1937 perror("[moldyn] realloc (runs)");
1941 schedule->runs[count-1]=runs;
1943 ptr=realloc(schedule->tau,count*sizeof(double));
1945 perror("[moldyn] realloc (tau)");
1949 schedule->tau[count-1]=tau;
1951 printf("[moldyn] schedule added:\n");
1952 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1958 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1960 moldyn->schedule.hook=hook;
1961 moldyn->schedule.hook_params=hook_params;
1968 * 'integration of newtons equation' - algorithms
1972 /* start the integration */
1974 int moldyn_integrate(t_moldyn *moldyn) {
1977 unsigned int e,m,s,v,p,t,a;
1979 t_moldyn_schedule *sched;
1984 double energy_scale;
1985 struct timeval t1,t2;
1988 #ifdef VISUAL_THREAD
1990 pthread_t io_thread;
1999 sched=&(moldyn->schedule);
2002 /* initialize linked cell method */
2003 link_cell_init(moldyn,VERBOSE);
2005 /* logging & visualization */
2014 /* sqaure of some variables */
2015 moldyn->tau_square=moldyn->tau*moldyn->tau;
2017 /* get current time */
2018 gettimeofday(&t1,NULL);
2020 /* calculate initial forces */
2021 potential_force_calc(moldyn);
2026 /* some stupid checks before we actually start calculating bullshit */
2027 if(moldyn->cutoff>0.5*moldyn->dim.x)
2028 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
2029 if(moldyn->cutoff>0.5*moldyn->dim.y)
2030 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
2031 if(moldyn->cutoff>0.5*moldyn->dim.z)
2032 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
2034 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
2035 if(ds>0.05*moldyn->nnd)
2036 printf("[moldyn] WARNING: forces too high / tau too small!\n");
2039 /* zero absolute time */
2040 // should have right values!
2042 //moldyn->total_steps=0;
2044 /* debugging, ignore */
2047 /* zero & init moldyn copy */
2048 #ifdef VISUAL_THREAD
2049 memset(&md_copy,0,sizeof(t_moldyn));
2050 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2051 if(atom_copy==NULL) {
2052 perror("[moldyn] malloc atom copy (init)");
2058 printf("##################\n");
2059 printf("# USING PTHREADS #\n");
2060 printf("##################\n");
2062 /* tell the world */
2063 printf("[moldyn] integration start, go get a coffee ...\n");
2065 /* executing the schedule */
2067 while(sched->count<sched->total_sched) {
2069 /* setting amount of runs and finite time step size */
2070 moldyn->tau=sched->tau[sched->count];
2071 moldyn->tau_square=moldyn->tau*moldyn->tau;
2072 moldyn->time_steps=sched->runs[sched->count];
2074 /* energy scaling factor (might change!) */
2075 energy_scale=moldyn->count*EV;
2077 /* integration according to schedule */
2079 for(i=0;i<moldyn->time_steps;i++) {
2081 /* integration step */
2082 moldyn->integrate(moldyn);
2084 /* calculate kinetic energy, temperature and pressure */
2086 temperature_calc(moldyn);
2088 pressure_calc(moldyn);
2090 thermodynamic_pressure_calc(moldyn);
2091 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
2095 /* calculate fluctuations + averages */
2096 average_and_fluctuation_calc(moldyn);
2099 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
2100 scale_velocity(moldyn,FALSE);
2101 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
2102 scale_volume(moldyn);
2104 /* check for log & visualization */
2106 if(!(moldyn->total_steps%e))
2107 dprintf(moldyn->efd,
2108 "%f %f %f %f %f %f\n",
2109 moldyn->time,moldyn->ekin/energy_scale,
2110 moldyn->energy/energy_scale,
2111 get_total_energy(moldyn)/energy_scale,
2112 moldyn->ekin/EV,moldyn->energy/EV);
2115 if(!(moldyn->total_steps%m)) {
2116 momentum=get_total_p(moldyn);
2117 dprintf(moldyn->mfd,
2118 "%f %f %f %f %f\n",moldyn->time,
2119 momentum.x,momentum.y,momentum.z,
2120 v3_norm(&momentum));
2124 if(!(moldyn->total_steps%p)) {
2125 dprintf(moldyn->pfd,
2126 "%f %f %f %f %f %f %f\n",moldyn->time,
2127 moldyn->p/BAR,moldyn->p_avg/BAR,
2128 moldyn->p/BAR,moldyn->p_avg/BAR,
2129 moldyn->tp/BAR,moldyn->tp_avg/BAR);
2133 if(!(moldyn->total_steps%t)) {
2134 dprintf(moldyn->tfd,
2136 moldyn->time,moldyn->t,moldyn->t_avg);
2140 if(!(moldyn->total_steps%v)) {
2141 dprintf(moldyn->vfd,
2142 "%f %f %f %f %f\n",moldyn->time,
2150 if(!(moldyn->total_steps%s)) {
2151 snprintf(dir,128,"%s/s-%08.f.save",
2152 moldyn->vlsdir,moldyn->time);
2153 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
2155 if(fd<0) perror("[moldyn] save fd open");
2157 write(fd,moldyn,sizeof(t_moldyn));
2158 write(fd,moldyn->atom,
2159 moldyn->count*sizeof(t_atom));
2165 if(!(moldyn->total_steps%a)) {
2166 #ifdef VISUAL_THREAD
2167 /* check whether thread has not terminated yet */
2169 ret=pthread_join(io_thread,NULL);
2172 /* prepare and start thread */
2173 if(moldyn->count!=md_copy.count) {
2177 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
2179 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2180 if(atom_copy==NULL) {
2181 perror("[moldyn] malloc atom copy (change)");
2185 md_copy.atom=atom_copy;
2186 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
2188 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
2190 perror("[moldyn] create visual atoms thread\n");
2194 visual_atoms(moldyn);
2199 /* display progress */
2203 /* get current time */
2204 gettimeofday(&t2,NULL);
2206 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
2207 sched->count,i,moldyn->total_steps,
2208 moldyn->t,moldyn->t_avg,
2210 moldyn->p/BAR,moldyn->p_avg/BAR,
2212 moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
2215 (int)(t2.tv_sec-t1.tv_sec));
2219 /* copy over time */
2225 /* increase absolute time */
2226 moldyn->time+=moldyn->tau;
2227 moldyn->total_steps+=1;
2231 /* check for hooks */
2233 printf("\n ## schedule hook %d start ##\n",
2235 sched->hook(moldyn,sched->hook_params);
2236 printf(" ## schedule hook end ##\n");
2239 /* increase the schedule counter */
2244 /* writing a final save file! */
2246 snprintf(dir,128,"%s/s-final.save",moldyn->vlsdir);
2247 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,S_IRUSR|S_IWUSR);
2248 if(fd<0) perror("[moldyn] save fd open");
2250 write(fd,moldyn,sizeof(t_moldyn));
2251 write(fd,moldyn->atom,
2252 moldyn->count*sizeof(t_atom));
2256 /* writing a final visual file! */
2258 visual_atoms(moldyn);
2263 /* velocity verlet */
2265 int velocity_verlet(t_moldyn *moldyn) {
2268 double tau,tau_square,h;
2271 #ifdef CONSTRAINT_11X_5832
2277 count=moldyn->count;
2279 tau_square=moldyn->tau_square;
2281 #ifdef CONSTRAINT_110_5832
2283 atom[5832].f.x=0.5*(atom[5832].f.x-atom[5832].f.y);
2284 atom[5832].f.y=-atom[5832].f.x;
2287 #ifdef CONSTRAINT_11X_5832
2291 yt=atom[5832].f.y*cos(-0.16935129)-atom[5832].f.z*sin(-0.16935129);
2292 zt=atom[5832].f.y*sin(-0.16935129)+atom[5832].f.z*cos(-0.16935129);
2294 xtt=xt*cos(M_PI/4.0)-yt*sin(M_PI/4.0);
2295 ytt=xt*sin(M_PI/4.0)+yt*sin(M_PI/4.0);
2297 // apply constraints
2299 // first trafo backwards
2300 xt=xtt*cos(M_PI/4.0)+ytt*sin(M_PI/4.0);
2301 yt=-xtt*sin(M_PI/4.0)+ytt*sin(M_PI/4.0);
2303 // second trafo backwards
2305 atom[5832].f.y=yt*cos(-0.16935129)+zt*sin(-0.16935129);
2306 atom[5832].f.z=-yt*sin(-0.16935129)+zt*cos(-0.16935129);
2309 for(i=0;i<count;i++) {
2310 /* check whether fixed atom */
2311 if(atom[i].attr&ATOM_ATTR_FP)
2315 v3_scale(&delta,&(atom[i].v),tau);
2316 #ifdef CONSTRAINT_110_5832
2322 #ifdef CONSTRAINT_11X_5832
2326 yt=delta.y*cos(-0.16935129)-delta.z*sin(-0.16935129);
2327 zt=delta.y*sin(-0.16935129)+delta.z*cos(-0.16935129);
2329 xtt=xt*cos(M_PI/4.0)-yt*sin(M_PI/4.0);
2330 ytt=xt*sin(M_PI/4.0)+yt*sin(M_PI/4.0);
2332 // apply constraints
2334 // first trafo backwards
2335 xt=xtt*cos(M_PI/4.0)+ytt*sin(M_PI/4.0);
2336 yt=-xtt*sin(M_PI/4.0)+ytt*sin(M_PI/4.0);
2338 // second trafo backwards
2340 delta.y=yt*cos(-0.16935129)+zt*sin(-0.16935129);
2341 delta.z=-yt*sin(-0.16935129)+zt*cos(-0.16935129);
2346 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2348 v3_scale(&delta,&(atom[i].f),h*tau_square);
2349 #ifdef CONSTRAINT_110_5832
2355 #ifdef CONSTRAINT_11X_5832
2359 yt=delta.y*cos(-0.16935129)-delta.z*sin(-0.16935129);
2360 zt=delta.y*sin(-0.16935129)+delta.z*cos(-0.16935129);
2362 xtt=xt*cos(M_PI/4.0)-yt*sin(M_PI/4.0);
2363 ytt=xt*sin(M_PI/4.0)+yt*sin(M_PI/4.0);
2365 // apply constraints
2367 // first trafo backwards
2368 xt=xtt*cos(M_PI/4.0)+ytt*sin(M_PI/4.0);
2369 yt=-xtt*sin(M_PI/4.0)+ytt*sin(M_PI/4.0);
2371 // second trafo backwards
2373 delta.y=yt*cos(-0.16935129)+zt*sin(-0.16935129);
2374 delta.z=-yt*sin(-0.16935129)+zt*cos(-0.16935129);
2378 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2379 //check_per_bound_and_care_for_pbc(moldyn,&(atom[i]));
2380 check_per_bound(moldyn,&(atom[i].r));
2382 /* velocities [actually v(t+tau/2)] */
2383 v3_scale(&delta,&(atom[i].f),h*tau);
2384 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2387 /* criticial check */
2388 moldyn_bc_check(moldyn);
2390 /* neighbour list update */
2391 link_cell_update(moldyn);
2393 /* forces depending on chosen potential */
2395 // if albe, use albe force calc routine
2396 //if(moldyn->func3b_j1==albe_mult_3bp_j1)
2397 // albe_potential_force_calc(moldyn);
2399 potential_force_calc(moldyn);
2401 albe_potential_force_calc(moldyn);
2404 #ifdef CONSTRAINT_110_5832
2406 atom[5832].f.x=0.5*(atom[5832].f.x-atom[5832].f.y);
2407 atom[5832].f.y=-atom[5832].f.x;
2410 #ifdef CONSTRAINT_11X_5832
2414 yt=atom[5832].f.y*cos(-0.16935129)-atom[5832].f.z*sin(-0.16935129);
2415 zt=atom[5832].f.y*sin(-0.16935129)+atom[5832].f.z*cos(-0.16935129);
2417 xtt=xt*cos(M_PI/4.0)-yt*sin(M_PI/4.0);
2418 ytt=xt*sin(M_PI/4.0)+yt*sin(M_PI/4.0);
2420 // apply constraints
2422 // first trafo backwards
2423 xt=xtt*cos(M_PI/4.0)+ytt*sin(M_PI/4.0);
2424 yt=-xtt*sin(M_PI/4.0)+ytt*sin(M_PI/4.0);
2426 // second trafo backwards
2428 atom[5832].f.y=yt*cos(-0.16935129)+zt*sin(-0.16935129);
2429 atom[5832].f.z=-yt*sin(-0.16935129)+zt*cos(-0.16935129);
2432 for(i=0;i<count;i++) {
2433 /* check whether fixed atom */
2434 if(atom[i].attr&ATOM_ATTR_FP)
2436 /* again velocities [actually v(t+tau)] */
2437 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2438 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2447 * potentials & corresponding forces & virial routine
2451 /* generic potential and force calculation */
2453 int potential_force_calc(t_moldyn *moldyn) {
2456 t_atom *itom,*jtom,*ktom;
2460 int *neighbour_i[27];
2464 int neighbour_i[27];
2467 t_list neighbour_i[27];
2468 t_list neighbour_i2[27];
2474 count=moldyn->count;
2484 /* reset global virial */
2485 memset(&(moldyn->gvir),0,sizeof(t_virial));
2487 /* reset force, site energy and virial of every atom */
2489 i=omp_get_thread_num();
2490 #pragma omp parallel for private(virial)
2492 for(i=0;i<count;i++) {
2495 v3_zero(&(itom[i].f));
2498 virial=(&(itom[i].virial));
2506 /* reset site energy */
2511 /* get energy, force and virial of every atom */
2513 /* first (and only) loop over atoms i */
2514 for(i=0;i<count;i++) {
2516 /* single particle potential/force */
2517 if(itom[i].attr&ATOM_ATTR_1BP)
2519 moldyn->func1b(moldyn,&(itom[i]));
2521 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2524 /* 2 body pair potential/force */
2526 link_cell_neighbour_index(moldyn,
2527 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2528 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2529 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2534 /* first loop over atoms j */
2535 if(moldyn->func2b) {
2542 while(neighbour_i[j][p]!=-1) {
2544 jtom=&(atom[neighbour_i[j][p]]);
2552 p=lc->subcell->list[p];
2554 this=&(neighbour_i[j]);
2557 if(this->start==NULL)
2561 jtom=this->current->data;
2564 if(jtom==&(itom[i]))
2567 if((jtom->attr&ATOM_ATTR_2BP)&
2568 (itom[i].attr&ATOM_ATTR_2BP)) {
2569 moldyn->func2b(moldyn,
2579 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2585 /* 3 body potential/force */
2587 if(!(itom[i].attr&ATOM_ATTR_3BP))
2590 /* copy the neighbour lists */
2592 /* no copy needed for static lists */
2594 /* no copy needed for lowmem lists */
2596 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2599 /* second loop over atoms j */
2606 while(neighbour_i[j][p]!=-1) {
2608 jtom=&(atom[neighbour_i[j][p]]);
2616 p=lc->subcell->list[p];
2618 this=&(neighbour_i[j]);
2621 if(this->start==NULL)
2626 jtom=this->current->data;
2629 if(jtom==&(itom[i]))
2632 if(!(jtom->attr&ATOM_ATTR_3BP))
2638 if(moldyn->func3b_j1)
2639 moldyn->func3b_j1(moldyn,
2644 /* in first j loop, 3bp run can be skipped */
2645 if(!(moldyn->run3bp))
2648 /* first loop over atoms k */
2649 if(moldyn->func3b_k1) {
2657 while(neighbour_i[k][q]!=-1) {
2659 ktom=&(atom[neighbour_i[k][q]]);
2667 q=lc->subcell->list[q];
2669 that=&(neighbour_i2[k]);
2672 if(that->start==NULL)
2676 ktom=that->current->data;
2679 if(!(ktom->attr&ATOM_ATTR_3BP))
2685 if(ktom==&(itom[i]))
2688 moldyn->func3b_k1(moldyn,
2699 } while(list_next_f(that)!=\
2707 if(moldyn->func3b_j2)
2708 moldyn->func3b_j2(moldyn,
2713 /* second loop over atoms k */
2714 if(moldyn->func3b_k2) {
2722 while(neighbour_i[k][q]!=-1) {
2724 ktom=&(atom[neighbour_i[k][q]]);
2732 q=lc->subcell->list[q];
2734 that=&(neighbour_i2[k]);
2737 if(that->start==NULL)
2741 ktom=that->current->data;
2744 if(!(ktom->attr&ATOM_ATTR_3BP))
2750 if(ktom==&(itom[i]))
2753 moldyn->func3b_k2(moldyn,
2764 } while(list_next_f(that)!=\
2772 /* 2bp post function */
2773 if(moldyn->func3b_j3) {
2774 moldyn->func3b_j3(moldyn,
2783 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2798 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2799 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2801 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2802 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2803 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2807 /* some postprocessing */
2809 #pragma omp parallel for
2811 for(i=0;i<count;i++) {
2812 /* calculate global virial */
2813 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2814 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2815 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2816 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2817 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2818 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2820 /* check forces regarding the given timestep */
2821 if(v3_norm(&(itom[i].f))>\
2822 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2823 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2831 * virial calculation
2834 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2835 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2837 a->virial.xx+=f->x*d->x;
2838 a->virial.yy+=f->y*d->y;
2839 a->virial.zz+=f->z*d->z;
2840 a->virial.xy+=f->x*d->y;
2841 a->virial.xz+=f->x*d->z;
2842 a->virial.yz+=f->y*d->z;
2848 * periodic boundary checking
2851 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2852 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2863 if(moldyn->status&MOLDYN_STAT_PBX) {
2864 if(a->x>=x) a->x-=dim->x;
2865 else if(-a->x>x) a->x+=dim->x;
2867 if(moldyn->status&MOLDYN_STAT_PBY) {
2868 if(a->y>=y) a->y-=dim->y;
2869 else if(-a->y>y) a->y+=dim->y;
2871 if(moldyn->status&MOLDYN_STAT_PBZ) {
2872 if(a->z>=z) a->z-=dim->z;
2873 else if(-a->z>z) a->z+=dim->z;
2879 int check_per_bound_and_care_for_pbc(t_moldyn *moldyn,t_atom *a) {
2890 if(moldyn->status&MOLDYN_STAT_PBX) {
2895 else if(-a->r.x>x) {
2900 if(moldyn->status&MOLDYN_STAT_PBY) {
2905 else if(-a->r.y>y) {
2910 if(moldyn->status&MOLDYN_STAT_PBZ) {
2915 else if(-a->r.z>z) {
2925 * debugging / critical check functions
2928 int moldyn_bc_check(t_moldyn *moldyn) {
2941 for(i=0;i<moldyn->count;i++) {
2942 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2943 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2944 i,atom[i].r.x,dim->x/2);
2945 printf("diagnostic:\n");
2946 printf("-----------\natom.r.x:\n");
2948 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2951 ((byte)&(1<<k))?1:0,
2954 printf("---------------\nx=dim.x/2:\n");
2956 memcpy(&byte,(u8 *)(&x)+j,1);
2959 ((byte)&(1<<k))?1:0,
2962 if(atom[i].r.x==x) printf("the same!\n");
2963 else printf("different!\n");
2965 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2966 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2967 i,atom[i].r.y,dim->y/2);
2968 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2969 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2970 i,atom[i].r.z,dim->z/2);
2980 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2987 fd=open(file,O_RDONLY);
2989 perror("[moldyn] load save file open");
2993 fsize=lseek(fd,0,SEEK_END);
2994 lseek(fd,0,SEEK_SET);
2996 size=sizeof(t_moldyn);
2999 cnt=read(fd,moldyn,size);
3001 perror("[moldyn] load save file read (moldyn)");
3007 size=moldyn->count*sizeof(t_atom);
3009 /* correcting possible atom data offset */
3011 if(fsize!=sizeof(t_moldyn)+size) {
3012 corr=fsize-sizeof(t_moldyn)-size;
3013 printf("[moldyn] WARNING: lsf (illegal file size)\n");
3014 printf(" modifying offset:\n");
3015 printf(" - current pos: %d\n",sizeof(t_moldyn));
3016 printf(" - atom size: %d\n",size);
3017 printf(" - file size: %d\n",fsize);
3018 printf(" => correction: %d\n",corr);
3019 lseek(fd,corr,SEEK_CUR);
3022 moldyn->atom=(t_atom *)malloc(size);
3023 if(moldyn->atom==NULL) {
3024 perror("[moldyn] load save file malloc (atoms)");
3029 cnt=read(fd,moldyn->atom,size);
3031 perror("[moldyn] load save file read (atoms)");
3038 amutex=malloc(moldyn->count*sizeof(pthread_mutex_t));
3040 perror("[moldyn] load save file (mutexes)");
3043 for(cnt=0;cnt<moldyn->count;cnt++)
3044 pthread_mutex_init(&(amutex[cnt]),NULL);
3052 int moldyn_free_save_file(t_moldyn *moldyn) {
3059 int moldyn_load(t_moldyn *moldyn) {
3067 * function to find/callback all combinations of 2 body bonds
3070 int process_2b_bonds(t_moldyn *moldyn,void *data,
3071 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3072 void *data,u8 bc)) {
3082 t_list neighbour[27];
3092 for(i=0;i<moldyn->count;i++) {
3093 /* neighbour indexing */
3094 link_cell_neighbour_index(moldyn,
3095 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
3096 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
3097 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
3102 bc=(j<lc->dnlc)?0:1;
3107 while(neighbour[j][p]!=-1) {
3109 jtom=&(moldyn->atom[neighbour[j][p]]);
3117 p=lc->subcell->list[p];
3119 this=&(neighbour[j]);
3122 if(this->start==NULL)
3127 jtom=this->current->data;
3131 process(moldyn,&(itom[i]),jtom,data,bc);
3138 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3148 * function to find neighboured atoms
3151 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
3152 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
3153 void *data,u8 bc)) {
3163 t_list neighbour[27];
3172 /* neighbour indexing */
3173 link_cell_neighbour_index(moldyn,
3174 (atom->r.x+moldyn->dim.x/2)/lc->x,
3175 (atom->r.y+moldyn->dim.y/2)/lc->x,
3176 (atom->r.z+moldyn->dim.z/2)/lc->x,
3181 bc=(j<lc->dnlc)?0:1;
3186 while(neighbour[j][p]!=-1) {
3188 natom=&(moldyn->atom[neighbour[j][p]]);
3195 natom=&(moldyn->atom[p]);
3196 p=lc->subcell->list[p];
3198 this=&(neighbour[j]);
3201 if(this->start==NULL)
3206 natom=this->current->data;
3210 process(moldyn,atom,natom,data,bc);
3217 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3226 * post processing functions
3229 int get_line(int fd,char *line,int max) {
3236 if(count==max) return count;
3237 ret=read(fd,line+count,1);
3238 if(ret<=0) return ret;
3239 if(line[count]=='\n') {
3240 memset(line+count,0,max-count-1);
3248 int pair_correlation_init(t_moldyn *moldyn,double dr) {
3254 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
3271 for(i=0;i<moldyn->count;i++) {
3273 /* care for pb crossing */
3274 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3275 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3276 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3277 /* calculate distance */
3278 v3_sub(&dist,&final_r,&(atom[i].r_0));
3279 d2=v3_absolute_square(&dist);
3293 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
3294 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
3295 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
3300 int calculate_msd(t_moldyn *moldyn,double *msd) {
3317 for(i=0;i<moldyn->count;i++) {
3319 /* care for pb crossing */
3320 if(atom[i].pbc[0]|atom[i].pbc[1]|atom[i].pbc[2]) {
3321 printf("[moldyn] msd pb crossings for atom %d\n",i);
3322 printf(" x: %d y: %d z: %d\n",
3323 atom[i].pbc[0],atom[i].pbc[1],atom[i].pbc[2]);
3325 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3326 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3327 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3328 /* calculate distance */
3329 v3_sub(&dist,&final_r,&(atom[i].r_0));
3330 d2=v3_absolute_square(&dist);
3346 msd[2]/=moldyn->count;
3351 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
3356 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
3357 t_atom *jtom,void *data,u8 bc) {
3364 /* only count pairs once,
3365 * skip same atoms */
3366 if(itom->tag>=jtom->tag)
3370 * pair correlation calc
3377 v3_sub(&dist,&(jtom->r),&(itom->r));
3378 if(bc) check_per_bound(moldyn,&dist);
3379 d=v3_absolute_square(&dist);
3381 /* ignore if greater cutoff */
3382 if(d>moldyn->cutoff_square)
3385 /* fill the slots */
3389 /* should never happen but it does 8) -
3390 * related to -ffloat-store problem! */
3392 printf("[moldyn] WARNING: pcc (%d/%d)",
3398 if(itom->brand!=jtom->brand) {
3403 /* type a - type a bonds */
3405 pcc->stat[s+pcc->o1]+=1;
3407 /* type b - type b bonds */
3408 pcc->stat[s+pcc->o2]+=1;
3414 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
3421 pcc.o1=moldyn->cutoff/dr;
3424 if(pcc.o1*dr<=moldyn->cutoff)
3425 printf("[moldyn] WARNING: pcc (low #slots)\n");
3427 printf("[moldyn] pair correlation calc info:\n");
3428 printf(" time: %f\n",moldyn->time);
3429 printf(" count: %d\n",moldyn->count);
3430 printf(" cutoff: %f\n",moldyn->cutoff);
3431 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
3434 pcc.stat=(double *)ptr;
3437 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
3438 if(pcc.stat==NULL) {
3439 perror("[moldyn] pair correlation malloc");
3444 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
3447 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
3450 for(i=1;i<pcc.o1;i++) {
3451 // normalization: 4 pi r^2 dr
3452 // here: not double counting pairs -> 2 pi r r dr
3453 // ... and actually it's a constant times r^2
3456 pcc.stat[pcc.o1+i]/=norm;
3457 pcc.stat[pcc.o2+i]/=norm;
3462 /* todo: store/print pair correlation function */
3469 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3476 if(itom->tag>=jtom->tag)
3480 v3_sub(&dist,&(jtom->r),&(itom->r));
3481 if(bc) check_per_bound(moldyn,&dist);
3482 d=v3_absolute_square(&dist);
3484 /* ignore if greater or equal cutoff */
3485 if(d>moldyn->cutoff_square)
3488 /* check for potential bond */
3489 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3492 /* now count this bonding ... */
3495 /* increase total bond counter
3500 ba->acnt[jtom->tag]+=1;
3502 ba->bcnt[jtom->tag]+=1;
3505 ba->acnt[itom->tag]+=1;
3507 ba->bcnt[itom->tag]+=1;
3512 int bond_analyze(t_moldyn *moldyn,double *quality) {
3523 ba.acnt=malloc(moldyn->count*sizeof(int));
3525 perror("[moldyn] bond analyze malloc (a)");
3528 memset(ba.acnt,0,moldyn->count*sizeof(int));
3530 ba.bcnt=malloc(moldyn->count*sizeof(int));
3532 perror("[moldyn] bond analyze malloc (b)");
3535 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3544 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3546 for(i=0;i<moldyn->count;i++) {
3547 if(atom[i].brand==0) {
3548 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3550 if((ba.acnt[i]==0)&(ba.bcnt[i]==3))
3554 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3558 if((ba.acnt[i]==3)&(ba.bcnt[i]==0)) {
3567 quality[0]=1.0*ccnt4/bcnt;
3568 quality[1]=1.0*ccnt3/bcnt;
3571 printf("[moldyn] bond analyze: %f %f\n",
3572 1.0*ccnt4/bcnt,1.0*ccnt3/bcnt);
3579 * visualization code
3582 int visual_init(t_moldyn *moldyn,char *filebase) {
3584 strncpy(moldyn->vis.fb,filebase,128);
3589 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3596 if(itom->tag>=jtom->tag)
3599 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3602 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3603 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3604 itom->r.x,itom->r.y,itom->r.z,
3605 jtom->r.x,jtom->r.y,jtom->r.z);
3610 #ifdef VISUAL_THREAD
3611 void *visual_atoms(void *ptr) {
3613 int visual_atoms(t_moldyn *moldyn) {
3624 #ifdef VISUAL_THREAD
3638 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3639 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3641 perror("open visual save file fd");
3642 #ifndef VISUAL_THREAD
3647 /* write the actual data file */
3650 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3651 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3653 // atomic configuration
3654 for(i=0;i<moldyn->count;i++) {
3655 v3_sub(&strain,&(atom[i].r),&(atom[i].r_0));
3656 check_per_bound(moldyn,&strain);
3657 // atom type, positions, color and kinetic energy
3658 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3662 pse_col[atom[i].element],
3664 sqrt(v3_absolute_square(&strain)));
3667 // bonds between atoms
3668 #ifndef VISUAL_THREAD
3669 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3674 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3675 -dim.x/2,-dim.y/2,-dim.z/2,
3676 dim.x/2,-dim.y/2,-dim.z/2);
3677 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3678 -dim.x/2,-dim.y/2,-dim.z/2,
3679 -dim.x/2,dim.y/2,-dim.z/2);
3680 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3681 dim.x/2,dim.y/2,-dim.z/2,
3682 dim.x/2,-dim.y/2,-dim.z/2);
3683 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3684 -dim.x/2,dim.y/2,-dim.z/2,
3685 dim.x/2,dim.y/2,-dim.z/2);
3687 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3688 -dim.x/2,-dim.y/2,dim.z/2,
3689 dim.x/2,-dim.y/2,dim.z/2);
3690 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3691 -dim.x/2,-dim.y/2,dim.z/2,
3692 -dim.x/2,dim.y/2,dim.z/2);
3693 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3694 dim.x/2,dim.y/2,dim.z/2,
3695 dim.x/2,-dim.y/2,dim.z/2);
3696 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3697 -dim.x/2,dim.y/2,dim.z/2,
3698 dim.x/2,dim.y/2,dim.z/2);
3700 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3701 -dim.x/2,-dim.y/2,dim.z/2,
3702 -dim.x/2,-dim.y/2,-dim.z/2);
3703 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3704 -dim.x/2,dim.y/2,dim.z/2,
3705 -dim.x/2,dim.y/2,-dim.z/2);
3706 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3707 dim.x/2,-dim.y/2,dim.z/2,
3708 dim.x/2,-dim.y/2,-dim.z/2);
3709 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3710 dim.x/2,dim.y/2,dim.z/2,
3711 dim.x/2,dim.y/2,-dim.z/2);
3716 #ifdef VISUAL_THREAD
3727 * fpu cntrol functions
3730 // set rounding to double (eliminates -ffloat-store!)
3731 int fpu_set_rtd(void) {
3737 ctrl&=~_FPU_EXTENDED;