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");
88 int moldyn_shutdown(t_moldyn *moldyn) {
94 printf("[moldyn] shutdown\n");
97 for(i=0;i<moldyn->count;i++)
98 pthread_mutex_destroy(&(amutex[i]));
100 pthread_mutex_destroy(&emutex);
103 moldyn_log_shutdown(moldyn);
104 link_cell_shutdown(moldyn);
105 rand_close(&(moldyn->random));
111 int set_int_alg(t_moldyn *moldyn,u8 algo) {
113 printf("[moldyn] integration algorithm: ");
116 case MOLDYN_INTEGRATE_VERLET:
117 moldyn->integrate=velocity_verlet;
118 printf("velocity verlet\n");
121 printf("unknown integration algorithm: %02x\n",algo);
129 int set_cutoff(t_moldyn *moldyn,double cutoff) {
131 moldyn->cutoff=cutoff;
132 moldyn->cutoff_square=cutoff*cutoff;
134 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
139 int set_temperature(t_moldyn *moldyn,double t_ref) {
143 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
148 int set_pressure(t_moldyn *moldyn,double p_ref) {
152 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
157 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
159 moldyn->pt_scale&=(~(P_SCALE_MASK));
160 moldyn->pt_scale|=ptype;
163 printf("[moldyn] p scaling:\n");
165 printf(" p: %s",ptype?"yes":"no ");
167 printf(" | type: %02x | factor: %f",ptype,ptc);
173 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
175 moldyn->pt_scale&=(~(T_SCALE_MASK));
176 moldyn->pt_scale|=ttype;
179 printf("[moldyn] t scaling:\n");
181 printf(" t: %s",ttype?"yes":"no ");
183 printf(" | type: %02x | factor: %f",ttype,ttc);
189 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
191 moldyn->pt_scale=(ptype|ttype);
195 printf("[moldyn] p/t scaling:\n");
197 printf(" p: %s",ptype?"yes":"no ");
199 printf(" | type: %02x | factor: %f",ptype,ptc);
202 printf(" t: %s",ttype?"yes":"no ");
204 printf(" | type: %02x | factor: %f",ttype,ttc);
210 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
216 moldyn->volume=x*y*z;
224 printf("[moldyn] dimensions in A and A^3 respectively:\n");
225 printf(" x: %f\n",moldyn->dim.x);
226 printf(" y: %f\n",moldyn->dim.y);
227 printf(" z: %f\n",moldyn->dim.z);
228 printf(" volume: %f\n",moldyn->volume);
229 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
234 int set_nn_dist(t_moldyn *moldyn,double dist) {
241 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
243 printf("[moldyn] periodic boundary conditions:\n");
246 moldyn->status|=MOLDYN_STAT_PBX;
249 moldyn->status|=MOLDYN_STAT_PBY;
252 moldyn->status|=MOLDYN_STAT_PBZ;
254 printf(" x: %s\n",x?"yes":"no");
255 printf(" y: %s\n",y?"yes":"no");
256 printf(" z: %s\n",z?"yes":"no");
261 int set_potential(t_moldyn *moldyn,u8 type) {
264 case MOLDYN_POTENTIAL_TM:
265 //moldyn->func1b=tersoff_mult_1bp;
266 moldyn->func3b_j1=tersoff_mult_3bp_j1;
267 moldyn->func3b_k1=tersoff_mult_3bp_k1;
268 moldyn->func3b_j2=tersoff_mult_3bp_j2;
269 moldyn->func3b_k2=tersoff_mult_3bp_k2;
270 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
272 case MOLDYN_POTENTIAL_AM:
273 moldyn->func3b_j1=albe_mult_3bp_j1;
274 moldyn->func3b_k1=albe_mult_3bp_k1;
275 moldyn->func3b_j2=albe_mult_3bp_j2;
276 moldyn->func3b_k2=albe_mult_3bp_k2;
277 moldyn->check_2b_bond=albe_mult_check_2b_bond;
279 case MOLDYN_POTENTIAL_HO:
280 moldyn->func2b=harmonic_oscillator;
281 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
283 case MOLDYN_POTENTIAL_LJ:
284 moldyn->func2b=lennard_jones;
285 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
288 printf("[moldyn] set potential: unknown type %02x\n",
296 int set_avg_skip(t_moldyn *moldyn,int skip) {
298 printf("[moldyn] skip %d steps before starting average calc\n",skip);
299 moldyn->avg_skip=skip;
304 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
306 strncpy(moldyn->vlsdir,dir,127);
311 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
313 strncpy(moldyn->rauthor,author,63);
314 strncpy(moldyn->rtitle,title,63);
319 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
324 printf("[moldyn] set log: ");
327 case LOG_TOTAL_ENERGY:
328 moldyn->ewrite=timer;
329 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
330 moldyn->efd=open(filename,
331 O_WRONLY|O_CREAT|O_EXCL,
334 perror("[moldyn] energy log fd open");
337 dprintf(moldyn->efd,"# total energy log file\n");
338 printf("total energy (%d)\n",timer);
340 case LOG_TOTAL_MOMENTUM:
341 moldyn->mwrite=timer;
342 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
343 moldyn->mfd=open(filename,
344 O_WRONLY|O_CREAT|O_EXCL,
347 perror("[moldyn] momentum log fd open");
350 dprintf(moldyn->efd,"# total momentum log file\n");
351 printf("total momentum (%d)\n",timer);
354 moldyn->pwrite=timer;
355 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
356 moldyn->pfd=open(filename,
357 O_WRONLY|O_CREAT|O_EXCL,
360 perror("[moldyn] pressure log file\n");
363 dprintf(moldyn->pfd,"# pressure log file\n");
364 printf("pressure (%d)\n",timer);
366 case LOG_TEMPERATURE:
367 moldyn->twrite=timer;
368 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
369 moldyn->tfd=open(filename,
370 O_WRONLY|O_CREAT|O_EXCL,
373 perror("[moldyn] temperature log file\n");
376 dprintf(moldyn->tfd,"# temperature log file\n");
377 printf("temperature (%d)\n",timer);
380 moldyn->vwrite=timer;
381 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
382 moldyn->vfd=open(filename,
383 O_WRONLY|O_CREAT|O_EXCL,
386 perror("[moldyn] volume log file\n");
389 dprintf(moldyn->vfd,"# volume log file\n");
390 printf("volume (%d)\n",timer);
393 moldyn->swrite=timer;
394 printf("save file (%d)\n",timer);
397 moldyn->awrite=timer;
398 ret=visual_init(moldyn,moldyn->vlsdir);
400 printf("[moldyn] visual init failure\n");
403 printf("visual file (%d)\n",timer);
406 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
407 moldyn->rfd=open(filename,
408 O_WRONLY|O_CREAT|O_EXCL,
411 perror("[moldyn] report fd open");
414 printf("report -> ");
416 snprintf(filename,127,"%s/e_plot.scr",
418 moldyn->epfd=open(filename,
419 O_WRONLY|O_CREAT|O_EXCL,
422 perror("[moldyn] energy plot fd open");
425 dprintf(moldyn->epfd,e_plot_script);
430 snprintf(filename,127,"%s/pressure_plot.scr",
432 moldyn->ppfd=open(filename,
433 O_WRONLY|O_CREAT|O_EXCL,
436 perror("[moldyn] p plot fd open");
439 dprintf(moldyn->ppfd,pressure_plot_script);
444 snprintf(filename,127,"%s/temperature_plot.scr",
446 moldyn->tpfd=open(filename,
447 O_WRONLY|O_CREAT|O_EXCL,
450 perror("[moldyn] t plot fd open");
453 dprintf(moldyn->tpfd,temperature_plot_script);
455 printf("temperature ");
457 dprintf(moldyn->rfd,report_start,
458 moldyn->rauthor,moldyn->rtitle);
462 printf("unknown log type: %02x\n",type);
469 int moldyn_log_shutdown(t_moldyn *moldyn) {
473 printf("[moldyn] log shutdown\n");
477 dprintf(moldyn->rfd,report_energy);
478 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
483 if(moldyn->mfd) close(moldyn->mfd);
487 dprintf(moldyn->rfd,report_pressure);
488 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
495 dprintf(moldyn->rfd,report_temperature);
496 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
501 dprintf(moldyn->rfd,report_end);
503 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
506 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
509 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
518 * creating lattice functions
521 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,
522 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin,
523 t_part_params *p_params,t_defect_params *d_params,
524 t_offset_params *o_params) {
533 pthread_mutex_t *mutex;
539 /* how many atoms do we expect */
542 printf("[moldyn] WARNING: create 'none' lattice called");
544 if(type==CUBIC) new*=1;
545 if(type==FCC) new*=4;
546 if(type==DIAMOND) new*=8;
550 switch(d_params->stype) {
551 case DEFECT_STYPE_DB_X:
552 case DEFECT_STYPE_DB_Y:
553 case DEFECT_STYPE_DB_Z:
554 case DEFECT_STYPE_DB_R:
558 printf("[moldyn] WARNING: cl unknown defect\n");
563 /* allocate space for atoms */
564 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
566 perror("[moldyn] realloc (create lattice)");
570 atom=&(moldyn->atom[count]);
573 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
575 perror("[moldyn] mutex realloc (add atom)");
579 mutex=&(amutex[count]);
582 /* no atoms on the boundaries (only reason: it looks better!) */
597 v3_add(&orig,&orig,&(o_params->o));
598 set_nn_dist(moldyn,lc);
599 ret=cubic_init(a,b,c,lc,atom,&orig,p_params,d_params);
600 strcpy(name,"cubic");
604 v3_scale(&orig,&orig,0.5);
606 v3_add(&orig,&orig,&(o_params->o));
607 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
608 ret=fcc_init(a,b,c,lc,atom,&orig,p_params,d_params);
613 v3_scale(&orig,&orig,0.25);
615 v3_add(&orig,&orig,&(o_params->o));
616 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
617 ret=diamond_init(a,b,c,lc,atom,&orig,p_params,d_params);
618 strcpy(name,"diamond");
621 printf("unknown lattice type (%02x)\n",type);
627 printf("[moldyn] creating %s lattice (lc=%f) incomplete\n",
629 printf(" (ignore for partial lattice creation)\n");
630 printf(" amount of atoms\n");
631 printf(" - expected: %d\n",new);
632 printf(" - created: %d\n",ret);
637 printf("[moldyn] created %s lattice with %d atoms\n",name,ret);
639 for(new=0;new<ret;new++) {
640 atom[new].element=element;
641 atom[new].mass=pse_mass[element];
643 atom[new].brand=brand;
644 atom[new].tag=count+new;
645 check_per_bound(moldyn,&(atom[new].r));
646 atom[new].r_0=atom[new].r;
648 pthread_mutex_init(&(mutex[new]),NULL);
652 atom[new].element=d_params->element;
653 atom[new].mass=pse_mass[d_params->element];
654 atom[new].attr=d_params->attr;
655 atom[new].brand=d_params->brand;
656 atom[new].tag=count+new;
657 check_per_bound(moldyn,&(atom[new].r));
658 atom[new].r_0=atom[new].r;
660 pthread_mutex_init(&(mutex[new]),NULL);
666 ptr=realloc(moldyn->atom,moldyn->count*sizeof(t_atom));
668 perror("[moldyn] realloc (create lattice - alloc fix)");
673 // WHAT ABOUT AMUTEX !!!!
676 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
678 perror("[moldyn] list realloc (create lattice)");
681 moldyn->lc.subcell->list=ptr;
684 /* update total system mass */
685 total_mass_calc(moldyn);
690 int add_atom(t_moldyn *moldyn,int element,u8 brand,u8 attr,
691 t_3dvec *r,t_3dvec *v) {
698 count=(moldyn->count)++; // asshole style!
700 ptr=realloc(atom,(count+1)*sizeof(t_atom));
702 perror("[moldyn] realloc (add atom)");
708 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
710 perror("[moldyn] list realloc (add atom)");
713 moldyn->lc.subcell->list=ptr;
717 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
719 perror("[moldyn] mutex realloc (add atom)");
723 pthread_mutex_init(&(amutex[count]),NULL);
728 /* initialize new atom */
729 memset(&(atom[count]),0,sizeof(t_atom));
732 atom[count].element=element;
733 atom[count].mass=pse_mass[element];
734 atom[count].brand=brand;
735 atom[count].tag=count;
736 atom[count].attr=attr;
737 check_per_bound(moldyn,&(atom[count].r));
738 atom[count].r_0=atom[count].r;
740 /* update total system mass */
741 total_mass_calc(moldyn);
746 int del_atom(t_moldyn *moldyn,int tag) {
750 #if defined LOWMEM_LISTS || defined PTHREADS
756 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
758 perror("[moldyn]malloc (del atom)");
762 for(cnt=0;cnt<tag;cnt++)
765 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
767 new[cnt-1].tag=cnt-1;
776 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
778 perror("[moldyn] list realloc (del atom)");
781 moldyn->lc.subcell->list=ptr;
783 link_cell_update(moldyn);
787 ptr=realloc(amutex,moldyn->count*sizeof(pthread_mutex_t));
789 perror("[moldyn] mutex realloc (add atom)");
793 pthread_mutex_destroy(&(amutex[moldyn->count+1]));
800 #define set_atom_positions(pos) \
801 if(d_params->type) {\
802 d_o.x=0; d_o.y=0; d_o.z=0;\
803 d_d.x=0; d_d.y=0; d_d.z=0;\
804 switch(d_params->stype) {\
805 case DEFECT_STYPE_DB_X:\
809 case DEFECT_STYPE_DB_Y:\
813 case DEFECT_STYPE_DB_Z:\
817 case DEFECT_STYPE_DB_R:\
820 printf("[moldyn] WARNING: unknown defect\n");\
823 v3_add(&dr,&pos,&d_o);\
824 v3_copy(&(atom[count].r),&dr);\
826 v3_add(&dr,&pos,&d_d);\
827 v3_copy(&(atom[count].r),&dr);\
831 v3_copy(&(atom[count].r),&pos);\
836 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
837 t_part_params *p_params,t_defect_params *d_params) {
857 /* shift partition values */
859 p.x=p_params->p.x+(a*lc)/2.0;
860 p.y=p_params->p.y+(b*lc)/2.0;
861 p.z=p_params->p.z+(c*lc)/2.0;
870 switch(p_params->type) {
873 if(v3_absolute_square(&dist)<
874 (p_params->r*p_params->r)) {
875 set_atom_positions(r);
880 if(v3_absolute_square(&dist)>=
881 (p_params->r*p_params->r)) {
882 set_atom_positions(r);
887 if((fabs(dist.x)<p_params->d.x)&&
888 (fabs(dist.y)<p_params->d.y)&&
889 (fabs(dist.z)<p_params->d.z)) {
890 set_atom_positions(r);
895 if((fabs(dist.x)>=p_params->d.x)||
896 (fabs(dist.y)>=p_params->d.y)||
897 (fabs(dist.z)>=p_params->d.z)) {
898 set_atom_positions(r);
902 set_atom_positions(r);
912 for(i=0;i<count;i++) {
913 atom[i].r.x-=(a*lc)/2.0;
914 atom[i].r.y-=(b*lc)/2.0;
915 atom[i].r.z-=(c*lc)/2.0;
921 /* fcc lattice init */
922 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
923 t_part_params *p_params,t_defect_params *d_params) {
941 /* construct the basis */
942 memset(basis,0,3*sizeof(t_3dvec));
950 /* shift partition values */
952 p.x=p_params->p.x+(a*lc)/2.0;
953 p.y=p_params->p.y+(b*lc)/2.0;
954 p.z=p_params->p.z+(c*lc)/2.0;
957 /* fill up the room */
965 switch(p_params->type) {
968 if(v3_absolute_square(&dist)<
969 (p_params->r*p_params->r)) {
970 set_atom_positions(r);
975 if(v3_absolute_square(&dist)>=
976 (p_params->r*p_params->r)) {
977 set_atom_positions(r);
982 if((fabs(dist.x)<p_params->d.x)&&
983 (fabs(dist.y)<p_params->d.y)&&
984 (fabs(dist.z)<p_params->d.z)) {
985 set_atom_positions(r);
990 if((fabs(dist.x)>=p_params->d.x)||
991 (fabs(dist.y)>=p_params->d.y)||
992 (fabs(dist.z)>=p_params->d.z)) {
993 set_atom_positions(r);
997 set_atom_positions(r);
1000 /* the three face centered atoms */
1002 v3_add(&n,&r,&basis[l]);
1003 switch(p_params->type) {
1005 v3_sub(&dist,&n,&p);
1006 if(v3_absolute_square(&dist)<
1007 (p_params->r*p_params->r)) {
1008 set_atom_positions(n);
1011 case PART_OUTSIDE_R:
1012 v3_sub(&dist,&n,&p);
1013 if(v3_absolute_square(&dist)>=
1014 (p_params->r*p_params->r)) {
1015 set_atom_positions(n);
1019 v3_sub(&dist,&n,&p);
1020 if((fabs(dist.x)<p_params->d.x)&&
1021 (fabs(dist.y)<p_params->d.y)&&
1022 (fabs(dist.z)<p_params->d.z)) {
1023 set_atom_positions(n);
1026 case PART_OUTSIDE_D:
1027 v3_sub(&dist,&n,&p);
1028 if((fabs(dist.x)>=p_params->d.x)||
1029 (fabs(dist.y)>=p_params->d.y)||
1030 (fabs(dist.z)>=p_params->d.z)) {
1031 set_atom_positions(n);
1035 set_atom_positions(n);
1046 /* coordinate transformation */
1047 for(i=0;i<count;i++) {
1048 atom[i].r.x-=(a*lc)/2.0;
1049 atom[i].r.y-=(b*lc)/2.0;
1050 atom[i].r.z-=(c*lc)/2.0;
1056 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
1057 t_part_params *p_params,t_defect_params *d_params) {
1062 count=fcc_init(a,b,c,lc,atom,origin,p_params,d_params);
1068 if(origin) v3_add(&o,&o,origin);
1070 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_params,d_params);
1075 int destroy_atoms(t_moldyn *moldyn) {
1077 if(moldyn->atom) free(moldyn->atom);
1082 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
1085 * - gaussian distribution of velocities
1086 * - zero total momentum
1087 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1092 t_3dvec p_total,delta;
1097 random=&(moldyn->random);
1099 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
1101 /* gaussian distribution of velocities */
1103 for(i=0;i<moldyn->count;i++) {
1104 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
1106 v=sigma*rand_get_gauss(random);
1108 p_total.x+=atom[i].mass*v;
1110 v=sigma*rand_get_gauss(random);
1112 p_total.y+=atom[i].mass*v;
1114 v=sigma*rand_get_gauss(random);
1116 p_total.z+=atom[i].mass*v;
1119 /* zero total momentum */
1120 v3_scale(&p_total,&p_total,1.0/moldyn->count);
1121 for(i=0;i<moldyn->count;i++) {
1122 v3_scale(&delta,&p_total,1.0/atom[i].mass);
1123 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
1126 /* velocity scaling */
1127 scale_velocity(moldyn,equi_init);
1132 double total_mass_calc(t_moldyn *moldyn) {
1138 for(i=0;i<moldyn->count;i++)
1139 moldyn->mass+=moldyn->atom[i].mass;
1141 return moldyn->mass;
1144 double temperature_calc(t_moldyn *moldyn) {
1146 /* assume up to date kinetic energy, which is 3/2 N k_B T */
1149 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
1155 double get_temperature(t_moldyn *moldyn) {
1160 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
1170 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1173 /* get kinetic energy / temperature & count involved atoms */
1176 for(i=0;i<moldyn->count;i++) {
1177 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
1178 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1183 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1184 else return 0; /* no atoms involved in scaling! */
1186 /* (temporary) hack for e,t = 0 */
1189 if(moldyn->t_ref!=0.0) {
1190 thermal_init(moldyn,equi_init);
1194 return 0; /* no scaling needed */
1198 /* get scaling factor */
1199 scale=moldyn->t_ref/moldyn->t;
1203 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1204 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1207 /* velocity scaling */
1208 for(i=0;i<moldyn->count;i++) {
1209 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1210 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1216 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1220 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1225 double virial_sum(t_moldyn *moldyn) {
1230 /* virial (sum over atom virials) */
1238 for(i=0;i<moldyn->count;i++) {
1239 virial=&(moldyn->atom[i].virial);
1240 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1241 moldyn->vir.xx+=virial->xx;
1242 moldyn->vir.yy+=virial->yy;
1243 moldyn->vir.zz+=virial->zz;
1244 moldyn->vir.xy+=virial->xy;
1245 moldyn->vir.xz+=virial->xz;
1246 moldyn->vir.yz+=virial->yz;
1249 /* global virial (absolute coordinates) */
1250 //virial=&(moldyn->gvir);
1251 //moldyn->gv=virial->xx+virial->yy+virial->zz;
1253 return moldyn->virial;
1256 double pressure_calc(t_moldyn *moldyn) {
1260 * with W = 1/3 sum_i f_i r_i (- skipped!)
1261 * virial = sum_i f_i r_i
1263 * => P = (2 Ekin + virial) / (3V)
1266 /* assume up to date virial & up to date kinetic energy */
1268 /* pressure (atom virials) */
1269 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1270 moldyn->p/=(3.0*moldyn->volume);
1272 //moldyn->px=2.0*moldyn->ekinx+moldyn->vir.xx;
1273 //moldyn->px/=moldyn->volume;
1274 //moldyn->py=2.0*moldyn->ekiny+moldyn->vir.yy;
1275 //moldyn->py/=moldyn->volume;
1276 //moldyn->pz=2.0*moldyn->ekinz+moldyn->vir.zz;
1277 //moldyn->pz/=moldyn->volume;
1279 /* pressure (absolute coordinates) */
1280 //moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1281 //moldyn->gp/=(3.0*moldyn->volume);
1286 int average_reset(t_moldyn *moldyn) {
1288 printf("[moldyn] average reset\n");
1290 /* update skip value */
1291 moldyn->avg_skip=moldyn->total_steps;
1293 /* kinetic energy */
1297 /* potential energy */
1305 moldyn->virial_sum=0.0;
1306 //moldyn->gv_sum=0.0;
1310 //moldyn->gp_sum=0.0;
1316 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1320 if(moldyn->total_steps<moldyn->avg_skip)
1323 denom=moldyn->total_steps+1-moldyn->avg_skip;
1325 /* assume up to date energies, temperature, pressure etc */
1327 /* kinetic energy */
1328 moldyn->k_sum+=moldyn->ekin;
1329 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1330 moldyn->k_avg=moldyn->k_sum/denom;
1331 moldyn->k2_avg=moldyn->k2_sum/denom;
1332 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1334 /* potential energy */
1335 moldyn->v_sum+=moldyn->energy;
1336 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1337 moldyn->v_avg=moldyn->v_sum/denom;
1338 moldyn->v2_avg=moldyn->v2_sum/denom;
1339 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1342 moldyn->t_sum+=moldyn->t;
1343 moldyn->t_avg=moldyn->t_sum/denom;
1346 moldyn->virial_sum+=moldyn->virial;
1347 moldyn->virial_avg=moldyn->virial_sum/denom;
1348 //moldyn->gv_sum+=moldyn->gv;
1349 //moldyn->gv_avg=moldyn->gv_sum/denom;
1352 moldyn->p_sum+=moldyn->p;
1353 moldyn->p_avg=moldyn->p_sum/denom;
1354 //moldyn->gp_sum+=moldyn->gp;
1355 //moldyn->gp_avg=moldyn->gp_sum/denom;
1356 moldyn->tp_sum+=moldyn->tp;
1357 moldyn->tp_avg=moldyn->tp_sum/denom;
1362 int get_heat_capacity(t_moldyn *moldyn) {
1366 /* averages needed for heat capacity calc */
1367 if(moldyn->total_steps<moldyn->avg_skip)
1370 /* (temperature average)^2 */
1371 temp2=moldyn->t_avg*moldyn->t_avg;
1372 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1375 /* ideal gas contribution */
1376 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1377 printf(" ideal gas contribution: %f\n",
1378 ighc/moldyn->mass*KILOGRAM/JOULE);
1380 /* specific heat for nvt ensemble */
1381 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1382 moldyn->c_v_nvt/=moldyn->mass;
1384 /* specific heat for nve ensemble */
1385 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1386 moldyn->c_v_nve/=moldyn->mass;
1388 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1389 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1390 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)));
1395 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1411 /* store atomic configuration + dimension */
1412 store=malloc(moldyn->count*sizeof(t_atom));
1414 printf("[moldyn] allocating store mem failed\n");
1417 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1422 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1423 su=pow(2.0-h,ONE_THIRD)-1.0;
1424 dv=(1.0-h)*moldyn->volume;
1426 /* scale up dimension and atom positions */
1427 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1428 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1429 link_cell_shutdown(moldyn);
1430 link_cell_init(moldyn,QUIET);
1431 potential_force_calc(moldyn);
1434 /* restore atomic configuration + dim */
1435 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1438 /* scale down dimension and atom positions */
1439 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1440 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1441 link_cell_shutdown(moldyn);
1442 link_cell_init(moldyn,QUIET);
1443 potential_force_calc(moldyn);
1446 /* calculate pressure */
1447 moldyn->tp=-(y1-y0)/(2.0*dv);
1449 /* restore atomic configuration */
1450 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1452 link_cell_shutdown(moldyn);
1453 link_cell_init(moldyn,QUIET);
1454 //potential_force_calc(moldyn);
1456 /* free store buffer */
1463 double get_pressure(t_moldyn *moldyn) {
1469 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1481 if(x) dim->x*=scale;
1482 if(y) dim->y*=scale;
1483 if(z) dim->z*=scale;
1488 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1499 for(i=0;i<moldyn->count;i++) {
1500 r=&(moldyn->atom[i].r);
1509 int scale_atoms_ind(t_moldyn *moldyn,double x,double y,double z) {
1514 for(i=0;i<moldyn->count;i++) {
1515 r=&(moldyn->atom[i].r);
1524 int scale_dim_ind(t_moldyn *moldyn,double x,double y,double z) {
1537 int scale_volume(t_moldyn *moldyn) {
1544 vdim=&(moldyn->vis.dim);
1548 /* scaling factor */
1549 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1550 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1551 scale=pow(scale,ONE_THIRD);
1554 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1559 sx=1.0-(moldyn->p_ref-moldyn->px)*moldyn->p_tc*moldyn->tau;
1560 sy=1.0-(moldyn->p_ref-moldyn->py)*moldyn->p_tc*moldyn->tau;
1561 sz=1.0-(moldyn->p_ref-moldyn->pz)*moldyn->p_tc*moldyn->tau;
1562 sx=pow(sx,ONE_THIRD);
1563 sy=pow(sy,ONE_THIRD);
1564 sz=pow(sz,ONE_THIRD);
1567 /* scale the atoms and dimensions */
1568 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1569 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1570 //scale_atoms_ind(moldyn,sx,sy,sz);
1571 //scale_dim_ind(moldyn,sx,sy,sz);
1573 /* visualize dimensions */
1580 /* recalculate scaled volume */
1581 moldyn->volume=dim->x*dim->y*dim->z;
1583 /* adjust/reinit linkcell */
1584 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1585 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1586 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1587 link_cell_shutdown(moldyn);
1588 link_cell_init(moldyn,QUIET);
1602 double e_kin_calc(t_moldyn *moldyn) {
1609 //moldyn->ekinx=0.0;
1610 //moldyn->ekiny=0.0;
1611 //moldyn->ekinz=0.0;
1613 for(i=0;i<moldyn->count;i++) {
1614 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1615 moldyn->ekin+=atom[i].ekin;
1616 //moldyn->ekinx+=0.5*atom[i].mass*atom[i].v.x*atom[i].v.x;
1617 //moldyn->ekiny+=0.5*atom[i].mass*atom[i].v.y*atom[i].v.y;
1618 //moldyn->ekinz+=0.5*atom[i].mass*atom[i].v.z*atom[i].v.z;
1621 return moldyn->ekin;
1624 double get_total_energy(t_moldyn *moldyn) {
1626 return(moldyn->ekin+moldyn->energy);
1629 t_3dvec get_total_p(t_moldyn *moldyn) {
1638 for(i=0;i<moldyn->count;i++) {
1639 v3_scale(&p,&(atom[i].v),atom[i].mass);
1640 v3_add(&p_total,&p_total,&p);
1646 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1650 /* nn_dist is the nearest neighbour distance */
1652 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1661 /* linked list / cell method */
1663 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1666 #ifndef LOWMEM_LISTS
1672 /* partitioning the md cell */
1673 lc->nx=moldyn->dim.x/moldyn->cutoff;
1674 lc->x=moldyn->dim.x/lc->nx;
1675 lc->ny=moldyn->dim.y/moldyn->cutoff;
1676 lc->y=moldyn->dim.y/lc->ny;
1677 lc->nz=moldyn->dim.z/moldyn->cutoff;
1678 lc->z=moldyn->dim.z/lc->nz;
1679 lc->cells=lc->nx*lc->ny*lc->nz;
1682 lc->subcell=malloc(lc->cells*sizeof(int*));
1684 lc->subcell=malloc(sizeof(t_lowmem_list));
1686 lc->subcell=malloc(lc->cells*sizeof(t_list));
1689 if(lc->subcell==NULL) {
1690 perror("[moldyn] cell init (malloc)");
1695 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1700 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1703 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1706 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1709 printf(" x: %d x %f A\n",lc->nx,lc->x);
1710 printf(" y: %d x %f A\n",lc->ny,lc->y);
1711 printf(" z: %d x %f A\n",lc->nz,lc->z);
1716 for(i=0;i<lc->cells;i++) {
1717 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1718 if(lc->subcell[i]==NULL) {
1719 perror("[moldyn] list init (malloc)");
1724 printf(" ---> %d malloc %p (%p)\n",
1725 i,lc->subcell[0],lc->subcell);
1729 lc->subcell->head=malloc(lc->cells*sizeof(int));
1730 if(lc->subcell->head==NULL) {
1731 perror("[moldyn] head init (malloc)");
1734 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1735 if(lc->subcell->list==NULL) {
1736 perror("[moldyn] list init (malloc)");
1740 for(i=0;i<lc->cells;i++)
1741 list_init_f(&(lc->subcell[i]));
1744 /* update the list */
1745 link_cell_update(moldyn);
1750 int link_cell_update(t_moldyn *moldyn) {
1768 for(i=0;i<lc->cells;i++)
1770 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1772 lc->subcell->head[i]=-1;
1774 list_destroy_f(&(lc->subcell[i]));
1777 for(count=0;count<moldyn->count;count++) {
1778 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1779 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1780 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1784 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1787 if(p>=MAX_ATOMS_PER_LIST) {
1788 printf("[moldyn] FATAL: amount of atoms too high!\n");
1792 lc->subcell[i+j*nx+k*nxy][p]=count;
1795 lc->subcell->list[count]=lc->subcell->head[p];
1796 lc->subcell->head[p]=count;
1798 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1802 printf(" ---> %d %d malloc %p (%p)\n",
1803 i,count,lc->subcell[i].current,lc->subcell);
1811 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1837 if(i>=nx||j>=ny||k>=nz)
1838 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1841 #ifndef LOWMEM_LISTS
1842 cell[0]=lc->subcell[i+j*nx+k*a];
1844 cell[0]=lc->subcell->head[i+j*nx+k*a];
1846 for(ci=-1;ci<=1;ci++) {
1849 if((x<0)||(x>=nx)) {
1853 for(cj=-1;cj<=1;cj++) {
1856 if((y<0)||(y>=ny)) {
1860 for(ck=-1;ck<=1;ck++) {
1863 if((z<0)||(z>=nz)) {
1867 if(!(ci|cj|ck)) continue;
1869 #ifndef LOWMEM_LISTS
1870 cell[--count2]=lc->subcell[x+y*nx+z*a];
1872 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1877 #ifndef LOWMEM_LISTS
1878 cell[count1++]=lc->subcell[x+y*nx+z*a];
1880 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1892 int link_cell_shutdown(t_moldyn *moldyn) {
1894 #ifndef LOWMEM_LISTS
1902 free(lc->subcell->head);
1903 free(lc->subcell->list);
1907 for(i=0;i<lc->cells;i++) {
1909 free(lc->subcell[i]);
1911 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1912 list_destroy_f(&(lc->subcell[i]));
1922 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1926 t_moldyn_schedule *schedule;
1928 schedule=&(moldyn->schedule);
1929 count=++(schedule->total_sched);
1931 ptr=realloc(schedule->runs,count*sizeof(int));
1933 perror("[moldyn] realloc (runs)");
1937 schedule->runs[count-1]=runs;
1939 ptr=realloc(schedule->tau,count*sizeof(double));
1941 perror("[moldyn] realloc (tau)");
1945 schedule->tau[count-1]=tau;
1947 printf("[moldyn] schedule added:\n");
1948 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1954 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1956 moldyn->schedule.hook=hook;
1957 moldyn->schedule.hook_params=hook_params;
1964 * 'integration of newtons equation' - algorithms
1968 /* start the integration */
1970 int moldyn_integrate(t_moldyn *moldyn) {
1973 unsigned int e,m,s,v,p,t,a;
1975 t_moldyn_schedule *sched;
1980 double energy_scale;
1981 struct timeval t1,t2;
1984 #ifdef VISUAL_THREAD
1986 pthread_t io_thread;
1995 sched=&(moldyn->schedule);
1998 /* initialize linked cell method */
1999 link_cell_init(moldyn,VERBOSE);
2001 /* logging & visualization */
2010 /* sqaure of some variables */
2011 moldyn->tau_square=moldyn->tau*moldyn->tau;
2013 /* get current time */
2014 gettimeofday(&t1,NULL);
2016 /* calculate initial forces */
2017 potential_force_calc(moldyn);
2022 /* some stupid checks before we actually start calculating bullshit */
2023 if(moldyn->cutoff>0.5*moldyn->dim.x)
2024 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
2025 if(moldyn->cutoff>0.5*moldyn->dim.y)
2026 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
2027 if(moldyn->cutoff>0.5*moldyn->dim.z)
2028 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
2030 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
2031 if(ds>0.05*moldyn->nnd)
2032 printf("[moldyn] WARNING: forces too high / tau too small!\n");
2035 /* zero absolute time */
2036 // should have right values!
2038 //moldyn->total_steps=0;
2040 /* debugging, ignore */
2043 /* zero & init moldyn copy */
2044 #ifdef VISUAL_THREAD
2045 memset(&md_copy,0,sizeof(t_moldyn));
2046 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2047 if(atom_copy==NULL) {
2048 perror("[moldyn] malloc atom copy (init)");
2054 printf("##################\n");
2055 printf("# USING PTHREADS #\n");
2056 printf("##################\n");
2058 /* tell the world */
2059 printf("[moldyn] integration start, go get a coffee ...\n");
2061 /* executing the schedule */
2063 while(sched->count<sched->total_sched) {
2065 /* setting amount of runs and finite time step size */
2066 moldyn->tau=sched->tau[sched->count];
2067 moldyn->tau_square=moldyn->tau*moldyn->tau;
2068 moldyn->time_steps=sched->runs[sched->count];
2070 /* energy scaling factor (might change!) */
2071 energy_scale=moldyn->count*EV;
2073 /* integration according to schedule */
2075 for(i=0;i<moldyn->time_steps;i++) {
2077 /* integration step */
2078 moldyn->integrate(moldyn);
2080 /* calculate kinetic energy, temperature and pressure */
2082 temperature_calc(moldyn);
2084 pressure_calc(moldyn);
2086 thermodynamic_pressure_calc(moldyn);
2087 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
2091 /* calculate fluctuations + averages */
2092 average_and_fluctuation_calc(moldyn);
2095 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
2096 scale_velocity(moldyn,FALSE);
2097 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
2098 scale_volume(moldyn);
2100 /* check for log & visualization */
2102 if(!(moldyn->total_steps%e))
2103 dprintf(moldyn->efd,
2104 "%f %f %f %f %f %f\n",
2105 moldyn->time,moldyn->ekin/energy_scale,
2106 moldyn->energy/energy_scale,
2107 get_total_energy(moldyn)/energy_scale,
2108 moldyn->ekin/EV,moldyn->energy/EV);
2111 if(!(moldyn->total_steps%m)) {
2112 momentum=get_total_p(moldyn);
2113 dprintf(moldyn->mfd,
2114 "%f %f %f %f %f\n",moldyn->time,
2115 momentum.x,momentum.y,momentum.z,
2116 v3_norm(&momentum));
2120 if(!(moldyn->total_steps%p)) {
2121 dprintf(moldyn->pfd,
2122 "%f %f %f %f %f %f %f\n",moldyn->time,
2123 moldyn->p/BAR,moldyn->p_avg/BAR,
2124 moldyn->p/BAR,moldyn->p_avg/BAR,
2125 moldyn->tp/BAR,moldyn->tp_avg/BAR);
2129 if(!(moldyn->total_steps%t)) {
2130 dprintf(moldyn->tfd,
2132 moldyn->time,moldyn->t,moldyn->t_avg);
2136 if(!(moldyn->total_steps%v)) {
2137 dprintf(moldyn->vfd,
2138 "%f %f %f %f %f\n",moldyn->time,
2146 if(!(moldyn->total_steps%s)) {
2147 snprintf(dir,128,"%s/s-%08.f.save",
2148 moldyn->vlsdir,moldyn->time);
2149 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
2151 if(fd<0) perror("[moldyn] save fd open");
2153 write(fd,moldyn,sizeof(t_moldyn));
2154 write(fd,moldyn->atom,
2155 moldyn->count*sizeof(t_atom));
2161 if(!(moldyn->total_steps%a)) {
2162 #ifdef VISUAL_THREAD
2163 /* check whether thread has not terminated yet */
2165 ret=pthread_join(io_thread,NULL);
2168 /* prepare and start thread */
2169 if(moldyn->count!=md_copy.count) {
2173 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
2175 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2176 if(atom_copy==NULL) {
2177 perror("[moldyn] malloc atom copy (change)");
2181 md_copy.atom=atom_copy;
2182 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
2184 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
2186 perror("[moldyn] create visual atoms thread\n");
2190 visual_atoms(moldyn);
2195 /* display progress */
2199 /* get current time */
2200 gettimeofday(&t2,NULL);
2202 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
2203 sched->count,i,moldyn->total_steps,
2204 moldyn->t,moldyn->t_avg,
2206 moldyn->p/BAR,moldyn->p_avg/BAR,
2208 moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
2211 (int)(t2.tv_sec-t1.tv_sec));
2215 /* copy over time */
2221 /* increase absolute time */
2222 moldyn->time+=moldyn->tau;
2223 moldyn->total_steps+=1;
2227 /* check for hooks */
2229 printf("\n ## schedule hook %d start ##\n",
2231 sched->hook(moldyn,sched->hook_params);
2232 printf(" ## schedule hook end ##\n");
2235 /* increase the schedule counter */
2240 /* writing a final save file! */
2242 snprintf(dir,128,"%s/s-final.save",moldyn->vlsdir);
2243 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,S_IRUSR|S_IWUSR);
2244 if(fd<0) perror("[moldyn] save fd open");
2246 write(fd,moldyn,sizeof(t_moldyn));
2247 write(fd,moldyn->atom,
2248 moldyn->count*sizeof(t_atom));
2252 /* writing a final visual file! */
2254 visual_atoms(moldyn);
2259 /* velocity verlet */
2261 int velocity_verlet(t_moldyn *moldyn) {
2264 double tau,tau_square,h;
2269 count=moldyn->count;
2271 tau_square=moldyn->tau_square;
2273 #ifdef CONSTRAINT_110_5832
2275 atom[5832].f.x=0.5*(atom[5832].f.x-atom[5832].f.y);
2276 atom[5832].f.y=-atom[5832].f.x;
2279 for(i=0;i<count;i++) {
2280 /* check whether fixed atom */
2281 if(atom[i].attr&ATOM_ATTR_FP)
2285 v3_scale(&delta,&(atom[i].v),tau);
2286 #ifdef CONSTRAINT_110_5832
2292 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2294 v3_scale(&delta,&(atom[i].f),h*tau_square);
2295 #ifdef CONSTRAINT_110_5832
2300 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2301 //check_per_bound_and_care_for_pbc(moldyn,&(atom[i]));
2302 check_per_bound(moldyn,&(atom[i].r));
2304 /* velocities [actually v(t+tau/2)] */
2305 v3_scale(&delta,&(atom[i].f),h*tau);
2306 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2309 /* criticial check */
2310 moldyn_bc_check(moldyn);
2312 /* neighbour list update */
2313 link_cell_update(moldyn);
2315 /* forces depending on chosen potential */
2317 // if albe, use albe force calc routine
2318 //if(moldyn->func3b_j1==albe_mult_3bp_j1)
2319 // albe_potential_force_calc(moldyn);
2321 potential_force_calc(moldyn);
2323 albe_potential_force_calc(moldyn);
2326 #ifdef CONSTRAINT_110_5832
2328 atom[5832].f.x=0.5*(atom[5832].f.x-atom[5832].f.y);
2329 atom[5832].f.y=-atom[5832].f.x;
2332 for(i=0;i<count;i++) {
2333 /* check whether fixed atom */
2334 if(atom[i].attr&ATOM_ATTR_FP)
2336 /* again velocities [actually v(t+tau)] */
2337 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2338 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2347 * potentials & corresponding forces & virial routine
2351 /* generic potential and force calculation */
2353 int potential_force_calc(t_moldyn *moldyn) {
2356 t_atom *itom,*jtom,*ktom;
2360 int *neighbour_i[27];
2364 int neighbour_i[27];
2367 t_list neighbour_i[27];
2368 t_list neighbour_i2[27];
2374 count=moldyn->count;
2384 /* reset global virial */
2385 memset(&(moldyn->gvir),0,sizeof(t_virial));
2387 /* reset force, site energy and virial of every atom */
2389 i=omp_get_thread_num();
2390 #pragma omp parallel for private(virial)
2392 for(i=0;i<count;i++) {
2395 v3_zero(&(itom[i].f));
2398 virial=(&(itom[i].virial));
2406 /* reset site energy */
2411 /* get energy, force and virial of every atom */
2413 /* first (and only) loop over atoms i */
2414 for(i=0;i<count;i++) {
2416 /* single particle potential/force */
2417 if(itom[i].attr&ATOM_ATTR_1BP)
2419 moldyn->func1b(moldyn,&(itom[i]));
2421 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2424 /* 2 body pair potential/force */
2426 link_cell_neighbour_index(moldyn,
2427 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2428 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2429 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2434 /* first loop over atoms j */
2435 if(moldyn->func2b) {
2442 while(neighbour_i[j][p]!=-1) {
2444 jtom=&(atom[neighbour_i[j][p]]);
2452 p=lc->subcell->list[p];
2454 this=&(neighbour_i[j]);
2457 if(this->start==NULL)
2461 jtom=this->current->data;
2464 if(jtom==&(itom[i]))
2467 if((jtom->attr&ATOM_ATTR_2BP)&
2468 (itom[i].attr&ATOM_ATTR_2BP)) {
2469 moldyn->func2b(moldyn,
2479 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2485 /* 3 body potential/force */
2487 if(!(itom[i].attr&ATOM_ATTR_3BP))
2490 /* copy the neighbour lists */
2492 /* no copy needed for static lists */
2494 /* no copy needed for lowmem lists */
2496 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2499 /* second loop over atoms j */
2506 while(neighbour_i[j][p]!=-1) {
2508 jtom=&(atom[neighbour_i[j][p]]);
2516 p=lc->subcell->list[p];
2518 this=&(neighbour_i[j]);
2521 if(this->start==NULL)
2526 jtom=this->current->data;
2529 if(jtom==&(itom[i]))
2532 if(!(jtom->attr&ATOM_ATTR_3BP))
2538 if(moldyn->func3b_j1)
2539 moldyn->func3b_j1(moldyn,
2544 /* in first j loop, 3bp run can be skipped */
2545 if(!(moldyn->run3bp))
2548 /* first loop over atoms k */
2549 if(moldyn->func3b_k1) {
2557 while(neighbour_i[k][q]!=-1) {
2559 ktom=&(atom[neighbour_i[k][q]]);
2567 q=lc->subcell->list[q];
2569 that=&(neighbour_i2[k]);
2572 if(that->start==NULL)
2576 ktom=that->current->data;
2579 if(!(ktom->attr&ATOM_ATTR_3BP))
2585 if(ktom==&(itom[i]))
2588 moldyn->func3b_k1(moldyn,
2599 } while(list_next_f(that)!=\
2607 if(moldyn->func3b_j2)
2608 moldyn->func3b_j2(moldyn,
2613 /* second loop over atoms k */
2614 if(moldyn->func3b_k2) {
2622 while(neighbour_i[k][q]!=-1) {
2624 ktom=&(atom[neighbour_i[k][q]]);
2632 q=lc->subcell->list[q];
2634 that=&(neighbour_i2[k]);
2637 if(that->start==NULL)
2641 ktom=that->current->data;
2644 if(!(ktom->attr&ATOM_ATTR_3BP))
2650 if(ktom==&(itom[i]))
2653 moldyn->func3b_k2(moldyn,
2664 } while(list_next_f(that)!=\
2672 /* 2bp post function */
2673 if(moldyn->func3b_j3) {
2674 moldyn->func3b_j3(moldyn,
2683 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2698 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2699 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2701 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2702 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2703 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2707 /* some postprocessing */
2709 #pragma omp parallel for
2711 for(i=0;i<count;i++) {
2712 /* calculate global virial */
2713 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2714 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2715 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2716 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2717 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2718 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2720 /* check forces regarding the given timestep */
2721 if(v3_norm(&(itom[i].f))>\
2722 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2723 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2731 * virial calculation
2734 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2735 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2737 a->virial.xx+=f->x*d->x;
2738 a->virial.yy+=f->y*d->y;
2739 a->virial.zz+=f->z*d->z;
2740 a->virial.xy+=f->x*d->y;
2741 a->virial.xz+=f->x*d->z;
2742 a->virial.yz+=f->y*d->z;
2748 * periodic boundary checking
2751 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2752 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2763 if(moldyn->status&MOLDYN_STAT_PBX) {
2764 if(a->x>=x) a->x-=dim->x;
2765 else if(-a->x>x) a->x+=dim->x;
2767 if(moldyn->status&MOLDYN_STAT_PBY) {
2768 if(a->y>=y) a->y-=dim->y;
2769 else if(-a->y>y) a->y+=dim->y;
2771 if(moldyn->status&MOLDYN_STAT_PBZ) {
2772 if(a->z>=z) a->z-=dim->z;
2773 else if(-a->z>z) a->z+=dim->z;
2779 int check_per_bound_and_care_for_pbc(t_moldyn *moldyn,t_atom *a) {
2790 if(moldyn->status&MOLDYN_STAT_PBX) {
2795 else if(-a->r.x>x) {
2800 if(moldyn->status&MOLDYN_STAT_PBY) {
2805 else if(-a->r.y>y) {
2810 if(moldyn->status&MOLDYN_STAT_PBZ) {
2815 else if(-a->r.z>z) {
2825 * debugging / critical check functions
2828 int moldyn_bc_check(t_moldyn *moldyn) {
2841 for(i=0;i<moldyn->count;i++) {
2842 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2843 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2844 i,atom[i].r.x,dim->x/2);
2845 printf("diagnostic:\n");
2846 printf("-----------\natom.r.x:\n");
2848 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2851 ((byte)&(1<<k))?1:0,
2854 printf("---------------\nx=dim.x/2:\n");
2856 memcpy(&byte,(u8 *)(&x)+j,1);
2859 ((byte)&(1<<k))?1:0,
2862 if(atom[i].r.x==x) printf("the same!\n");
2863 else printf("different!\n");
2865 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2866 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2867 i,atom[i].r.y,dim->y/2);
2868 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2869 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2870 i,atom[i].r.z,dim->z/2);
2880 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2887 fd=open(file,O_RDONLY);
2889 perror("[moldyn] load save file open");
2893 fsize=lseek(fd,0,SEEK_END);
2894 lseek(fd,0,SEEK_SET);
2896 size=sizeof(t_moldyn);
2899 cnt=read(fd,moldyn,size);
2901 perror("[moldyn] load save file read (moldyn)");
2907 size=moldyn->count*sizeof(t_atom);
2909 /* correcting possible atom data offset */
2911 if(fsize!=sizeof(t_moldyn)+size) {
2912 corr=fsize-sizeof(t_moldyn)-size;
2913 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2914 printf(" modifying offset:\n");
2915 printf(" - current pos: %d\n",sizeof(t_moldyn));
2916 printf(" - atom size: %d\n",size);
2917 printf(" - file size: %d\n",fsize);
2918 printf(" => correction: %d\n",corr);
2919 lseek(fd,corr,SEEK_CUR);
2922 moldyn->atom=(t_atom *)malloc(size);
2923 if(moldyn->atom==NULL) {
2924 perror("[moldyn] load save file malloc (atoms)");
2929 cnt=read(fd,moldyn->atom,size);
2931 perror("[moldyn] load save file read (atoms)");
2938 amutex=malloc(moldyn->count*sizeof(pthread_mutex_t));
2940 perror("[moldyn] load save file (mutexes)");
2943 for(cnt=0;cnt<moldyn->count;cnt++)
2944 pthread_mutex_init(&(amutex[cnt]),NULL);
2952 int moldyn_free_save_file(t_moldyn *moldyn) {
2959 int moldyn_load(t_moldyn *moldyn) {
2967 * function to find/callback all combinations of 2 body bonds
2970 int process_2b_bonds(t_moldyn *moldyn,void *data,
2971 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2972 void *data,u8 bc)) {
2982 t_list neighbour[27];
2992 for(i=0;i<moldyn->count;i++) {
2993 /* neighbour indexing */
2994 link_cell_neighbour_index(moldyn,
2995 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2996 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2997 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
3002 bc=(j<lc->dnlc)?0:1;
3007 while(neighbour[j][p]!=-1) {
3009 jtom=&(moldyn->atom[neighbour[j][p]]);
3017 p=lc->subcell->list[p];
3019 this=&(neighbour[j]);
3022 if(this->start==NULL)
3027 jtom=this->current->data;
3031 process(moldyn,&(itom[i]),jtom,data,bc);
3038 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3048 * function to find neighboured atoms
3051 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
3052 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
3053 void *data,u8 bc)) {
3063 t_list neighbour[27];
3072 /* neighbour indexing */
3073 link_cell_neighbour_index(moldyn,
3074 (atom->r.x+moldyn->dim.x/2)/lc->x,
3075 (atom->r.y+moldyn->dim.y/2)/lc->x,
3076 (atom->r.z+moldyn->dim.z/2)/lc->x,
3081 bc=(j<lc->dnlc)?0:1;
3086 while(neighbour[j][p]!=-1) {
3088 natom=&(moldyn->atom[neighbour[j][p]]);
3095 natom=&(moldyn->atom[p]);
3096 p=lc->subcell->list[p];
3098 this=&(neighbour[j]);
3101 if(this->start==NULL)
3106 natom=this->current->data;
3110 process(moldyn,atom,natom,data,bc);
3117 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3126 * post processing functions
3129 int get_line(int fd,char *line,int max) {
3136 if(count==max) return count;
3137 ret=read(fd,line+count,1);
3138 if(ret<=0) return ret;
3139 if(line[count]=='\n') {
3140 memset(line+count,0,max-count-1);
3148 int pair_correlation_init(t_moldyn *moldyn,double dr) {
3154 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
3171 for(i=0;i<moldyn->count;i++) {
3173 /* care for pb crossing */
3174 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3175 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3176 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3177 /* calculate distance */
3178 v3_sub(&dist,&final_r,&(atom[i].r_0));
3179 d2=v3_absolute_square(&dist);
3193 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
3194 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
3195 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
3200 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
3205 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
3206 t_atom *jtom,void *data,u8 bc) {
3213 /* only count pairs once,
3214 * skip same atoms */
3215 if(itom->tag>=jtom->tag)
3219 * pair correlation calc
3226 v3_sub(&dist,&(jtom->r),&(itom->r));
3227 if(bc) check_per_bound(moldyn,&dist);
3228 d=v3_absolute_square(&dist);
3230 /* ignore if greater cutoff */
3231 if(d>moldyn->cutoff_square)
3234 /* fill the slots */
3238 /* should never happen but it does 8) -
3239 * related to -ffloat-store problem! */
3241 printf("[moldyn] WARNING: pcc (%d/%d)",
3247 if(itom->brand!=jtom->brand) {
3252 /* type a - type a bonds */
3254 pcc->stat[s+pcc->o1]+=1;
3256 /* type b - type b bonds */
3257 pcc->stat[s+pcc->o2]+=1;
3263 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
3270 pcc.o1=moldyn->cutoff/dr;
3273 if(pcc.o1*dr<=moldyn->cutoff)
3274 printf("[moldyn] WARNING: pcc (low #slots)\n");
3276 printf("[moldyn] pair correlation calc info:\n");
3277 printf(" time: %f\n",moldyn->time);
3278 printf(" count: %d\n",moldyn->count);
3279 printf(" cutoff: %f\n",moldyn->cutoff);
3280 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
3283 pcc.stat=(double *)ptr;
3286 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
3287 if(pcc.stat==NULL) {
3288 perror("[moldyn] pair correlation malloc");
3293 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
3296 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
3299 for(i=1;i<pcc.o1;i++) {
3300 // normalization: 4 pi r^2 dr
3301 // here: not double counting pairs -> 2 pi r r dr
3302 // ... and actually it's a constant times r^2
3305 pcc.stat[pcc.o1+i]/=norm;
3306 pcc.stat[pcc.o2+i]/=norm;
3311 /* todo: store/print pair correlation function */
3318 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3325 if(itom->tag>=jtom->tag)
3329 v3_sub(&dist,&(jtom->r),&(itom->r));
3330 if(bc) check_per_bound(moldyn,&dist);
3331 d=v3_absolute_square(&dist);
3333 /* ignore if greater or equal cutoff */
3334 if(d>moldyn->cutoff_square)
3337 /* check for potential bond */
3338 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3341 /* now count this bonding ... */
3344 /* increase total bond counter
3349 ba->acnt[jtom->tag]+=1;
3351 ba->bcnt[jtom->tag]+=1;
3354 ba->acnt[itom->tag]+=1;
3356 ba->bcnt[itom->tag]+=1;
3361 int bond_analyze(t_moldyn *moldyn,double *quality) {
3372 ba.acnt=malloc(moldyn->count*sizeof(int));
3374 perror("[moldyn] bond analyze malloc (a)");
3377 memset(ba.acnt,0,moldyn->count*sizeof(int));
3379 ba.bcnt=malloc(moldyn->count*sizeof(int));
3381 perror("[moldyn] bond analyze malloc (b)");
3384 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3393 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3395 for(i=0;i<moldyn->count;i++) {
3396 if(atom[i].brand==0) {
3397 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3399 if((ba.acnt[i]==0)&(ba.bcnt[i]==3))
3403 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3407 if((ba.acnt[i]==3)&(ba.bcnt[i]==0)) {
3416 quality[0]=1.0*ccnt4/bcnt;
3417 quality[1]=1.0*ccnt3/bcnt;
3420 printf("[moldyn] bond analyze: %f %f\n",
3421 1.0*ccnt4/bcnt,1.0*ccnt3/bcnt);
3428 * visualization code
3431 int visual_init(t_moldyn *moldyn,char *filebase) {
3433 strncpy(moldyn->vis.fb,filebase,128);
3438 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3445 if(itom->tag>=jtom->tag)
3448 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3451 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3452 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3453 itom->r.x,itom->r.y,itom->r.z,
3454 jtom->r.x,jtom->r.y,jtom->r.z);
3459 #ifdef VISUAL_THREAD
3460 void *visual_atoms(void *ptr) {
3462 int visual_atoms(t_moldyn *moldyn) {
3473 #ifdef VISUAL_THREAD
3487 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3488 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3490 perror("open visual save file fd");
3491 #ifndef VISUAL_THREAD
3496 /* write the actual data file */
3499 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3500 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3502 // atomic configuration
3503 for(i=0;i<moldyn->count;i++) {
3504 v3_sub(&strain,&(atom[i].r),&(atom[i].r_0));
3505 check_per_bound(moldyn,&strain);
3506 // atom type, positions, color and kinetic energy
3507 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3511 pse_col[atom[i].element],
3513 sqrt(v3_absolute_square(&strain)));
3516 // bonds between atoms
3517 #ifndef VISUAL_THREAD
3518 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3523 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3524 -dim.x/2,-dim.y/2,-dim.z/2,
3525 dim.x/2,-dim.y/2,-dim.z/2);
3526 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3527 -dim.x/2,-dim.y/2,-dim.z/2,
3528 -dim.x/2,dim.y/2,-dim.z/2);
3529 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3530 dim.x/2,dim.y/2,-dim.z/2,
3531 dim.x/2,-dim.y/2,-dim.z/2);
3532 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3533 -dim.x/2,dim.y/2,-dim.z/2,
3534 dim.x/2,dim.y/2,-dim.z/2);
3536 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3537 -dim.x/2,-dim.y/2,dim.z/2,
3538 dim.x/2,-dim.y/2,dim.z/2);
3539 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3540 -dim.x/2,-dim.y/2,dim.z/2,
3541 -dim.x/2,dim.y/2,dim.z/2);
3542 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3543 dim.x/2,dim.y/2,dim.z/2,
3544 dim.x/2,-dim.y/2,dim.z/2);
3545 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3546 -dim.x/2,dim.y/2,dim.z/2,
3547 dim.x/2,dim.y/2,dim.z/2);
3549 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3550 -dim.x/2,-dim.y/2,dim.z/2,
3551 -dim.x/2,-dim.y/2,-dim.z/2);
3552 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3553 -dim.x/2,dim.y/2,dim.z/2,
3554 -dim.x/2,dim.y/2,-dim.z/2);
3555 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3556 dim.x/2,-dim.y/2,dim.z/2,
3557 dim.x/2,-dim.y/2,-dim.z/2);
3558 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3559 dim.x/2,dim.y/2,dim.z/2,
3560 dim.x/2,dim.y/2,-dim.z/2);
3565 #ifdef VISUAL_THREAD
3576 * fpu cntrol functions
3579 // set rounding to double (eliminates -ffloat-store!)
3580 int fpu_set_rtd(void) {
3586 ctrl&=~_FPU_EXTENDED;