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"
37 #include "potentials/albe_orig.h"
39 #include "potentials/tersoff_orig.h"
41 #include "potentials/tersoff.h"
55 pthread_mutex_t *amutex;
56 pthread_mutex_t emutex;
59 /* fully constrained relaxation technique - global pointers */
65 * the moldyn functions
68 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
70 printf("[moldyn] init\n");
72 /* only needed if compiled without -msse2 (float-store prob!) */
75 memset(moldyn,0,sizeof(t_moldyn));
80 rand_init(&(moldyn->random),NULL,1);
81 moldyn->random.status|=RAND_STAT_VERBOSE;
84 pthread_mutex_init(&emutex,NULL);
90 int moldyn_shutdown(t_moldyn *moldyn) {
96 printf("[moldyn] shutdown\n");
99 for(i=0;i<moldyn->count;i++)
100 pthread_mutex_destroy(&(amutex[i]));
102 pthread_mutex_destroy(&emutex);
105 moldyn_log_shutdown(moldyn);
106 link_cell_shutdown(moldyn);
107 rand_close(&(moldyn->random));
113 int set_int_alg(t_moldyn *moldyn,u8 algo) {
115 printf("[moldyn] integration algorithm: ");
118 case MOLDYN_INTEGRATE_VERLET:
119 moldyn->integrate=velocity_verlet;
120 printf("velocity verlet\n");
123 printf("unknown integration algorithm: %02x\n",algo);
131 int set_cutoff(t_moldyn *moldyn,double cutoff) {
133 moldyn->cutoff=cutoff;
134 moldyn->cutoff_square=cutoff*cutoff;
136 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
141 int set_temperature(t_moldyn *moldyn,double t_ref) {
145 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
150 int set_pressure(t_moldyn *moldyn,double p_ref) {
154 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
159 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
161 moldyn->pt_scale&=(~(P_SCALE_MASK));
162 moldyn->pt_scale|=ptype;
165 printf("[moldyn] p scaling:\n");
167 printf(" p: %s",ptype?"yes":"no ");
169 printf(" | type: %02x | factor: %f",ptype,ptc);
175 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
177 moldyn->pt_scale&=(~(T_SCALE_MASK));
178 moldyn->pt_scale|=ttype;
181 printf("[moldyn] t scaling:\n");
183 printf(" t: %s",ttype?"yes":"no ");
185 printf(" | type: %02x | factor: %f",ttype,ttc);
191 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
193 moldyn->pt_scale=(ptype|ttype);
197 printf("[moldyn] p/t scaling:\n");
199 printf(" p: %s",ptype?"yes":"no ");
201 printf(" | type: %02x | factor: %f",ptype,ptc);
204 printf(" t: %s",ttype?"yes":"no ");
206 printf(" | type: %02x | factor: %f",ttype,ttc);
212 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
218 moldyn->volume=x*y*z;
226 printf("[moldyn] dimensions in A and A^3 respectively:\n");
227 printf(" x: %f\n",moldyn->dim.x);
228 printf(" y: %f\n",moldyn->dim.y);
229 printf(" z: %f\n",moldyn->dim.z);
230 printf(" volume: %f\n",moldyn->volume);
231 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
236 int set_nn_dist(t_moldyn *moldyn,double dist) {
243 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
245 printf("[moldyn] periodic boundary conditions:\n");
248 moldyn->status|=MOLDYN_STAT_PBX;
251 moldyn->status|=MOLDYN_STAT_PBY;
254 moldyn->status|=MOLDYN_STAT_PBZ;
256 printf(" x: %s\n",x?"yes":"no");
257 printf(" y: %s\n",y?"yes":"no");
258 printf(" z: %s\n",z?"yes":"no");
263 int set_potential(t_moldyn *moldyn,u8 type) {
266 case MOLDYN_POTENTIAL_TM:
267 //moldyn->func1b=tersoff_mult_1bp;
268 moldyn->func3b_j1=tersoff_mult_3bp_j1;
269 moldyn->func3b_k1=tersoff_mult_3bp_k1;
270 moldyn->func3b_j2=tersoff_mult_3bp_j2;
271 moldyn->func3b_k2=tersoff_mult_3bp_k2;
272 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
274 case MOLDYN_POTENTIAL_AO:
275 moldyn->func_j1=albe_orig_mult_3bp_j1;
276 moldyn->func_j1_k0=albe_orig_mult_3bp_k1;
277 moldyn->func_j1c=albe_orig_mult_3bp_j2;
278 moldyn->func_j1_k1=albe_orig_mult_3bp_k2;
279 moldyn->check_2b_bond=albe_orig_mult_check_2b_bond;
281 case MOLDYN_POTENTIAL_AM:
282 moldyn->func_i0=albe_mult_i0;
283 moldyn->func_j0=albe_mult_i0_j0;
284 moldyn->func_j0_k0=albe_mult_i0_j0_k0;
285 moldyn->func_j0e=albe_mult_i0_j1;
286 moldyn->func_j1=albe_mult_i0_j2;
287 moldyn->func_j1_k0=albe_mult_i0_j2_k0;
288 moldyn->func_j1c=albe_mult_i0_j3;
289 moldyn->check_2b_bond=albe_mult_check_2b_bond;
291 case MOLDYN_POTENTIAL_HO:
292 moldyn->func_j0=harmonic_oscillator;
293 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
295 case MOLDYN_POTENTIAL_LJ:
296 moldyn->func_j0=lennard_jones;
297 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
300 printf("[moldyn] set potential: unknown type %02x\n",
308 int set_avg_skip(t_moldyn *moldyn,int skip) {
310 printf("[moldyn] skip %d steps before starting average calc\n",skip);
311 moldyn->avg_skip=skip;
316 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
318 strncpy(moldyn->vlsdir,dir,127);
323 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
325 strncpy(moldyn->rauthor,author,63);
326 strncpy(moldyn->rtitle,title,63);
331 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
336 printf("[moldyn] set log: ");
339 case LOG_TOTAL_ENERGY:
340 moldyn->ewrite=timer;
341 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
342 moldyn->efd=open(filename,
343 O_WRONLY|O_CREAT|O_EXCL,
346 perror("[moldyn] energy log fd open");
349 dprintf(moldyn->efd,"# total energy log file\n");
350 printf("total energy (%d)\n",timer);
352 case LOG_TOTAL_MOMENTUM:
353 moldyn->mwrite=timer;
354 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
355 moldyn->mfd=open(filename,
356 O_WRONLY|O_CREAT|O_EXCL,
359 perror("[moldyn] momentum log fd open");
362 dprintf(moldyn->efd,"# total momentum log file\n");
363 printf("total momentum (%d)\n",timer);
366 moldyn->pwrite=timer;
367 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
368 moldyn->pfd=open(filename,
369 O_WRONLY|O_CREAT|O_EXCL,
372 perror("[moldyn] pressure log file\n");
375 dprintf(moldyn->pfd,"# pressure log file\n");
376 printf("pressure (%d)\n",timer);
378 case LOG_TEMPERATURE:
379 moldyn->twrite=timer;
380 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
381 moldyn->tfd=open(filename,
382 O_WRONLY|O_CREAT|O_EXCL,
385 perror("[moldyn] temperature log file\n");
388 dprintf(moldyn->tfd,"# temperature log file\n");
389 printf("temperature (%d)\n",timer);
392 moldyn->vwrite=timer;
393 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
394 moldyn->vfd=open(filename,
395 O_WRONLY|O_CREAT|O_EXCL,
398 perror("[moldyn] volume log file\n");
401 dprintf(moldyn->vfd,"# volume log file\n");
402 printf("volume (%d)\n",timer);
405 moldyn->swrite=timer;
406 printf("save file (%d)\n",timer);
409 moldyn->awrite=timer;
410 ret=visual_init(moldyn,moldyn->vlsdir);
412 printf("[moldyn] visual init failure\n");
415 printf("visual file (%d)\n",timer);
418 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
419 moldyn->rfd=open(filename,
420 O_WRONLY|O_CREAT|O_EXCL,
423 perror("[moldyn] report fd open");
426 printf("report -> ");
428 snprintf(filename,127,"%s/e_plot.scr",
430 moldyn->epfd=open(filename,
431 O_WRONLY|O_CREAT|O_EXCL,
434 perror("[moldyn] energy plot fd open");
437 dprintf(moldyn->epfd,e_plot_script);
442 snprintf(filename,127,"%s/pressure_plot.scr",
444 moldyn->ppfd=open(filename,
445 O_WRONLY|O_CREAT|O_EXCL,
448 perror("[moldyn] p plot fd open");
451 dprintf(moldyn->ppfd,pressure_plot_script);
456 snprintf(filename,127,"%s/temperature_plot.scr",
458 moldyn->tpfd=open(filename,
459 O_WRONLY|O_CREAT|O_EXCL,
462 perror("[moldyn] t plot fd open");
465 dprintf(moldyn->tpfd,temperature_plot_script);
467 printf("temperature ");
469 dprintf(moldyn->rfd,report_start,
470 moldyn->rauthor,moldyn->rtitle);
474 printf("unknown log type: %02x\n",type);
481 int moldyn_log_shutdown(t_moldyn *moldyn) {
485 printf("[moldyn] log shutdown\n");
489 dprintf(moldyn->rfd,report_energy);
490 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
495 if(moldyn->mfd) close(moldyn->mfd);
499 dprintf(moldyn->rfd,report_pressure);
500 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
507 dprintf(moldyn->rfd,report_temperature);
508 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
513 dprintf(moldyn->rfd,report_end);
515 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
518 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
521 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
530 * creating lattice functions
533 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,
534 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin,
535 t_part_params *p_params,t_defect_params *d_params,
536 t_offset_params *o_params) {
545 pthread_mutex_t *mutex;
551 /* how many atoms do we expect */
554 printf("[moldyn] WARNING: create 'none' lattice called");
556 if(type==CUBIC) new*=1;
557 if(type==FCC) new*=4;
558 if(type==DIAMOND) new*=8;
562 switch(d_params->stype) {
563 case DEFECT_STYPE_DB_X:
564 case DEFECT_STYPE_DB_Y:
565 case DEFECT_STYPE_DB_Z:
566 case DEFECT_STYPE_DB_R:
570 printf("[moldyn] WARNING: cl unknown defect\n");
575 /* allocate space for atoms */
576 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
578 perror("[moldyn] realloc (create lattice)");
582 atom=&(moldyn->atom[count]);
585 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
587 perror("[moldyn] mutex realloc (add atom)");
591 mutex=&(amutex[count]);
594 /* no atoms on the boundaries (only reason: it looks better!) */
609 v3_add(&orig,&orig,&(o_params->o));
610 set_nn_dist(moldyn,lc);
611 ret=cubic_init(a,b,c,lc,atom,&orig,p_params,d_params);
612 strcpy(name,"cubic");
616 v3_scale(&orig,&orig,0.5);
618 v3_add(&orig,&orig,&(o_params->o));
619 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
620 ret=fcc_init(a,b,c,lc,atom,&orig,p_params,d_params);
625 v3_scale(&orig,&orig,0.25);
627 v3_add(&orig,&orig,&(o_params->o));
628 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
629 ret=diamond_init(a,b,c,lc,atom,&orig,p_params,d_params);
630 strcpy(name,"diamond");
633 printf("unknown lattice type (%02x)\n",type);
639 printf("[moldyn] creating %s lattice (lc=%f) incomplete\n",
641 printf(" (ignore for partial lattice creation)\n");
642 printf(" amount of atoms\n");
643 printf(" - expected: %d\n",new);
644 printf(" - created: %d\n",ret);
649 printf("[moldyn] created %s lattice with %d atoms\n",name,ret);
651 for(new=0;new<ret;new++) {
652 atom[new].element=element;
653 atom[new].mass=pse_mass[element];
655 atom[new].brand=brand;
656 atom[new].tag=count+new;
657 check_per_bound(moldyn,&(atom[new].r));
658 atom[new].r_0=atom[new].r;
660 pthread_mutex_init(&(mutex[new]),NULL);
664 atom[new].element=d_params->element;
665 atom[new].mass=pse_mass[d_params->element];
666 atom[new].attr=d_params->attr;
667 atom[new].brand=d_params->brand;
668 atom[new].tag=count+new;
669 check_per_bound(moldyn,&(atom[new].r));
670 atom[new].r_0=atom[new].r;
672 pthread_mutex_init(&(mutex[new]),NULL);
678 ptr=realloc(moldyn->atom,moldyn->count*sizeof(t_atom));
680 perror("[moldyn] realloc (create lattice - alloc fix)");
685 // WHAT ABOUT AMUTEX !!!!
688 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
690 perror("[moldyn] list realloc (create lattice)");
693 moldyn->lc.subcell->list=ptr;
696 /* update total system mass */
697 total_mass_calc(moldyn);
702 int add_atom(t_moldyn *moldyn,int element,u8 brand,u8 attr,
703 t_3dvec *r,t_3dvec *v) {
710 count=(moldyn->count)++; // asshole style!
712 ptr=realloc(atom,(count+1)*sizeof(t_atom));
714 perror("[moldyn] realloc (add atom)");
720 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
722 perror("[moldyn] list realloc (add atom)");
725 moldyn->lc.subcell->list=ptr;
729 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
731 perror("[moldyn] mutex realloc (add atom)");
735 pthread_mutex_init(&(amutex[count]),NULL);
740 /* initialize new atom */
741 memset(&(atom[count]),0,sizeof(t_atom));
744 atom[count].element=element;
745 atom[count].mass=pse_mass[element];
746 atom[count].brand=brand;
747 atom[count].tag=count;
748 atom[count].attr=attr;
749 check_per_bound(moldyn,&(atom[count].r));
750 atom[count].r_0=atom[count].r;
752 /* update total system mass */
753 total_mass_calc(moldyn);
758 int del_atom(t_moldyn *moldyn,int tag) {
762 #if defined LOWMEM_LISTS || defined PTHREADS
768 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
770 perror("[moldyn]malloc (del atom)");
774 for(cnt=0;cnt<tag;cnt++)
777 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
779 new[cnt-1].tag=cnt-1;
788 ptr=realloc(moldyn->lc.subcell->list,moldyn->count*sizeof(int));
790 perror("[moldyn] list realloc (del atom)");
793 moldyn->lc.subcell->list=ptr;
795 link_cell_update(moldyn);
799 ptr=realloc(amutex,moldyn->count*sizeof(pthread_mutex_t));
801 perror("[moldyn] mutex realloc (add atom)");
805 pthread_mutex_destroy(&(amutex[moldyn->count+1]));
812 #define set_atom_positions(pos) \
813 if(d_params->type) {\
814 d_o.x=0; d_o.y=0; d_o.z=0;\
815 d_d.x=0; d_d.y=0; d_d.z=0;\
816 switch(d_params->stype) {\
817 case DEFECT_STYPE_DB_X:\
821 case DEFECT_STYPE_DB_Y:\
825 case DEFECT_STYPE_DB_Z:\
829 case DEFECT_STYPE_DB_R:\
832 printf("[moldyn] WARNING: unknown defect\n");\
835 v3_add(&dr,&pos,&d_o);\
836 v3_copy(&(atom[count].r),&dr);\
838 v3_add(&dr,&pos,&d_d);\
839 v3_copy(&(atom[count].r),&dr);\
843 v3_copy(&(atom[count].r),&pos);\
848 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
849 t_part_params *p_params,t_defect_params *d_params) {
869 /* shift partition values */
871 p.x=p_params->p.x+(a*lc)/2.0;
872 p.y=p_params->p.y+(b*lc)/2.0;
873 p.z=p_params->p.z+(c*lc)/2.0;
882 switch(p_params->type) {
885 if(v3_absolute_square(&dist)<
886 (p_params->r*p_params->r)) {
887 set_atom_positions(r);
892 if(v3_absolute_square(&dist)>=
893 (p_params->r*p_params->r)) {
894 set_atom_positions(r);
899 if((fabs(dist.x)<p_params->d.x)&&
900 (fabs(dist.y)<p_params->d.y)&&
901 (fabs(dist.z)<p_params->d.z)) {
902 set_atom_positions(r);
907 if((fabs(dist.x)>=p_params->d.x)||
908 (fabs(dist.y)>=p_params->d.y)||
909 (fabs(dist.z)>=p_params->d.z)) {
910 set_atom_positions(r);
914 set_atom_positions(r);
924 for(i=0;i<count;i++) {
925 atom[i].r.x-=(a*lc)/2.0;
926 atom[i].r.y-=(b*lc)/2.0;
927 atom[i].r.z-=(c*lc)/2.0;
933 /* fcc lattice init */
934 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
935 t_part_params *p_params,t_defect_params *d_params) {
953 /* construct the basis */
954 memset(basis,0,3*sizeof(t_3dvec));
962 /* shift partition values */
964 p.x=p_params->p.x+(a*lc)/2.0;
965 p.y=p_params->p.y+(b*lc)/2.0;
966 p.z=p_params->p.z+(c*lc)/2.0;
969 /* fill up the room */
977 switch(p_params->type) {
980 if(v3_absolute_square(&dist)<
981 (p_params->r*p_params->r)) {
982 set_atom_positions(r);
987 if(v3_absolute_square(&dist)>=
988 (p_params->r*p_params->r)) {
989 set_atom_positions(r);
994 if((fabs(dist.x)<p_params->d.x)&&
995 (fabs(dist.y)<p_params->d.y)&&
996 (fabs(dist.z)<p_params->d.z)) {
997 set_atom_positions(r);
1000 case PART_OUTSIDE_D:
1001 v3_sub(&dist,&r,&p);
1002 if((fabs(dist.x)>=p_params->d.x)||
1003 (fabs(dist.y)>=p_params->d.y)||
1004 (fabs(dist.z)>=p_params->d.z)) {
1005 set_atom_positions(r);
1009 set_atom_positions(r);
1012 /* the three face centered atoms */
1014 v3_add(&n,&r,&basis[l]);
1015 switch(p_params->type) {
1017 v3_sub(&dist,&n,&p);
1018 if(v3_absolute_square(&dist)<
1019 (p_params->r*p_params->r)) {
1020 set_atom_positions(n);
1023 case PART_OUTSIDE_R:
1024 v3_sub(&dist,&n,&p);
1025 if(v3_absolute_square(&dist)>=
1026 (p_params->r*p_params->r)) {
1027 set_atom_positions(n);
1031 v3_sub(&dist,&n,&p);
1032 if((fabs(dist.x)<p_params->d.x)&&
1033 (fabs(dist.y)<p_params->d.y)&&
1034 (fabs(dist.z)<p_params->d.z)) {
1035 set_atom_positions(n);
1038 case PART_OUTSIDE_D:
1039 v3_sub(&dist,&n,&p);
1040 if((fabs(dist.x)>=p_params->d.x)||
1041 (fabs(dist.y)>=p_params->d.y)||
1042 (fabs(dist.z)>=p_params->d.z)) {
1043 set_atom_positions(n);
1047 set_atom_positions(n);
1058 /* coordinate transformation */
1059 for(i=0;i<count;i++) {
1060 atom[i].r.x-=(a*lc)/2.0;
1061 atom[i].r.y-=(b*lc)/2.0;
1062 atom[i].r.z-=(c*lc)/2.0;
1068 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
1069 t_part_params *p_params,t_defect_params *d_params) {
1074 count=fcc_init(a,b,c,lc,atom,origin,p_params,d_params);
1080 if(origin) v3_add(&o,&o,origin);
1082 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_params,d_params);
1087 int destroy_atoms(t_moldyn *moldyn) {
1089 if(moldyn->atom) free(moldyn->atom);
1094 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
1097 * - gaussian distribution of velocities
1098 * - zero total momentum
1099 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1104 t_3dvec p_total,delta;
1109 random=&(moldyn->random);
1111 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
1113 /* gaussian distribution of velocities */
1115 for(i=0;i<moldyn->count;i++) {
1116 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
1118 v=sigma*rand_get_gauss(random);
1120 p_total.x+=atom[i].mass*v;
1122 v=sigma*rand_get_gauss(random);
1124 p_total.y+=atom[i].mass*v;
1126 v=sigma*rand_get_gauss(random);
1128 p_total.z+=atom[i].mass*v;
1131 /* zero total momentum */
1132 v3_scale(&p_total,&p_total,1.0/moldyn->count);
1133 for(i=0;i<moldyn->count;i++) {
1134 v3_scale(&delta,&p_total,1.0/atom[i].mass);
1135 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
1138 /* velocity scaling */
1139 scale_velocity(moldyn,equi_init);
1144 double total_mass_calc(t_moldyn *moldyn) {
1150 for(i=0;i<moldyn->count;i++)
1151 moldyn->mass+=moldyn->atom[i].mass;
1153 return moldyn->mass;
1156 double temperature_calc(t_moldyn *moldyn) {
1158 /* assume up to date kinetic energy, which is 3/2 N k_B T */
1161 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
1167 double get_temperature(t_moldyn *moldyn) {
1172 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
1182 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
1185 /* get kinetic energy / temperature & count involved atoms */
1188 for(i=0;i<moldyn->count;i++) {
1189 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
1190 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1195 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1196 else return 0; /* no atoms involved in scaling! */
1198 /* (temporary) hack for e,t = 0 */
1201 if(moldyn->t_ref!=0.0) {
1202 thermal_init(moldyn,equi_init);
1206 return 0; /* no scaling needed */
1210 /* get scaling factor */
1211 scale=moldyn->t_ref/moldyn->t;
1215 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1216 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1219 /* velocity scaling */
1220 for(i=0;i<moldyn->count;i++) {
1221 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1222 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1228 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1232 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1237 double virial_sum(t_moldyn *moldyn) {
1242 /* virial (sum over atom virials) */
1250 for(i=0;i<moldyn->count;i++) {
1251 virial=&(moldyn->atom[i].virial);
1252 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1253 moldyn->vir.xx+=virial->xx;
1254 moldyn->vir.yy+=virial->yy;
1255 moldyn->vir.zz+=virial->zz;
1256 moldyn->vir.xy+=virial->xy;
1257 moldyn->vir.xz+=virial->xz;
1258 moldyn->vir.yz+=virial->yz;
1261 /* global virial (absolute coordinates) */
1262 //virial=&(moldyn->gvir);
1263 //moldyn->gv=virial->xx+virial->yy+virial->zz;
1265 return moldyn->virial;
1268 double pressure_calc(t_moldyn *moldyn) {
1272 * with W = 1/3 sum_i f_i r_i (- skipped!)
1273 * virial = sum_i f_i r_i
1275 * => P = (2 Ekin + virial) / (3V)
1278 /* assume up to date virial & up to date kinetic energy */
1280 /* pressure (atom virials) */
1281 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1282 moldyn->p/=(3.0*moldyn->volume);
1284 //moldyn->px=2.0*moldyn->ekinx+moldyn->vir.xx;
1285 //moldyn->px/=moldyn->volume;
1286 //moldyn->py=2.0*moldyn->ekiny+moldyn->vir.yy;
1287 //moldyn->py/=moldyn->volume;
1288 //moldyn->pz=2.0*moldyn->ekinz+moldyn->vir.zz;
1289 //moldyn->pz/=moldyn->volume;
1291 /* pressure (absolute coordinates) */
1292 //moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1293 //moldyn->gp/=(3.0*moldyn->volume);
1298 int average_reset(t_moldyn *moldyn) {
1300 printf("[moldyn] average reset\n");
1302 /* update skip value */
1303 moldyn->avg_skip=moldyn->total_steps;
1305 /* kinetic energy */
1309 /* potential energy */
1317 moldyn->virial_sum=0.0;
1318 //moldyn->gv_sum=0.0;
1322 //moldyn->gp_sum=0.0;
1328 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1332 if(moldyn->total_steps<moldyn->avg_skip)
1335 denom=moldyn->total_steps+1-moldyn->avg_skip;
1337 /* assume up to date energies, temperature, pressure etc */
1339 /* kinetic energy */
1340 moldyn->k_sum+=moldyn->ekin;
1341 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1342 moldyn->k_avg=moldyn->k_sum/denom;
1343 moldyn->k2_avg=moldyn->k2_sum/denom;
1344 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1346 /* potential energy */
1347 moldyn->v_sum+=moldyn->energy;
1348 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1349 moldyn->v_avg=moldyn->v_sum/denom;
1350 moldyn->v2_avg=moldyn->v2_sum/denom;
1351 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1354 moldyn->t_sum+=moldyn->t;
1355 moldyn->t_avg=moldyn->t_sum/denom;
1358 moldyn->virial_sum+=moldyn->virial;
1359 moldyn->virial_avg=moldyn->virial_sum/denom;
1360 //moldyn->gv_sum+=moldyn->gv;
1361 //moldyn->gv_avg=moldyn->gv_sum/denom;
1364 moldyn->p_sum+=moldyn->p;
1365 moldyn->p_avg=moldyn->p_sum/denom;
1366 //moldyn->gp_sum+=moldyn->gp;
1367 //moldyn->gp_avg=moldyn->gp_sum/denom;
1368 moldyn->tp_sum+=moldyn->tp;
1369 moldyn->tp_avg=moldyn->tp_sum/denom;
1374 int get_heat_capacity(t_moldyn *moldyn) {
1378 /* averages needed for heat capacity calc */
1379 if(moldyn->total_steps<moldyn->avg_skip)
1382 /* (temperature average)^2 */
1383 temp2=moldyn->t_avg*moldyn->t_avg;
1384 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1387 /* ideal gas contribution */
1388 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1389 printf(" ideal gas contribution: %f\n",
1390 ighc/moldyn->mass*KILOGRAM/JOULE);
1392 /* specific heat for nvt ensemble */
1393 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1394 moldyn->c_v_nvt/=moldyn->mass;
1396 /* specific heat for nve ensemble */
1397 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1398 moldyn->c_v_nve/=moldyn->mass;
1400 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1401 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1402 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)));
1407 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1423 /* store atomic configuration + dimension */
1424 store=malloc(moldyn->count*sizeof(t_atom));
1426 printf("[moldyn] allocating store mem failed\n");
1429 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1434 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1435 su=pow(2.0-h,ONE_THIRD)-1.0;
1436 dv=(1.0-h)*moldyn->volume;
1438 /* scale up dimension and atom positions */
1439 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1440 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1441 link_cell_shutdown(moldyn);
1442 link_cell_init(moldyn,QUIET);
1443 potential_force_calc(moldyn);
1446 /* restore atomic configuration + dim */
1447 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1450 /* scale down dimension and atom positions */
1451 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1452 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1453 link_cell_shutdown(moldyn);
1454 link_cell_init(moldyn,QUIET);
1455 potential_force_calc(moldyn);
1458 /* calculate pressure */
1459 moldyn->tp=-(y1-y0)/(2.0*dv);
1461 /* restore atomic configuration */
1462 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1464 link_cell_shutdown(moldyn);
1465 link_cell_init(moldyn,QUIET);
1466 //potential_force_calc(moldyn);
1468 /* free store buffer */
1475 double get_pressure(t_moldyn *moldyn) {
1481 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1493 if(x) dim->x*=scale;
1494 if(y) dim->y*=scale;
1495 if(z) dim->z*=scale;
1500 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1511 for(i=0;i<moldyn->count;i++) {
1512 r=&(moldyn->atom[i].r);
1521 int scale_atoms_ind(t_moldyn *moldyn,double x,double y,double z) {
1526 for(i=0;i<moldyn->count;i++) {
1527 r=&(moldyn->atom[i].r);
1536 int scale_dim_ind(t_moldyn *moldyn,double x,double y,double z) {
1549 int scale_volume(t_moldyn *moldyn) {
1556 vdim=&(moldyn->vis.dim);
1560 /* scaling factor */
1561 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1562 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1563 scale=pow(scale,ONE_THIRD);
1566 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1571 sx=1.0-(moldyn->p_ref-moldyn->px)*moldyn->p_tc*moldyn->tau;
1572 sy=1.0-(moldyn->p_ref-moldyn->py)*moldyn->p_tc*moldyn->tau;
1573 sz=1.0-(moldyn->p_ref-moldyn->pz)*moldyn->p_tc*moldyn->tau;
1574 sx=pow(sx,ONE_THIRD);
1575 sy=pow(sy,ONE_THIRD);
1576 sz=pow(sz,ONE_THIRD);
1579 /* scale the atoms and dimensions */
1580 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1581 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1582 //scale_atoms_ind(moldyn,sx,sy,sz);
1583 //scale_dim_ind(moldyn,sx,sy,sz);
1585 /* visualize dimensions */
1592 /* recalculate scaled volume */
1593 moldyn->volume=dim->x*dim->y*dim->z;
1595 /* adjust/reinit linkcell */
1596 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1597 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1598 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1599 link_cell_shutdown(moldyn);
1600 link_cell_init(moldyn,QUIET);
1614 double e_kin_calc(t_moldyn *moldyn) {
1621 //moldyn->ekinx=0.0;
1622 //moldyn->ekiny=0.0;
1623 //moldyn->ekinz=0.0;
1625 for(i=0;i<moldyn->count;i++) {
1626 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1627 moldyn->ekin+=atom[i].ekin;
1628 //moldyn->ekinx+=0.5*atom[i].mass*atom[i].v.x*atom[i].v.x;
1629 //moldyn->ekiny+=0.5*atom[i].mass*atom[i].v.y*atom[i].v.y;
1630 //moldyn->ekinz+=0.5*atom[i].mass*atom[i].v.z*atom[i].v.z;
1633 return moldyn->ekin;
1636 double get_total_energy(t_moldyn *moldyn) {
1638 return(moldyn->ekin+moldyn->energy);
1641 t_3dvec get_total_p(t_moldyn *moldyn) {
1650 for(i=0;i<moldyn->count;i++) {
1651 v3_scale(&p,&(atom[i].v),atom[i].mass);
1652 v3_add(&p_total,&p_total,&p);
1658 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1662 /* nn_dist is the nearest neighbour distance */
1664 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1673 /* linked list / cell method */
1675 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1678 #ifndef LOWMEM_LISTS
1684 /* partitioning the md cell */
1685 lc->nx=moldyn->dim.x/moldyn->cutoff;
1686 lc->x=moldyn->dim.x/lc->nx;
1687 lc->ny=moldyn->dim.y/moldyn->cutoff;
1688 lc->y=moldyn->dim.y/lc->ny;
1689 lc->nz=moldyn->dim.z/moldyn->cutoff;
1690 lc->z=moldyn->dim.z/lc->nz;
1691 lc->cells=lc->nx*lc->ny*lc->nz;
1694 lc->subcell=malloc(lc->cells*sizeof(int*));
1696 lc->subcell=malloc(sizeof(t_lowmem_list));
1698 lc->subcell=malloc(lc->cells*sizeof(t_list));
1701 if(lc->subcell==NULL) {
1702 perror("[moldyn] cell init (malloc)");
1707 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1712 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1715 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1718 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1721 printf(" x: %d x %f A\n",lc->nx,lc->x);
1722 printf(" y: %d x %f A\n",lc->ny,lc->y);
1723 printf(" z: %d x %f A\n",lc->nz,lc->z);
1728 for(i=0;i<lc->cells;i++) {
1729 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1730 if(lc->subcell[i]==NULL) {
1731 perror("[moldyn] list init (malloc)");
1736 printf(" ---> %d malloc %p (%p)\n",
1737 i,lc->subcell[0],lc->subcell);
1741 lc->subcell->head=malloc(lc->cells*sizeof(int));
1742 if(lc->subcell->head==NULL) {
1743 perror("[moldyn] head init (malloc)");
1746 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1747 if(lc->subcell->list==NULL) {
1748 perror("[moldyn] list init (malloc)");
1752 for(i=0;i<lc->cells;i++)
1753 list_init_f(&(lc->subcell[i]));
1756 /* update the list */
1757 link_cell_update(moldyn);
1762 int link_cell_update(t_moldyn *moldyn) {
1780 for(i=0;i<lc->cells;i++)
1782 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1784 lc->subcell->head[i]=-1;
1786 list_destroy_f(&(lc->subcell[i]));
1789 for(count=0;count<moldyn->count;count++) {
1790 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1791 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1792 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1796 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1799 if(p>=MAX_ATOMS_PER_LIST) {
1800 printf("[moldyn] FATAL: amount of atoms too high!\n");
1804 lc->subcell[i+j*nx+k*nxy][p]=count;
1807 lc->subcell->list[count]=lc->subcell->head[p];
1808 lc->subcell->head[p]=count;
1810 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1814 printf(" ---> %d %d malloc %p (%p)\n",
1815 i,count,lc->subcell[i].current,lc->subcell);
1823 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1849 if(i>=nx||j>=ny||k>=nz)
1850 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1853 #ifndef LOWMEM_LISTS
1854 cell[0]=lc->subcell[i+j*nx+k*a];
1856 cell[0]=lc->subcell->head[i+j*nx+k*a];
1858 for(ci=-1;ci<=1;ci++) {
1861 if((x<0)||(x>=nx)) {
1865 for(cj=-1;cj<=1;cj++) {
1868 if((y<0)||(y>=ny)) {
1872 for(ck=-1;ck<=1;ck++) {
1875 if((z<0)||(z>=nz)) {
1879 if(!(ci|cj|ck)) continue;
1881 #ifndef LOWMEM_LISTS
1882 cell[--count2]=lc->subcell[x+y*nx+z*a];
1884 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1889 #ifndef LOWMEM_LISTS
1890 cell[count1++]=lc->subcell[x+y*nx+z*a];
1892 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1904 int link_cell_shutdown(t_moldyn *moldyn) {
1906 #ifndef LOWMEM_LISTS
1914 free(lc->subcell->head);
1915 free(lc->subcell->list);
1919 for(i=0;i<lc->cells;i++) {
1921 free(lc->subcell[i]);
1923 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1924 list_destroy_f(&(lc->subcell[i]));
1934 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1938 t_moldyn_schedule *schedule;
1940 schedule=&(moldyn->schedule);
1941 count=++(schedule->total_sched);
1943 ptr=realloc(schedule->runs,count*sizeof(int));
1945 perror("[moldyn] realloc (runs)");
1949 schedule->runs[count-1]=runs;
1951 ptr=realloc(schedule->tau,count*sizeof(double));
1953 perror("[moldyn] realloc (tau)");
1957 schedule->tau[count-1]=tau;
1959 printf("[moldyn] schedule added:\n");
1960 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1966 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1968 moldyn->schedule.hook=hook;
1969 moldyn->schedule.hook_params=hook_params;
1976 * 'integration of newtons equation' - algorithms
1980 /* start the integration */
1982 int moldyn_integrate(t_moldyn *moldyn) {
1985 unsigned int e,m,s,v,p,t,a;
1987 t_moldyn_schedule *sched;
1992 double energy_scale;
1993 struct timeval t1,t2;
1996 #ifdef VISUAL_THREAD
1998 pthread_t io_thread;
2007 sched=&(moldyn->schedule);
2010 /* initialize linked cell method */
2011 link_cell_init(moldyn,VERBOSE);
2013 /* logging & visualization */
2022 /* sqaure of some variables */
2023 moldyn->tau_square=moldyn->tau*moldyn->tau;
2025 /* get current time */
2026 gettimeofday(&t1,NULL);
2028 /* calculate initial forces */
2029 potential_force_calc(moldyn);
2034 /* some stupid checks before we actually start calculating bullshit */
2035 if(moldyn->cutoff>0.5*moldyn->dim.x)
2036 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
2037 if(moldyn->cutoff>0.5*moldyn->dim.y)
2038 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
2039 if(moldyn->cutoff>0.5*moldyn->dim.z)
2040 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
2042 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
2043 if(ds>0.05*moldyn->nnd)
2044 printf("[moldyn] WARNING: forces too high / tau too small!\n");
2047 /* zero absolute time */
2048 // should have right values!
2050 //moldyn->total_steps=0;
2052 /* debugging, ignore */
2055 /* zero & init moldyn copy */
2056 #ifdef VISUAL_THREAD
2057 memset(&md_copy,0,sizeof(t_moldyn));
2058 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2059 if(atom_copy==NULL) {
2060 perror("[moldyn] malloc atom copy (init)");
2066 printf("##################\n");
2067 printf("# USING PTHREADS #\n");
2068 printf("##################\n");
2070 /* tell the world */
2071 printf("[moldyn] integration start, go get a coffee ...\n");
2073 /* executing the schedule */
2075 while(sched->count<sched->total_sched) {
2077 /* setting amount of runs and finite time step size */
2078 moldyn->tau=sched->tau[sched->count];
2079 moldyn->tau_square=moldyn->tau*moldyn->tau;
2080 moldyn->time_steps=sched->runs[sched->count];
2082 /* energy scaling factor (might change!) */
2083 energy_scale=moldyn->count*EV;
2085 /* integration according to schedule */
2087 for(i=0;i<moldyn->time_steps;i++) {
2089 /* integration step */
2090 moldyn->integrate(moldyn);
2092 /* calculate kinetic energy, temperature and pressure */
2094 temperature_calc(moldyn);
2096 pressure_calc(moldyn);
2098 thermodynamic_pressure_calc(moldyn);
2099 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
2103 /* calculate fluctuations + averages */
2104 average_and_fluctuation_calc(moldyn);
2107 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
2108 scale_velocity(moldyn,FALSE);
2109 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
2110 scale_volume(moldyn);
2112 /* check for log & visualization */
2114 if(!(moldyn->total_steps%e))
2115 dprintf(moldyn->efd,
2116 "%f %f %f %f %f %f\n",
2117 moldyn->time,moldyn->ekin/energy_scale,
2118 moldyn->energy/energy_scale,
2119 get_total_energy(moldyn)/energy_scale,
2120 moldyn->ekin/EV,moldyn->energy/EV);
2123 if(!(moldyn->total_steps%m)) {
2124 momentum=get_total_p(moldyn);
2125 dprintf(moldyn->mfd,
2126 "%f %f %f %f %f\n",moldyn->time,
2127 momentum.x,momentum.y,momentum.z,
2128 v3_norm(&momentum));
2132 if(!(moldyn->total_steps%p)) {
2133 dprintf(moldyn->pfd,
2134 "%f %f %f %f %f %f %f\n",moldyn->time,
2135 moldyn->p/BAR,moldyn->p_avg/BAR,
2136 moldyn->p/BAR,moldyn->p_avg/BAR,
2137 moldyn->tp/BAR,moldyn->tp_avg/BAR);
2141 if(!(moldyn->total_steps%t)) {
2142 dprintf(moldyn->tfd,
2144 moldyn->time,moldyn->t,moldyn->t_avg);
2148 if(!(moldyn->total_steps%v)) {
2149 dprintf(moldyn->vfd,
2150 "%f %f %f %f %f\n",moldyn->time,
2158 if(!(moldyn->total_steps%s)) {
2159 snprintf(dir,128,"%s/s-%08.f.save",
2160 moldyn->vlsdir,moldyn->time);
2161 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
2163 if(fd<0) perror("[moldyn] save fd open");
2165 write(fd,moldyn,sizeof(t_moldyn));
2166 write(fd,moldyn->atom,
2167 moldyn->count*sizeof(t_atom));
2173 if(!(moldyn->total_steps%a)) {
2174 #ifdef VISUAL_THREAD
2175 /* check whether thread has not terminated yet */
2177 ret=pthread_join(io_thread,NULL);
2180 /* prepare and start thread */
2181 if(moldyn->count!=md_copy.count) {
2185 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
2187 atom_copy=malloc(moldyn->count*sizeof(t_atom));
2188 if(atom_copy==NULL) {
2189 perror("[moldyn] malloc atom copy (change)");
2193 md_copy.atom=atom_copy;
2194 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
2196 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
2198 perror("[moldyn] create visual atoms thread\n");
2202 visual_atoms(moldyn);
2207 /* display progress */
2211 /* get current time */
2212 gettimeofday(&t2,NULL);
2214 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
2215 sched->count,i,moldyn->total_steps,
2216 moldyn->t,moldyn->t_avg,
2218 moldyn->p/BAR,moldyn->p_avg/BAR,
2220 moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
2223 (int)(t2.tv_sec-t1.tv_sec));
2227 /* copy over time */
2233 /* increase absolute time */
2234 moldyn->time+=moldyn->tau;
2235 moldyn->total_steps+=1;
2239 /* check for hooks */
2241 printf("\n ## schedule hook %d start ##\n",
2243 sched->hook(moldyn,sched->hook_params);
2244 printf(" ## schedule hook end ##\n");
2247 /* increase the schedule counter */
2252 /* writing a final save file! */
2254 snprintf(dir,128,"%s/s-final.save",moldyn->vlsdir);
2255 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,S_IRUSR|S_IWUSR);
2256 if(fd<0) perror("[moldyn] save fd open");
2258 write(fd,moldyn,sizeof(t_moldyn));
2259 write(fd,moldyn->atom,
2260 moldyn->count*sizeof(t_atom));
2264 /* writing a final visual file! */
2266 visual_atoms(moldyn);
2276 int basis_trafo(t_3dvec *r,u8 dir,double z,double x) {
2283 r->x=cos(z)*tmp.x-sin(z)*tmp.y;
2284 r->y=sin(z)*tmp.x+cos(z)*tmp.y;
2288 r->y=cos(x)*tmp.y-sin(x)*tmp.z;
2289 r->z=sin(x)*tmp.y+cos(x)*tmp.z;
2295 r->y=cos(-x)*tmp.y-sin(-x)*tmp.z;
2296 r->z=sin(-x)*tmp.y+cos(-x)*tmp.z;
2300 r->x=cos(-z)*tmp.x-sin(-z)*tmp.y;
2301 r->y=sin(-z)*tmp.x+cos(-z)*tmp.y;
2308 /* velocity verlet */
2310 int velocity_verlet(t_moldyn *moldyn) {
2313 double tau,tau_square,h;
2318 count=moldyn->count;
2320 tau_square=moldyn->tau_square;
2322 for(i=0;i<count;i++) {
2324 /* check whether fixed atom */
2325 if(atom[i].attr&ATOM_ATTR_FP)
2331 /* constraint relaxation */
2334 basis_trafo(&(atom[i].f),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].f),BACKWARD,
2343 trafo_angle[2*i],trafo_angle[2*i+1]);
2345 basis_trafo(&(atom[i].v),FORWARD,
2346 trafo_angle[2*i],trafo_angle[2*i+1]);
2347 if(constraints[3*i])
2349 if(constraints[3*i+1])
2351 if(constraints[3*i+2])
2353 basis_trafo(&(atom[i].v),BACKWARD,
2354 trafo_angle[2*i],trafo_angle[2*i+1]);
2358 v3_scale(&delta,&(atom[i].v),tau);
2359 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2361 v3_scale(&delta,&(atom[i].f),h*tau_square);
2362 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2363 //check_per_bound_and_care_for_pbc(moldyn,&(atom[i]));
2364 check_per_bound(moldyn,&(atom[i].r));
2366 /* velocities [actually v(t+tau/2)] */
2367 v3_scale(&delta,&(atom[i].f),h*tau);
2368 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2371 /* criticial check */
2372 moldyn_bc_check(moldyn);
2374 /* neighbour list update */
2375 link_cell_update(moldyn);
2377 /* forces depending on chosen potential */
2379 // if albe, use albe force calc routine
2380 //if(moldyn->func3b_j1==albe_mult_3bp_j1)
2381 // albe_potential_force_calc(moldyn);
2383 potential_force_calc(moldyn);
2385 albe_potential_force_calc(moldyn);
2388 for(i=0;i<count;i++) {
2389 /* check whether fixed atom */
2390 if(atom[i].attr&ATOM_ATTR_FP)
2393 /* constraint relaxation */
2396 basis_trafo(&(atom[i].f),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].f),BACKWARD,
2405 trafo_angle[2*i],trafo_angle[2*i+1]);
2407 basis_trafo(&(atom[i].v),FORWARD,
2408 trafo_angle[2*i],trafo_angle[2*i+1]);
2409 if(constraints[3*i])
2411 if(constraints[3*i+1])
2413 if(constraints[3*i+2])
2415 basis_trafo(&(atom[i].v),BACKWARD,
2416 trafo_angle[2*i],trafo_angle[2*i+1]);
2419 /* again velocities [actually v(t+tau)] */
2420 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2421 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2430 * potentials & corresponding forces & virial routine
2434 /* generic potential and force calculation */
2436 int potential_force_calc(t_moldyn *moldyn) {
2439 t_atom *itom,*jtom,*ktom;
2443 int *neighbour_i[27];
2447 int neighbour_i[27];
2450 t_list neighbour_i[27];
2451 t_list neighbour_i2[27];
2457 count=moldyn->count;
2467 /* reset global virial */
2468 memset(&(moldyn->gvir),0,sizeof(t_virial));
2470 /* reset force, site energy and virial of every atom */
2472 i=omp_get_thread_num();
2473 #pragma omp parallel for private(virial)
2475 for(i=0;i<count;i++) {
2478 v3_zero(&(itom[i].f));
2481 virial=(&(itom[i].virial));
2489 /* reset site energy */
2494 /* get energy, force and virial of every atom */
2496 /* first (and only) loop over atoms i */
2497 for(i=0;i<count;i++) {
2499 /* single particle potential/force */
2500 if(itom[i].attr&ATOM_ATTR_1BP)
2502 moldyn->func_i0(moldyn,&(itom[i]));
2504 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2507 /* 2 body pair potential/force */
2509 link_cell_neighbour_index(moldyn,
2510 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2511 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2512 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2517 #ifndef STATIC_LISTS
2518 /* check for later 3 body interaction */
2519 if(itom[i].attr&ATOM_ATTR_3BP)
2520 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2523 /* first loop over atoms j */
2524 if(moldyn->func_j0) {
2531 while(neighbour_i[j][p]!=-1) {
2533 jtom=&(atom[neighbour_i[j][p]]);
2541 p=lc->subcell->list[p];
2543 this=&(neighbour_i[j]);
2546 if(this->start==NULL)
2550 jtom=this->current->data;
2553 if(jtom==&(itom[i]))
2559 if((jtom->attr&ATOM_ATTR_2BP)&
2560 (itom[i].attr&ATOM_ATTR_2BP)) {
2561 moldyn->func_j0(moldyn,
2567 /* 3 body potential/force */
2569 /* in j loop, 3bp run can be skipped */
2570 if(!(moldyn->run3bp))
2573 if(!(itom[i].attr&ATOM_ATTR_3BP))
2576 if(!(jtom->attr&ATOM_ATTR_3BP))
2579 if(moldyn->func_j0_k0==NULL)
2582 /* first loop over atoms k in first j loop */
2589 while(neighbour_i[j][q]!=0) {
2591 ktom=&(atom[neighbour_i[k][q]]);
2594 that=&(neighbour_i2[k]);
2597 if(that->start==NULL)
2601 ktom=that->current->data;
2604 if(!(ktom->attr&ATOM_ATTR_3BP))
2610 if(ktom==&(itom[i]))
2613 moldyn->func_j0_k0(moldyn,
2625 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2631 /* continued 3 body potential/force */
2633 if(!(itom[i].attr&ATOM_ATTR_3BP))
2636 /* copy the neighbour lists */
2638 /* no copy needed for static lists */
2640 /* no copy needed for lowmem lists */
2642 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2645 /* second loop over atoms j */
2652 while(neighbour_i[j][p]!=-1) {
2654 jtom=&(atom[neighbour_i[j][p]]);
2662 p=lc->subcell->list[p];
2664 this=&(neighbour_i[j]);
2667 if(this->start==NULL)
2672 jtom=this->current->data;
2675 if(jtom==&(itom[i]))
2678 if(!(jtom->attr&ATOM_ATTR_3BP))
2685 moldyn->func_j1(moldyn,
2690 /* in j loop, 3bp run can be skipped */
2691 if(!(moldyn->run3bp))
2694 /* first loop over atoms k in second j loop */
2695 if(moldyn->func_j1_k0) {
2703 while(neighbour_i[k][q]!=-1) {
2705 ktom=&(atom[neighbour_i[k][q]]);
2713 q=lc->subcell->list[q];
2715 that=&(neighbour_i2[k]);
2718 if(that->start==NULL)
2722 ktom=that->current->data;
2725 if(!(ktom->attr&ATOM_ATTR_3BP))
2731 if(ktom==&(itom[i]))
2734 moldyn->func3b_k1(moldyn,
2745 } while(list_next_f(that)!=\
2753 /* continued j loop after first k loop */
2754 if(moldyn->func_j1c)
2755 moldyn->func_j1c(moldyn,
2760 /* second loop over atoms k */
2761 if(moldyn->func_j1_k1) {
2769 while(neighbour_i[k][q]!=-1) {
2771 ktom=&(atom[neighbour_i[k][q]]);
2779 q=lc->subcell->list[q];
2781 that=&(neighbour_i2[k]);
2784 if(that->start==NULL)
2788 ktom=that->current->data;
2791 if(!(ktom->attr&ATOM_ATTR_3BP))
2797 if(ktom==&(itom[i]))
2800 moldyn->func_j1_k1(moldyn,
2811 } while(list_next_f(that)!=\
2819 /* finish of j loop after second k loop */
2820 if(moldyn->func_j1e) {
2821 moldyn->func_j1e(moldyn,
2830 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2845 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2846 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2848 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2849 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2850 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2854 /* some postprocessing */
2856 #pragma omp parallel for
2858 for(i=0;i<count;i++) {
2859 /* calculate global virial */
2860 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2861 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2862 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2863 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2864 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2865 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2867 /* check forces regarding the given timestep */
2868 if(v3_norm(&(itom[i].f))>\
2869 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2870 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2878 * virial calculation
2881 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2882 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2884 a->virial.xx+=f->x*d->x;
2885 a->virial.yy+=f->y*d->y;
2886 a->virial.zz+=f->z*d->z;
2887 a->virial.xy+=f->x*d->y;
2888 a->virial.xz+=f->x*d->z;
2889 a->virial.yz+=f->y*d->z;
2895 * periodic boundary checking
2898 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2899 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2910 if(moldyn->status&MOLDYN_STAT_PBX) {
2911 if(a->x>=x) a->x-=dim->x;
2912 else if(-a->x>x) a->x+=dim->x;
2914 if(moldyn->status&MOLDYN_STAT_PBY) {
2915 if(a->y>=y) a->y-=dim->y;
2916 else if(-a->y>y) a->y+=dim->y;
2918 if(moldyn->status&MOLDYN_STAT_PBZ) {
2919 if(a->z>=z) a->z-=dim->z;
2920 else if(-a->z>z) a->z+=dim->z;
2926 int check_per_bound_and_care_for_pbc(t_moldyn *moldyn,t_atom *a) {
2937 if(moldyn->status&MOLDYN_STAT_PBX) {
2942 else if(-a->r.x>x) {
2947 if(moldyn->status&MOLDYN_STAT_PBY) {
2952 else if(-a->r.y>y) {
2957 if(moldyn->status&MOLDYN_STAT_PBZ) {
2962 else if(-a->r.z>z) {
2972 * debugging / critical check functions
2975 int moldyn_bc_check(t_moldyn *moldyn) {
2988 for(i=0;i<moldyn->count;i++) {
2989 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2990 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2991 i,atom[i].r.x,dim->x/2);
2992 printf("diagnostic:\n");
2993 printf("-----------\natom.r.x:\n");
2995 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2998 ((byte)&(1<<k))?1:0,
3001 printf("---------------\nx=dim.x/2:\n");
3003 memcpy(&byte,(u8 *)(&x)+j,1);
3006 ((byte)&(1<<k))?1:0,
3009 if(atom[i].r.x==x) printf("the same!\n");
3010 else printf("different!\n");
3012 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
3013 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
3014 i,atom[i].r.y,dim->y/2);
3015 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
3016 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
3017 i,atom[i].r.z,dim->z/2);
3027 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
3034 fd=open(file,O_RDONLY);
3036 perror("[moldyn] load save file open");
3040 fsize=lseek(fd,0,SEEK_END);
3041 lseek(fd,0,SEEK_SET);
3043 size=sizeof(t_moldyn);
3046 cnt=read(fd,moldyn,size);
3048 perror("[moldyn] load save file read (moldyn)");
3054 size=moldyn->count*sizeof(t_atom);
3056 /* correcting possible atom data offset */
3058 if(fsize!=sizeof(t_moldyn)+size) {
3059 corr=fsize-sizeof(t_moldyn)-size;
3060 printf("[moldyn] WARNING: lsf (illegal file size)\n");
3061 printf(" modifying offset:\n");
3062 printf(" - current pos: %d\n",sizeof(t_moldyn));
3063 printf(" - atom size: %d\n",size);
3064 printf(" - file size: %d\n",fsize);
3065 printf(" => correction: %d\n",corr);
3066 lseek(fd,corr,SEEK_CUR);
3069 moldyn->atom=(t_atom *)malloc(size);
3070 if(moldyn->atom==NULL) {
3071 perror("[moldyn] load save file malloc (atoms)");
3076 cnt=read(fd,moldyn->atom,size);
3078 perror("[moldyn] load save file read (atoms)");
3085 amutex=malloc(moldyn->count*sizeof(pthread_mutex_t));
3087 perror("[moldyn] load save file (mutexes)");
3090 for(cnt=0;cnt<moldyn->count;cnt++)
3091 pthread_mutex_init(&(amutex[cnt]),NULL);
3099 int moldyn_free_save_file(t_moldyn *moldyn) {
3106 int moldyn_load(t_moldyn *moldyn) {
3114 * function to find/callback all combinations of 2 body bonds
3117 int process_2b_bonds(t_moldyn *moldyn,void *data,
3118 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3119 void *data,u8 bc)) {
3129 t_list neighbour[27];
3139 for(i=0;i<moldyn->count;i++) {
3140 /* neighbour indexing */
3141 link_cell_neighbour_index(moldyn,
3142 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
3143 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
3144 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
3149 bc=(j<lc->dnlc)?0:1;
3154 while(neighbour[j][p]!=-1) {
3156 jtom=&(moldyn->atom[neighbour[j][p]]);
3164 p=lc->subcell->list[p];
3166 this=&(neighbour[j]);
3169 if(this->start==NULL)
3174 jtom=this->current->data;
3178 process(moldyn,&(itom[i]),jtom,data,bc);
3185 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3195 * function to find neighboured atoms
3198 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
3199 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
3200 void *data,u8 bc)) {
3210 t_list neighbour[27];
3219 /* neighbour indexing */
3220 link_cell_neighbour_index(moldyn,
3221 (atom->r.x+moldyn->dim.x/2)/lc->x,
3222 (atom->r.y+moldyn->dim.y/2)/lc->x,
3223 (atom->r.z+moldyn->dim.z/2)/lc->x,
3228 bc=(j<lc->dnlc)?0:1;
3233 while(neighbour[j][p]!=-1) {
3235 natom=&(moldyn->atom[neighbour[j][p]]);
3242 natom=&(moldyn->atom[p]);
3243 p=lc->subcell->list[p];
3245 this=&(neighbour[j]);
3248 if(this->start==NULL)
3253 natom=this->current->data;
3257 process(moldyn,atom,natom,data,bc);
3264 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
3273 * post processing functions
3276 int get_line(int fd,char *line,int max) {
3283 if(count==max) return count;
3284 ret=read(fd,line+count,1);
3285 if(ret<=0) return ret;
3286 if(line[count]=='\n') {
3287 memset(line+count,0,max-count-1);
3295 int pair_correlation_init(t_moldyn *moldyn,double dr) {
3301 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
3318 for(i=0;i<moldyn->count;i++) {
3320 /* care for pb crossing */
3321 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3322 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3323 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3324 /* calculate distance */
3325 v3_sub(&dist,&final_r,&(atom[i].r_0));
3326 d2=v3_absolute_square(&dist);
3340 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
3341 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
3342 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
3347 int calculate_msd(t_moldyn *moldyn,double *msd) {
3364 for(i=0;i<moldyn->count;i++) {
3366 /* care for pb crossing */
3367 if(atom[i].pbc[0]|atom[i].pbc[1]|atom[i].pbc[2]) {
3368 printf("[moldyn] msd pb crossings for atom %d\n",i);
3369 printf(" x: %d y: %d z: %d\n",
3370 atom[i].pbc[0],atom[i].pbc[1],atom[i].pbc[2]);
3372 final_r.x=atom[i].r.x+atom[i].pbc[0]*moldyn->dim.x;
3373 final_r.y=atom[i].r.y+atom[i].pbc[1]*moldyn->dim.y;
3374 final_r.z=atom[i].r.z+atom[i].pbc[2]*moldyn->dim.z;
3375 /* calculate distance */
3376 v3_sub(&dist,&final_r,&(atom[i].r_0));
3377 d2=v3_absolute_square(&dist);
3393 msd[2]/=moldyn->count;
3398 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
3403 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
3404 t_atom *jtom,void *data,u8 bc) {
3411 /* only count pairs once,
3412 * skip same atoms */
3413 if(itom->tag>=jtom->tag)
3417 * pair correlation calc
3424 v3_sub(&dist,&(jtom->r),&(itom->r));
3425 if(bc) check_per_bound(moldyn,&dist);
3426 d=v3_absolute_square(&dist);
3428 /* ignore if greater cutoff */
3429 if(d>moldyn->cutoff_square)
3432 /* fill the slots */
3436 /* should never happen but it does 8) -
3437 * related to -ffloat-store problem! */
3439 printf("[moldyn] WARNING: pcc (%d/%d)",
3445 if(itom->brand!=jtom->brand) {
3450 /* type a - type a bonds */
3452 pcc->stat[s+pcc->o1]+=1;
3454 /* type b - type b bonds */
3455 pcc->stat[s+pcc->o2]+=1;
3461 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
3468 pcc.o1=moldyn->cutoff/dr;
3471 if(pcc.o1*dr<=moldyn->cutoff)
3472 printf("[moldyn] WARNING: pcc (low #slots)\n");
3474 printf("[moldyn] pair correlation calc info:\n");
3475 printf(" time: %f\n",moldyn->time);
3476 printf(" count: %d\n",moldyn->count);
3477 printf(" cutoff: %f\n",moldyn->cutoff);
3478 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
3481 pcc.stat=(double *)ptr;
3484 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
3485 if(pcc.stat==NULL) {
3486 perror("[moldyn] pair correlation malloc");
3491 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
3494 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
3497 for(i=1;i<pcc.o1;i++) {
3498 // normalization: 4 pi r^2 dr
3499 // here: not double counting pairs -> 2 pi r r dr
3500 // ... and actually it's a constant times r^2
3503 pcc.stat[pcc.o1+i]/=norm;
3504 pcc.stat[pcc.o2+i]/=norm;
3509 /* todo: store/print pair correlation function */
3516 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3523 if(itom->tag>=jtom->tag)
3527 v3_sub(&dist,&(jtom->r),&(itom->r));
3528 if(bc) check_per_bound(moldyn,&dist);
3529 d=v3_absolute_square(&dist);
3531 /* ignore if greater or equal cutoff */
3532 if(d>moldyn->cutoff_square)
3535 /* check for potential bond */
3536 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3539 /* now count this bonding ... */
3542 /* increase total bond counter
3547 ba->acnt[jtom->tag]+=1;
3549 ba->bcnt[jtom->tag]+=1;
3552 ba->acnt[itom->tag]+=1;
3554 ba->bcnt[itom->tag]+=1;
3559 int bond_analyze(t_moldyn *moldyn,double *quality) {
3570 ba.acnt=malloc(moldyn->count*sizeof(int));
3572 perror("[moldyn] bond analyze malloc (a)");
3575 memset(ba.acnt,0,moldyn->count*sizeof(int));
3577 ba.bcnt=malloc(moldyn->count*sizeof(int));
3579 perror("[moldyn] bond analyze malloc (b)");
3582 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3591 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3593 for(i=0;i<moldyn->count;i++) {
3594 if(atom[i].brand==0) {
3595 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3597 if((ba.acnt[i]==0)&(ba.bcnt[i]==3))
3601 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3605 if((ba.acnt[i]==3)&(ba.bcnt[i]==0)) {
3614 quality[0]=1.0*ccnt4/bcnt;
3615 quality[1]=1.0*ccnt3/bcnt;
3618 printf("[moldyn] bond analyze: %f %f\n",
3619 1.0*ccnt4/bcnt,1.0*ccnt3/bcnt);
3626 * visualization code
3629 int visual_init(t_moldyn *moldyn,char *filebase) {
3631 strncpy(moldyn->vis.fb,filebase,128);
3636 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3643 if(itom->tag>=jtom->tag)
3646 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3649 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3650 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3651 itom->r.x,itom->r.y,itom->r.z,
3652 jtom->r.x,jtom->r.y,jtom->r.z);
3657 #ifdef VISUAL_THREAD
3658 void *visual_atoms(void *ptr) {
3660 int visual_atoms(t_moldyn *moldyn) {
3671 #ifdef VISUAL_THREAD
3685 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3686 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3688 perror("open visual save file fd");
3689 #ifndef VISUAL_THREAD
3694 /* write the actual data file */
3697 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3698 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3700 // atomic configuration
3701 for(i=0;i<moldyn->count;i++) {
3702 v3_sub(&strain,&(atom[i].r),&(atom[i].r_0));
3703 check_per_bound(moldyn,&strain);
3704 // atom type, positions, color and kinetic energy
3705 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3709 pse_col[atom[i].element],
3711 sqrt(v3_absolute_square(&strain)));
3714 // bonds between atoms
3715 #ifndef VISUAL_THREAD
3716 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3721 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3722 -dim.x/2,-dim.y/2,-dim.z/2,
3723 dim.x/2,-dim.y/2,-dim.z/2);
3724 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3725 -dim.x/2,-dim.y/2,-dim.z/2,
3726 -dim.x/2,dim.y/2,-dim.z/2);
3727 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3728 dim.x/2,dim.y/2,-dim.z/2,
3729 dim.x/2,-dim.y/2,-dim.z/2);
3730 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3731 -dim.x/2,dim.y/2,-dim.z/2,
3732 dim.x/2,dim.y/2,-dim.z/2);
3734 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3735 -dim.x/2,-dim.y/2,dim.z/2,
3736 dim.x/2,-dim.y/2,dim.z/2);
3737 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3738 -dim.x/2,-dim.y/2,dim.z/2,
3739 -dim.x/2,dim.y/2,dim.z/2);
3740 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3741 dim.x/2,dim.y/2,dim.z/2,
3742 dim.x/2,-dim.y/2,dim.z/2);
3743 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3744 -dim.x/2,dim.y/2,dim.z/2,
3745 dim.x/2,dim.y/2,dim.z/2);
3747 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3748 -dim.x/2,-dim.y/2,dim.z/2,
3749 -dim.x/2,-dim.y/2,-dim.z/2);
3750 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3751 -dim.x/2,dim.y/2,dim.z/2,
3752 -dim.x/2,dim.y/2,-dim.z/2);
3753 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3754 dim.x/2,-dim.y/2,dim.z/2,
3755 dim.x/2,-dim.y/2,-dim.z/2);
3756 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3757 dim.x/2,dim.y/2,dim.z/2,
3758 dim.x/2,dim.y/2,-dim.z/2);
3763 #ifdef VISUAL_THREAD
3774 * fpu cntrol functions
3777 // set rounding to double (eliminates -ffloat-store!)
3778 int fpu_set_rtd(void) {
3784 ctrl&=~_FPU_EXTENDED;