implemented basic p ctrl stuff + video with 13 fps
[physik/posic.git] / moldyn.c
index edda007..91fead7 100644 (file)
--- a/moldyn.c
+++ b/moldyn.c
@@ -191,15 +191,12 @@ int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
                moldyn->vis.dim.z=z;
        }
 
-       moldyn->dv=0.000001*moldyn->volume;
-
        printf("[moldyn] dimensions in A and A^3 respectively:\n");
        printf("  x: %f\n",moldyn->dim.x);
        printf("  y: %f\n",moldyn->dim.y);
        printf("  z: %f\n",moldyn->dim.z);
        printf("  volume: %f\n",moldyn->volume);
        printf("  visualize simulation box: %s\n",visualize?"yes":"no");
-       printf("  delta volume (pressure calc): %f\n",moldyn->dv);
 
        return 0;
 }
@@ -370,12 +367,25 @@ int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer) {
                        dprintf(moldyn->tfd,"# temperature log file\n");
                        printf("temperature (%d)\n",timer);
                        break;
+               case LOG_VOLUME:
+                       moldyn->vwrite=timer;
+                       snprintf(filename,127,"%s/volume",moldyn->vlsdir);
+                       moldyn->vfd=open(filename,
+                                        O_WRONLY|O_CREAT|O_EXCL,
+                                        S_IRUSR|S_IWUSR);
+                       if(moldyn->vfd<0) {
+                               perror("[moldyn] volume log file\n");
+                               return moldyn->vfd;
+                       }
+                       dprintf(moldyn->vfd,"# volume log file\n");
+                       printf("volume (%d)\n",timer);
+                       break;
                case SAVE_STEP:
                        moldyn->swrite=timer;
                        printf("save file (%d)\n",timer);
                        break;
                case VISUAL_STEP:
-                       moldyn->vwrite=timer;
+                       moldyn->awrite=timer;
                        ret=visual_init(moldyn,moldyn->vlsdir);
                        if(ret<0) {
                                printf("[moldyn] visual init failure\n");
@@ -595,7 +605,7 @@ int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
        int count;
        
        atom=moldyn->atom;
-       count=(moldyn->count)++;
+       count=(moldyn->count)++;        // asshole style!
 
        ptr=realloc(atom,(count+1)*sizeof(t_atom));
        if(!ptr) {
@@ -605,6 +615,9 @@ int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
        moldyn->atom=ptr;
 
        atom=moldyn->atom;
+
+       /* initialize new atom */
+       memset(&(atom[count]),0,sizeof(t_atom));
        atom[count].r=*r;
        atom[count].v=*v;
        atom[count].element=element;
@@ -913,16 +926,26 @@ double ideal_gas_law_pressure(t_moldyn *moldyn) {
 double virial_sum(t_moldyn *moldyn) {
 
        int i;
-       double v;
        t_virial *virial;
 
        /* virial (sum over atom virials) */
-       v=0.0;
+       moldyn->virial=0.0;
+       moldyn->vir.xx=0.0;
+       moldyn->vir.yy=0.0;
+       moldyn->vir.zz=0.0;
+       moldyn->vir.xy=0.0;
+       moldyn->vir.xz=0.0;
+       moldyn->vir.yz=0.0;
        for(i=0;i<moldyn->count;i++) {
                virial=&(moldyn->atom[i].virial);
-               v+=(virial->xx+virial->yy+virial->zz);
+               moldyn->virial+=(virial->xx+virial->yy+virial->zz);
+               moldyn->vir.xx+=virial->xx;
+               moldyn->vir.yy+=virial->yy;
+               moldyn->vir.zz+=virial->zz;
+               moldyn->vir.xy+=virial->xy;
+               moldyn->vir.xz+=virial->xz;
+               moldyn->vir.yz+=virial->yz;
        }
-       moldyn->virial=v;
 
        /* global virial (absolute coordinates) */
        virial=&(moldyn->gvir);
@@ -954,6 +977,36 @@ double pressure_calc(t_moldyn *moldyn) {
        return moldyn->p;
 }
 
+int average_reset(t_moldyn *moldyn) {
+
+       printf("[moldyn] average reset\n");
+
+       /* update skip value */
+       moldyn->avg_skip=moldyn->total_steps;
+
+       /* kinetic energy */
+       moldyn->k_sum=0.0;
+       moldyn->k2_sum=0.0;
+       
+       /* potential energy */
+       moldyn->v_sum=0.0;
+       moldyn->v2_sum=0.0;
+
+       /* temperature */
+       moldyn->t_sum=0.0;
+
+       /* virial */
+       moldyn->virial_sum=0.0;
+       moldyn->gv_sum=0.0;
+
+       /* pressure */
+       moldyn->p_sum=0.0;
+       moldyn->gp_sum=0.0;
+       moldyn->tp_sum=0.0;
+
+       return 0;
+}
+
 int average_and_fluctuation_calc(t_moldyn *moldyn) {
 
        int denom;
@@ -994,6 +1047,8 @@ int average_and_fluctuation_calc(t_moldyn *moldyn) {
        moldyn->p_avg=moldyn->p_sum/denom;
        moldyn->gp_sum+=moldyn->gp;
        moldyn->gp_avg=moldyn->gp_sum/denom;
+       moldyn->tp_sum+=moldyn->tp;
+       moldyn->tp_avg=moldyn->tp_sum/denom;
 
        return 0;
 }
@@ -1035,8 +1090,9 @@ double thermodynamic_pressure_calc(t_moldyn *moldyn) {
 
        t_3dvec dim;
        //t_3dvec *tp;
-       double u_up,u_down,dv;
-       double scale,p;
+       double h,dv;
+       double y0,y1;
+       double su,sd;
        t_atom *store;
 
        /*
@@ -1046,54 +1102,56 @@ double thermodynamic_pressure_calc(t_moldyn *moldyn) {
         *
         */
 
-       scale=0.00001;
-       dv=8*scale*scale*scale*moldyn->volume;
-
+       /* store atomic configuration + dimension */
        store=malloc(moldyn->count*sizeof(t_atom));
        if(store==NULL) {
                printf("[moldyn] allocating store mem failed\n");
                return -1;
        }
-
-       /* save unscaled potential energy + atom/dim configuration */
        memcpy(store,moldyn->atom,moldyn->count*sizeof(t_atom));
        dim=moldyn->dim;
 
+       /* x1, y1 */
+       sd=0.00001;
+       h=(1.0-sd)*(1.0-sd)*(1.0-sd);
+       su=pow(2.0-h,ONE_THIRD)-1.0;
+       dv=(1.0-h)*moldyn->volume;
+
        /* scale up dimension and atom positions */
-       scale_dim(moldyn,SCALE_UP,scale,TRUE,TRUE,TRUE);
-       scale_atoms(moldyn,SCALE_UP,scale,TRUE,TRUE,TRUE);
+       scale_dim(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
+       scale_atoms(moldyn,SCALE_UP,su,TRUE,TRUE,TRUE);
        link_cell_shutdown(moldyn);
        link_cell_init(moldyn,QUIET);
        potential_force_calc(moldyn);
-       u_up=moldyn->energy;
+       y1=moldyn->energy;
 
        /* restore atomic configuration + dim */
        memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
        moldyn->dim=dim;
 
        /* scale down dimension and atom positions */
-       scale_dim(moldyn,SCALE_DOWN,scale,TRUE,TRUE,TRUE);
-       scale_atoms(moldyn,SCALE_DOWN,scale,TRUE,TRUE,TRUE);
+       scale_dim(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
+       scale_atoms(moldyn,SCALE_DOWN,sd,TRUE,TRUE,TRUE);
        link_cell_shutdown(moldyn);
        link_cell_init(moldyn,QUIET);
        potential_force_calc(moldyn);
-       u_down=moldyn->energy;
+       y0=moldyn->energy;
        
        /* calculate pressure */
-       p=-(u_up-u_down)/dv;
-printf("-------> %.10f %.10f %f\n",u_up/EV/moldyn->count,u_down/EV/moldyn->count,p/BAR);
+       moldyn->tp=-(y1-y0)/(2.0*dv);
 
-       /* restore atomic configuration + dim */
+       /* restore atomic configuration */
        memcpy(moldyn->atom,store,moldyn->count*sizeof(t_atom));
        moldyn->dim=dim;
-
-       /* restore energy */
-       potential_force_calc(moldyn);
-
        link_cell_shutdown(moldyn);
        link_cell_init(moldyn,QUIET);
+       //potential_force_calc(moldyn);
 
-       return p;
+       /* free store buffer */
+       if(store)
+               free(store);
+
+       return moldyn->tp;
 }
 
 double get_pressure(t_moldyn *moldyn) {
@@ -1154,13 +1212,12 @@ int scale_volume(t_moldyn *moldyn) {
 
        /* scaling factor */
        if(moldyn->pt_scale&P_SCALE_BERENDSEN) {
-               scale=1.0-(moldyn->p_ref-moldyn->p)/moldyn->p_tc;
+               scale=1.0-(moldyn->p_ref-moldyn->p)*moldyn->p_tc;
                scale=pow(scale,ONE_THIRD);
        }
        else {
                scale=pow(moldyn->p/moldyn->p_ref,ONE_THIRD);
        }
-moldyn->debug=scale;
 
        /* scale the atoms and dimensions */
        scale_atoms(moldyn,SCALE_DIRECT,scale,TRUE,TRUE,TRUE);
@@ -1342,7 +1399,7 @@ int link_cell_update(t_moldyn *moldyn) {
                i=((atom[count].r.x+(moldyn->dim.x/2))/lc->x);
                j=((atom[count].r.y+(moldyn->dim.y/2))/lc->y);
                k=((atom[count].r.z+(moldyn->dim.z/2))/lc->z);
-
+       
 #ifdef STATIC_LISTS
                p=0;
                while(lc->subcell[i+j*nx+k*nx*ny][p]!=0)
@@ -1392,6 +1449,10 @@ int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,
        count2=27;
        a=nx*ny;
 
+       if(i>=nx||j>=ny||k>=nz)
+               printf("[moldyn] WARNING: lcni %d/%d %d/%d %d/%d\n",
+                      i,nx,j,ny,k,nz);
+
        cell[0]=lc->subcell[i+j*nx+k*a];
        for(ci=-1;ci<=1;ci++) {
                bx=0;
@@ -1502,7 +1563,7 @@ int moldyn_set_schedule_hook(t_moldyn *moldyn,set_hook hook,void *hook_params) {
 int moldyn_integrate(t_moldyn *moldyn) {
 
        int i;
-       unsigned int e,m,s,v,p,t;
+       unsigned int e,m,s,v,p,t,a;
        t_3dvec momentum;
        t_moldyn_schedule *sched;
        t_atom *atom;
@@ -1524,6 +1585,7 @@ int moldyn_integrate(t_moldyn *moldyn) {
        m=moldyn->mwrite;
        s=moldyn->swrite;
        v=moldyn->vwrite;
+       a=moldyn->awrite;
        p=moldyn->pwrite;
        t=moldyn->twrite;
 
@@ -1542,14 +1604,14 @@ int moldyn_integrate(t_moldyn *moldyn) {
 
        /* some stupid checks before we actually start calculating bullshit */
        if(moldyn->cutoff>0.5*moldyn->dim.x)
-               printf("[moldyn] warning: cutoff > 0.5 x dim.x\n");
+               printf("[moldyn] WARNING: cutoff > 0.5 x dim.x\n");
        if(moldyn->cutoff>0.5*moldyn->dim.y)
-               printf("[moldyn] warning: cutoff > 0.5 x dim.y\n");
+               printf("[moldyn] WARNING: cutoff > 0.5 x dim.y\n");
        if(moldyn->cutoff>0.5*moldyn->dim.z)
-               printf("[moldyn] warning: cutoff > 0.5 x dim.z\n");
+               printf("[moldyn] WARNING: cutoff > 0.5 x dim.z\n");
        ds=0.5*atom[0].f.x*moldyn->tau_square/atom[0].mass;
        if(ds>0.05*moldyn->nnd)
-               printf("[moldyn] warning: forces too high / tau too small!\n");
+               printf("[moldyn] WARNING: forces too high / tau too small!\n");
 
        /* zero absolute time */
        moldyn->time=0.0;
@@ -1585,6 +1647,9 @@ int moldyn_integrate(t_moldyn *moldyn) {
                temperature_calc(moldyn);
                virial_sum(moldyn);
                pressure_calc(moldyn);
+               //thermodynamic_pressure_calc(moldyn);
+
+               /* calculate fluctuations + averages */
                average_and_fluctuation_calc(moldyn);
 
                /* p/t scaling */
@@ -1614,9 +1679,10 @@ int moldyn_integrate(t_moldyn *moldyn) {
                if(p) {
                        if(!(moldyn->total_steps%p)) {
                                dprintf(moldyn->pfd,
-                                       "%f %f %f %f %f\n",moldyn->time,
+                                       "%f %f %f %f %f %f %f\n",moldyn->time,
                                         moldyn->p/BAR,moldyn->p_avg/BAR,
-                                        moldyn->gp/BAR,moldyn->gp_avg/BAR);
+                                        moldyn->gp/BAR,moldyn->gp_avg/BAR,
+                                        moldyn->tp/BAR,moldyn->tp_avg/BAR);
                        }
                }
                if(t) {
@@ -1626,6 +1692,12 @@ int moldyn_integrate(t_moldyn *moldyn) {
                                        moldyn->time,moldyn->t,moldyn->t_avg);
                        }
                }
+               if(v) {
+                       if(!(moldyn->total_steps%v)) {
+                               dprintf(moldyn->vfd,
+                                       "%f %f\n",moldyn->time,moldyn->volume);
+                       }
+               }
                if(s) {
                        if(!(moldyn->total_steps%s)) {
                                snprintf(dir,128,"%s/s-%07.f.save",
@@ -1641,8 +1713,8 @@ int moldyn_integrate(t_moldyn *moldyn) {
                                close(fd);
                        }       
                }
-               if(v) {
-                       if(!(moldyn->total_steps%v)) {
+               if(a) {
+                       if(!(moldyn->total_steps%a)) {
                                visual_atoms(moldyn);
                        }
                }
@@ -1652,10 +1724,10 @@ int moldyn_integrate(t_moldyn *moldyn) {
                        /* get current time */
                        gettimeofday(&t2,NULL);
 
-printf("\rsched:%d, steps:%d/%d, T:%3.1f/%3.1f P:%4.1f/%4.1f V:%6.1f (%d)",
+printf("\rsched:%d, steps:%d/%d, T:%4.1f/%4.1f P:%4.1f/%4.1f V:%6.1f (%d)",
        sched->count,i,moldyn->total_steps,
        moldyn->t,moldyn->t_avg,
-       moldyn->p_avg/BAR,moldyn->gp_avg/BAR,
+       moldyn->p/BAR,moldyn->p_avg/BAR,
        moldyn->volume,
        (int)(t2.tv_sec-t1.tv_sec));
 
@@ -1702,6 +1774,9 @@ int velocity_verlet(t_moldyn *moldyn) {
        tau_square=moldyn->tau_square;
 
        for(i=0;i<count;i++) {
+               /* check whether fixed atom */
+               if(atom[i].attr&ATOM_ATTR_FP)
+                       continue;
                /* new positions */
                h=0.5/atom[i].mass;
                v3_scale(&delta,&(atom[i].v),tau);
@@ -1715,6 +1790,9 @@ int velocity_verlet(t_moldyn *moldyn) {
                v3_add(&(atom[i].v),&(atom[i].v),&delta);
        }
 
+       /* criticial check */
+       moldyn_bc_check(moldyn);
+
        /* neighbour list update */
        link_cell_update(moldyn);
 
@@ -1722,6 +1800,9 @@ int velocity_verlet(t_moldyn *moldyn) {
        potential_force_calc(moldyn);
 
        for(i=0;i<count;i++) {
+               /* check whether fixed atom */
+               if(atom[i].attr&ATOM_ATTR_FP)
+                       continue;
                /* again velocities [actually v(t+tau)] */
                v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
                v3_add(&(atom[i].v),&(atom[i].v),&delta);
@@ -2049,20 +2130,27 @@ int potential_force_calc(t_moldyn *moldyn) {
        //printf("\nATOM 0: %f %f %f\n\n",itom->f.x,itom->f.y,itom->f.z);
        if(moldyn->time>DSTART&&moldyn->time<DEND) {
                printf("force:\n");
-               printf("  x: %0.40f\n",moldyn->atom[5832].f.x);
-               printf("  y: %0.40f\n",moldyn->atom[5832].f.y);
-               printf("  z: %0.40f\n",moldyn->atom[5832].f.z);
+               printf("  x: %0.40f\n",moldyn->atom[DATOM].f.x);
+               printf("  y: %0.40f\n",moldyn->atom[DATOM].f.y);
+               printf("  z: %0.40f\n",moldyn->atom[DATOM].f.z);
        }
 #endif
 
-       /* calculate global virial */
+       /* some postprocessing */
        for(i=0;i<count;i++) {
-               moldyn->gvir.xx+=moldyn->atom[i].r.x*moldyn->atom[i].f.x;
-               moldyn->gvir.yy+=moldyn->atom[i].r.y*moldyn->atom[i].f.y;
-               moldyn->gvir.zz+=moldyn->atom[i].r.z*moldyn->atom[i].f.z;
-               moldyn->gvir.xy+=moldyn->atom[i].r.y*moldyn->atom[i].f.x;
-               moldyn->gvir.xz+=moldyn->atom[i].r.z*moldyn->atom[i].f.x;
-               moldyn->gvir.yz+=moldyn->atom[i].r.z*moldyn->atom[i].f.y;
+               /* calculate global virial */
+               moldyn->gvir.xx+=itom[i].r.x*itom[i].f.x;
+               moldyn->gvir.yy+=itom[i].r.y*itom[i].f.y;
+               moldyn->gvir.zz+=itom[i].r.z*itom[i].f.z;
+               moldyn->gvir.xy+=itom[i].r.y*itom[i].f.x;
+               moldyn->gvir.xz+=itom[i].r.z*itom[i].f.x;
+               moldyn->gvir.yz+=itom[i].r.z*itom[i].f.y;
+
+               /* check forces regarding the given timestep */
+               if(v3_norm(&(itom[i].f))>\
+                  0.1*moldyn->nnd*itom[i].mass/moldyn->tau_square)
+                       printf("[moldyn] WARNING: pfc (high force: atom %d)\n",
+                              i);
        }
 
        return 0;
@@ -2177,6 +2265,8 @@ int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
 
        int fd;
        int cnt,size;
+       int fsize;
+       int corr;
 
        fd=open(file,O_RDONLY);
        if(fd<0) {
@@ -2184,6 +2274,9 @@ int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
                return fd;
        }
 
+       fsize=lseek(fd,0,SEEK_END);
+       lseek(fd,0,SEEK_SET);
+
        size=sizeof(t_moldyn);
 
        while(size) {
@@ -2197,6 +2290,19 @@ int moldyn_read_save_file(t_moldyn *moldyn,char *file) {
 
        size=moldyn->count*sizeof(t_atom);
 
+       /* correcting possible atom data offset */
+       corr=0;
+       if(fsize!=sizeof(t_moldyn)+size) {
+               corr=fsize-sizeof(t_moldyn)-size;
+               printf("[moldyn] WARNING: lsf (illegal file size)\n");
+               printf("  moifying offset:\n");
+               printf("  - current pos: %d\n",sizeof(t_moldyn));
+               printf("  - atom size: %d\n",size);
+               printf("  - file size: %d\n",fsize);
+               printf("  => correction: %d\n",corr);
+               lseek(fd,corr,SEEK_CUR);
+       }
+
        moldyn->atom=(t_atom *)malloc(size);
        if(moldyn->atom==NULL) {
                perror("[moldyn] load save file malloc (atoms)");
@@ -2252,6 +2358,47 @@ int pair_correlation_init(t_moldyn *moldyn,double dr) {
        return 0;
 }
 
+int calculate_diffusion_coefficient(t_moldyn *moldyn,double *dc) {
+
+       int i;
+       t_atom *atom;
+       t_3dvec dist;
+       double d2;
+       int a_cnt;
+       int b_cnt;
+
+       atom=moldyn->atom;
+       dc[0]=0;
+       dc[1]=0;
+       dc[2]=0;
+       a_cnt=0;
+       b_cnt=0;
+
+       for(i=0;i<moldyn->count;i++) {
+
+               v3_sub(&dist,&(atom[i].r),&(atom[i].r_0));
+               check_per_bound(moldyn,&dist);
+               d2=v3_absolute_square(&dist);
+
+               if(atom[i].brand) {
+                       b_cnt+=1;
+                       dc[1]+=d2;
+               }
+               else {
+                       a_cnt+=1;
+                       dc[0]+=d2;
+               }
+
+               dc[2]+=d2;
+       }
+
+       dc[0]*=(1.0/(6.0*moldyn->time*a_cnt));
+       dc[1]*=(1.0/(6.0*moldyn->time*b_cnt));
+       dc[2]*=(1.0/(6.0*moldyn->time*moldyn->count));
+               
+       return 0;
+}
+
 int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
 
        int slots;
@@ -2269,15 +2416,18 @@ int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
        unsigned char bc;
        t_3dvec dist;
        double d;
-       //double norm;
+       double norm;
        int o,s;
        unsigned char ibrand;
 
        lc=&(moldyn->lc);
 
-       slots=moldyn->cutoff/dr;
+       slots=2.0*moldyn->cutoff/dr;
        o=2*slots;
 
+       if(slots*dr<=moldyn->cutoff)
+               printf("[moldyn] WARNING: pcc (low #slots)\n");
+
        printf("[moldyn] pair correlation calc info:\n");
        printf("  time: %f\n",moldyn->time);
        printf("  count: %d\n",moldyn->count);
@@ -2349,7 +2499,7 @@ int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
                                d=v3_absolute_square(&dist);
 
                                /* ignore if greater or equal cutoff */
-                               if(d>=moldyn->cutoff_square)
+                               if(d>=4.0*moldyn->cutoff_square)
                                        continue;
 
                                /* fill the slots */
@@ -2358,7 +2508,12 @@ int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
 
                                /* should never happen but it does 8) -
                                 * related to -ffloat-store problem! */
-                               if(s>=slots) s=slots-1;
+                               if(s>=slots) {
+                                       printf("[moldyn] WARNING: pcc (%d/%d)",
+                                              s,slots);
+                                       printf("\n");
+                                       s=slots-1;
+                               }
 
                                if(ibrand!=jtom->brand) {
                                        /* mixed */
@@ -2380,16 +2535,17 @@ int calculate_pair_correlation(t_moldyn *moldyn,double dr,void *ptr) {
                }
        }
 
-       /* normalization 
+       /* normalization */
        for(i=1;i<slots;i++) {
-                // normalization: 4 pi r r dr
+                // normalization: 4 pi r^2 dr
                 // here: not double counting pairs -> 2 pi r r dr
-               norm=2*M_PI*moldyn->count*(i*dr*i*dr)*dr;
+                // ... and actually it's a constant times r^2
+               norm=i*i*dr*dr;
                stat[i]/=norm;
                stat[slots+i]/=norm;
                stat[o+i]/=norm;
        }
-       */
+       /* */
 
        if(ptr==NULL) {
                /* todo: store/print pair correlation function */