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) {
532 pthread_mutex_t *mutex;
538 /* how many atoms do we expect */
541 printf("[moldyn] WARNING: create 'none' lattice called");
543 if(type==CUBIC) new*=1;
544 if(type==FCC) new*=4;
545 if(type==DIAMOND) new*=8;
549 switch(d_params->stype) {
550 case DEFECT_STYPE_DB_X:
551 case DEFECT_STYPE_DB_Y:
552 case DEFECT_STYPE_DB_Z:
553 case DEFECT_STYPE_DB_R:
557 printf("[moldyn] WARNING: cl unknown defect\n");
562 /* allocate space for atoms */
563 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
565 perror("[moldyn] realloc (create lattice)");
569 atom=&(moldyn->atom[count]);
572 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
574 perror("[moldyn] mutex realloc (add atom)");
578 mutex=&(amutex[count]);
581 /* no atoms on the boundaries (only reason: it looks better!) */
595 set_nn_dist(moldyn,lc);
596 ret=cubic_init(a,b,c,lc,atom,&orig,p_params,d_params);
597 strcpy(name,"cubic");
601 v3_scale(&orig,&orig,0.5);
602 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
603 ret=fcc_init(a,b,c,lc,atom,&orig,p_params,d_params);
608 v3_scale(&orig,&orig,0.25);
609 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
610 ret=diamond_init(a,b,c,lc,atom,&orig,p_params,d_params);
611 strcpy(name,"diamond");
614 printf("unknown lattice type (%02x)\n",type);
620 printf("[moldyn] creating %s lattice (lc=%f) incomplete\n",
622 printf(" (ignore for partial lattice creation)\n");
623 printf(" amount of atoms\n");
624 printf(" - expected: %d\n",new);
625 printf(" - created: %d\n",ret);
630 printf("[moldyn] created %s lattice with %d atoms\n",name,ret);
632 for(new=0;new<ret;new++) {
633 atom[new].element=element;
634 atom[new].mass=pse_mass[element];
636 atom[new].brand=brand;
637 atom[new].tag=count+new;
638 check_per_bound(moldyn,&(atom[new].r));
639 atom[new].r_0=atom[new].r;
641 pthread_mutex_init(&(mutex[new]),NULL);
645 atom[new].element=d_params->element;
646 atom[new].mass=pse_mass[d_params->element];
647 atom[new].attr=d_params->attr;
648 atom[new].brand=d_params->brand;
649 atom[new].tag=count+new;
650 check_per_bound(moldyn,&(atom[new].r));
651 atom[new].r_0=atom[new].r;
653 pthread_mutex_init(&(mutex[new]),NULL);
659 ptr=realloc(moldyn->atom,moldyn->count*sizeof(t_atom));
661 perror("[moldyn] realloc (create lattice - alloc fix)");
666 // WHAT ABOUT AMUTEX !!!!
669 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
671 perror("[moldyn] list realloc (create lattice)");
674 moldyn->lc.subcell->list=ptr;
677 /* update total system mass */
678 total_mass_calc(moldyn);
683 int add_atom(t_moldyn *moldyn,int element,u8 brand,u8 attr,
684 t_3dvec *r,t_3dvec *v) {
691 count=(moldyn->count)++; // asshole style!
693 ptr=realloc(atom,(count+1)*sizeof(t_atom));
695 perror("[moldyn] realloc (add atom)");
701 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
703 perror("[moldyn] list realloc (add atom)");
706 moldyn->lc.subcell->list=ptr;
710 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
712 perror("[moldyn] mutex realloc (add atom)");
716 pthread_mutex_init(&(amutex[count]),NULL);
721 /* initialize new atom */
722 memset(&(atom[count]),0,sizeof(t_atom));
725 atom[count].element=element;
726 atom[count].mass=pse_mass[element];
727 atom[count].brand=brand;
728 atom[count].tag=count;
729 atom[count].attr=attr;
730 check_per_bound(moldyn,&(atom[count].r));
731 atom[count].r_0=atom[count].r;
733 /* update total system mass */
734 total_mass_calc(moldyn);
739 int del_atom(t_moldyn *moldyn,int tag) {
743 #if defined LOWMEM_LISTS || defined PTHREADS
749 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
751 perror("[moldyn]malloc (del atom)");
755 for(cnt=0;cnt<tag;cnt++)
758 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
760 new[cnt-1].tag=cnt-1;
769 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
771 perror("[moldyn] list realloc (del atom)");
774 moldyn->lc.subcell->list=ptr;
776 link_cell_update(moldyn);
780 ptr=realloc(amutex,moldyn->count*sizeof(pthread_mutex_t));
782 perror("[moldyn] mutex realloc (add atom)");
786 pthread_mutex_destroy(&(amutex[moldyn->count+1]));
793 #define set_atom_positions(pos) \
794 if(d_params->type) {\
795 d_o.x=0; d_o.y=0; d_o.z=0;\
796 d_d.x=0; d_d.y=0; d_d.z=0;\
797 switch(d_params->stype) {\
798 case DEFECT_STYPE_DB_X:\
802 case DEFECT_STYPE_DB_Y:\
806 case DEFECT_STYPE_DB_Z:\
812 case DEFECT_STYPE_DB_R:\
815 printf("[moldyn] WARNING: unknown defect\n");\
818 v3_add(&dr,&pos,&d_o);\
819 v3_copy(&(atom[count].r),&dr);\
821 v3_add(&dr,&pos,&d_d);\
822 v3_copy(&(atom[count].r),&dr);\
826 v3_copy(&(atom[count].r),&pos);\
831 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
832 t_part_params *p_params,t_defect_params *d_params) {
852 /* shift partition values */
854 p.x=p_params->p.x+(a*lc)/2.0;
855 p.y=p_params->p.y+(b*lc)/2.0;
856 p.z=p_params->p.z+(c*lc)/2.0;
865 switch(p_params->type) {
868 if(v3_absolute_square(&dist)<
869 (p_params->r*p_params->r)) {
870 set_atom_positions(r);
875 if(v3_absolute_square(&dist)>=
876 (p_params->r*p_params->r)) {
877 set_atom_positions(r);
882 if((fabs(dist.x)<p_params->d.x)&&
883 (fabs(dist.y)<p_params->d.y)&&
884 (fabs(dist.z)<p_params->d.z)) {
885 set_atom_positions(r);
890 if((fabs(dist.x)>=p_params->d.x)||
891 (fabs(dist.y)>=p_params->d.y)||
892 (fabs(dist.z)>=p_params->d.z)) {
893 set_atom_positions(r);
897 set_atom_positions(r);
907 for(i=0;i<count;i++) {
908 atom[i].r.x-=(a*lc)/2.0;
909 atom[i].r.y-=(b*lc)/2.0;
910 atom[i].r.z-=(c*lc)/2.0;
916 /* fcc lattice init */
917 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
918 t_part_params *p_params,t_defect_params *d_params) {
936 /* construct the basis */
937 memset(basis,0,3*sizeof(t_3dvec));
945 /* shift partition values */
947 p.x=p_params->p.x+(a*lc)/2.0;
948 p.y=p_params->p.y+(b*lc)/2.0;
949 p.z=p_params->p.z+(c*lc)/2.0;
952 /* fill up the room */
960 switch(p_params->type) {
963 if(v3_absolute_square(&dist)<
964 (p_params->r*p_params->r)) {
965 set_atom_positions(r);
970 if(v3_absolute_square(&dist)>=
971 (p_params->r*p_params->r)) {
972 set_atom_positions(r);
977 if((fabs(dist.x)<p_params->d.x)&&
978 (fabs(dist.y)<p_params->d.y)&&
979 (fabs(dist.z)<p_params->d.z)) {
980 set_atom_positions(r);
985 if((fabs(dist.x)>=p_params->d.x)||
986 (fabs(dist.y)>=p_params->d.y)||
987 (fabs(dist.z)>=p_params->d.z)) {
988 set_atom_positions(r);
992 set_atom_positions(r);
995 /* the three face centered atoms */
997 v3_add(&n,&r,&basis[l]);
998 switch(p_params->type) {
1000 v3_sub(&dist,&n,&p);
1001 if(v3_absolute_square(&dist)<
1002 (p_params->r*p_params->r)) {
1003 set_atom_positions(n);
1006 case PART_OUTSIDE_R:
1007 v3_sub(&dist,&n,&p);
1008 if(v3_absolute_square(&dist)>=
1009 (p_params->r*p_params->r)) {
1010 set_atom_positions(n);
1014 v3_sub(&dist,&n,&p);
1015 if((fabs(dist.x)<p_params->d.x)&&
1016 (fabs(dist.y)<p_params->d.y)&&
1017 (fabs(dist.z)<p_params->d.z)) {
1018 set_atom_positions(n);
1021 case PART_OUTSIDE_D:
1022 v3_sub(&dist,&n,&p);
1023 if((fabs(dist.x)>=p_params->d.x)||
1024 (fabs(dist.y)>=p_params->d.y)||
1025 (fabs(dist.z)>=p_params->d.z)) {
1026 set_atom_positions(n);
1030 set_atom_positions(n);
1041 /* coordinate transformation */
1042 for(i=0;i<count;i++) {
1043 atom[i].r.x-=(a*lc)/2.0;
1044 atom[i].r.y-=(b*lc)/2.0;
1045 atom[i].r.z-=(c*lc)/2.0;
1051 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
1052 t_part_params *p_params,t_defect_params *d_params) {
1057 count=fcc_init(a,b,c,lc,atom,origin,p_params,d_params);
1063 if(origin) v3_add(&o,&o,origin);
1065 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_params,d_params);
1070 int destroy_atoms(t_moldyn *moldyn) {
1072 if(moldyn->atom) free(moldyn->atom);
1077 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
1080 * - gaussian distribution of velocities
1081 * - zero total momentum
1082 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1087 t_3dvec p_total,delta;
1092 random=&(moldyn->random);
1094 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
1096 /* gaussian distribution of velocities */
1098 for(i=0;i<moldyn->count;i++) {
1099 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
1101 v=sigma*rand_get_gauss(random);
1103 p_total.x+=atom[i].mass*v;
1105 v=sigma*rand_get_gauss(random);
1107 p_total.y+=atom[i].mass*v;
1109 v=sigma*rand_get_gauss(random);
1111 p_total.z+=atom[i].mass*v;
1114 /* zero total momentum */
1115 v3_scale(&p_total,&p_total,1.0/moldyn->count);
1116 for(i=0;i<moldyn->count;i++) {
1117 v3_scale(&delta,&p_total,1.0/atom[i].mass);
1118 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
1121 /* velocity scaling */
1122 scale_velocity(moldyn,equi_init);
1127 double total_mass_calc(t_moldyn *moldyn) {
1133 for(i=0;i<moldyn->count;i++)
1134 moldyn->mass+=moldyn->atom[i].mass;
1136 return moldyn->mass;
1139 double temperature_calc(t_moldyn *moldyn) {
1141 /* assume up to date kinetic energy, which is 3/2 N k_B T */
1144 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
1150 double get_temperature(t_moldyn *moldyn) {
1155 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
1165 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1168 /* get kinetic energy / temperature & count involved atoms */
1171 for(i=0;i<moldyn->count;i++) {
1172 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
1173 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1178 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1179 else return 0; /* no atoms involved in scaling! */
1181 /* (temporary) hack for e,t = 0 */
1184 if(moldyn->t_ref!=0.0) {
1185 thermal_init(moldyn,equi_init);
1189 return 0; /* no scaling needed */
1193 /* get scaling factor */
1194 scale=moldyn->t_ref/moldyn->t;
1198 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1199 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1202 /* velocity scaling */
1203 for(i=0;i<moldyn->count;i++) {
1204 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1205 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1211 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1215 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1220 double virial_sum(t_moldyn *moldyn) {
1225 /* virial (sum over atom virials) */
1233 for(i=0;i<moldyn->count;i++) {
1234 virial=&(moldyn->atom[i].virial);
1235 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1236 moldyn->vir.xx+=virial->xx;
1237 moldyn->vir.yy+=virial->yy;
1238 moldyn->vir.zz+=virial->zz;
1239 moldyn->vir.xy+=virial->xy;
1240 moldyn->vir.xz+=virial->xz;
1241 moldyn->vir.yz+=virial->yz;
1244 /* global virial (absolute coordinates) */
1245 //virial=&(moldyn->gvir);
1246 //moldyn->gv=virial->xx+virial->yy+virial->zz;
1248 return moldyn->virial;
1251 double pressure_calc(t_moldyn *moldyn) {
1255 * with W = 1/3 sum_i f_i r_i (- skipped!)
1256 * virial = sum_i f_i r_i
1258 * => P = (2 Ekin + virial) / (3V)
1261 /* assume up to date virial & up to date kinetic energy */
1263 /* pressure (atom virials) */
1264 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1265 moldyn->p/=(3.0*moldyn->volume);
1267 //moldyn->px=2.0*moldyn->ekinx+moldyn->vir.xx;
1268 //moldyn->px/=moldyn->volume;
1269 //moldyn->py=2.0*moldyn->ekiny+moldyn->vir.yy;
1270 //moldyn->py/=moldyn->volume;
1271 //moldyn->pz=2.0*moldyn->ekinz+moldyn->vir.zz;
1272 //moldyn->pz/=moldyn->volume;
1274 /* pressure (absolute coordinates) */
1275 //moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1276 //moldyn->gp/=(3.0*moldyn->volume);
1281 int average_reset(t_moldyn *moldyn) {
1283 printf("[moldyn] average reset\n");
1285 /* update skip value */
1286 moldyn->avg_skip=moldyn->total_steps;
1288 /* kinetic energy */
1292 /* potential energy */
1300 moldyn->virial_sum=0.0;
1301 //moldyn->gv_sum=0.0;
1305 //moldyn->gp_sum=0.0;
1311 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1315 if(moldyn->total_steps<moldyn->avg_skip)
1318 denom=moldyn->total_steps+1-moldyn->avg_skip;
1320 /* assume up to date energies, temperature, pressure etc */
1322 /* kinetic energy */
1323 moldyn->k_sum+=moldyn->ekin;
1324 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1325 moldyn->k_avg=moldyn->k_sum/denom;
1326 moldyn->k2_avg=moldyn->k2_sum/denom;
1327 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1329 /* potential energy */
1330 moldyn->v_sum+=moldyn->energy;
1331 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1332 moldyn->v_avg=moldyn->v_sum/denom;
1333 moldyn->v2_avg=moldyn->v2_sum/denom;
1334 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1337 moldyn->t_sum+=moldyn->t;
1338 moldyn->t_avg=moldyn->t_sum/denom;
1341 moldyn->virial_sum+=moldyn->virial;
1342 moldyn->virial_avg=moldyn->virial_sum/denom;
1343 //moldyn->gv_sum+=moldyn->gv;
1344 //moldyn->gv_avg=moldyn->gv_sum/denom;
1347 moldyn->p_sum+=moldyn->p;
1348 moldyn->p_avg=moldyn->p_sum/denom;
1349 //moldyn->gp_sum+=moldyn->gp;
1350 //moldyn->gp_avg=moldyn->gp_sum/denom;
1351 moldyn->tp_sum+=moldyn->tp;
1352 moldyn->tp_avg=moldyn->tp_sum/denom;
1357 int get_heat_capacity(t_moldyn *moldyn) {
1361 /* averages needed for heat capacity calc */
1362 if(moldyn->total_steps<moldyn->avg_skip)
1365 /* (temperature average)^2 */
1366 temp2=moldyn->t_avg*moldyn->t_avg;
1367 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1370 /* ideal gas contribution */
1371 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1372 printf(" ideal gas contribution: %f\n",
1373 ighc/moldyn->mass*KILOGRAM/JOULE);
1375 /* specific heat for nvt ensemble */
1376 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1377 moldyn->c_v_nvt/=moldyn->mass;
1379 /* specific heat for nve ensemble */
1380 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1381 moldyn->c_v_nve/=moldyn->mass;
1383 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1384 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1385 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)));
1390 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1406 /* store atomic configuration + dimension */
1407 store=malloc(moldyn->count*sizeof(t_atom));
1409 printf("[moldyn] allocating store mem failed\n");
1412 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1417 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1418 su=pow(2.0-h,ONE_THIRD)-1.0;
1419 dv=(1.0-h)*moldyn->volume;
1421 /* scale up dimension and atom positions */
1422 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1423 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1424 link_cell_shutdown(moldyn);
1425 link_cell_init(moldyn,QUIET);
1426 potential_force_calc(moldyn);
1429 /* restore atomic configuration + dim */
1430 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1433 /* scale down dimension and atom positions */
1434 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1435 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1436 link_cell_shutdown(moldyn);
1437 link_cell_init(moldyn,QUIET);
1438 potential_force_calc(moldyn);
1441 /* calculate pressure */
1442 moldyn->tp=-(y1-y0)/(2.0*dv);
1444 /* restore atomic configuration */
1445 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1447 link_cell_shutdown(moldyn);
1448 link_cell_init(moldyn,QUIET);
1449 //potential_force_calc(moldyn);
1451 /* free store buffer */
1458 double get_pressure(t_moldyn *moldyn) {
1464 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1476 if(x) dim->x*=scale;
1477 if(y) dim->y*=scale;
1478 if(z) dim->z*=scale;
1483 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1494 for(i=0;i<moldyn->count;i++) {
1495 r=&(moldyn->atom[i].r);
1504 int scale_atoms_ind(t_moldyn *moldyn,double x,double y,double z) {
1509 for(i=0;i<moldyn->count;i++) {
1510 r=&(moldyn->atom[i].r);
1519 int scale_dim_ind(t_moldyn *moldyn,double x,double y,double z) {
1532 int scale_volume(t_moldyn *moldyn) {
1539 vdim=&(moldyn->vis.dim);
1543 /* scaling factor */
1544 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1545 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1546 scale=pow(scale,ONE_THIRD);
1549 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1554 sx=1.0-(moldyn->p_ref-moldyn->px)*moldyn->p_tc*moldyn->tau;
1555 sy=1.0-(moldyn->p_ref-moldyn->py)*moldyn->p_tc*moldyn->tau;
1556 sz=1.0-(moldyn->p_ref-moldyn->pz)*moldyn->p_tc*moldyn->tau;
1557 sx=pow(sx,ONE_THIRD);
1558 sy=pow(sy,ONE_THIRD);
1559 sz=pow(sz,ONE_THIRD);
1562 /* scale the atoms and dimensions */
1563 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1564 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1565 //scale_atoms_ind(moldyn,sx,sy,sz);
1566 //scale_dim_ind(moldyn,sx,sy,sz);
1568 /* visualize dimensions */
1575 /* recalculate scaled volume */
1576 moldyn->volume=dim->x*dim->y*dim->z;
1578 /* adjust/reinit linkcell */
1579 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1580 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1581 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1582 link_cell_shutdown(moldyn);
1583 link_cell_init(moldyn,QUIET);
1597 double e_kin_calc(t_moldyn *moldyn) {
1604 //moldyn->ekinx=0.0;
1605 //moldyn->ekiny=0.0;
1606 //moldyn->ekinz=0.0;
1608 for(i=0;i<moldyn->count;i++) {
1609 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1610 moldyn->ekin+=atom[i].ekin;
1611 //moldyn->ekinx+=0.5*atom[i].mass*atom[i].v.x*atom[i].v.x;
1612 //moldyn->ekiny+=0.5*atom[i].mass*atom[i].v.y*atom[i].v.y;
1613 //moldyn->ekinz+=0.5*atom[i].mass*atom[i].v.z*atom[i].v.z;
1616 return moldyn->ekin;
1619 double get_total_energy(t_moldyn *moldyn) {
1621 return(moldyn->ekin+moldyn->energy);
1624 t_3dvec get_total_p(t_moldyn *moldyn) {
1633 for(i=0;i<moldyn->count;i++) {
1634 v3_scale(&p,&(atom[i].v),atom[i].mass);
1635 v3_add(&p_total,&p_total,&p);
1641 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1645 /* nn_dist is the nearest neighbour distance */
1647 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1656 /* linked list / cell method */
1658 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1661 #ifndef LOWMEM_LISTS
1667 /* partitioning the md cell */
1668 lc->nx=moldyn->dim.x/moldyn->cutoff;
1669 lc->x=moldyn->dim.x/lc->nx;
1670 lc->ny=moldyn->dim.y/moldyn->cutoff;
1671 lc->y=moldyn->dim.y/lc->ny;
1672 lc->nz=moldyn->dim.z/moldyn->cutoff;
1673 lc->z=moldyn->dim.z/lc->nz;
1674 lc->cells=lc->nx*lc->ny*lc->nz;
1677 lc->subcell=malloc(lc->cells*sizeof(int*));
1679 lc->subcell=malloc(sizeof(t_lowmem_list));
1681 lc->subcell=malloc(lc->cells*sizeof(t_list));
1684 if(lc->subcell==NULL) {
1685 perror("[moldyn] cell init (malloc)");
1690 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1695 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1698 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1701 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1704 printf(" x: %d x %f A\n",lc->nx,lc->x);
1705 printf(" y: %d x %f A\n",lc->ny,lc->y);
1706 printf(" z: %d x %f A\n",lc->nz,lc->z);
1711 for(i=0;i<lc->cells;i++) {
1712 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1713 if(lc->subcell[i]==NULL) {
1714 perror("[moldyn] list init (malloc)");
1719 printf(" ---> %d malloc %p (%p)\n",
1720 i,lc->subcell[0],lc->subcell);
1724 lc->subcell->head=malloc(lc->cells*sizeof(int));
1725 if(lc->subcell->head==NULL) {
1726 perror("[moldyn] head init (malloc)");
1729 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1730 if(lc->subcell->list==NULL) {
1731 perror("[moldyn] list init (malloc)");
1735 for(i=0;i<lc->cells;i++)
1736 list_init_f(&(lc->subcell[i]));
1739 /* update the list */
1740 link_cell_update(moldyn);
1745 int link_cell_update(t_moldyn *moldyn) {
1763 for(i=0;i<lc->cells;i++)
1765 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1767 lc->subcell->head[i]=-1;
1769 list_destroy_f(&(lc->subcell[i]));
1772 for(count=0;count<moldyn->count;count++) {
1773 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1774 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1775 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1779 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1782 if(p>=MAX_ATOMS_PER_LIST) {
1783 printf("[moldyn] FATAL: amount of atoms too high!\n");
1787 lc->subcell[i+j*nx+k*nxy][p]=count;
1790 lc->subcell->list[count]=lc->subcell->head[p];
1791 lc->subcell->head[p]=count;
1793 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1797 printf(" ---> %d %d malloc %p (%p)\n",
1798 i,count,lc->subcell[i].current,lc->subcell);
1806 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1832 if(i>=nx||j>=ny||k>=nz)
1833 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1836 #ifndef LOWMEM_LISTS
1837 cell[0]=lc->subcell[i+j*nx+k*a];
1839 cell[0]=lc->subcell->head[i+j*nx+k*a];
1841 for(ci=-1;ci<=1;ci++) {
1844 if((x<0)||(x>=nx)) {
1848 for(cj=-1;cj<=1;cj++) {
1851 if((y<0)||(y>=ny)) {
1855 for(ck=-1;ck<=1;ck++) {
1858 if((z<0)||(z>=nz)) {
1862 if(!(ci|cj|ck)) continue;
1864 #ifndef LOWMEM_LISTS
1865 cell[--count2]=lc->subcell[x+y*nx+z*a];
1867 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1872 #ifndef LOWMEM_LISTS
1873 cell[count1++]=lc->subcell[x+y*nx+z*a];
1875 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1887 int link_cell_shutdown(t_moldyn *moldyn) {
1889 #ifndef LOWMEM_LISTS
1897 free(lc->subcell->head);
1898 free(lc->subcell->list);
1902 for(i=0;i<lc->cells;i++) {
1904 free(lc->subcell[i]);
1906 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1907 list_destroy_f(&(lc->subcell[i]));
1917 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1921 t_moldyn_schedule *schedule;
1923 schedule=&(moldyn->schedule);
1924 count=++(schedule->total_sched);
1926 ptr=realloc(schedule->runs,count*sizeof(int));
1928 perror("[moldyn] realloc (runs)");
1932 schedule->runs[count-1]=runs;
1934 ptr=realloc(schedule->tau,count*sizeof(double));
1936 perror("[moldyn] realloc (tau)");
1940 schedule->tau[count-1]=tau;
1942 printf("[moldyn] schedule added:\n");
1943 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1949 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1951 moldyn->schedule.hook=hook;
1952 moldyn->schedule.hook_params=hook_params;
1959 * 'integration of newtons equation' - algorithms
1963 /* start the integration */
1965 int moldyn_integrate(t_moldyn *moldyn) {
1968 unsigned int e,m,s,v,p,t,a;
1970 t_moldyn_schedule *sched;
1975 double energy_scale;
1976 struct timeval t1,t2;
1979 #ifdef VISUAL_THREAD
1981 pthread_t io_thread;
1990 sched=&(moldyn->schedule);
1993 /* initialize linked cell method */
1994 link_cell_init(moldyn,VERBOSE);
1996 /* logging & visualization */
2005 /* sqaure of some variables */
2006 moldyn->tau_square=moldyn->tau*moldyn->tau;
2008 /* get current time */
2009 gettimeofday(&t1,NULL);
2011 /* calculate initial forces */
2012 potential_force_calc(moldyn);
2017 /* some stupid checks before we actually start calculating bullshit */
2018 if(moldyn->cutoff>0.5*moldyn->dim.x)
2019 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
2020 if(moldyn->cutoff>0.5*moldyn->dim.y)
2021 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
2022 if(moldyn->cutoff>0.5*moldyn->dim.z)
2023 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
2025 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
2026 if(ds>0.05*moldyn->nnd)
2027 printf("[moldyn] WARNING: forces too high / tau too small!\n");
2030 /* zero absolute time */
2031 // should have right values!
2033 //moldyn->total_steps=0;
2035 /* debugging, ignore */
2038 /* zero & init moldyn copy */
2039 #ifdef VISUAL_THREAD
2040 memset(&md_copy,0,sizeof(t_moldyn));
2041 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2042 if(atom_copy==NULL) {
2043 perror("[moldyn] malloc atom copy (init)");
2049 printf("##################\n");
2050 printf("# USING PTHREADS #\n");
2051 printf("##################\n");
2053 /* tell the world */
2054 printf("[moldyn] integration start, go get a coffee ...\n");
2056 /* executing the schedule */
2058 while(sched->count<sched->total_sched) {
2060 /* setting amount of runs and finite time step size */
2061 moldyn->tau=sched->tau[sched->count];
2062 moldyn->tau_square=moldyn->tau*moldyn->tau;
2063 moldyn->time_steps=sched->runs[sched->count];
2065 /* energy scaling factor (might change!) */
2066 energy_scale=moldyn->count*EV;
2068 /* integration according to schedule */
2070 for(i=0;i<moldyn->time_steps;i++) {
2072 /* integration step */
2073 moldyn->integrate(moldyn);
2075 /* calculate kinetic energy, temperature and pressure */
2077 temperature_calc(moldyn);
2079 pressure_calc(moldyn);
2081 thermodynamic_pressure_calc(moldyn);
2082 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
2086 /* calculate fluctuations + averages */
2087 average_and_fluctuation_calc(moldyn);
2090 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
2091 scale_velocity(moldyn,FALSE);
2092 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
2093 scale_volume(moldyn);
2095 /* check for log & visualization */
2097 if(!(moldyn->total_steps%e))
2098 dprintf(moldyn->efd,
2099 "%f %f %f %f %f %f\n",
2100 moldyn->time,moldyn->ekin/energy_scale,
2101 moldyn->energy/energy_scale,
2102 get_total_energy(moldyn)/energy_scale,
2103 moldyn->ekin/EV,moldyn->energy/EV);
2106 if(!(moldyn->total_steps%m)) {
2107 momentum=get_total_p(moldyn);
2108 dprintf(moldyn->mfd,
2109 "%f %f %f %f %f\n",moldyn->time,
2110 momentum.x,momentum.y,momentum.z,
2111 v3_norm(&momentum));
2115 if(!(moldyn->total_steps%p)) {
2116 dprintf(moldyn->pfd,
2117 "%f %f %f %f %f %f %f\n",moldyn->time,
2118 moldyn->p/BAR,moldyn->p_avg/BAR,
2119 moldyn->p/BAR,moldyn->p_avg/BAR,
2120 moldyn->tp/BAR,moldyn->tp_avg/BAR);
2124 if(!(moldyn->total_steps%t)) {
2125 dprintf(moldyn->tfd,
2127 moldyn->time,moldyn->t,moldyn->t_avg);
2131 if(!(moldyn->total_steps%v)) {
2132 dprintf(moldyn->vfd,
2133 "%f %f %f %f %f\n",moldyn->time,
2141 if(!(moldyn->total_steps%s)) {
2142 snprintf(dir,128,"%s/s-%08.f.save",
2143 moldyn->vlsdir,moldyn->time);
2144 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
2146 if(fd<0) perror("[moldyn] save fd open");
2148 write(fd,moldyn,sizeof(t_moldyn));
2149 write(fd,moldyn->atom,
2150 moldyn->count*sizeof(t_atom));
2156 if(!(moldyn->total_steps%a)) {
2157 #ifdef VISUAL_THREAD
2158 /* check whether thread has not terminated yet */
2160 ret=pthread_join(io_thread,NULL);
2163 /* prepare and start thread */
2164 if(moldyn->count!=md_copy.count) {
2168 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
2170 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2171 if(atom_copy==NULL) {
2172 perror("[moldyn] malloc atom copy (change)");
2176 md_copy.atom=atom_copy;
2177 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
2179 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
2181 perror("[moldyn] create visual atoms thread\n");
2185 visual_atoms(moldyn);
2190 /* display progress */
2194 /* get current time */
2195 gettimeofday(&t2,NULL);
2197 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
2198 sched->count,i,moldyn->total_steps,
2199 moldyn->t,moldyn->t_avg,
2201 moldyn->p/BAR,moldyn->p_avg/BAR,
2203 moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
2206 (int)(t2.tv_sec-t1.tv_sec));
2210 /* copy over time */
2216 /* increase absolute time */
2217 moldyn->time+=moldyn->tau;
2218 moldyn->total_steps+=1;
2222 /* check for hooks */
2224 printf("\n ## schedule hook %d start ##\n",
2226 sched->hook(moldyn,sched->hook_params);
2227 printf(" ## schedule hook end ##\n");
2230 /* increase the schedule counter */
2235 /* writing a final save file! */
2237 snprintf(dir,128,"%s/s-final.save",moldyn->vlsdir);
2238 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,S_IRUSR|S_IWUSR);
2239 if(fd<0) perror("[moldyn] save fd open");
2241 write(fd,moldyn,sizeof(t_moldyn));
2242 write(fd,moldyn->atom,
2243 moldyn->count*sizeof(t_atom));
2247 /* writing a final visual file! */
2249 visual_atoms(moldyn);
2254 /* velocity verlet */
2256 int velocity_verlet(t_moldyn *moldyn) {
2259 double tau,tau_square,h;
2264 count=moldyn->count;
2266 tau_square=moldyn->tau_square;
2268 #ifdef CONSTRAINT_110_5832
2270 atom[5832].f.x=0.5*(atom[5832].f.x-atom[5832].f.y);
2271 atom[5832].f.y=-atom[5832].f.x;
2274 for(i=0;i<count;i++) {
2275 /* check whether fixed atom */
2276 if(atom[i].attr&ATOM_ATTR_FP)
2280 v3_scale(&delta,&(atom[i].v),tau);
2281 #ifdef CONSTRAINT_110_5832
2287 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2289 v3_scale(&delta,&(atom[i].f),h*tau_square);
2290 #ifdef CONSTRAINT_110_5832
2295 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2296 //check_per_bound_and_care_for_pbc(moldyn,&(atom[i]));
2297 check_per_bound(moldyn,&(atom[i].r));
2299 /* velocities [actually v(t+tau/2)] */
2300 v3_scale(&delta,&(atom[i].f),h*tau);
2301 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2304 /* criticial check */
2305 moldyn_bc_check(moldyn);
2307 /* neighbour list update */
2308 link_cell_update(moldyn);
2310 /* forces depending on chosen potential */
2312 // if albe, use albe force calc routine
2313 //if(moldyn->func3b_j1==albe_mult_3bp_j1)
2314 // albe_potential_force_calc(moldyn);
2316 potential_force_calc(moldyn);
2318 albe_potential_force_calc(moldyn);
2321 #ifdef CONSTRAINT_110_5832
2323 atom[5832].f.x=0.5*(atom[5832].f.x-atom[5832].f.y);
2324 atom[5832].f.y=-atom[5832].f.x;
2327 for(i=0;i<count;i++) {
2328 /* check whether fixed atom */
2329 if(atom[i].attr&ATOM_ATTR_FP)
2331 /* again velocities [actually v(t+tau)] */
2332 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2333 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2342 * potentials & corresponding forces & virial routine
2346 /* generic potential and force calculation */
2348 int potential_force_calc(t_moldyn *moldyn) {
2351 t_atom *itom,*jtom,*ktom;
2355 int *neighbour_i[27];
2359 int neighbour_i[27];
2362 t_list neighbour_i[27];
2363 t_list neighbour_i2[27];
2369 count=moldyn->count;
2379 /* reset global virial */
2380 memset(&(moldyn->gvir),0,sizeof(t_virial));
2382 /* reset force, site energy and virial of every atom */
2384 i=omp_get_thread_num();
2385 #pragma omp parallel for private(virial)
2387 for(i=0;i<count;i++) {
2390 v3_zero(&(itom[i].f));
2393 virial=(&(itom[i].virial));
2401 /* reset site energy */
2406 /* get energy, force and virial of every atom */
2408 /* first (and only) loop over atoms i */
2409 for(i=0;i<count;i++) {
2411 /* single particle potential/force */
2412 if(itom[i].attr&ATOM_ATTR_1BP)
2414 moldyn->func1b(moldyn,&(itom[i]));
2416 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2419 /* 2 body pair potential/force */
2421 link_cell_neighbour_index(moldyn,
2422 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2423 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2424 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2429 /* first loop over atoms j */
2430 if(moldyn->func2b) {
2437 while(neighbour_i[j][p]!=-1) {
2439 jtom=&(atom[neighbour_i[j][p]]);
2447 p=lc->subcell->list[p];
2449 this=&(neighbour_i[j]);
2452 if(this->start==NULL)
2456 jtom=this->current->data;
2459 if(jtom==&(itom[i]))
2462 if((jtom->attr&ATOM_ATTR_2BP)&
2463 (itom[i].attr&ATOM_ATTR_2BP)) {
2464 moldyn->func2b(moldyn,
2474 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2480 /* 3 body potential/force */
2482 if(!(itom[i].attr&ATOM_ATTR_3BP))
2485 /* copy the neighbour lists */
2487 /* no copy needed for static lists */
2489 /* no copy needed for lowmem lists */
2491 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2494 /* second loop over atoms j */
2501 while(neighbour_i[j][p]!=-1) {
2503 jtom=&(atom[neighbour_i[j][p]]);
2511 p=lc->subcell->list[p];
2513 this=&(neighbour_i[j]);
2516 if(this->start==NULL)
2521 jtom=this->current->data;
2524 if(jtom==&(itom[i]))
2527 if(!(jtom->attr&ATOM_ATTR_3BP))
2533 if(moldyn->func3b_j1)
2534 moldyn->func3b_j1(moldyn,
2539 /* in first j loop, 3bp run can be skipped */
2540 if(!(moldyn->run3bp))
2543 /* first loop over atoms k */
2544 if(moldyn->func3b_k1) {
2552 while(neighbour_i[k][q]!=-1) {
2554 ktom=&(atom[neighbour_i[k][q]]);
2562 q=lc->subcell->list[q];
2564 that=&(neighbour_i2[k]);
2567 if(that->start==NULL)
2571 ktom=that->current->data;
2574 if(!(ktom->attr&ATOM_ATTR_3BP))
2580 if(ktom==&(itom[i]))
2583 moldyn->func3b_k1(moldyn,
2594 } while(list_next_f(that)!=\
2602 if(moldyn->func3b_j2)
2603 moldyn->func3b_j2(moldyn,
2608 /* second loop over atoms k */
2609 if(moldyn->func3b_k2) {
2617 while(neighbour_i[k][q]!=-1) {
2619 ktom=&(atom[neighbour_i[k][q]]);
2627 q=lc->subcell->list[q];
2629 that=&(neighbour_i2[k]);
2632 if(that->start==NULL)
2636 ktom=that->current->data;
2639 if(!(ktom->attr&ATOM_ATTR_3BP))
2645 if(ktom==&(itom[i]))
2648 moldyn->func3b_k2(moldyn,
2659 } while(list_next_f(that)!=\
2667 /* 2bp post function */
2668 if(moldyn->func3b_j3) {
2669 moldyn->func3b_j3(moldyn,
2678 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2693 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2694 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2696 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2697 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2698 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2702 /* some postprocessing */
2704 #pragma omp parallel for
2706 for(i=0;i<count;i++) {
2707 /* calculate global virial */
2708 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2709 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2710 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2711 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2712 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2713 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2715 /* check forces regarding the given timestep */
2716 if(v3_norm(&(itom[i].f))>\
2717 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2718 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2726 * virial calculation
2729 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2730 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2732 a->virial.xx+=f->x*d->x;
2733 a->virial.yy+=f->y*d->y;
2734 a->virial.zz+=f->z*d->z;
2735 a->virial.xy+=f->x*d->y;
2736 a->virial.xz+=f->x*d->z;
2737 a->virial.yz+=f->y*d->z;
2743 * periodic boundary checking
2746 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2747 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2758 if(moldyn->status&MOLDYN_STAT_PBX) {
2759 if(a->x>=x) a->x-=dim->x;
2760 else if(-a->x>x) a->x+=dim->x;
2762 if(moldyn->status&MOLDYN_STAT_PBY) {
2763 if(a->y>=y) a->y-=dim->y;
2764 else if(-a->y>y) a->y+=dim->y;
2766 if(moldyn->status&MOLDYN_STAT_PBZ) {
2767 if(a->z>=z) a->z-=dim->z;
2768 else if(-a->z>z) a->z+=dim->z;
2774 int check_per_bound_and_care_for_pbc(t_moldyn *moldyn,t_atom *a) {
2785 if(moldyn->status&MOLDYN_STAT_PBX) {
2790 else if(-a->r.x>x) {
2795 if(moldyn->status&MOLDYN_STAT_PBY) {
2800 else if(-a->r.y>y) {
2805 if(moldyn->status&MOLDYN_STAT_PBZ) {
2810 else if(-a->r.z>z) {
2820 * debugging / critical check functions
2823 int moldyn_bc_check(t_moldyn *moldyn) {
2836 for(i=0;i<moldyn->count;i++) {
2837 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2838 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2839 i,atom[i].r.x,dim->x/2);
2840 printf("diagnostic:\n");
2841 printf("-----------\natom.r.x:\n");
2843 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2846 ((byte)&(1<<k))?1:0,
2849 printf("---------------\nx=dim.x/2:\n");
2851 memcpy(&byte,(u8 *)(&x)+j,1);
2854 ((byte)&(1<<k))?1:0,
2857 if(atom[i].r.x==x) printf("the same!\n");
2858 else printf("different!\n");
2860 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2861 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2862 i,atom[i].r.y,dim->y/2);
2863 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2864 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2865 i,atom[i].r.z,dim->z/2);
2875 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2882 fd=open(file,O_RDONLY);
2884 perror("[moldyn] load save file open");
2888 fsize=lseek(fd,0,SEEK_END);
2889 lseek(fd,0,SEEK_SET);
2891 size=sizeof(t_moldyn);
2894 cnt=read(fd,moldyn,size);
2896 perror("[moldyn] load save file read (moldyn)");
2902 size=moldyn->count*sizeof(t_atom);
2904 /* correcting possible atom data offset */
2906 if(fsize!=sizeof(t_moldyn)+size) {
2907 corr=fsize-sizeof(t_moldyn)-size;
2908 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2909 printf(" modifying offset:\n");
2910 printf(" - current pos: %d\n",sizeof(t_moldyn));
2911 printf(" - atom size: %d\n",size);
2912 printf(" - file size: %d\n",fsize);
2913 printf(" => correction: %d\n",corr);
2914 lseek(fd,corr,SEEK_CUR);
2917 moldyn->atom=(t_atom *)malloc(size);
2918 if(moldyn->atom==NULL) {
2919 perror("[moldyn] load save file malloc (atoms)");
2924 cnt=read(fd,moldyn->atom,size);
2926 perror("[moldyn] load save file read (atoms)");
2933 amutex=malloc(moldyn->count*sizeof(pthread_mutex_t));
2935 perror("[moldyn] load save file (mutexes)");
2938 for(cnt=0;cnt<moldyn->count;cnt++)
2939 pthread_mutex_init(&(amutex[cnt]),NULL);
2947 int moldyn_free_save_file(t_moldyn *moldyn) {
2954 int moldyn_load(t_moldyn *moldyn) {
2962 * function to find/callback all combinations of 2 body bonds
2965 int process_2b_bonds(t_moldyn *moldyn,void *data,
2966 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2967 void *data,u8 bc)) {
2977 t_list neighbour[27];
2987 for(i=0;i<moldyn->count;i++) {
2988 /* neighbour indexing */
2989 link_cell_neighbour_index(moldyn,
2990 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2991 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2992 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
2997 bc=(j<lc->dnlc)?0:1;
3002 while(neighbour[j][p]!=-1) {
3004 jtom=&(moldyn->atom[neighbour[j][p]]);
3012 p=lc->subcell->list[p];
3014 this=&(neighbour[j]);
3017 if(this->start==NULL)
3022 jtom=this->current->data;
3026 process(moldyn,&(itom[i]),jtom,data,bc);
3033 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3043 * function to find neighboured atoms
3046 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
3047 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
3048 void *data,u8 bc)) {
3058 t_list neighbour[27];
3067 /* neighbour indexing */
3068 link_cell_neighbour_index(moldyn,
3069 (atom->r.x+moldyn->dim.x/2)/lc->x,
3070 (atom->r.y+moldyn->dim.y/2)/lc->x,
3071 (atom->r.z+moldyn->dim.z/2)/lc->x,
3076 bc=(j<lc->dnlc)?0:1;
3081 while(neighbour[j][p]!=-1) {
3083 natom=&(moldyn->atom[neighbour[j][p]]);
3090 natom=&(moldyn->atom[p]);
3091 p=lc->subcell->list[p];
3093 this=&(neighbour[j]);
3096 if(this->start==NULL)
3101 natom=this->current->data;
3105 process(moldyn,atom,natom,data,bc);
3112 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3121 * post processing functions
3124 int get_line(int fd,char *line,int max) {
3131 if(count==max) return count;
3132 ret=read(fd,line+count,1);
3133 if(ret<=0) return ret;
3134 if(line[count]=='\n') {
3135 memset(line+count,0,max-count-1);
3143 int pair_correlation_init(t_moldyn *moldyn,double dr) {
3149 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
3166 for(i=0;i<moldyn->count;i++) {
3168 /* care for pb crossing */
3169 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3170 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3171 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3172 /* calculate distance */
3173 v3_sub(&dist,&final_r,&(atom[i].r_0));
3174 d2=v3_absolute_square(&dist);
3188 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
3189 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
3190 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
3195 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
3200 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
3201 t_atom *jtom,void *data,u8 bc) {
3208 /* only count pairs once,
3209 * skip same atoms */
3210 if(itom->tag>=jtom->tag)
3214 * pair correlation calc
3221 v3_sub(&dist,&(jtom->r),&(itom->r));
3222 if(bc) check_per_bound(moldyn,&dist);
3223 d=v3_absolute_square(&dist);
3225 /* ignore if greater cutoff */
3226 if(d>moldyn->cutoff_square)
3229 /* fill the slots */
3233 /* should never happen but it does 8) -
3234 * related to -ffloat-store problem! */
3236 printf("[moldyn] WARNING: pcc (%d/%d)",
3242 if(itom->brand!=jtom->brand) {
3247 /* type a - type a bonds */
3249 pcc->stat[s+pcc->o1]+=1;
3251 /* type b - type b bonds */
3252 pcc->stat[s+pcc->o2]+=1;
3258 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
3265 pcc.o1=moldyn->cutoff/dr;
3268 if(pcc.o1*dr<=moldyn->cutoff)
3269 printf("[moldyn] WARNING: pcc (low #slots)\n");
3271 printf("[moldyn] pair correlation calc info:\n");
3272 printf(" time: %f\n",moldyn->time);
3273 printf(" count: %d\n",moldyn->count);
3274 printf(" cutoff: %f\n",moldyn->cutoff);
3275 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
3278 pcc.stat=(double *)ptr;
3281 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
3282 if(pcc.stat==NULL) {
3283 perror("[moldyn] pair correlation malloc");
3288 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
3291 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
3294 for(i=1;i<pcc.o1;i++) {
3295 // normalization: 4 pi r^2 dr
3296 // here: not double counting pairs -> 2 pi r r dr
3297 // ... and actually it's a constant times r^2
3300 pcc.stat[pcc.o1+i]/=norm;
3301 pcc.stat[pcc.o2+i]/=norm;
3306 /* todo: store/print pair correlation function */
3313 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3320 if(itom->tag>=jtom->tag)
3324 v3_sub(&dist,&(jtom->r),&(itom->r));
3325 if(bc) check_per_bound(moldyn,&dist);
3326 d=v3_absolute_square(&dist);
3328 /* ignore if greater or equal cutoff */
3329 if(d>moldyn->cutoff_square)
3332 /* check for potential bond */
3333 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3336 /* now count this bonding ... */
3339 /* increase total bond counter
3344 ba->acnt[jtom->tag]+=1;
3346 ba->bcnt[jtom->tag]+=1;
3349 ba->acnt[itom->tag]+=1;
3351 ba->bcnt[itom->tag]+=1;
3356 int bond_analyze(t_moldyn *moldyn,double *quality) {
3367 ba.acnt=malloc(moldyn->count*sizeof(int));
3369 perror("[moldyn] bond analyze malloc (a)");
3372 memset(ba.acnt,0,moldyn->count*sizeof(int));
3374 ba.bcnt=malloc(moldyn->count*sizeof(int));
3376 perror("[moldyn] bond analyze malloc (b)");
3379 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3388 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3390 for(i=0;i<moldyn->count;i++) {
3391 if(atom[i].brand==0) {
3392 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3394 if((ba.acnt[i]==0)&(ba.bcnt[i]==3))
3398 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3402 if((ba.acnt[i]==3)&(ba.bcnt[i]==0)) {
3411 quality[0]=1.0*ccnt4/bcnt;
3412 quality[1]=1.0*ccnt3/bcnt;
3415 printf("[moldyn] bond analyze: %f %f\n",
3416 1.0*ccnt4/bcnt,1.0*ccnt3/bcnt);
3423 * visualization code
3426 int visual_init(t_moldyn *moldyn,char *filebase) {
3428 strncpy(moldyn->vis.fb,filebase,128);
3433 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3440 if(itom->tag>=jtom->tag)
3443 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3446 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3447 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3448 itom->r.x,itom->r.y,itom->r.z,
3449 jtom->r.x,jtom->r.y,jtom->r.z);
3454 #ifdef VISUAL_THREAD
3455 void *visual_atoms(void *ptr) {
3457 int visual_atoms(t_moldyn *moldyn) {
3468 #ifdef VISUAL_THREAD
3482 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3483 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3485 perror("open visual save file fd");
3486 #ifndef VISUAL_THREAD
3491 /* write the actual data file */
3494 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3495 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3497 // atomic configuration
3498 for(i=0;i<moldyn->count;i++) {
3499 v3_sub(&strain,&(atom[i].r),&(atom[i].r_0));
3500 check_per_bound(moldyn,&strain);
3501 // atom type, positions, color and kinetic energy
3502 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3506 pse_col[atom[i].element],
3508 sqrt(v3_absolute_square(&strain)));
3511 // bonds between atoms
3512 #ifndef VISUAL_THREAD
3513 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3518 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3519 -dim.x/2,-dim.y/2,-dim.z/2,
3520 dim.x/2,-dim.y/2,-dim.z/2);
3521 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3522 -dim.x/2,-dim.y/2,-dim.z/2,
3523 -dim.x/2,dim.y/2,-dim.z/2);
3524 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3525 dim.x/2,dim.y/2,-dim.z/2,
3526 dim.x/2,-dim.y/2,-dim.z/2);
3527 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3528 -dim.x/2,dim.y/2,-dim.z/2,
3529 dim.x/2,dim.y/2,-dim.z/2);
3531 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3532 -dim.x/2,-dim.y/2,dim.z/2,
3533 dim.x/2,-dim.y/2,dim.z/2);
3534 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3535 -dim.x/2,-dim.y/2,dim.z/2,
3536 -dim.x/2,dim.y/2,dim.z/2);
3537 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3538 dim.x/2,dim.y/2,dim.z/2,
3539 dim.x/2,-dim.y/2,dim.z/2);
3540 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3541 -dim.x/2,dim.y/2,dim.z/2,
3542 dim.x/2,dim.y/2,dim.z/2);
3544 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3545 -dim.x/2,-dim.y/2,dim.z/2,
3546 -dim.x/2,-dim.y/2,-dim.z/2);
3547 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3548 -dim.x/2,dim.y/2,dim.z/2,
3549 -dim.x/2,dim.y/2,-dim.z/2);
3550 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3551 dim.x/2,-dim.y/2,dim.z/2,
3552 dim.x/2,-dim.y/2,-dim.z/2);
3553 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3554 dim.x/2,dim.y/2,dim.z/2,
3555 dim.x/2,dim.y/2,-dim.z/2);
3560 #ifdef VISUAL_THREAD
3571 * fpu cntrol functions
3574 // set rounding to double (eliminates -ffloat-store!)
3575 int fpu_set_rtd(void) {
3581 ctrl&=~_FPU_EXTENDED;