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"
52 pthread_mutex_t *amutex;
53 pthread_mutex_t emutex;
57 * the moldyn functions
60 int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
62 printf("[moldyn] init\n");
64 /* only needed if compiled without -msse2 (float-store prob!) */
67 memset(moldyn,0,sizeof(t_moldyn));
72 rand_init(&(moldyn->random),NULL,1);
73 moldyn->random.status|=RAND_STAT_VERBOSE;
76 pthread_mutex_init(&emutex,NULL);
82 int moldyn_shutdown(t_moldyn *moldyn) {
88 printf("[moldyn] shutdown\n");
91 for(i=0;i<moldyn->count;i++)
92 pthread_mutex_destroy(&(amutex[i]));
94 pthread_mutex_destroy(&emutex);
97 moldyn_log_shutdown(moldyn);
98 link_cell_shutdown(moldyn);
99 rand_close(&(moldyn->random));
105 int set_int_alg(t_moldyn *moldyn,u8 algo) {
107 printf("[moldyn] integration algorithm: ");
110 case MOLDYN_INTEGRATE_VERLET:
111 moldyn->integrate=velocity_verlet;
112 printf("velocity verlet\n");
115 printf("unknown integration algorithm: %02x\n",algo);
123 int set_cutoff(t_moldyn *moldyn,double cutoff) {
125 moldyn->cutoff=cutoff;
126 moldyn->cutoff_square=cutoff*cutoff;
128 printf("[moldyn] cutoff [A]: %f\n",moldyn->cutoff);
133 int set_temperature(t_moldyn *moldyn,double t_ref) {
137 printf("[moldyn] temperature [K]: %f\n",moldyn->t_ref);
142 int set_pressure(t_moldyn *moldyn,double p_ref) {
146 printf("[moldyn] pressure [bar]: %f\n",moldyn->p_ref/BAR);
151 int set_p_scale(t_moldyn *moldyn,u8 ptype,double ptc) {
153 moldyn->pt_scale&=(~(P_SCALE_MASK));
154 moldyn->pt_scale|=ptype;
157 printf("[moldyn] p scaling:\n");
159 printf(" p: %s",ptype?"yes":"no ");
161 printf(" | type: %02x | factor: %f",ptype,ptc);
167 int set_t_scale(t_moldyn *moldyn,u8 ttype,double ttc) {
169 moldyn->pt_scale&=(~(T_SCALE_MASK));
170 moldyn->pt_scale|=ttype;
173 printf("[moldyn] t scaling:\n");
175 printf(" t: %s",ttype?"yes":"no ");
177 printf(" | type: %02x | factor: %f",ttype,ttc);
183 int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
185 moldyn->pt_scale=(ptype|ttype);
189 printf("[moldyn] p/t scaling:\n");
191 printf(" p: %s",ptype?"yes":"no ");
193 printf(" | type: %02x | factor: %f",ptype,ptc);
196 printf(" t: %s",ttype?"yes":"no ");
198 printf(" | type: %02x | factor: %f",ttype,ttc);
204 int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
210 moldyn->volume=x*y*z;
218 printf("[moldyn] dimensions in A and A^3 respectively:\n");
219 printf(" x: %f\n",moldyn->dim.x);
220 printf(" y: %f\n",moldyn->dim.y);
221 printf(" z: %f\n",moldyn->dim.z);
222 printf(" volume: %f\n",moldyn->volume);
223 printf(" visualize simulation box: %s\n",visualize?"yes":"no");
228 int set_nn_dist(t_moldyn *moldyn,double dist) {
235 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
237 printf("[moldyn] periodic boundary conditions:\n");
240 moldyn->status|=MOLDYN_STAT_PBX;
243 moldyn->status|=MOLDYN_STAT_PBY;
246 moldyn->status|=MOLDYN_STAT_PBZ;
248 printf(" x: %s\n",x?"yes":"no");
249 printf(" y: %s\n",y?"yes":"no");
250 printf(" z: %s\n",z?"yes":"no");
255 int set_potential(t_moldyn *moldyn,u8 type) {
258 case MOLDYN_POTENTIAL_TM:
259 moldyn->func1b=tersoff_mult_1bp;
260 moldyn->func3b_j1=tersoff_mult_3bp_j1;
261 moldyn->func3b_k1=tersoff_mult_3bp_k1;
262 moldyn->func3b_j2=tersoff_mult_3bp_j2;
263 moldyn->func3b_k2=tersoff_mult_3bp_k2;
264 moldyn->check_2b_bond=tersoff_mult_check_2b_bond;
266 case MOLDYN_POTENTIAL_AM:
267 moldyn->func3b_j1=albe_mult_3bp_j1;
268 moldyn->func3b_k1=albe_mult_3bp_k1;
269 moldyn->func3b_j2=albe_mult_3bp_j2;
270 moldyn->func3b_k2=albe_mult_3bp_k2;
271 moldyn->check_2b_bond=albe_mult_check_2b_bond;
273 case MOLDYN_POTENTIAL_HO:
274 moldyn->func2b=harmonic_oscillator;
275 moldyn->check_2b_bond=harmonic_oscillator_check_2b_bond;
277 case MOLDYN_POTENTIAL_LJ:
278 moldyn->func2b=lennard_jones;
279 moldyn->check_2b_bond=lennard_jones_check_2b_bond;
282 printf("[moldyn] set potential: unknown type %02x\n",
290 int set_avg_skip(t_moldyn *moldyn,int skip) {
292 printf("[moldyn] skip %d steps before starting average calc\n",skip);
293 moldyn->avg_skip=skip;
298 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir) {
300 strncpy(moldyn->vlsdir,dir,127);
305 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title) {
307 strncpy(moldyn->rauthor,author,63);
308 strncpy(moldyn->rtitle,title,63);
313 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
318 printf("[moldyn] set log: ");
321 case LOG_TOTAL_ENERGY:
322 moldyn->ewrite=timer;
323 snprintf(filename,127,"%s/energy",moldyn->vlsdir);
324 moldyn->efd=open(filename,
325 O_WRONLY|O_CREAT|O_EXCL,
328 perror("[moldyn] energy log fd open");
331 dprintf(moldyn->efd,"# total energy log file\n");
332 printf("total energy (%d)\n",timer);
334 case LOG_TOTAL_MOMENTUM:
335 moldyn->mwrite=timer;
336 snprintf(filename,127,"%s/momentum",moldyn->vlsdir);
337 moldyn->mfd=open(filename,
338 O_WRONLY|O_CREAT|O_EXCL,
341 perror("[moldyn] momentum log fd open");
344 dprintf(moldyn->efd,"# total momentum log file\n");
345 printf("total momentum (%d)\n",timer);
348 moldyn->pwrite=timer;
349 snprintf(filename,127,"%s/pressure",moldyn->vlsdir);
350 moldyn->pfd=open(filename,
351 O_WRONLY|O_CREAT|O_EXCL,
354 perror("[moldyn] pressure log file\n");
357 dprintf(moldyn->pfd,"# pressure log file\n");
358 printf("pressure (%d)\n",timer);
360 case LOG_TEMPERATURE:
361 moldyn->twrite=timer;
362 snprintf(filename,127,"%s/temperature",moldyn->vlsdir);
363 moldyn->tfd=open(filename,
364 O_WRONLY|O_CREAT|O_EXCL,
367 perror("[moldyn] temperature log file\n");
370 dprintf(moldyn->tfd,"# temperature log file\n");
371 printf("temperature (%d)\n",timer);
374 moldyn->vwrite=timer;
375 snprintf(filename,127,"%s/volume",moldyn->vlsdir);
376 moldyn->vfd=open(filename,
377 O_WRONLY|O_CREAT|O_EXCL,
380 perror("[moldyn] volume log file\n");
383 dprintf(moldyn->vfd,"# volume log file\n");
384 printf("volume (%d)\n",timer);
387 moldyn->swrite=timer;
388 printf("save file (%d)\n",timer);
391 moldyn->awrite=timer;
392 ret=visual_init(moldyn,moldyn->vlsdir);
394 printf("[moldyn] visual init failure\n");
397 printf("visual file (%d)\n",timer);
400 snprintf(filename,127,"%s/report.tex",moldyn->vlsdir);
401 moldyn->rfd=open(filename,
402 O_WRONLY|O_CREAT|O_EXCL,
405 perror("[moldyn] report fd open");
408 printf("report -> ");
410 snprintf(filename,127,"%s/e_plot.scr",
412 moldyn->epfd=open(filename,
413 O_WRONLY|O_CREAT|O_EXCL,
416 perror("[moldyn] energy plot fd open");
419 dprintf(moldyn->epfd,e_plot_script);
424 snprintf(filename,127,"%s/pressure_plot.scr",
426 moldyn->ppfd=open(filename,
427 O_WRONLY|O_CREAT|O_EXCL,
430 perror("[moldyn] p plot fd open");
433 dprintf(moldyn->ppfd,pressure_plot_script);
438 snprintf(filename,127,"%s/temperature_plot.scr",
440 moldyn->tpfd=open(filename,
441 O_WRONLY|O_CREAT|O_EXCL,
444 perror("[moldyn] t plot fd open");
447 dprintf(moldyn->tpfd,temperature_plot_script);
449 printf("temperature ");
451 dprintf(moldyn->rfd,report_start,
452 moldyn->rauthor,moldyn->rtitle);
456 printf("unknown log type: %02x\n",type);
463 int moldyn_log_shutdown(t_moldyn *moldyn) {
467 printf("[moldyn] log shutdown\n");
471 dprintf(moldyn->rfd,report_energy);
472 snprintf(sc,255,"cd %s && gnuplot e_plot.scr",
477 if(moldyn->mfd) close(moldyn->mfd);
481 dprintf(moldyn->rfd,report_pressure);
482 snprintf(sc,255,"cd %s && gnuplot pressure_plot.scr",
489 dprintf(moldyn->rfd,report_temperature);
490 snprintf(sc,255,"cd %s && gnuplot temperature_plot.scr",
495 dprintf(moldyn->rfd,report_end);
497 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
500 snprintf(sc,255,"cd %s && pdflatex report >/dev/null 2>&1",
503 snprintf(sc,255,"cd %s && dvipdf report >/dev/null 2>&1",
512 * creating lattice functions
515 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
516 u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin,
517 u8 p_type,t_part_vals *p_vals) {
526 pthread_mutex_t *mutex;
532 /* how many atoms do we expect */
535 printf("[moldyn] WARNING: create 'none' lattice called");
537 if(type==CUBIC) new*=1;
538 if(type==FCC) new*=4;
539 if(type==DIAMOND) new*=8;
541 /* allocate space for atoms */
542 ptr=realloc(moldyn->atom,(count+new)*sizeof(t_atom));
544 perror("[moldyn] realloc (create lattice)");
548 atom=&(moldyn->atom[count]);
551 ptr=realloc(amutex,(count+new)*sizeof(pthread_mutex_t));
553 perror("[moldyn] mutex realloc (add atom)");
557 mutex=&(amutex[count]);
560 /* no atoms on the boundaries (only reason: it looks better!) */
574 set_nn_dist(moldyn,lc);
575 ret=cubic_init(a,b,c,lc,atom,&orig,p_type,p_vals);
576 strcpy(name,"cubic");
580 v3_scale(&orig,&orig,0.5);
581 set_nn_dist(moldyn,0.5*sqrt(2.0)*lc);
582 ret=fcc_init(a,b,c,lc,atom,&orig,p_type,p_vals);
587 v3_scale(&orig,&orig,0.25);
588 set_nn_dist(moldyn,0.25*sqrt(3.0)*lc);
589 ret=diamond_init(a,b,c,lc,atom,&orig,p_type,p_vals);
590 strcpy(name,"diamond");
593 printf("unknown lattice type (%02x)\n",type);
599 printf("[moldyn] creating lattice failed\n");
600 printf(" amount of atoms\n");
601 printf(" - expected: %d\n",new);
602 printf(" - created: %d\n",ret);
607 printf("[moldyn] created %s lattice with %d atoms\n",name,new);
609 for(ret=0;ret<new;ret++) {
610 atom[ret].element=element;
613 atom[ret].brand=brand;
614 atom[ret].tag=count+ret;
615 check_per_bound(moldyn,&(atom[ret].r));
616 atom[ret].r_0=atom[ret].r;
618 pthread_mutex_init(&(mutex[ret]),NULL);
622 /* update total system mass */
623 total_mass_calc(moldyn);
628 int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
629 t_3dvec *r,t_3dvec *v) {
636 count=(moldyn->count)++; // asshole style!
638 ptr=realloc(atom,(count+1)*sizeof(t_atom));
640 perror("[moldyn] realloc (add atom)");
646 ptr=realloc(moldyn->lc.subcell->list,(count+1)*sizeof(int));
648 perror("[moldyn] list realloc (add atom)");
651 moldyn->lc.subcell->list=ptr;
655 ptr=realloc(amutex,(count+1)*sizeof(pthread_mutex_t));
657 perror("[moldyn] mutex realloc (add atom)");
661 pthread_mutex_init(&(amutex[count]),NULL);
666 /* initialize new atom */
667 memset(&(atom[count]),0,sizeof(t_atom));
670 atom[count].element=element;
671 atom[count].mass=mass;
672 atom[count].brand=brand;
673 atom[count].tag=count;
674 atom[count].attr=attr;
675 check_per_bound(moldyn,&(atom[count].r));
676 atom[count].r_0=atom[count].r;
678 /* update total system mass */
679 total_mass_calc(moldyn);
684 int del_atom(t_moldyn *moldyn,int tag) {
691 new=(t_atom *)malloc((moldyn->count-1)*sizeof(t_atom));
693 perror("[moldyn]malloc (del atom)");
697 for(cnt=0;cnt<tag;cnt++)
700 for(cnt=tag+1;cnt<moldyn->count;cnt++) {
702 new[cnt-1].tag=cnt-1;
714 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
715 u8 p_type,t_part_vals *p_vals) {
729 /* shift partition values */
731 p_vals->p.x+=(a*lc)/2.0;
732 p_vals->p.y+=(b*lc)/2.0;
733 p_vals->p.z+=(c*lc)/2.0;
744 v3_sub(&dist,&r,&(p_vals->p));
745 if(v3_absolute_square(&dist)<(p_vals->r*p_vals->r)) {
746 v3_copy(&(atom[count].r),&r);
751 v3_sub(&dist,&r,&(p_vals->p));
752 if(v3_absolute_square(&dist)>=(p_vals->r*p_vals->r)) {
753 v3_copy(&(atom[count].r),&r);
758 v3_copy(&(atom[count].r),&r);
769 for(i=0;i<count;i++) {
770 atom[i].r.x-=(a*lc)/2.0;
771 atom[i].r.y-=(b*lc)/2.0;
772 atom[i].r.z-=(c*lc)/2.0;
778 /* fcc lattice init */
779 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
780 u8 p_type,t_part_vals *p_vals) {
794 /* construct the basis */
795 memset(basis,0,3*sizeof(t_3dvec));
803 /* shift partition values */
805 p_vals->p.x+=(a*lc)/2.0;
806 p_vals->p.y+=(b*lc)/2.0;
807 p_vals->p.z+=(c*lc)/2.0;
810 /* fill up the room */
820 v3_sub(&dist,&r,&(p_vals->p));
821 if(v3_absolute_square(&dist)<(p_vals->r*p_vals->r)) {
822 v3_copy(&(atom[count].r),&r);
827 v3_sub(&dist,&r,&(p_vals->p));
828 if(v3_absolute_square(&dist)>=(p_vals->r*p_vals->r)) {
829 v3_copy(&(atom[count].r),&r);
834 v3_copy(&(atom[count].r),&r);
838 /* the three face centered atoms */
840 v3_add(&n,&r,&basis[l]);
843 v3_sub(&dist,&n,&(p_vals->p));
844 if(v3_absolute_square(&dist)<(p_vals->r*p_vals->r)) {
845 v3_copy(&(atom[count].r),&r);
850 v3_sub(&dist,&n,&(p_vals->p));
851 if(v3_absolute_square(&dist)>=(p_vals->r*p_vals->r)) {
852 v3_copy(&(atom[count].r),&r);
857 v3_copy(&(atom[count].r),&n);
869 /* coordinate transformation */
870 for(i=0;i<count;i++) {
871 atom[i].r.x-=(a*lc)/2.0;
872 atom[i].r.y-=(b*lc)/2.0;
873 atom[i].r.z-=(c*lc)/2.0;
879 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin,
880 u8 p_type,t_part_vals *p_vals) {
885 count=fcc_init(a,b,c,lc,atom,origin,p_type,p_vals);
891 if(origin) v3_add(&o,&o,origin);
893 count+=fcc_init(a,b,c,lc,&atom[count],&o,p_type,p_vals);
898 int destroy_atoms(t_moldyn *moldyn) {
900 if(moldyn->atom) free(moldyn->atom);
905 int thermal_init(t_moldyn *moldyn,u8 equi_init) {
908 * - gaussian distribution of velocities
909 * - zero total momentum
910 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
915 t_3dvec p_total,delta;
920 random=&(moldyn->random);
922 printf("[moldyn] thermal init (equi init: %s)\n",equi_init?"yes":"no");
924 /* gaussian distribution of velocities */
926 for(i=0;i<moldyn->count;i++) {
927 sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
929 v=sigma*rand_get_gauss(random);
931 p_total.x+=atom[i].mass*v;
933 v=sigma*rand_get_gauss(random);
935 p_total.y+=atom[i].mass*v;
937 v=sigma*rand_get_gauss(random);
939 p_total.z+=atom[i].mass*v;
942 /* zero total momentum */
943 v3_scale(&p_total,&p_total,1.0/moldyn->count);
944 for(i=0;i<moldyn->count;i++) {
945 v3_scale(&delta,&p_total,1.0/atom[i].mass);
946 v3_sub(&(atom[i].v),&(atom[i].v),&delta);
949 /* velocity scaling */
950 scale_velocity(moldyn,equi_init);
955 double total_mass_calc(t_moldyn *moldyn) {
961 for(i=0;i<moldyn->count;i++)
962 moldyn->mass+=moldyn->atom[i].mass;
967 double temperature_calc(t_moldyn *moldyn) {
969 /* assume up to date kinetic energy, which is 3/2 N k_B T */
971 moldyn->t=(2.0*moldyn->ekin)/(3.0*K_BOLTZMANN*moldyn->count);
976 double get_temperature(t_moldyn *moldyn) {
981 int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
991 * - velocity scaling (E = 3/2 N k T), E: kinetic energy
994 /* get kinetic energy / temperature & count involved atoms */
997 for(i=0;i<moldyn->count;i++) {
998 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
999 e+=atom[i].mass*v3_absolute_square(&(atom[i].v));
1004 if(count!=0) moldyn->t=e/(1.5*count*K_BOLTZMANN);
1005 else return 0; /* no atoms involved in scaling! */
1007 /* (temporary) hack for e,t = 0 */
1010 if(moldyn->t_ref!=0.0) {
1011 thermal_init(moldyn,equi_init);
1015 return 0; /* no scaling needed */
1019 /* get scaling factor */
1020 scale=moldyn->t_ref/moldyn->t;
1024 if(moldyn->pt_scale&T_SCALE_BERENDSEN)
1025 scale=1.0+(scale-1.0)*moldyn->tau/moldyn->t_tc;
1028 /* velocity scaling */
1029 for(i=0;i<moldyn->count;i++) {
1030 if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
1031 v3_scale(&(atom[i].v),&(atom[i].v),scale);
1037 double ideal_gas_law_pressure(t_moldyn *moldyn) {
1041 p=moldyn->count*moldyn->t*K_BOLTZMANN/moldyn->volume;
1046 double virial_sum(t_moldyn *moldyn) {
1051 /* virial (sum over atom virials) */
1059 for(i=0;i<moldyn->count;i++) {
1060 virial=&(moldyn->atom[i].virial);
1061 moldyn->virial+=(virial->xx+virial->yy+virial->zz);
1062 moldyn->vir.xx+=virial->xx;
1063 moldyn->vir.yy+=virial->yy;
1064 moldyn->vir.zz+=virial->zz;
1065 moldyn->vir.xy+=virial->xy;
1066 moldyn->vir.xz+=virial->xz;
1067 moldyn->vir.yz+=virial->yz;
1070 /* global virial (absolute coordinates) */
1071 virial=&(moldyn->gvir);
1072 moldyn->gv=virial->xx+virial->yy+virial->zz;
1074 return moldyn->virial;
1077 double pressure_calc(t_moldyn *moldyn) {
1081 * with W = 1/3 sum_i f_i r_i (- skipped!)
1082 * virial = sum_i f_i r_i
1084 * => P = (2 Ekin + virial) / (3V)
1087 /* assume up to date virial & up to date kinetic energy */
1089 /* pressure (atom virials) */
1090 moldyn->p=2.0*moldyn->ekin+moldyn->virial;
1091 moldyn->p/=(3.0*moldyn->volume);
1093 /* pressure (absolute coordinates) */
1094 moldyn->gp=2.0*moldyn->ekin+moldyn->gv;
1095 moldyn->gp/=(3.0*moldyn->volume);
1100 int average_reset(t_moldyn *moldyn) {
1102 printf("[moldyn] average reset\n");
1104 /* update skip value */
1105 moldyn->avg_skip=moldyn->total_steps;
1107 /* kinetic energy */
1111 /* potential energy */
1119 moldyn->virial_sum=0.0;
1130 int average_and_fluctuation_calc(t_moldyn *moldyn) {
1134 if(moldyn->total_steps<moldyn->avg_skip)
1137 denom=moldyn->total_steps+1-moldyn->avg_skip;
1139 /* assume up to date energies, temperature, pressure etc */
1141 /* kinetic energy */
1142 moldyn->k_sum+=moldyn->ekin;
1143 moldyn->k2_sum+=(moldyn->ekin*moldyn->ekin);
1144 moldyn->k_avg=moldyn->k_sum/denom;
1145 moldyn->k2_avg=moldyn->k2_sum/denom;
1146 moldyn->dk2_avg=moldyn->k2_avg-(moldyn->k_avg*moldyn->k_avg);
1148 /* potential energy */
1149 moldyn->v_sum+=moldyn->energy;
1150 moldyn->v2_sum+=(moldyn->energy*moldyn->energy);
1151 moldyn->v_avg=moldyn->v_sum/denom;
1152 moldyn->v2_avg=moldyn->v2_sum/denom;
1153 moldyn->dv2_avg=moldyn->v2_avg-(moldyn->v_avg*moldyn->v_avg);
1156 moldyn->t_sum+=moldyn->t;
1157 moldyn->t_avg=moldyn->t_sum/denom;
1160 moldyn->virial_sum+=moldyn->virial;
1161 moldyn->virial_avg=moldyn->virial_sum/denom;
1162 moldyn->gv_sum+=moldyn->gv;
1163 moldyn->gv_avg=moldyn->gv_sum/denom;
1166 moldyn->p_sum+=moldyn->p;
1167 moldyn->p_avg=moldyn->p_sum/denom;
1168 moldyn->gp_sum+=moldyn->gp;
1169 moldyn->gp_avg=moldyn->gp_sum/denom;
1170 moldyn->tp_sum+=moldyn->tp;
1171 moldyn->tp_avg=moldyn->tp_sum/denom;
1176 int get_heat_capacity(t_moldyn *moldyn) {
1180 /* averages needed for heat capacity calc */
1181 if(moldyn->total_steps<moldyn->avg_skip)
1184 /* (temperature average)^2 */
1185 temp2=moldyn->t_avg*moldyn->t_avg;
1186 printf("[moldyn] specific heat capacity for T=%f K [J/(kg K)]\n",
1189 /* ideal gas contribution */
1190 ighc=3.0*moldyn->count*K_BOLTZMANN/2.0;
1191 printf(" ideal gas contribution: %f\n",
1192 ighc/moldyn->mass*KILOGRAM/JOULE);
1194 /* specific heat for nvt ensemble */
1195 moldyn->c_v_nvt=moldyn->dv2_avg/(K_BOLTZMANN*temp2)+ighc;
1196 moldyn->c_v_nvt/=moldyn->mass;
1198 /* specific heat for nve ensemble */
1199 moldyn->c_v_nve=ighc/(1.0-(moldyn->dv2_avg/(ighc*K_BOLTZMANN*temp2)));
1200 moldyn->c_v_nve/=moldyn->mass;
1202 printf(" NVE: %f\n",moldyn->c_v_nve*KILOGRAM/JOULE);
1203 printf(" NVT: %f\n",moldyn->c_v_nvt*KILOGRAM/JOULE);
1204 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)));
1209 double thermodynamic_pressure_calc(t_moldyn *moldyn) {
1225 /* store atomic configuration + dimension */
1226 store=malloc(moldyn->count*sizeof(t_atom));
1228 printf("[moldyn] allocating store mem failed\n");
1231 memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
1236 h=(1.0-sd)*(1.0-sd)*(1.0-sd);
1237 su=pow(2.0-h,ONE_THIRD)-1.0;
1238 dv=(1.0-h)*moldyn->volume;
1240 /* scale up dimension and atom positions */
1241 scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1242 scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
1243 link_cell_shutdown(moldyn);
1244 link_cell_init(moldyn,QUIET);
1245 potential_force_calc(moldyn);
1248 /* restore atomic configuration + dim */
1249 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1252 /* scale down dimension and atom positions */
1253 scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1254 scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
1255 link_cell_shutdown(moldyn);
1256 link_cell_init(moldyn,QUIET);
1257 potential_force_calc(moldyn);
1260 /* calculate pressure */
1261 moldyn->tp=-(y1-y0)/(2.0*dv);
1263 /* restore atomic configuration */
1264 memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
1266 link_cell_shutdown(moldyn);
1267 link_cell_init(moldyn,QUIET);
1268 //potential_force_calc(moldyn);
1270 /* free store buffer */
1277 double get_pressure(t_moldyn *moldyn) {
1283 int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1295 if(x) dim->x*=scale;
1296 if(y) dim->y*=scale;
1297 if(z) dim->z*=scale;
1302 int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z) {
1313 for(i=0;i<moldyn->count;i++) {
1314 r=&(moldyn->atom[i].r);
1323 int scale_volume(t_moldyn *moldyn) {
1329 vdim=&(moldyn->vis.dim);
1333 /* scaling factor */
1334 if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
1335 scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc*moldyn->tau;
1336 scale=pow(scale,ONE_THIRD);
1339 scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
1342 /* scale the atoms and dimensions */
1343 scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1344 scale_dim(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
1346 /* visualize dimensions */
1353 /* recalculate scaled volume */
1354 moldyn->volume=dim->x*dim->y*dim->z;
1356 /* adjust/reinit linkcell */
1357 if(((int)(dim->x/moldyn->cutoff)!=lc->nx)||
1358 ((int)(dim->y/moldyn->cutoff)!=lc->ny)||
1359 ((int)(dim->z/moldyn->cutoff)!=lc->nx)) {
1360 link_cell_shutdown(moldyn);
1361 link_cell_init(moldyn,QUIET);
1372 double e_kin_calc(t_moldyn *moldyn) {
1380 for(i=0;i<moldyn->count;i++) {
1381 atom[i].ekin=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
1382 moldyn->ekin+=atom[i].ekin;
1385 return moldyn->ekin;
1388 double get_total_energy(t_moldyn *moldyn) {
1390 return(moldyn->ekin+moldyn->energy);
1393 t_3dvec get_total_p(t_moldyn *moldyn) {
1402 for(i=0;i<moldyn->count;i++) {
1403 v3_scale(&p,&(atom[i].v),atom[i].mass);
1404 v3_add(&p_total,&p_total,&p);
1410 double estimate_time_step(t_moldyn *moldyn,double nn_dist) {
1414 /* nn_dist is the nearest neighbour distance */
1416 tau=(0.05*nn_dist*moldyn->atom[0].mass)/sqrt(3.0*K_BOLTZMANN*moldyn->t);
1425 /* linked list / cell method */
1427 int link_cell_init(t_moldyn *moldyn,u8 vol) {
1430 #ifndef LOWMEM_LISTS
1436 /* partitioning the md cell */
1437 lc->nx=moldyn->dim.x/moldyn->cutoff;
1438 lc->x=moldyn->dim.x/lc->nx;
1439 lc->ny=moldyn->dim.y/moldyn->cutoff;
1440 lc->y=moldyn->dim.y/lc->ny;
1441 lc->nz=moldyn->dim.z/moldyn->cutoff;
1442 lc->z=moldyn->dim.z/lc->nz;
1443 lc->cells=lc->nx*lc->ny*lc->nz;
1446 lc->subcell=malloc(lc->cells*sizeof(int*));
1448 lc->subcell=malloc(sizeof(t_lowmem_list));
1450 lc->subcell=malloc(lc->cells*sizeof(t_list));
1453 if(lc->subcell==NULL) {
1454 perror("[moldyn] cell init (malloc)");
1459 printf("[moldyn] FATAL: less then 27 subcells! (%d)\n",
1464 printf("[moldyn] initializing 'static' linked cells (%d)\n",
1467 printf("[moldyn] initializing 'lowmem' linked cells (%d)\n",
1470 printf("[moldyn] initializing 'dynamic' linked cells (%d)\n",
1473 printf(" x: %d x %f A\n",lc->nx,lc->x);
1474 printf(" y: %d x %f A\n",lc->ny,lc->y);
1475 printf(" z: %d x %f A\n",lc->nz,lc->z);
1480 for(i=0;i<lc->cells;i++) {
1481 lc->subcell[i]=malloc((MAX_ATOMS_PER_LIST+1)*sizeof(int));
1482 if(lc->subcell[i]==NULL) {
1483 perror("[moldyn] list init (malloc)");
1488 printf(" ---> %d malloc %p (%p)\n",
1489 i,lc->subcell[0],lc->subcell);
1493 lc->subcell->head=malloc(lc->cells*sizeof(int));
1494 if(lc->subcell->head==NULL) {
1495 perror("[moldyn] head init (malloc)");
1498 lc->subcell->list=malloc(moldyn->count*sizeof(int));
1499 if(lc->subcell->list==NULL) {
1500 perror("[moldyn] list init (malloc)");
1504 for(i=0;i<lc->cells;i++)
1505 list_init_f(&(lc->subcell[i]));
1508 /* update the list */
1509 link_cell_update(moldyn);
1514 int link_cell_update(t_moldyn *moldyn) {
1532 for(i=0;i<lc->cells;i++)
1534 memset(lc->subcell[i],-1,(MAX_ATOMS_PER_LIST+1)*sizeof(int));
1536 lc->subcell->head[i]=-1;
1538 list_destroy_f(&(lc->subcell[i]));
1541 for(count=0;count<moldyn->count;count++) {
1542 i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
1543 j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
1544 k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
1548 while(lc->subcell[i+j*nx+k*nxy][p]!=-1)
1551 if(p>=MAX_ATOMS_PER_LIST) {
1552 printf("[moldyn] FATAL: amount of atoms too high!\n");
1556 lc->subcell[i+j*nx+k*nxy][p]=count;
1559 lc->subcell->list[count]=lc->subcell->head[p];
1560 lc->subcell->head[p]=count;
1562 list_add_immediate_f(&(lc->subcell[i+j*nx+k*nxy]),
1566 printf(" ---> %d %d malloc %p (%p)\n",
1567 i,count,lc->subcell[i].current,lc->subcell);
1575 int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
1601 if(i>=nx||j>=ny||k>=nz)
1602 printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
1605 #ifndef LOWMEM_LISTS
1606 cell[0]=lc->subcell[i+j*nx+k*a];
1608 cell[0]=lc->subcell->head[i+j*nx+k*a];
1610 for(ci=-1;ci<=1;ci++) {
1613 if((x<0)||(x>=nx)) {
1617 for(cj=-1;cj<=1;cj++) {
1620 if((y<0)||(y>=ny)) {
1624 for(ck=-1;ck<=1;ck++) {
1627 if((z<0)||(z>=nz)) {
1631 if(!(ci|cj|ck)) continue;
1633 #ifndef LOWMEM_LISTS
1634 cell[--count2]=lc->subcell[x+y*nx+z*a];
1636 cell[--count2]=lc->subcell->head[x+y*nx+z*a];
1641 #ifndef LOWMEM_LISTS
1642 cell[count1++]=lc->subcell[x+y*nx+z*a];
1644 cell[count1++]=lc->subcell->head[x+y*nx+z*a];
1656 int link_cell_shutdown(t_moldyn *moldyn) {
1658 #ifndef LOWMEM_LISTS
1666 free(lc->subcell->head);
1667 free(lc->subcell->list);
1671 for(i=0;i<lc->cells;i++) {
1673 free(lc->subcell[i]);
1675 //printf(" ---> %d free %p\n",i,lc->subcell[i].start);
1676 list_destroy_f(&(lc->subcell[i]));
1686 int moldyn_add_schedule(t_moldyn *moldyn,int runs,double tau) {
1690 t_moldyn_schedule *schedule;
1692 schedule=&(moldyn->schedule);
1693 count=++(schedule->total_sched);
1695 ptr=realloc(schedule->runs,count*sizeof(int));
1697 perror("[moldyn] realloc (runs)");
1701 schedule->runs[count-1]=runs;
1703 ptr=realloc(schedule->tau,count*sizeof(double));
1705 perror("[moldyn] realloc (tau)");
1709 schedule->tau[count-1]=tau;
1711 printf("[moldyn] schedule added:\n");
1712 printf(" number: %d | runs: %d | tau: %f\n",count-1,runs,tau);
1718 int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
1720 moldyn->schedule.hook=hook;
1721 moldyn->schedule.hook_params=hook_params;
1728 * 'integration of newtons equation' - algorithms
1732 /* start the integration */
1734 int moldyn_integrate(t_moldyn *moldyn) {
1737 unsigned int e,m,s,v,p,t,a;
1739 t_moldyn_schedule *sched;
1744 double energy_scale;
1745 struct timeval t1,t2;
1748 #ifdef VISUAL_THREAD
1750 pthread_t io_thread;
1759 sched=&(moldyn->schedule);
1762 /* initialize linked cell method */
1763 link_cell_init(moldyn,VERBOSE);
1765 /* logging & visualization */
1774 /* sqaure of some variables */
1775 moldyn->tau_square=moldyn->tau*moldyn->tau;
1777 /* get current time */
1778 gettimeofday(&t1,NULL);
1780 /* calculate initial forces */
1781 potential_force_calc(moldyn);
1786 /* some stupid checks before we actually start calculating bullshit */
1787 if(moldyn->cutoff>0.5*moldyn->dim.x)
1788 printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
1789 if(moldyn->cutoff>0.5*moldyn->dim.y)
1790 printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
1791 if(moldyn->cutoff>0.5*moldyn->dim.z)
1792 printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
1794 ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
1795 if(ds>0.05*moldyn->nnd)
1796 printf("[moldyn] WARNING: forces too high / tau too small!\n");
1799 /* zero absolute time */
1800 // should have right values!
1802 //moldyn->total_steps=0;
1804 /* debugging, ignore */
1807 /* zero & init moldyn copy */
1808 #ifdef VISUAL_THREAD
1809 memset(&md_copy,0,sizeof(t_moldyn));
1810 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1811 if(atom_copy==NULL) {
1812 perror("[moldyn] malloc atom copy (init)");
1818 printf("##################\n");
1819 printf("# USING PTHREADS #\n");
1820 printf("##################\n");
1822 /* tell the world */
1823 printf("[moldyn] integration start, go get a coffee ...\n");
1825 /* executing the schedule */
1827 while(sched->count<sched->total_sched) {
1829 /* setting amount of runs and finite time step size */
1830 moldyn->tau=sched->tau[sched->count];
1831 moldyn->tau_square=moldyn->tau*moldyn->tau;
1832 moldyn->time_steps=sched->runs[sched->count];
1834 /* energy scaling factor (might change!) */
1835 energy_scale=moldyn->count*EV;
1837 /* integration according to schedule */
1839 for(i=0;i<moldyn->time_steps;i++) {
1841 /* integration step */
1842 moldyn->integrate(moldyn);
1844 /* calculate kinetic energy, temperature and pressure */
1846 temperature_calc(moldyn);
1848 pressure_calc(moldyn);
1850 thermodynamic_pressure_calc(moldyn);
1851 printf("\n\nDEBUG: numeric pressure calc: %f\n\n",
1855 /* calculate fluctuations + averages */
1856 average_and_fluctuation_calc(moldyn);
1859 if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
1860 scale_velocity(moldyn,FALSE);
1861 if(moldyn->pt_scale&(P_SCALE_BERENDSEN|P_SCALE_DIRECT))
1862 scale_volume(moldyn);
1864 /* check for log & visualization */
1866 if(!(moldyn->total_steps%e))
1867 dprintf(moldyn->efd,
1869 moldyn->time,moldyn->ekin/energy_scale,
1870 moldyn->energy/energy_scale,
1871 get_total_energy(moldyn)/energy_scale);
1874 if(!(moldyn->total_steps%m)) {
1875 momentum=get_total_p(moldyn);
1876 dprintf(moldyn->mfd,
1877 "%f %f %f %f %f\n",moldyn->time,
1878 momentum.x,momentum.y,momentum.z,
1879 v3_norm(&momentum));
1883 if(!(moldyn->total_steps%p)) {
1884 dprintf(moldyn->pfd,
1885 "%f %f %f %f %f %f %f\n",moldyn->time,
1886 moldyn->p/BAR,moldyn->p_avg/BAR,
1887 moldyn->gp/BAR,moldyn->gp_avg/BAR,
1888 moldyn->tp/BAR,moldyn->tp_avg/BAR);
1892 if(!(moldyn->total_steps%t)) {
1893 dprintf(moldyn->tfd,
1895 moldyn->time,moldyn->t,moldyn->t_avg);
1899 if(!(moldyn->total_steps%v)) {
1900 dprintf(moldyn->vfd,
1901 "%f %f\n",moldyn->time,moldyn->volume);
1905 if(!(moldyn->total_steps%s)) {
1906 snprintf(dir,128,"%s/s-%08.f.save",
1907 moldyn->vlsdir,moldyn->time);
1908 fd=open(dir,O_WRONLY|O_TRUNC|O_CREAT,
1910 if(fd<0) perror("[moldyn] save fd open");
1912 write(fd,moldyn,sizeof(t_moldyn));
1913 write(fd,moldyn->atom,
1914 moldyn->count*sizeof(t_atom));
1920 if(!(moldyn->total_steps%a)) {
1921 #ifdef VISUAL_THREAD
1922 /* check whether thread has not terminated yet */
1924 ret=pthread_join(io_thread,NULL);
1927 /* prepare and start thread */
1928 if(moldyn->count!=md_copy.count) {
1932 memcpy(&md_copy,moldyn,sizeof(t_moldyn));
1934 atom_copy=malloc(moldyn->count*sizeof(t_atom));
1935 if(atom_copy==NULL) {
1936 perror("[moldyn] malloc atom copy (change)");
1940 md_copy.atom=atom_copy;
1941 memcpy(atom_copy,moldyn->atom,moldyn->count*sizeof(t_atom));
1943 ret=pthread_create(&io_thread,NULL,visual_atoms,&md_copy);
1945 perror("[moldyn] create visual atoms thread\n");
1949 visual_atoms(moldyn);
1954 /* display progress */
1956 /* get current time */
1957 gettimeofday(&t2,NULL);
1959 printf("sched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)\n",
1960 sched->count,i,moldyn->total_steps,
1961 moldyn->t,moldyn->t_avg,
1962 moldyn->p/BAR,moldyn->p_avg/BAR,
1963 //moldyn->p/BAR,(moldyn->p-2.0*moldyn->ekin/(3.0*moldyn->volume))/BAR,
1965 (int)(t2.tv_sec-t1.tv_sec));
1969 /* copy over time */
1973 /* increase absolute time */
1974 moldyn->time+=moldyn->tau;
1975 moldyn->total_steps+=1;
1979 /* check for hooks */
1981 printf("\n ## schedule hook %d start ##\n",
1983 sched->hook(moldyn,sched->hook_params);
1984 printf(" ## schedule hook end ##\n");
1987 /* increase the schedule counter */
1995 /* velocity verlet */
1997 int velocity_verlet(t_moldyn *moldyn) {
2000 double tau,tau_square,h;
2005 count=moldyn->count;
2007 tau_square=moldyn->tau_square;
2009 for(i=0;i<count;i++) {
2010 /* check whether fixed atom */
2011 if(atom[i].attr&ATOM_ATTR_FP)
2015 v3_scale(&delta,&(atom[i].v),tau);
2016 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2017 v3_scale(&delta,&(atom[i].f),h*tau_square);
2018 v3_add(&(atom[i].r),&(atom[i].r),&delta);
2019 check_per_bound(moldyn,&(atom[i].r));
2021 /* velocities [actually v(t+tau/2)] */
2022 v3_scale(&delta,&(atom[i].f),h*tau);
2023 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2026 /* criticial check */
2027 moldyn_bc_check(moldyn);
2029 /* neighbour list update */
2030 link_cell_update(moldyn);
2032 /* forces depending on chosen potential */
2034 potential_force_calc(moldyn);
2036 albe_potential_force_calc(moldyn);
2039 for(i=0;i<count;i++) {
2040 /* check whether fixed atom */
2041 if(atom[i].attr&ATOM_ATTR_FP)
2043 /* again velocities [actually v(t+tau)] */
2044 v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
2045 v3_add(&(atom[i].v),&(atom[i].v),&delta);
2054 * potentials & corresponding forces & virial routine
2058 /* generic potential and force calculation */
2060 int potential_force_calc(t_moldyn *moldyn) {
2063 t_atom *itom,*jtom,*ktom;
2067 int *neighbour_i[27];
2071 int neighbour_i[27];
2074 t_list neighbour_i[27];
2075 t_list neighbour_i2[27];
2081 count=moldyn->count;
2091 /* reset global virial */
2092 memset(&(moldyn->gvir),0,sizeof(t_virial));
2094 /* reset force, site energy and virial of every atom */
2096 i=omp_get_thread_num();
2097 #pragma omp parallel for private(virial)
2099 for(i=0;i<count;i++) {
2102 v3_zero(&(itom[i].f));
2105 virial=(&(itom[i].virial));
2113 /* reset site energy */
2118 /* get energy, force and virial of every atom */
2120 /* first (and only) loop over atoms i */
2121 for(i=0;i<count;i++) {
2123 /* single particle potential/force */
2124 if(itom[i].attr&ATOM_ATTR_1BP)
2126 moldyn->func1b(moldyn,&(itom[i]));
2128 if(!(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)))
2131 /* 2 body pair potential/force */
2133 link_cell_neighbour_index(moldyn,
2134 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2135 (itom[i].r.y+moldyn->dim.y/2)/lc->y,
2136 (itom[i].r.z+moldyn->dim.z/2)/lc->z,
2141 /* first loop over atoms j */
2142 if(moldyn->func2b) {
2149 while(neighbour_i[j][p]!=-1) {
2151 jtom=&(atom[neighbour_i[j][p]]);
2159 p=lc->subcell->list[p];
2161 this=&(neighbour_i[j]);
2164 if(this->start==NULL)
2168 jtom=this->current->data;
2171 if(jtom==&(itom[i]))
2174 if((jtom->attr&ATOM_ATTR_2BP)&
2175 (itom[i].attr&ATOM_ATTR_2BP)) {
2176 moldyn->func2b(moldyn,
2186 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2192 /* 3 body potential/force */
2194 if(!(itom[i].attr&ATOM_ATTR_3BP))
2197 /* copy the neighbour lists */
2199 /* no copy needed for static lists */
2201 /* no copy needed for lowmem lists */
2203 memcpy(neighbour_i2,neighbour_i,27*sizeof(t_list));
2206 /* second loop over atoms j */
2213 while(neighbour_i[j][p]!=-1) {
2215 jtom=&(atom[neighbour_i[j][p]]);
2223 p=lc->subcell->list[p];
2225 this=&(neighbour_i[j]);
2228 if(this->start==NULL)
2233 jtom=this->current->data;
2236 if(jtom==&(itom[i]))
2239 if(!(jtom->attr&ATOM_ATTR_3BP))
2245 if(moldyn->func3b_j1)
2246 moldyn->func3b_j1(moldyn,
2251 /* in first j loop, 3bp run can be skipped */
2252 if(!(moldyn->run3bp))
2255 /* first loop over atoms k */
2256 if(moldyn->func3b_k1) {
2264 while(neighbour_i[k][q]!=-1) {
2266 ktom=&(atom[neighbour_i[k][q]]);
2274 q=lc->subcell->list[q];
2276 that=&(neighbour_i2[k]);
2279 if(that->start==NULL)
2283 ktom=that->current->data;
2286 if(!(ktom->attr&ATOM_ATTR_3BP))
2292 if(ktom==&(itom[i]))
2295 moldyn->func3b_k1(moldyn,
2306 } while(list_next_f(that)!=\
2314 if(moldyn->func3b_j2)
2315 moldyn->func3b_j2(moldyn,
2320 /* second loop over atoms k */
2321 if(moldyn->func3b_k2) {
2329 while(neighbour_i[k][q]!=-1) {
2331 ktom=&(atom[neighbour_i[k][q]]);
2339 q=lc->subcell->list[q];
2341 that=&(neighbour_i2[k]);
2344 if(that->start==NULL)
2348 ktom=that->current->data;
2351 if(!(ktom->attr&ATOM_ATTR_3BP))
2357 if(ktom==&(itom[i]))
2360 moldyn->func3b_k2(moldyn,
2371 } while(list_next_f(that)!=\
2379 /* 2bp post function */
2380 if(moldyn->func3b_j3) {
2381 moldyn->func3b_j3(moldyn,
2390 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2405 //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
2406 if(moldyn->time>DSTART&&moldyn->time<DEND) {
2408 printf(" x: %0.40f\n",moldyn->atom[DATOM].f.x);
2409 printf(" y: %0.40f\n",moldyn->atom[DATOM].f.y);
2410 printf(" z: %0.40f\n",moldyn->atom[DATOM].f.z);
2414 /* some postprocessing */
2416 #pragma omp parallel for
2418 for(i=0;i<count;i++) {
2419 /* calculate global virial */
2420 moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
2421 moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
2422 moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
2423 moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
2424 moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
2425 moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
2427 /* check forces regarding the given timestep */
2428 if(v3_norm(&(itom[i].f))>\
2429 0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
2430 printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
2438 * virial calculation
2441 //inline int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2442 int virial_calc(t_atom *a,t_3dvec *f,t_3dvec *d) {
2444 a->virial.xx+=f->x*d->x;
2445 a->virial.yy+=f->y*d->y;
2446 a->virial.zz+=f->z*d->z;
2447 a->virial.xy+=f->x*d->y;
2448 a->virial.xz+=f->x*d->z;
2449 a->virial.yz+=f->y*d->z;
2455 * periodic boundary checking
2458 //inline int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2459 int check_per_bound(t_moldyn *moldyn,t_3dvec *a) {
2470 if(moldyn->status&MOLDYN_STAT_PBX) {
2471 if(a->x>=x) a->x-=dim->x;
2472 else if(-a->x>x) a->x+=dim->x;
2474 if(moldyn->status&MOLDYN_STAT_PBY) {
2475 if(a->y>=y) a->y-=dim->y;
2476 else if(-a->y>y) a->y+=dim->y;
2478 if(moldyn->status&MOLDYN_STAT_PBZ) {
2479 if(a->z>=z) a->z-=dim->z;
2480 else if(-a->z>z) a->z+=dim->z;
2487 * debugging / critical check functions
2490 int moldyn_bc_check(t_moldyn *moldyn) {
2503 for(i=0;i<moldyn->count;i++) {
2504 if(atom[i].r.x>=dim->x/2||-atom[i].r.x>dim->x/2) {
2505 printf("FATAL: atom %d: x: %.20f (%.20f)\n",
2506 i,atom[i].r.x,dim->x/2);
2507 printf("diagnostic:\n");
2508 printf("-----------\natom.r.x:\n");
2510 memcpy(&byte,(u8 *)(&(atom[i].r.x))+j,1);
2513 ((byte)&(1<<k))?1:0,
2516 printf("---------------\nx=dim.x/2:\n");
2518 memcpy(&byte,(u8 *)(&x)+j,1);
2521 ((byte)&(1<<k))?1:0,
2524 if(atom[i].r.x==x) printf("the same!\n");
2525 else printf("different!\n");
2527 if(atom[i].r.y>=dim->y/2||-atom[i].r.y>dim->y/2)
2528 printf("FATAL: atom %d: y: %.20f (%.20f)\n",
2529 i,atom[i].r.y,dim->y/2);
2530 if(atom[i].r.z>=dim->z/2||-atom[i].r.z>dim->z/2)
2531 printf("FATAL: atom %d: z: %.20f (%.20f)\n",
2532 i,atom[i].r.z,dim->z/2);
2542 int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
2549 fd=open(file,O_RDONLY);
2551 perror("[moldyn] load save file open");
2555 fsize=lseek(fd,0,SEEK_END);
2556 lseek(fd,0,SEEK_SET);
2558 size=sizeof(t_moldyn);
2561 cnt=read(fd,moldyn,size);
2563 perror("[moldyn] load save file read (moldyn)");
2569 size=moldyn->count*sizeof(t_atom);
2571 /* correcting possible atom data offset */
2573 if(fsize!=sizeof(t_moldyn)+size) {
2574 corr=fsize-sizeof(t_moldyn)-size;
2575 printf("[moldyn] WARNING: lsf (illegal file size)\n");
2576 printf(" moifying offset:\n");
2577 printf(" - current pos: %d\n",sizeof(t_moldyn));
2578 printf(" - atom size: %d\n",size);
2579 printf(" - file size: %d\n",fsize);
2580 printf(" => correction: %d\n",corr);
2581 lseek(fd,corr,SEEK_CUR);
2584 moldyn->atom=(t_atom *)malloc(size);
2585 if(moldyn->atom==NULL) {
2586 perror("[moldyn] load save file malloc (atoms)");
2591 cnt=read(fd,moldyn->atom,size);
2593 perror("[moldyn] load save file read (atoms)");
2604 int moldyn_free_save_file(t_moldyn *moldyn) {
2611 int moldyn_load(t_moldyn *moldyn) {
2619 * function to find/callback all combinations of 2 body bonds
2622 int process_2b_bonds(t_moldyn *moldyn,void *data,
2623 int (*process)(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2624 void *data,u8 bc)) {
2634 t_list neighbour[27];
2644 for(i=0;i<moldyn->count;i++) {
2645 /* neighbour indexing */
2646 link_cell_neighbour_index(moldyn,
2647 (itom[i].r.x+moldyn->dim.x/2)/lc->x,
2648 (itom[i].r.y+moldyn->dim.y/2)/lc->x,
2649 (itom[i].r.z+moldyn->dim.z/2)/lc->x,
2654 bc=(j<lc->dnlc)?0:1;
2659 while(neighbour[j][p]!=-1) {
2661 jtom=&(moldyn->atom[neighbour[j][p]]);
2669 p=lc->subcell->list[p];
2671 this=&(neighbour[j]);
2674 if(this->start==NULL)
2679 jtom=this->current->data;
2683 process(moldyn,&(itom[i]),jtom,data,bc);
2690 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2700 * function to find neighboured atoms
2703 int process_neighbours(t_moldyn *moldyn,void *data,t_atom *atom,
2704 int (*process)(t_moldyn *moldyn,t_atom *atom,t_atom *natom,
2705 void *data,u8 bc)) {
2715 t_list neighbour[27];
2724 /* neighbour indexing */
2725 link_cell_neighbour_index(moldyn,
2726 (atom->r.x+moldyn->dim.x/2)/lc->x,
2727 (atom->r.y+moldyn->dim.y/2)/lc->x,
2728 (atom->r.z+moldyn->dim.z/2)/lc->x,
2733 bc=(j<lc->dnlc)?0:1;
2738 while(neighbour[j][p]!=-1) {
2740 natom=&(moldyn->atom[neighbour[j][p]]);
2747 natom=&(moldyn->atom[p]);
2748 p=lc->subcell->list[p];
2750 this=&(neighbour[j]);
2753 if(this->start==NULL)
2758 natom=this->current->data;
2762 process(moldyn,atom,natom,data,bc);
2769 } while(list_next_f(this)!=L_NO_NEXT_ELEMENT);
2778 * post processing functions
2781 int get_line(int fd,char *line,int max) {
2788 if(count==max) return count;
2789 ret=read(fd,line+count,1);
2790 if(ret<=0) return ret;
2791 if(line[count]=='\n') {
2792 memset(line+count,0,max-count-1);
2800 int pair_correlation_init(t_moldyn *moldyn,double dr) {
2806 int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
2822 for(i=0;i<moldyn->count;i++) {
2824 v3_sub(&dist,&(atom[i].r),&(atom[i].r_0));
2825 check_per_bound(moldyn,&dist);
2826 d2=v3_absolute_square(&dist);
2840 dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
2841 dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
2842 dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
2847 int bonding_analyze(t_moldyn *moldyn,double *cnt) {
2852 int calculate_pair_correlation_process(t_moldyn *moldyn,t_atom *itom,
2853 t_atom *jtom,void *data,u8 bc) {
2860 /* only count pairs once,
2861 * skip same atoms */
2862 if(itom->tag>=jtom->tag)
2866 * pair correlation calc
2873 v3_sub(&dist,&(jtom->r),&(itom->r));
2874 if(bc) check_per_bound(moldyn,&dist);
2875 d=v3_absolute_square(&dist);
2877 /* ignore if greater cutoff */
2878 if(d>moldyn->cutoff_square)
2881 /* fill the slots */
2885 /* should never happen but it does 8) -
2886 * related to -ffloat-store problem! */
2888 printf("[moldyn] WARNING: pcc (%d/%d)",
2894 if(itom->brand!=jtom->brand) {
2899 /* type a - type a bonds */
2901 pcc->stat[s+pcc->o1]+=1;
2903 /* type b - type b bonds */
2904 pcc->stat[s+pcc->o2]+=1;
2910 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
2917 pcc.o1=moldyn->cutoff/dr;
2920 if(pcc.o1*dr<=moldyn->cutoff)
2921 printf("[moldyn] WARNING: pcc (low #slots)\n");
2923 printf("[moldyn] pair correlation calc info:\n");
2924 printf(" time: %f\n",moldyn->time);
2925 printf(" count: %d\n",moldyn->count);
2926 printf(" cutoff: %f\n",moldyn->cutoff);
2927 printf(" temperature: cur=%f avg=%f\n",moldyn->t,moldyn->t_avg);
2930 pcc.stat=(double *)ptr;
2933 pcc.stat=(double *)malloc(3*pcc.o1*sizeof(double));
2934 if(pcc.stat==NULL) {
2935 perror("[moldyn] pair correlation malloc");
2940 memset(pcc.stat,0,3*pcc.o1*sizeof(double));
2943 process_2b_bonds(moldyn,&pcc,calculate_pair_correlation_process);
2946 for(i=1;i<pcc.o1;i++) {
2947 // normalization: 4 pi r^2 dr
2948 // here: not double counting pairs -> 2 pi r r dr
2949 // ... and actually it's a constant times r^2
2952 pcc.stat[pcc.o1+i]/=norm;
2953 pcc.stat[pcc.o2+i]/=norm;
2958 /* todo: store/print pair correlation function */
2965 int bond_analyze_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
2972 if(itom->tag>=jtom->tag)
2976 v3_sub(&dist,&(jtom->r),&(itom->r));
2977 if(bc) check_per_bound(moldyn,&dist);
2978 d=v3_absolute_square(&dist);
2980 /* ignore if greater or equal cutoff */
2981 if(d>moldyn->cutoff_square)
2984 /* check for potential bond */
2985 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
2988 /* now count this bonding ... */
2991 /* increase total bond counter
2992 * ... double counting!
2997 ba->acnt[jtom->tag]+=1;
2999 ba->bcnt[jtom->tag]+=1;
3002 ba->acnt[itom->tag]+=1;
3004 ba->bcnt[itom->tag]+=1;
3009 int bond_analyze(t_moldyn *moldyn,double *quality) {
3011 // by now: # bonds of type 'a-4b' and 'b-4a' / # bonds total
3019 ba.acnt=malloc(moldyn->count*sizeof(int));
3021 perror("[moldyn] bond analyze malloc (a)");
3024 memset(ba.acnt,0,moldyn->count*sizeof(int));
3026 ba.bcnt=malloc(moldyn->count*sizeof(int));
3028 perror("[moldyn] bond analyze malloc (b)");
3031 memset(ba.bcnt,0,moldyn->count*sizeof(int));
3040 process_2b_bonds(moldyn,&ba,bond_analyze_process);
3042 for(i=0;i<moldyn->count;i++) {
3043 if(atom[i].brand==0) {
3044 if((ba.acnt[i]==0)&(ba.bcnt[i]==4))
3048 if((ba.acnt[i]==4)&(ba.bcnt[i]==0)) {
3056 printf("[moldyn] bond analyze: c_cnt=%d | set=%d\n",ccnt,cset);
3057 printf("[moldyn] bond analyze: q_cnt=%d | tot=%d\n",qcnt,ba.tcnt);
3060 quality[0]=1.0*ccnt/cset;
3061 quality[1]=1.0*qcnt/ba.tcnt;
3064 printf("[moldyn] bond analyze: c_bnd_q=%f\n",1.0*qcnt/ba.tcnt);
3065 printf("[moldyn] bond analyze: tot_q=%f\n",1.0*qcnt/ba.tcnt);
3072 * visualization code
3075 int visual_init(t_moldyn *moldyn,char *filebase) {
3077 strncpy(moldyn->vis.fb,filebase,128);
3082 int visual_bonds_process(t_moldyn *moldyn,t_atom *itom,t_atom *jtom,
3089 if(itom->tag>=jtom->tag)
3092 if(moldyn->check_2b_bond(moldyn,itom,jtom,bc)==FALSE)
3095 if((itom->attr&ATOM_ATTR_VB)|(jtom->attr&ATOM_ATTR_VB))
3096 dprintf(vb->fd,"# [B] %f %f %f %f %f %f\n",
3097 itom->r.x,itom->r.y,itom->r.z,
3098 jtom->r.x,jtom->r.y,jtom->r.z);
3103 #ifdef VISUAL_THREAD
3104 void *visual_atoms(void *ptr) {
3106 int visual_atoms(t_moldyn *moldyn) {
3116 #ifdef VISUAL_THREAD
3130 sprintf(file,"%s/atomic_conf_%08.f.xyz",v->fb,moldyn->time);
3131 vb.fd=open(file,O_WRONLY|O_CREAT|O_TRUNC,S_IRUSR|S_IWUSR);
3133 perror("open visual save file fd");
3134 #ifndef VISUAL_THREAD
3139 /* write the actual data file */
3142 dprintf(vb.fd,"# [P] %d %08.f <%f,%f,%f>\n",
3143 moldyn->count,moldyn->time,help/40.0,help/40.0,-0.8*help);
3145 // atomic configuration
3146 for(i=0;i<moldyn->count;i++)
3147 // atom type, positions, color and kinetic energy
3148 dprintf(vb.fd,"%s %f %f %f %s %f\n",pse_name[atom[i].element],
3152 pse_col[atom[i].element],
3155 // bonds between atoms
3156 #ifndef VISUAL_THREAD
3157 process_2b_bonds(moldyn,&vb,visual_bonds_process);
3162 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3163 -dim.x/2,-dim.y/2,-dim.z/2,
3164 dim.x/2,-dim.y/2,-dim.z/2);
3165 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3166 -dim.x/2,-dim.y/2,-dim.z/2,
3167 -dim.x/2,dim.y/2,-dim.z/2);
3168 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3169 dim.x/2,dim.y/2,-dim.z/2,
3170 dim.x/2,-dim.y/2,-dim.z/2);
3171 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3172 -dim.x/2,dim.y/2,-dim.z/2,
3173 dim.x/2,dim.y/2,-dim.z/2);
3175 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3176 -dim.x/2,-dim.y/2,dim.z/2,
3177 dim.x/2,-dim.y/2,dim.z/2);
3178 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3179 -dim.x/2,-dim.y/2,dim.z/2,
3180 -dim.x/2,dim.y/2,dim.z/2);
3181 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3182 dim.x/2,dim.y/2,dim.z/2,
3183 dim.x/2,-dim.y/2,dim.z/2);
3184 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3185 -dim.x/2,dim.y/2,dim.z/2,
3186 dim.x/2,dim.y/2,dim.z/2);
3188 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3189 -dim.x/2,-dim.y/2,dim.z/2,
3190 -dim.x/2,-dim.y/2,-dim.z/2);
3191 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3192 -dim.x/2,dim.y/2,dim.z/2,
3193 -dim.x/2,dim.y/2,-dim.z/2);
3194 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3195 dim.x/2,-dim.y/2,dim.z/2,
3196 dim.x/2,-dim.y/2,-dim.z/2);
3197 dprintf(vb.fd,"# [D] %f %f %f %f %f %f\n",
3198 dim.x/2,dim.y/2,dim.z/2,
3199 dim.x/2,dim.y/2,-dim.z/2);
3204 #ifdef VISUAL_THREAD
3215 * fpu cntrol functions
3218 // set rounding to double (eliminates -ffloat-store!)
3219 int fpu_set_rtd(void) {
3225 ctrl&=~_FPU_EXTENDED;