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;
58 /* fully constrained relaxation technique - global pointers */
64 * the moldyn functions
67 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
69 printf("[moldyn] init\n");
71 /* only needed if compiled without -msse2 (float-store prob!) */
74 memset(moldyn,0,sizeof(t_moldyn));
79 rand_init(&(moldyn->random),NULL,1);
80 moldyn->random.status|=RAND_STAT_VERBOSE;
83 pthread_mutex_init(&emutex,NULL);
89 int moldyn_shutdown(t_moldyn *moldyn) {
95 printf("[moldyn] shutdown\n");
98 for(i=0;i<moldyn->count;i++)
99 pthread_mutex_destroy(&(amutex[i]));
101 pthread_mutex_destroy(&emutex);
104 moldyn_log_shutdown(moldyn);
105 link_cell_shutdown(moldyn);
106 rand_close(&(moldyn->random));
112 int set_int_alg(t_moldyn *moldyn,u8 algo) {
114 printf("[moldyn] integration algorithm: ");
117 case MOLDYN_INTEGRATE_VERLET:
118 moldyn->integrate=velocity_verlet;
119 printf("velocity verlet\n");
122 printf("unknown integration algorithm: %02x\n",algo);
130 int set_cutoff(t_moldyn *moldyn,double cutoff) {
132 moldyn->cutoff=cutoff;
133 moldyn->cutoff_square=cutoff*cutoff;
135 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
140 int set_temperature(t_moldyn *moldyn,double t_ref) {
144 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
149 int set_pressure(t_moldyn *moldyn,double p_ref) {
153 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
158 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
160 moldyn->pt_scale&=(~(P_SCALE_MASK));
161 moldyn->pt_scale|=ptype;
164 printf("[moldyn] p scaling:\n");
166 printf(" p: %s",ptype?"yes":"no ");
168 printf(" | type: %02x | factor: %f",ptype,ptc);
174 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
176 moldyn->pt_scale&=(~(T_SCALE_MASK));
177 moldyn->pt_scale|=ttype;
180 printf("[moldyn] t scaling:\n");
182 printf(" t: %s",ttype?"yes":"no ");
184 printf(" | type: %02x | factor: %f",ttype,ttc);
190 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
192 moldyn->pt_scale=(ptype|ttype);
196 printf("[moldyn] p/t scaling:\n");
198 printf(" p: %s",ptype?"yes":"no ");
200 printf(" | type: %02x | factor: %f",ptype,ptc);
203 printf(" t: %s",ttype?"yes":"no ");
205 printf(" | type: %02x | factor: %f",ttype,ttc);
211 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
217 moldyn->volume=x*y*z;
225 printf("[moldyn] dimensions in A and A^3 respectively:\n");
226 printf(" x: %f\n",moldyn->dim.x);
227 printf(" y: %f\n",moldyn->dim.y);
228 printf(" z: %f\n",moldyn->dim.z);
229 printf(" volume: %f\n",moldyn->volume);
230 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
235 int set_nn_dist(t_moldyn *moldyn,double dist) {
242 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
244 printf("[moldyn] periodic boundary conditions:\n");
247 moldyn->status|=MOLDYN_STAT_PBX;
250 moldyn->status|=MOLDYN_STAT_PBY;
253 moldyn->status|=MOLDYN_STAT_PBZ;
255 printf(" x: %s\n",x?"yes":"no");
256 printf(" y: %s\n",y?"yes":"no");
257 printf(" z: %s\n",z?"yes":"no");
262 int set_potential(t_moldyn *moldyn,u8 type) {
265 case MOLDYN_POTENTIAL_TM:
266 //moldyn->func1b=tersoff_mult_1bp;
267 moldyn->func3b_j1=tersoff_mult_3bp_j1;
268 moldyn->func3b_k1=tersoff_mult_3bp_k1;
269 moldyn->func3b_j2=tersoff_mult_3bp_j2;
270 moldyn->func3b_k2=tersoff_mult_3bp_k2;
271 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
273 case MOLDYN_POTENTIAL_AM:
274 moldyn->func3b_j1=albe_mult_3bp_j1;
275 moldyn->func3b_k1=albe_mult_3bp_k1;
276 moldyn->func3b_j2=albe_mult_3bp_j2;
277 moldyn->func3b_k2=albe_mult_3bp_k2;
278 moldyn->check_2b_bond=albe_mult_check_2b_bond;
280 case MOLDYN_POTENTIAL_HO:
281 moldyn->func2b=harmonic_oscillator;
282 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
284 case MOLDYN_POTENTIAL_LJ:
285 moldyn->func2b=lennard_jones;
286 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
289 printf("[moldyn] set potential: unknown type %02x\n",
297 int set_avg_skip(t_moldyn *moldyn,int skip) {
299 printf("[moldyn] skip %d steps before starting average calc\n",skip);
300 moldyn->avg_skip=skip;
305 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
307 strncpy(moldyn->vlsdir,dir,127);
312 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
314 strncpy(moldyn->rauthor,author,63);
315 strncpy(moldyn->rtitle,title,63);
320 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
325 printf("[moldyn] set log: ");
328 case LOG_TOTAL_ENERGY:
329 moldyn->ewrite=timer;
330 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
331 moldyn->efd=open(filename,
332 O_WRONLY|O_CREAT|O_EXCL,
335 perror("[moldyn] energy log fd open");
338 dprintf(moldyn->efd,"# total energy log file\n");
339 printf("total energy (%d)\n",timer);
341 case LOG_TOTAL_MOMENTUM:
342 moldyn->mwrite=timer;
343 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
344 moldyn->mfd=open(filename,
345 O_WRONLY|O_CREAT|O_EXCL,
348 perror("[moldyn] momentum log fd open");
351 dprintf(moldyn->efd,"# total momentum log file\n");
352 printf("total momentum (%d)\n",timer);
355 moldyn->pwrite=timer;
356 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
357 moldyn->pfd=open(filename,
358 O_WRONLY|O_CREAT|O_EXCL,
361 perror("[moldyn] pressure log file\n");
364 dprintf(moldyn->pfd,"# pressure log file\n");
365 printf("pressure (%d)\n",timer);
367 case LOG_TEMPERATURE:
368 moldyn->twrite=timer;
369 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
370 moldyn->tfd=open(filename,
371 O_WRONLY|O_CREAT|O_EXCL,
374 perror("[moldyn] temperature log file\n");
377 dprintf(moldyn->tfd,"# temperature log file\n");
378 printf("temperature (%d)\n",timer);
381 moldyn->vwrite=timer;
382 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
383 moldyn->vfd=open(filename,
384 O_WRONLY|O_CREAT|O_EXCL,
387 perror("[moldyn] volume log file\n");
390 dprintf(moldyn->vfd,"# volume log file\n");
391 printf("volume (%d)\n",timer);
394 moldyn->swrite=timer;
395 printf("save file (%d)\n",timer);
398 moldyn->awrite=timer;
399 ret=visual_init(moldyn,moldyn->vlsdir);
401 printf("[moldyn] visual init failure\n");
404 printf("visual file (%d)\n",timer);
407 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
408 moldyn->rfd=open(filename,
409 O_WRONLY|O_CREAT|O_EXCL,
412 perror("[moldyn] report fd open");
415 printf("report -> ");
417 snprintf(filename,127,"%s/e_plot.scr",
419 moldyn->epfd=open(filename,
420 O_WRONLY|O_CREAT|O_EXCL,
423 perror("[moldyn] energy plot fd open");
426 dprintf(moldyn->epfd,e_plot_script);
431 snprintf(filename,127,"%s/pressure_plot.scr",
433 moldyn->ppfd=open(filename,
434 O_WRONLY|O_CREAT|O_EXCL,
437 perror("[moldyn] p plot fd open");
440 dprintf(moldyn->ppfd,pressure_plot_script);
445 snprintf(filename,127,"%s/temperature_plot.scr",
447 moldyn->tpfd=open(filename,
448 O_WRONLY|O_CREAT|O_EXCL,
451 perror("[moldyn] t plot fd open");
454 dprintf(moldyn->tpfd,temperature_plot_script);
456 printf("temperature ");
458 dprintf(moldyn->rfd,report_start,
459 moldyn->rauthor,moldyn->rtitle);
463 printf("unknown log type: %02x\n",type);
470 int moldyn_log_shutdown(t_moldyn *moldyn) {
474 printf("[moldyn] log shutdown\n");
478 dprintf(moldyn->rfd,report_energy);
479 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
484 if(moldyn->mfd) close(moldyn->mfd);
488 dprintf(moldyn->rfd,report_pressure);
489 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
496 dprintf(moldyn->rfd,report_temperature);
497 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
502 dprintf(moldyn->rfd,report_end);
504 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
507 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
510 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
519 * creating lattice functions
522 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,
523 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin,
524 t_part_params *p_params,t_defect_params *d_params,
525 t_offset_params *o_params) {
534 pthread_mutex_t *mutex;
540 /* how many atoms do we expect */
543 printf("[moldyn] WARNING: create 'none' lattice called");
545 if(type==CUBIC) new*=1;
546 if(type==FCC) new*=4;
547 if(type==DIAMOND) new*=8;
551 switch(d_params->stype) {
552 case DEFECT_STYPE_DB_X:
553 case DEFECT_STYPE_DB_Y:
554 case DEFECT_STYPE_DB_Z:
555 case DEFECT_STYPE_DB_R:
559 printf("[moldyn] WARNING: cl unknown defect\n");
564 /* allocate space for atoms */
565 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
567 perror("[moldyn] realloc (create lattice)");
571 atom=&(moldyn->atom[count]);
574 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
576 perror("[moldyn] mutex realloc (add atom)");
580 mutex=&(amutex[count]);
583 /* no atoms on the boundaries (only reason: it looks better!) */
598 v3_add(&orig,&orig,&(o_params->o));
599 set_nn_dist(moldyn,lc);
600 ret=cubic_init(a,b,c,lc,atom,&orig,p_params,d_params);
601 strcpy(name,"cubic");
605 v3_scale(&orig,&orig,0.5);
607 v3_add(&orig,&orig,&(o_params->o));
608 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
609 ret=fcc_init(a,b,c,lc,atom,&orig,p_params,d_params);
614 v3_scale(&orig,&orig,0.25);
616 v3_add(&orig,&orig,&(o_params->o));
617 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
618 ret=diamond_init(a,b,c,lc,atom,&orig,p_params,d_params);
619 strcpy(name,"diamond");
622 printf("unknown lattice type (%02x)\n",type);
628 printf("[moldyn] creating %s lattice (lc=%f) incomplete\n",
630 printf(" (ignore for partial lattice creation)\n");
631 printf(" amount of atoms\n");
632 printf(" - expected: %d\n",new);
633 printf(" - created: %d\n",ret);
638 printf("[moldyn] created %s lattice with %d atoms\n",name,ret);
640 for(new=0;new<ret;new++) {
641 atom[new].element=element;
642 atom[new].mass=pse_mass[element];
644 atom[new].brand=brand;
645 atom[new].tag=count+new;
646 check_per_bound(moldyn,&(atom[new].r));
647 atom[new].r_0=atom[new].r;
649 pthread_mutex_init(&(mutex[new]),NULL);
653 atom[new].element=d_params->element;
654 atom[new].mass=pse_mass[d_params->element];
655 atom[new].attr=d_params->attr;
656 atom[new].brand=d_params->brand;
657 atom[new].tag=count+new;
658 check_per_bound(moldyn,&(atom[new].r));
659 atom[new].r_0=atom[new].r;
661 pthread_mutex_init(&(mutex[new]),NULL);
667 ptr=realloc(moldyn->atom,moldyn->count*sizeof(t_atom));
669 perror("[moldyn] realloc (create lattice - alloc fix)");
674 // WHAT ABOUT AMUTEX !!!!
677 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
679 perror("[moldyn] list realloc (create lattice)");
682 moldyn->lc.subcell->list=ptr;
685 /* update total system mass */
686 total_mass_calc(moldyn);
691 int add_atom(t_moldyn *moldyn,int element,u8 brand,u8 attr,
692 t_3dvec *r,t_3dvec *v) {
699 count=(moldyn->count)++; // asshole style!
701 ptr=realloc(atom,(count+1)*sizeof(t_atom));
703 perror("[moldyn] realloc (add atom)");
709 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
711 perror("[moldyn] list realloc (add atom)");
714 moldyn->lc.subcell->list=ptr;
718 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
720 perror("[moldyn] mutex realloc (add atom)");
724 pthread_mutex_init(&(amutex[count]),NULL);
729 /* initialize new atom */
730 memset(&(atom[count]),0,sizeof(t_atom));
733 atom[count].element=element;
734 atom[count].mass=pse_mass[element];
735 atom[count].brand=brand;
736 atom[count].tag=count;
737 atom[count].attr=attr;
738 check_per_bound(moldyn,&(atom[count].r));
739 atom[count].r_0=atom[count].r;
741 /* update total system mass */
742 total_mass_calc(moldyn);
747 int del_atom(t_moldyn *moldyn,int tag) {
751 #if defined LOWMEM_LISTS || defined PTHREADS
757 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
759 perror("[moldyn]malloc (del atom)");
763 for(cnt=0;cnt<tag;cnt++)
766 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
768 new[cnt-1].tag=cnt-1;
777 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
779 perror("[moldyn] list realloc (del atom)");
782 moldyn->lc.subcell->list=ptr;
784 link_cell_update(moldyn);
788 ptr=realloc(amutex,moldyn->count*sizeof(pthread_mutex_t));
790 perror("[moldyn] mutex realloc (add atom)");
794 pthread_mutex_destroy(&(amutex[moldyn->count+1]));
801 #define set_atom_positions(pos) \
802 if(d_params->type) {\
803 d_o.x=0; d_o.y=0; d_o.z=0;\
804 d_d.x=0; d_d.y=0; d_d.z=0;\
805 switch(d_params->stype) {\
806 case DEFECT_STYPE_DB_X:\
810 case DEFECT_STYPE_DB_Y:\
814 case DEFECT_STYPE_DB_Z:\
818 case DEFECT_STYPE_DB_R:\
821 printf("[moldyn] WARNING: unknown defect\n");\
824 v3_add(&dr,&pos,&d_o);\
825 v3_copy(&(atom[count].r),&dr);\
827 v3_add(&dr,&pos,&d_d);\
828 v3_copy(&(atom[count].r),&dr);\
832 v3_copy(&(atom[count].r),&pos);\
837 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
838 t_part_params *p_params,t_defect_params *d_params) {
858 /* shift partition values */
860 p.x=p_params->p.x+(a*lc)/2.0;
861 p.y=p_params->p.y+(b*lc)/2.0;
862 p.z=p_params->p.z+(c*lc)/2.0;
871 switch(p_params->type) {
874 if(v3_absolute_square(&dist)<
875 (p_params->r*p_params->r)) {
876 set_atom_positions(r);
881 if(v3_absolute_square(&dist)>=
882 (p_params->r*p_params->r)) {
883 set_atom_positions(r);
888 if((fabs(dist.x)<p_params->d.x)&&
889 (fabs(dist.y)<p_params->d.y)&&
890 (fabs(dist.z)<p_params->d.z)) {
891 set_atom_positions(r);
896 if((fabs(dist.x)>=p_params->d.x)||
897 (fabs(dist.y)>=p_params->d.y)||
898 (fabs(dist.z)>=p_params->d.z)) {
899 set_atom_positions(r);
903 set_atom_positions(r);
913 for(i=0;i<count;i++) {
914 atom[i].r.x-=(a*lc)/2.0;
915 atom[i].r.y-=(b*lc)/2.0;
916 atom[i].r.z-=(c*lc)/2.0;
922 /* fcc lattice init */
923 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
924 t_part_params *p_params,t_defect_params *d_params) {
942 /* construct the basis */
943 memset(basis,0,3*sizeof(t_3dvec));
951 /* shift partition values */
953 p.x=p_params->p.x+(a*lc)/2.0;
954 p.y=p_params->p.y+(b*lc)/2.0;
955 p.z=p_params->p.z+(c*lc)/2.0;
958 /* fill up the room */
966 switch(p_params->type) {
969 if(v3_absolute_square(&dist)<
970 (p_params->r*p_params->r)) {
971 set_atom_positions(r);
976 if(v3_absolute_square(&dist)>=
977 (p_params->r*p_params->r)) {
978 set_atom_positions(r);
983 if((fabs(dist.x)<p_params->d.x)&&
984 (fabs(dist.y)<p_params->d.y)&&
985 (fabs(dist.z)<p_params->d.z)) {
986 set_atom_positions(r);
991 if((fabs(dist.x)>=p_params->d.x)||
992 (fabs(dist.y)>=p_params->d.y)||
993 (fabs(dist.z)>=p_params->d.z)) {
994 set_atom_positions(r);
998 set_atom_positions(r);
1001 /* the three face centered atoms */
1003 v3_add(&n,&r,&basis[l]);
1004 switch(p_params->type) {
1006 v3_sub(&dist,&n,&p);
1007 if(v3_absolute_square(&dist)<
1008 (p_params->r*p_params->r)) {
1009 set_atom_positions(n);
1012 case PART_OUTSIDE_R:
1013 v3_sub(&dist,&n,&p);
1014 if(v3_absolute_square(&dist)>=
1015 (p_params->r*p_params->r)) {
1016 set_atom_positions(n);
1020 v3_sub(&dist,&n,&p);
1021 if((fabs(dist.x)<p_params->d.x)&&
1022 (fabs(dist.y)<p_params->d.y)&&
1023 (fabs(dist.z)<p_params->d.z)) {
1024 set_atom_positions(n);
1027 case PART_OUTSIDE_D:
1028 v3_sub(&dist,&n,&p);
1029 if((fabs(dist.x)>=p_params->d.x)||
1030 (fabs(dist.y)>=p_params->d.y)||
1031 (fabs(dist.z)>=p_params->d.z)) {
1032 set_atom_positions(n);
1036 set_atom_positions(n);
1047 /* coordinate transformation */
1048 for(i=0;i<count;i++) {
1049 atom[i].r.x-=(a*lc)/2.0;
1050 atom[i].r.y-=(b*lc)/2.0;
1051 atom[i].r.z-=(c*lc)/2.0;
1057 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
1058 t_part_params *p_params,t_defect_params *d_params) {
1063 count=fcc_init(a,b,c,lc,atom,origin,p_params,d_params);
1069 if(origin) v3_add(&o,&o,origin);
1071 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_params,d_params);
1076 int destroy_atoms(t_moldyn *moldyn) {
1078 if(moldyn->atom) free(moldyn->atom);
1083 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
1086 * - gaussian distribution of velocities
1087 * - zero total momentum
1088 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1093 t_3dvec p_total,delta;
1098 random=&(moldyn->random);
1100 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
1102 /* gaussian distribution of velocities */
1104 for(i=0;i<moldyn->count;i++) {
1105 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
1107 v=sigma*rand_get_gauss(random);
1109 p_total.x+=atom[i].mass*v;
1111 v=sigma*rand_get_gauss(random);
1113 p_total.y+=atom[i].mass*v;
1115 v=sigma*rand_get_gauss(random);
1117 p_total.z+=atom[i].mass*v;
1120 /* zero total momentum */
1121 v3_scale(&p_total,&p_total,1.0/moldyn->count);
1122 for(i=0;i<moldyn->count;i++) {
1123 v3_scale(&delta,&p_total,1.0/atom[i].mass);
1124 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
1127 /* velocity scaling */
1128 scale_velocity(moldyn,equi_init);
1133 double total_mass_calc(t_moldyn *moldyn) {
1139 for(i=0;i<moldyn->count;i++)
1140 moldyn->mass+=moldyn->atom[i].mass;
1142 return moldyn->mass;
1145 double temperature_calc(t_moldyn *moldyn) {
1147 /* assume up to date kinetic energy, which is 3/2 N k_B T */
1150 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
1156 double get_temperature(t_moldyn *moldyn) {
1161 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
1171 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1174 /* get kinetic energy / temperature & count involved atoms */
1177 for(i=0;i<moldyn->count;i++) {
1178 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
1179 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1184 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1185 else return 0; /* no atoms involved in scaling! */
1187 /* (temporary) hack for e,t = 0 */
1190 if(moldyn->t_ref!=0.0) {
1191 thermal_init(moldyn,equi_init);
1195 return 0; /* no scaling needed */
1199 /* get scaling factor */
1200 scale=moldyn->t_ref/moldyn->t;
1204 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1205 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1208 /* velocity scaling */
1209 for(i=0;i<moldyn->count;i++) {
1210 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1211 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1217 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1221 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1226 double virial_sum(t_moldyn *moldyn) {
1231 /* virial (sum over atom virials) */
1239 for(i=0;i<moldyn->count;i++) {
1240 virial=&(moldyn->atom[i].virial);
1241 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1242 moldyn->vir.xx+=virial->xx;
1243 moldyn->vir.yy+=virial->yy;
1244 moldyn->vir.zz+=virial->zz;
1245 moldyn->vir.xy+=virial->xy;
1246 moldyn->vir.xz+=virial->xz;
1247 moldyn->vir.yz+=virial->yz;
1250 /* global virial (absolute coordinates) */
1251 //virial=&(moldyn->gvir);
1252 //moldyn->gv=virial->xx+virial->yy+virial->zz;
1254 return moldyn->virial;
1257 double pressure_calc(t_moldyn *moldyn) {
1261 * with W = 1/3 sum_i f_i r_i (- skipped!)
1262 * virial = sum_i f_i r_i
1264 * => P = (2 Ekin + virial) / (3V)
1267 /* assume up to date virial & up to date kinetic energy */
1269 /* pressure (atom virials) */
1270 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1271 moldyn->p/=(3.0*moldyn->volume);
1273 //moldyn->px=2.0*moldyn->ekinx+moldyn->vir.xx;
1274 //moldyn->px/=moldyn->volume;
1275 //moldyn->py=2.0*moldyn->ekiny+moldyn->vir.yy;
1276 //moldyn->py/=moldyn->volume;
1277 //moldyn->pz=2.0*moldyn->ekinz+moldyn->vir.zz;
1278 //moldyn->pz/=moldyn->volume;
1280 /* pressure (absolute coordinates) */
1281 //moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1282 //moldyn->gp/=(3.0*moldyn->volume);
1287 int average_reset(t_moldyn *moldyn) {
1289 printf("[moldyn] average reset\n");
1291 /* update skip value */
1292 moldyn->avg_skip=moldyn->total_steps;
1294 /* kinetic energy */
1298 /* potential energy */
1306 moldyn->virial_sum=0.0;
1307 //moldyn->gv_sum=0.0;
1311 //moldyn->gp_sum=0.0;
1317 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1321 if(moldyn->total_steps<moldyn->avg_skip)
1324 denom=moldyn->total_steps+1-moldyn->avg_skip;
1326 /* assume up to date energies, temperature, pressure etc */
1328 /* kinetic energy */
1329 moldyn->k_sum+=moldyn->ekin;
1330 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1331 moldyn->k_avg=moldyn->k_sum/denom;
1332 moldyn->k2_avg=moldyn->k2_sum/denom;
1333 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1335 /* potential energy */
1336 moldyn->v_sum+=moldyn->energy;
1337 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1338 moldyn->v_avg=moldyn->v_sum/denom;
1339 moldyn->v2_avg=moldyn->v2_sum/denom;
1340 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1343 moldyn->t_sum+=moldyn->t;
1344 moldyn->t_avg=moldyn->t_sum/denom;
1347 moldyn->virial_sum+=moldyn->virial;
1348 moldyn->virial_avg=moldyn->virial_sum/denom;
1349 //moldyn->gv_sum+=moldyn->gv;
1350 //moldyn->gv_avg=moldyn->gv_sum/denom;
1353 moldyn->p_sum+=moldyn->p;
1354 moldyn->p_avg=moldyn->p_sum/denom;
1355 //moldyn->gp_sum+=moldyn->gp;
1356 //moldyn->gp_avg=moldyn->gp_sum/denom;
1357 moldyn->tp_sum+=moldyn->tp;
1358 moldyn->tp_avg=moldyn->tp_sum/denom;
1363 int get_heat_capacity(t_moldyn *moldyn) {
1367 /* averages needed for heat capacity calc */
1368 if(moldyn->total_steps<moldyn->avg_skip)
1371 /* (temperature average)^2 */
1372 temp2=moldyn->t_avg*moldyn->t_avg;
1373 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1376 /* ideal gas contribution */
1377 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1378 printf(" ideal gas contribution: %f\n",
1379 ighc/moldyn->mass*KILOGRAM/JOULE);
1381 /* specific heat for nvt ensemble */
1382 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1383 moldyn->c_v_nvt/=moldyn->mass;
1385 /* specific heat for nve ensemble */
1386 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1387 moldyn->c_v_nve/=moldyn->mass;
1389 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1390 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1391 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)));
1396 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1412 /* store atomic configuration + dimension */
1413 store=malloc(moldyn->count*sizeof(t_atom));
1415 printf("[moldyn] allocating store mem failed\n");
1418 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1423 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1424 su=pow(2.0-h,ONE_THIRD)-1.0;
1425 dv=(1.0-h)*moldyn->volume;
1427 /* scale up dimension and atom positions */
1428 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1429 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1430 link_cell_shutdown(moldyn);
1431 link_cell_init(moldyn,QUIET);
1432 potential_force_calc(moldyn);
1435 /* restore atomic configuration + dim */
1436 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1439 /* scale down dimension and atom positions */
1440 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1441 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1442 link_cell_shutdown(moldyn);
1443 link_cell_init(moldyn,QUIET);
1444 potential_force_calc(moldyn);
1447 /* calculate pressure */
1448 moldyn->tp=-(y1-y0)/(2.0*dv);
1450 /* restore atomic configuration */
1451 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1453 link_cell_shutdown(moldyn);
1454 link_cell_init(moldyn,QUIET);
1455 //potential_force_calc(moldyn);
1457 /* free store buffer */
1464 double get_pressure(t_moldyn *moldyn) {
1470 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1482 if(x) dim->x*=scale;
1483 if(y) dim->y*=scale;
1484 if(z) dim->z*=scale;
1489 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1500 for(i=0;i<moldyn->count;i++) {
1501 r=&(moldyn->atom[i].r);
1510 int scale_atoms_ind(t_moldyn *moldyn,double x,double y,double z) {
1515 for(i=0;i<moldyn->count;i++) {
1516 r=&(moldyn->atom[i].r);
1525 int scale_dim_ind(t_moldyn *moldyn,double x,double y,double z) {
1538 int scale_volume(t_moldyn *moldyn) {
1545 vdim=&(moldyn->vis.dim);
1549 /* scaling factor */
1550 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1551 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1552 scale=pow(scale,ONE_THIRD);
1555 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1560 sx=1.0-(moldyn->p_ref-moldyn->px)*moldyn->p_tc*moldyn->tau;
1561 sy=1.0-(moldyn->p_ref-moldyn->py)*moldyn->p_tc*moldyn->tau;
1562 sz=1.0-(moldyn->p_ref-moldyn->pz)*moldyn->p_tc*moldyn->tau;
1563 sx=pow(sx,ONE_THIRD);
1564 sy=pow(sy,ONE_THIRD);
1565 sz=pow(sz,ONE_THIRD);
1568 /* scale the atoms and dimensions */
1569 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1570 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1571 //scale_atoms_ind(moldyn,sx,sy,sz);
1572 //scale_dim_ind(moldyn,sx,sy,sz);
1574 /* visualize dimensions */
1581 /* recalculate scaled volume */
1582 moldyn->volume=dim->x*dim->y*dim->z;
1584 /* adjust/reinit linkcell */
1585 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1586 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1587 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1588 link_cell_shutdown(moldyn);
1589 link_cell_init(moldyn,QUIET);
1603 double e_kin_calc(t_moldyn *moldyn) {
1610 //moldyn->ekinx=0.0;
1611 //moldyn->ekiny=0.0;
1612 //moldyn->ekinz=0.0;
1614 for(i=0;i<moldyn->count;i++) {
1615 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1616 moldyn->ekin+=atom[i].ekin;
1617 //moldyn->ekinx+=0.5*atom[i].mass*atom[i].v.x*atom[i].v.x;
1618 //moldyn->ekiny+=0.5*atom[i].mass*atom[i].v.y*atom[i].v.y;
1619 //moldyn->ekinz+=0.5*atom[i].mass*atom[i].v.z*atom[i].v.z;
1622 return moldyn->ekin;
1625 double get_total_energy(t_moldyn *moldyn) {
1627 return(moldyn->ekin+moldyn->energy);
1630 t_3dvec get_total_p(t_moldyn *moldyn) {
1639 for(i=0;i<moldyn->count;i++) {
1640 v3_scale(&p,&(atom[i].v),atom[i].mass);
1641 v3_add(&p_total,&p_total,&p);
1647 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1651 /* nn_dist is the nearest neighbour distance */
1653 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1662 /* linked list / cell method */
1664 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1667 #ifndef LOWMEM_LISTS
1673 /* partitioning the md cell */
1674 lc->nx=moldyn->dim.x/moldyn->cutoff;
1675 lc->x=moldyn->dim.x/lc->nx;
1676 lc->ny=moldyn->dim.y/moldyn->cutoff;
1677 lc->y=moldyn->dim.y/lc->ny;
1678 lc->nz=moldyn->dim.z/moldyn->cutoff;
1679 lc->z=moldyn->dim.z/lc->nz;
1680 lc->cells=lc->nx*lc->ny*lc->nz;
1683 lc->subcell=malloc(lc->cells*sizeof(int*));
1685 lc->subcell=malloc(sizeof(t_lowmem_list));
1687 lc->subcell=malloc(lc->cells*sizeof(t_list));
1690 if(lc->subcell==NULL) {
1691 perror("[moldyn] cell init (malloc)");
1696 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1701 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1704 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1707 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1710 printf(" x: %d x %f A\n",lc->nx,lc->x);
1711 printf(" y: %d x %f A\n",lc->ny,lc->y);
1712 printf(" z: %d x %f A\n",lc->nz,lc->z);
1717 for(i=0;i<lc->cells;i++) {
1718 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1719 if(lc->subcell[i]==NULL) {
1720 perror("[moldyn] list init (malloc)");
1725 printf(" ---> %d malloc %p (%p)\n",
1726 i,lc->subcell[0],lc->subcell);
1730 lc->subcell->head=malloc(lc->cells*sizeof(int));
1731 if(lc->subcell->head==NULL) {
1732 perror("[moldyn] head init (malloc)");
1735 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1736 if(lc->subcell->list==NULL) {
1737 perror("[moldyn] list init (malloc)");
1741 for(i=0;i<lc->cells;i++)
1742 list_init_f(&(lc->subcell[i]));
1745 /* update the list */
1746 link_cell_update(moldyn);
1751 int link_cell_update(t_moldyn *moldyn) {
1769 for(i=0;i<lc->cells;i++)
1771 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1773 lc->subcell->head[i]=-1;
1775 list_destroy_f(&(lc->subcell[i]));
1778 for(count=0;count<moldyn->count;count++) {
1779 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1780 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1781 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1785 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1788 if(p>=MAX_ATOMS_PER_LIST) {
1789 printf("[moldyn] FATAL: amount of atoms too high!\n");
1793 lc->subcell[i+j*nx+k*nxy][p]=count;
1796 lc->subcell->list[count]=lc->subcell->head[p];
1797 lc->subcell->head[p]=count;
1799 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1803 printf(" ---> %d %d malloc %p (%p)\n",
1804 i,count,lc->subcell[i].current,lc->subcell);
1812 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1838 if(i>=nx||j>=ny||k>=nz)
1839 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1842 #ifndef LOWMEM_LISTS
1843 cell[0]=lc->subcell[i+j*nx+k*a];
1845 cell[0]=lc->subcell->head[i+j*nx+k*a];
1847 for(ci=-1;ci<=1;ci++) {
1850 if((x<0)||(x>=nx)) {
1854 for(cj=-1;cj<=1;cj++) {
1857 if((y<0)||(y>=ny)) {
1861 for(ck=-1;ck<=1;ck++) {
1864 if((z<0)||(z>=nz)) {
1868 if(!(ci|cj|ck)) continue;
1870 #ifndef LOWMEM_LISTS
1871 cell[--count2]=lc->subcell[x+y*nx+z*a];
1873 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1878 #ifndef LOWMEM_LISTS
1879 cell[count1++]=lc->subcell[x+y*nx+z*a];
1881 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1893 int link_cell_shutdown(t_moldyn *moldyn) {
1895 #ifndef LOWMEM_LISTS
1903 free(lc->subcell->head);
1904 free(lc->subcell->list);
1908 for(i=0;i<lc->cells;i++) {
1910 free(lc->subcell[i]);
1912 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1913 list_destroy_f(&(lc->subcell[i]));
1923 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1927 t_moldyn_schedule *schedule;
1929 schedule=&(moldyn->schedule);
1930 count=++(schedule->total_sched);
1932 ptr=realloc(schedule->runs,count*sizeof(int));
1934 perror("[moldyn] realloc (runs)");
1938 schedule->runs[count-1]=runs;
1940 ptr=realloc(schedule->tau,count*sizeof(double));
1942 perror("[moldyn] realloc (tau)");
1946 schedule->tau[count-1]=tau;
1948 printf("[moldyn] schedule added:\n");
1949 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1955 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1957 moldyn->schedule.hook=hook;
1958 moldyn->schedule.hook_params=hook_params;
1965 * 'integration of newtons equation' - algorithms
1969 /* start the integration */
1971 int moldyn_integrate(t_moldyn *moldyn) {
1974 unsigned int e,m,s,v,p,t,a;
1976 t_moldyn_schedule *sched;
1981 double energy_scale;
1982 struct timeval t1,t2;
1985 #ifdef VISUAL_THREAD
1987 pthread_t io_thread;
1996 sched=&(moldyn->schedule);
1999 /* initialize linked cell method */
2000 link_cell_init(moldyn,VERBOSE);
2002 /* logging & visualization */
2011 /* sqaure of some variables */
2012 moldyn->tau_square=moldyn->tau*moldyn->tau;
2014 /* get current time */
2015 gettimeofday(&t1,NULL);
2017 /* calculate initial forces */
2018 potential_force_calc(moldyn);
2023 /* some stupid checks before we actually start calculating bullshit */
2024 if(moldyn->cutoff>0.5*moldyn->dim.x)
2025 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
2026 if(moldyn->cutoff>0.5*moldyn->dim.y)
2027 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
2028 if(moldyn->cutoff>0.5*moldyn->dim.z)
2029 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
2031 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
2032 if(ds>0.05*moldyn->nnd)
2033 printf("[moldyn] WARNING: forces too high / tau too small!\n");
2036 /* zero absolute time */
2037 // should have right values!
2039 //moldyn->total_steps=0;
2041 /* debugging, ignore */
2044 /* zero & init moldyn copy */
2045 #ifdef VISUAL_THREAD
2046 memset(&md_copy,0,sizeof(t_moldyn));
2047 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2048 if(atom_copy==NULL) {
2049 perror("[moldyn] malloc atom copy (init)");
2055 printf("##################\n");
2056 printf("# USING PTHREADS #\n");
2057 printf("##################\n");
2059 /* tell the world */
2060 printf("[moldyn] integration start, go get a coffee ...\n");
2062 /* executing the schedule */
2064 while(sched->count<sched->total_sched) {
2066 /* setting amount of runs and finite time step size */
2067 moldyn->tau=sched->tau[sched->count];
2068 moldyn->tau_square=moldyn->tau*moldyn->tau;
2069 moldyn->time_steps=sched->runs[sched->count];
2071 /* energy scaling factor (might change!) */
2072 energy_scale=moldyn->count*EV;
2074 /* integration according to schedule */
2076 for(i=0;i<moldyn->time_steps;i++) {
2078 /* integration step */
2079 moldyn->integrate(moldyn);
2081 /* calculate kinetic energy, temperature and pressure */
2083 temperature_calc(moldyn);
2085 pressure_calc(moldyn);
2087 thermodynamic_pressure_calc(moldyn);
2088 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
2092 /* calculate fluctuations + averages */
2093 average_and_fluctuation_calc(moldyn);
2096 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
2097 scale_velocity(moldyn,FALSE);
2098 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
2099 scale_volume(moldyn);
2101 /* check for log & visualization */
2103 if(!(moldyn->total_steps%e))
2104 dprintf(moldyn->efd,
2105 "%f %f %f %f %f %f\n",
2106 moldyn->time,moldyn->ekin/energy_scale,
2107 moldyn->energy/energy_scale,
2108 get_total_energy(moldyn)/energy_scale,
2109 moldyn->ekin/EV,moldyn->energy/EV);
2112 if(!(moldyn->total_steps%m)) {
2113 momentum=get_total_p(moldyn);
2114 dprintf(moldyn->mfd,
2115 "%f %f %f %f %f\n",moldyn->time,
2116 momentum.x,momentum.y,momentum.z,
2117 v3_norm(&momentum));
2121 if(!(moldyn->total_steps%p)) {
2122 dprintf(moldyn->pfd,
2123 "%f %f %f %f %f %f %f\n",moldyn->time,
2124 moldyn->p/BAR,moldyn->p_avg/BAR,
2125 moldyn->p/BAR,moldyn->p_avg/BAR,
2126 moldyn->tp/BAR,moldyn->tp_avg/BAR);
2130 if(!(moldyn->total_steps%t)) {
2131 dprintf(moldyn->tfd,
2133 moldyn->time,moldyn->t,moldyn->t_avg);
2137 if(!(moldyn->total_steps%v)) {
2138 dprintf(moldyn->vfd,
2139 "%f %f %f %f %f\n",moldyn->time,
2147 if(!(moldyn->total_steps%s)) {
2148 snprintf(dir,128,"%s/s-%08.f.save",
2149 moldyn->vlsdir,moldyn->time);
2150 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
2152 if(fd<0) perror("[moldyn] save fd open");
2154 write(fd,moldyn,sizeof(t_moldyn));
2155 write(fd,moldyn->atom,
2156 moldyn->count*sizeof(t_atom));
2162 if(!(moldyn->total_steps%a)) {
2163 #ifdef VISUAL_THREAD
2164 /* check whether thread has not terminated yet */
2166 ret=pthread_join(io_thread,NULL);
2169 /* prepare and start thread */
2170 if(moldyn->count!=md_copy.count) {
2174 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
2176 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2177 if(atom_copy==NULL) {
2178 perror("[moldyn] malloc atom copy (change)");
2182 md_copy.atom=atom_copy;
2183 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
2185 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
2187 perror("[moldyn] create visual atoms thread\n");
2191 visual_atoms(moldyn);
2196 /* display progress */
2200 /* get current time */
2201 gettimeofday(&t2,NULL);
2203 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
2204 sched->count,i,moldyn->total_steps,
2205 moldyn->t,moldyn->t_avg,
2207 moldyn->p/BAR,moldyn->p_avg/BAR,
2209 moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
2212 (int)(t2.tv_sec-t1.tv_sec));
2216 /* copy over time */
2222 /* increase absolute time */
2223 moldyn->time+=moldyn->tau;
2224 moldyn->total_steps+=1;
2228 /* check for hooks */
2230 printf("\n ## schedule hook %d start ##\n",
2232 sched->hook(moldyn,sched->hook_params);
2233 printf(" ## schedule hook end ##\n");
2236 /* increase the schedule counter */
2241 /* writing a final save file! */
2243 snprintf(dir,128,"%s/s-final.save",moldyn->vlsdir);
2244 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,S_IRUSR|S_IWUSR);
2245 if(fd<0) perror("[moldyn] save fd open");
2247 write(fd,moldyn,sizeof(t_moldyn));
2248 write(fd,moldyn->atom,
2249 moldyn->count*sizeof(t_atom));
2253 /* writing a final visual file! */
2255 visual_atoms(moldyn);
2265 int basis_trafo(t_3dvec *r,u8 dir,double z,double x) {
2272 r->x=cos(z)*tmp.x-sin(z)*tmp.y;
2273 r->y=sin(z)*tmp.x+cos(z)*tmp.y;
2277 r->y=cos(x)*tmp.y-sin(x)*tmp.z;
2278 r->z=sin(x)*tmp.y+cos(x)*tmp.z;
2284 r->y=cos(-x)*tmp.y-sin(-x)*tmp.z;
2285 r->z=sin(-x)*tmp.y+cos(-x)*tmp.z;
2289 r->x=cos(-z)*tmp.x-sin(-z)*tmp.y;
2290 r->y=sin(-z)*tmp.x+cos(-z)*tmp.y;
2297 /* velocity verlet */
2299 int velocity_verlet(t_moldyn *moldyn) {
2302 double tau,tau_square,h;
2307 count=moldyn->count;
2309 tau_square=moldyn->tau_square;
2311 for(i=0;i<count;i++) {
2313 /* check whether fixed atom */
2314 if(atom[i].attr&ATOM_ATTR_FP)
2320 /* constraint relaxation */
2323 basis_trafo(&(atom[i].f),FORWARD,
2324 trafo_angle[2*i],trafo_angle[2*i+1]);
2325 if(constraints[3*i])
2327 if(constraints[3*i+1])
2329 if(constraints[3*i+2])
2331 basis_trafo(&(atom[i].f),BACKWARD,
2332 trafo_angle[2*i],trafo_angle[2*i+1]);
2334 basis_trafo(&(atom[i].v),FORWARD,
2335 trafo_angle[2*i],trafo_angle[2*i+1]);
2336 if(constraints[3*i])
2338 if(constraints[3*i+1])
2340 if(constraints[3*i+2])
2342 basis_trafo(&(atom[i].v),BACKWARD,
2343 trafo_angle[2*i],trafo_angle[2*i+1]);
2347 v3_scale(&delta,&(atom[i].v),tau);
2348 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2350 v3_scale(&delta,&(atom[i].f),h*tau_square);
2351 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2352 //check_per_bound_and_care_for_pbc(moldyn,&(atom[i]));
2353 check_per_bound(moldyn,&(atom[i].r));
2355 /* velocities [actually v(t+tau/2)] */
2356 v3_scale(&delta,&(atom[i].f),h*tau);
2357 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2360 /* criticial check */
2361 moldyn_bc_check(moldyn);
2363 /* neighbour list update */
2364 link_cell_update(moldyn);
2366 /* forces depending on chosen potential */
2368 // if albe, use albe force calc routine
2369 //if(moldyn->func3b_j1==albe_mult_3bp_j1)
2370 // albe_potential_force_calc(moldyn);
2372 potential_force_calc(moldyn);
2374 albe_potential_force_calc(moldyn);
2377 for(i=0;i<count;i++) {
2378 /* check whether fixed atom */
2379 if(atom[i].attr&ATOM_ATTR_FP)
2382 /* constraint relaxation */
2385 basis_trafo(&(atom[i].f),FORWARD,
2386 trafo_angle[2*i],trafo_angle[2*i+1]);
2387 if(constraints[3*i])
2389 if(constraints[3*i+1])
2391 if(constraints[3*i+2])
2393 basis_trafo(&(atom[i].f),BACKWARD,
2394 trafo_angle[2*i],trafo_angle[2*i+1]);
2396 basis_trafo(&(atom[i].v),FORWARD,
2397 trafo_angle[2*i],trafo_angle[2*i+1]);
2398 if(constraints[3*i])
2400 if(constraints[3*i+1])
2402 if(constraints[3*i+2])
2404 basis_trafo(&(atom[i].v),BACKWARD,
2405 trafo_angle[2*i],trafo_angle[2*i+1]);
2408 /* again velocities [actually v(t+tau)] */
2409 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2410 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2419 * potentials & corresponding forces & virial routine
2423 /* generic potential and force calculation */
2425 int potential_force_calc(t_moldyn *moldyn) {
2428 t_atom *itom,*jtom,*ktom;
2432 int *neighbour_i[27];
2436 int neighbour_i[27];
2439 t_list neighbour_i[27];
2440 t_list neighbour_i2[27];
2446 count=moldyn->count;
2456 /* reset global virial */
2457 memset(&(moldyn->gvir),0,sizeof(t_virial));
2459 /* reset force, site energy and virial of every atom */
2461 i=omp_get_thread_num();
2462 #pragma omp parallel for private(virial)
2464 for(i=0;i<count;i++) {
2467 v3_zero(&(itom[i].f));
2470 virial=(&(itom[i].virial));
2478 /* reset site energy */
2483 /* get energy, force and virial of every atom */
2485 /* first (and only) loop over atoms i */
2486 for(i=0;i<count;i++) {
2488 /* single particle potential/force */
2489 if(itom[i].attr&ATOM_ATTR_1BP)
2491 moldyn->func1b(moldyn,&(itom[i]));
2493 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2496 /* 2 body pair potential/force */
2498 link_cell_neighbour_index(moldyn,
2499 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2500 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2501 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2506 /* first loop over atoms j */
2507 if(moldyn->func2b) {
2514 while(neighbour_i[j][p]!=-1) {
2516 jtom=&(atom[neighbour_i[j][p]]);
2524 p=lc->subcell->list[p];
2526 this=&(neighbour_i[j]);
2529 if(this->start==NULL)
2533 jtom=this->current->data;
2536 if(jtom==&(itom[i]))
2539 if((jtom->attr&ATOM_ATTR_2BP)&
2540 (itom[i].attr&ATOM_ATTR_2BP)) {
2541 moldyn->func2b(moldyn,
2551 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2557 /* 3 body potential/force */
2559 if(!(itom[i].attr&ATOM_ATTR_3BP))
2562 /* copy the neighbour lists */
2564 /* no copy needed for static lists */
2566 /* no copy needed for lowmem lists */
2568 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2571 /* second loop over atoms j */
2578 while(neighbour_i[j][p]!=-1) {
2580 jtom=&(atom[neighbour_i[j][p]]);
2588 p=lc->subcell->list[p];
2590 this=&(neighbour_i[j]);
2593 if(this->start==NULL)
2598 jtom=this->current->data;
2601 if(jtom==&(itom[i]))
2604 if(!(jtom->attr&ATOM_ATTR_3BP))
2610 if(moldyn->func3b_j1)
2611 moldyn->func3b_j1(moldyn,
2616 /* in first j loop, 3bp run can be skipped */
2617 if(!(moldyn->run3bp))
2620 /* first loop over atoms k */
2621 if(moldyn->func3b_k1) {
2629 while(neighbour_i[k][q]!=-1) {
2631 ktom=&(atom[neighbour_i[k][q]]);
2639 q=lc->subcell->list[q];
2641 that=&(neighbour_i2[k]);
2644 if(that->start==NULL)
2648 ktom=that->current->data;
2651 if(!(ktom->attr&ATOM_ATTR_3BP))
2657 if(ktom==&(itom[i]))
2660 moldyn->func3b_k1(moldyn,
2671 } while(list_next_f(that)!=\
2679 if(moldyn->func3b_j2)
2680 moldyn->func3b_j2(moldyn,
2685 /* second loop over atoms k */
2686 if(moldyn->func3b_k2) {
2694 while(neighbour_i[k][q]!=-1) {
2696 ktom=&(atom[neighbour_i[k][q]]);
2704 q=lc->subcell->list[q];
2706 that=&(neighbour_i2[k]);
2709 if(that->start==NULL)
2713 ktom=that->current->data;
2716 if(!(ktom->attr&ATOM_ATTR_3BP))
2722 if(ktom==&(itom[i]))
2725 moldyn->func3b_k2(moldyn,
2736 } while(list_next_f(that)!=\
2744 /* 2bp post function */
2745 if(moldyn->func3b_j3) {
2746 moldyn->func3b_j3(moldyn,
2755 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2770 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2771 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2773 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2774 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2775 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2779 /* some postprocessing */
2781 #pragma omp parallel for
2783 for(i=0;i<count;i++) {
2784 /* calculate global virial */
2785 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2786 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2787 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2788 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2789 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2790 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2792 /* check forces regarding the given timestep */
2793 if(v3_norm(&(itom[i].f))>\
2794 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2795 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2803 * virial calculation
2806 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2807 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2809 a->virial.xx+=f->x*d->x;
2810 a->virial.yy+=f->y*d->y;
2811 a->virial.zz+=f->z*d->z;
2812 a->virial.xy+=f->x*d->y;
2813 a->virial.xz+=f->x*d->z;
2814 a->virial.yz+=f->y*d->z;
2820 * periodic boundary checking
2823 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2824 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2835 if(moldyn->status&MOLDYN_STAT_PBX) {
2836 if(a->x>=x) a->x-=dim->x;
2837 else if(-a->x>x) a->x+=dim->x;
2839 if(moldyn->status&MOLDYN_STAT_PBY) {
2840 if(a->y>=y) a->y-=dim->y;
2841 else if(-a->y>y) a->y+=dim->y;
2843 if(moldyn->status&MOLDYN_STAT_PBZ) {
2844 if(a->z>=z) a->z-=dim->z;
2845 else if(-a->z>z) a->z+=dim->z;
2851 int check_per_bound_and_care_for_pbc(t_moldyn *moldyn,t_atom *a) {
2862 if(moldyn->status&MOLDYN_STAT_PBX) {
2867 else if(-a->r.x>x) {
2872 if(moldyn->status&MOLDYN_STAT_PBY) {
2877 else if(-a->r.y>y) {
2882 if(moldyn->status&MOLDYN_STAT_PBZ) {
2887 else if(-a->r.z>z) {
2897 * debugging / critical check functions
2900 int moldyn_bc_check(t_moldyn *moldyn) {
2913 for(i=0;i<moldyn->count;i++) {
2914 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2915 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2916 i,atom[i].r.x,dim->x/2);
2917 printf("diagnostic:\n");
2918 printf("-----------\natom.r.x:\n");
2920 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2923 ((byte)&(1<<k))?1:0,
2926 printf("---------------\nx=dim.x/2:\n");
2928 memcpy(&byte,(u8 *)(&x)+j,1);
2931 ((byte)&(1<<k))?1:0,
2934 if(atom[i].r.x==x) printf("the same!\n");
2935 else printf("different!\n");
2937 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2938 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2939 i,atom[i].r.y,dim->y/2);
2940 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2941 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2942 i,atom[i].r.z,dim->z/2);
2952 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2959 fd=open(file,O_RDONLY);
2961 perror("[moldyn] load save file open");
2965 fsize=lseek(fd,0,SEEK_END);
2966 lseek(fd,0,SEEK_SET);
2968 size=sizeof(t_moldyn);
2971 cnt=read(fd,moldyn,size);
2973 perror("[moldyn] load save file read (moldyn)");
2979 size=moldyn->count*sizeof(t_atom);
2981 /* correcting possible atom data offset */
2983 if(fsize!=sizeof(t_moldyn)+size) {
2984 corr=fsize-sizeof(t_moldyn)-size;
2985 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2986 printf(" modifying offset:\n");
2987 printf(" - current pos: %d\n",sizeof(t_moldyn));
2988 printf(" - atom size: %d\n",size);
2989 printf(" - file size: %d\n",fsize);
2990 printf(" => correction: %d\n",corr);
2991 lseek(fd,corr,SEEK_CUR);
2994 moldyn->atom=(t_atom *)malloc(size);
2995 if(moldyn->atom==NULL) {
2996 perror("[moldyn] load save file malloc (atoms)");
3001 cnt=read(fd,moldyn->atom,size);
3003 perror("[moldyn] load save file read (atoms)");
3010 amutex=malloc(moldyn->count*sizeof(pthread_mutex_t));
3012 perror("[moldyn] load save file (mutexes)");
3015 for(cnt=0;cnt<moldyn->count;cnt++)
3016 pthread_mutex_init(&(amutex[cnt]),NULL);
3024 int moldyn_free_save_file(t_moldyn *moldyn) {
3031 int moldyn_load(t_moldyn *moldyn) {
3039 * function to find/callback all combinations of 2 body bonds
3042 int process_2b_bonds(t_moldyn *moldyn,void *data,
3043 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3044 void *data,u8 bc)) {
3054 t_list neighbour[27];
3064 for(i=0;i<moldyn->count;i++) {
3065 /* neighbour indexing */
3066 link_cell_neighbour_index(moldyn,
3067 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
3068 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
3069 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
3074 bc=(j<lc->dnlc)?0:1;
3079 while(neighbour[j][p]!=-1) {
3081 jtom=&(moldyn->atom[neighbour[j][p]]);
3089 p=lc->subcell->list[p];
3091 this=&(neighbour[j]);
3094 if(this->start==NULL)
3099 jtom=this->current->data;
3103 process(moldyn,&(itom[i]),jtom,data,bc);
3110 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3120 * function to find neighboured atoms
3123 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
3124 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
3125 void *data,u8 bc)) {
3135 t_list neighbour[27];
3144 /* neighbour indexing */
3145 link_cell_neighbour_index(moldyn,
3146 (atom->r.x+moldyn->dim.x/2)/lc->x,
3147 (atom->r.y+moldyn->dim.y/2)/lc->x,
3148 (atom->r.z+moldyn->dim.z/2)/lc->x,
3153 bc=(j<lc->dnlc)?0:1;
3158 while(neighbour[j][p]!=-1) {
3160 natom=&(moldyn->atom[neighbour[j][p]]);
3167 natom=&(moldyn->atom[p]);
3168 p=lc->subcell->list[p];
3170 this=&(neighbour[j]);
3173 if(this->start==NULL)
3178 natom=this->current->data;
3182 process(moldyn,atom,natom,data,bc);
3189 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3198 * post processing functions
3201 int get_line(int fd,char *line,int max) {
3208 if(count==max) return count;
3209 ret=read(fd,line+count,1);
3210 if(ret<=0) return ret;
3211 if(line[count]=='\n') {
3212 memset(line+count,0,max-count-1);
3220 int pair_correlation_init(t_moldyn *moldyn,double dr) {
3226 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
3243 for(i=0;i<moldyn->count;i++) {
3245 /* care for pb crossing */
3246 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3247 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3248 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3249 /* calculate distance */
3250 v3_sub(&dist,&final_r,&(atom[i].r_0));
3251 d2=v3_absolute_square(&dist);
3265 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
3266 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
3267 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
3272 int calculate_msd(t_moldyn *moldyn,double *msd) {
3289 for(i=0;i<moldyn->count;i++) {
3291 /* care for pb crossing */
3292 if(atom[i].pbc[0]|atom[i].pbc[1]|atom[i].pbc[2]) {
3293 printf("[moldyn] msd pb crossings for atom %d\n",i);
3294 printf(" x: %d y: %d z: %d\n",
3295 atom[i].pbc[0],atom[i].pbc[1],atom[i].pbc[2]);
3297 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3298 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3299 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3300 /* calculate distance */
3301 v3_sub(&dist,&final_r,&(atom[i].r_0));
3302 d2=v3_absolute_square(&dist);
3318 msd[2]/=moldyn->count;
3323 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
3328 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
3329 t_atom *jtom,void *data,u8 bc) {
3336 /* only count pairs once,
3337 * skip same atoms */
3338 if(itom->tag>=jtom->tag)
3342 * pair correlation calc
3349 v3_sub(&dist,&(jtom->r),&(itom->r));
3350 if(bc) check_per_bound(moldyn,&dist);
3351 d=v3_absolute_square(&dist);
3353 /* ignore if greater cutoff */
3354 if(d>moldyn->cutoff_square)
3357 /* fill the slots */
3361 /* should never happen but it does 8) -
3362 * related to -ffloat-store problem! */
3364 printf("[moldyn] WARNING: pcc (%d/%d)",
3370 if(itom->brand!=jtom->brand) {
3375 /* type a - type a bonds */
3377 pcc->stat[s+pcc->o1]+=1;
3379 /* type b - type b bonds */
3380 pcc->stat[s+pcc->o2]+=1;
3386 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
3393 pcc.o1=moldyn->cutoff/dr;
3396 if(pcc.o1*dr<=moldyn->cutoff)
3397 printf("[moldyn] WARNING: pcc (low #slots)\n");
3399 printf("[moldyn] pair correlation calc info:\n");
3400 printf(" time: %f\n",moldyn->time);
3401 printf(" count: %d\n",moldyn->count);
3402 printf(" cutoff: %f\n",moldyn->cutoff);
3403 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
3406 pcc.stat=(double *)ptr;
3409 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
3410 if(pcc.stat==NULL) {
3411 perror("[moldyn] pair correlation malloc");
3416 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
3419 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
3422 for(i=1;i<pcc.o1;i++) {
3423 // normalization: 4 pi r^2 dr
3424 // here: not double counting pairs -> 2 pi r r dr
3425 // ... and actually it's a constant times r^2
3428 pcc.stat[pcc.o1+i]/=norm;
3429 pcc.stat[pcc.o2+i]/=norm;
3434 /* todo: store/print pair correlation function */
3441 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3448 if(itom->tag>=jtom->tag)
3452 v3_sub(&dist,&(jtom->r),&(itom->r));
3453 if(bc) check_per_bound(moldyn,&dist);
3454 d=v3_absolute_square(&dist);
3456 /* ignore if greater or equal cutoff */
3457 if(d>moldyn->cutoff_square)
3460 /* check for potential bond */
3461 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3464 /* now count this bonding ... */
3467 /* increase total bond counter
3472 ba->acnt[jtom->tag]+=1;
3474 ba->bcnt[jtom->tag]+=1;
3477 ba->acnt[itom->tag]+=1;
3479 ba->bcnt[itom->tag]+=1;
3484 int bond_analyze(t_moldyn *moldyn,double *quality) {
3495 ba.acnt=malloc(moldyn->count*sizeof(int));
3497 perror("[moldyn] bond analyze malloc (a)");
3500 memset(ba.acnt,0,moldyn->count*sizeof(int));
3502 ba.bcnt=malloc(moldyn->count*sizeof(int));
3504 perror("[moldyn] bond analyze malloc (b)");
3507 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3516 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3518 for(i=0;i<moldyn->count;i++) {
3519 if(atom[i].brand==0) {
3520 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3522 if((ba.acnt[i]==0)&(ba.bcnt[i]==3))
3526 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3530 if((ba.acnt[i]==3)&(ba.bcnt[i]==0)) {
3539 quality[0]=1.0*ccnt4/bcnt;
3540 quality[1]=1.0*ccnt3/bcnt;
3543 printf("[moldyn] bond analyze: %f %f\n",
3544 1.0*ccnt4/bcnt,1.0*ccnt3/bcnt);
3551 * visualization code
3554 int visual_init(t_moldyn *moldyn,char *filebase) {
3556 strncpy(moldyn->vis.fb,filebase,128);
3561 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3568 if(itom->tag>=jtom->tag)
3571 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3574 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3575 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3576 itom->r.x,itom->r.y,itom->r.z,
3577 jtom->r.x,jtom->r.y,jtom->r.z);
3582 #ifdef VISUAL_THREAD
3583 void *visual_atoms(void *ptr) {
3585 int visual_atoms(t_moldyn *moldyn) {
3596 #ifdef VISUAL_THREAD
3610 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3611 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3613 perror("open visual save file fd");
3614 #ifndef VISUAL_THREAD
3619 /* write the actual data file */
3622 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3623 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3625 // atomic configuration
3626 for(i=0;i<moldyn->count;i++) {
3627 v3_sub(&strain,&(atom[i].r),&(atom[i].r_0));
3628 check_per_bound(moldyn,&strain);
3629 // atom type, positions, color and kinetic energy
3630 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3634 pse_col[atom[i].element],
3636 sqrt(v3_absolute_square(&strain)));
3639 // bonds between atoms
3640 #ifndef VISUAL_THREAD
3641 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3646 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3647 -dim.x/2,-dim.y/2,-dim.z/2,
3648 dim.x/2,-dim.y/2,-dim.z/2);
3649 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3650 -dim.x/2,-dim.y/2,-dim.z/2,
3651 -dim.x/2,dim.y/2,-dim.z/2);
3652 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3653 dim.x/2,dim.y/2,-dim.z/2,
3654 dim.x/2,-dim.y/2,-dim.z/2);
3655 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3656 -dim.x/2,dim.y/2,-dim.z/2,
3657 dim.x/2,dim.y/2,-dim.z/2);
3659 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3660 -dim.x/2,-dim.y/2,dim.z/2,
3661 dim.x/2,-dim.y/2,dim.z/2);
3662 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3663 -dim.x/2,-dim.y/2,dim.z/2,
3664 -dim.x/2,dim.y/2,dim.z/2);
3665 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3666 dim.x/2,dim.y/2,dim.z/2,
3667 dim.x/2,-dim.y/2,dim.z/2);
3668 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3669 -dim.x/2,dim.y/2,dim.z/2,
3670 dim.x/2,dim.y/2,dim.z/2);
3672 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3673 -dim.x/2,-dim.y/2,dim.z/2,
3674 -dim.x/2,-dim.y/2,-dim.z/2);
3675 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3676 -dim.x/2,dim.y/2,dim.z/2,
3677 -dim.x/2,dim.y/2,-dim.z/2);
3678 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3679 dim.x/2,-dim.y/2,dim.z/2,
3680 dim.x/2,-dim.y/2,-dim.z/2);
3681 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3682 dim.x/2,dim.y/2,dim.z/2,
3683 dim.x/2,dim.y/2,-dim.z/2);
3688 #ifdef VISUAL_THREAD
3699 * fpu cntrol functions
3702 // set rounding to double (eliminates -ffloat-store!)
3703 int fpu_set_rtd(void) {
3709 ctrl&=~_FPU_EXTENDED;