ready for tersoff 3bp debugging
[physik/posic.git] / moldyn.c
index 2f274c2..af0132d 100644 (file)
--- a/moldyn.c
+++ b/moldyn.c
@@ -41,9 +41,9 @@ int moldyn_init(t_moldyn *moldyn,int argc,char **argv) {
 
 int moldyn_shutdown(t_moldyn *moldyn) {
 
+       printf("[moldyn] shutdown\n");
        moldyn_log_shutdown(moldyn);
        link_cell_shutdown(moldyn);
-       moldyn_log_shutdown(moldyn);
        rand_close(&(moldyn->random));
        free(moldyn->atom);
 
@@ -71,9 +71,18 @@ int set_cutoff(t_moldyn *moldyn,double cutoff) {
        return 0;
 }
 
-int set_temperature(t_moldyn *moldyn,double t) {
-       
-       moldyn->t=t;
+int set_temperature(t_moldyn *moldyn,double t_ref) {
+
+       moldyn->t_ref=t_ref;
+
+       return 0;
+}
+
+int set_pt_scale(t_moldyn *moldyn,u8 ptype,double ptc,u8 ttype,double ttc) {
+
+       moldyn->pt_scale=(ptype|ttype);
+       moldyn->t_tc=ttc;
+       moldyn->p_tc=ptc;
 
        return 0;
 }
@@ -93,6 +102,13 @@ int set_dim(t_moldyn *moldyn,double x,double y,double z,u8 visualize) {
        return 0;
 }
 
+int set_nn_dist(t_moldyn *moldyn,double dist) {
+
+       moldyn->nnd=dist;
+
+       return 0;
+}
+
 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z) {
 
        if(x)
@@ -171,9 +187,10 @@ int moldyn_set_log(t_moldyn *moldyn,u8 type,char *fb,int timer) {
 
 int moldyn_log_shutdown(t_moldyn *moldyn) {
 
+       printf("[moldyn] log shutdown\n");
        if(moldyn->efd) close(moldyn->efd);
        if(moldyn->mfd) close(moldyn->mfd);
-       if(moldyn->visual) visual_tini(moldyn->visual);
+       if(&(moldyn->vis)) visual_tini(&(moldyn->vis));
 
        return 0;
 }
@@ -184,17 +201,15 @@ int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
        int count;
        int ret;
        t_3dvec origin;
-       t_atom *atom;
 
        count=a*b*c;
-       atom=moldyn->atom;
 
        if(type==FCC) count*=4;
 
        if(type==DIAMOND) count*=8;
 
-       atom=malloc(count*sizeof(t_atom));
-       if(atom==NULL) {
+       moldyn->atom=malloc(count*sizeof(t_atom));
+       if(moldyn->atom==NULL) {
                perror("malloc (atoms)");
                return -1;
        }
@@ -203,10 +218,10 @@ int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
 
        switch(type) {
                case FCC:
-                       ret=fcc_init(a,b,c,lc,atom,&origin);
+                       ret=fcc_init(a,b,c,lc,moldyn->atom,&origin);
                        break;
                case DIAMOND:
-                       ret=diamond_init(a,b,c,lc,atom,&origin);
+                       ret=diamond_init(a,b,c,lc,moldyn->atom,&origin);
                        break;
                default:
                        printf("unknown lattice type (%02x)\n",type);
@@ -222,15 +237,18 @@ int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
        }
 
        moldyn->count=count;
+       printf("[moldyn] created lattice with %d atoms\n",count);
 
        while(count) {
-               atom[count-1].element=element;
-               atom[count-1].mass=mass;
-               atom[count-1].attr=attr;
-               atom[count-1].bnum=bnum;
                count-=1;
+               moldyn->atom[count].element=element;
+               moldyn->atom[count].mass=mass;
+               moldyn->atom[count].attr=attr;
+               moldyn->atom[count].bnum=bnum;
+               check_per_bound(moldyn,&(moldyn->atom[count].r));
        }
 
+
        return ret;
 }
 
@@ -252,11 +270,12 @@ int add_atom(t_moldyn *moldyn,int element,double mass,u8 bnum,u8 attr,
        moldyn->atom=ptr;
 
        atom=moldyn->atom;
-       atom->r=*r;
-       atom->v=*v;
-       atom->element=element;
-       atom->bnum=bnum;
-       atom->attr=attr;
+       atom[count-1].r=*r;
+       atom[count-1].v=*v;
+       atom[count-1].element=element;
+       atom[count-1].mass=mass;
+       atom[count-1].bnum=bnum;
+       atom[count-1].attr=attr;
 
        return 0;
 }
@@ -268,7 +287,7 @@ int destroy_atoms(t_moldyn *moldyn) {
        return 0;
 }
 
-int thermal_init(t_moldyn *moldyn) {
+int thermal_init(t_moldyn *moldyn,u8 equi_init) {
 
        /*
         * - gaussian distribution of velocities
@@ -288,7 +307,7 @@ int thermal_init(t_moldyn *moldyn) {
        /* gaussian distribution of velocities */
        v3_zero(&p_total);
        for(i=0;i<moldyn->count;i++) {
-               sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t/atom[i].mass);
+               sigma=sqrt(2.0*K_BOLTZMANN*moldyn->t_ref/atom[i].mass);
                /* x direction */
                v=sigma*rand_get_gauss(random);
                atom[i].v.x=v;
@@ -311,28 +330,59 @@ int thermal_init(t_moldyn *moldyn) {
        }
 
        /* velocity scaling */
-       scale_velocity(moldyn);
+       scale_velocity(moldyn,equi_init);
 
        return 0;
 }
 
-int scale_velocity(t_moldyn *moldyn) {
+int scale_velocity(t_moldyn *moldyn,u8 equi_init) {
 
        int i;
-       double e,c;
+       double e,scale;
        t_atom *atom;
+       int count;
 
        atom=moldyn->atom;
 
        /*
         * - velocity scaling (E = 3/2 N k T), E: kinetic energy
         */
+
+       /* get kinetic energy / temperature & count involved atoms */
        e=0.0;
+       count=0;
+       for(i=0;i<moldyn->count;i++) {
+               if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB)) {
+                       e+=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
+                       count+=1;
+               }
+       }
+       if(count!=0) moldyn->t=(2.0*e)/(3.0*count*K_BOLTZMANN);
+       else return 0;  /* no atoms involved in scaling! */
+       
+       /* (temporary) hack for e,t = 0 */
+       if(e==0.0) {
+       moldyn->t=0.0;
+               if(moldyn->t_ref!=0.0)
+                       thermal_init(moldyn,equi_init);
+               else
+                       return 0; /* no scaling needed */
+       }
+
+
+       /* get scaling factor */
+       scale=moldyn->t_ref/moldyn->t;
+       if(equi_init&TRUE)
+               scale*=2.0;
+       else
+               if(moldyn->pt_scale&T_SCALE_BERENDSEN)
+                       scale=1.0+moldyn->tau*(scale-1.0)/moldyn->t_tc;
+       scale=sqrt(scale);
+
+       /* velocity scaling */
        for(i=0;i<moldyn->count;i++)
-               e+=0.5*atom[i].mass*v3_absolute_square(&(atom[i].v));
-       c=sqrt((2.0*e)/(3.0*moldyn->count*K_BOLTZMANN*moldyn->t));
-       for(i=0;i<moldyn->count;i++)
-               v3_scale(&(atom[i].v),&(atom[i].v),(1.0/c));
+               if((equi_init&TRUE)||(atom[i].attr&ATOM_ATTR_HB))
+                       v3_scale(&(atom[i].v),&(atom[i].v),scale);
 
        return 0;
 }
@@ -484,7 +534,6 @@ int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,t_list *cell) {
        count2=27;
        a=nx*ny;
 
-
        cell[0]=lc->subcell[i+j*nx+k*a];
        for(ci=-1;ci<=1;ci++) {
                bx=0;
@@ -518,7 +567,7 @@ int link_cell_neighbour_index(t_moldyn *moldyn,int i,int j,int k,t_list *cell) {
                }
        }
 
-       lc->dnlc=count2;
+       lc->dnlc=count1;
        lc->countn=27;
 
        return count2;
@@ -587,11 +636,13 @@ int moldyn_integrate(t_moldyn *moldyn) {
        unsigned int e,m,s,v;
        t_3dvec p;
        t_moldyn_schedule *schedule;
-
+       t_atom *atom;
        int fd;
        char fb[128];
+       double ds;
 
        schedule=&(moldyn->schedule);
+       atom=moldyn->atom;
 
        /* initialize linked cell method */
        link_cell_init(moldyn);
@@ -605,13 +656,22 @@ int moldyn_integrate(t_moldyn *moldyn) {
        /* sqaure of some variables */
        moldyn->tau_square=moldyn->tau*moldyn->tau;
        moldyn->cutoff_square=moldyn->cutoff*moldyn->cutoff;
-
        /* calculate initial forces */
        potential_force_calc(moldyn);
 
+       /* do some 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");
+       if(moldyn->cutoff>0.5*moldyn->dim.y)
+               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");
+       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");
+
        /* zero absolute time */
        moldyn->time=0.0;
-
        for(sched=0;sched<moldyn->schedule.content_count;sched++) {
 
                /* setting amount of runs and finite time step size */
@@ -626,6 +686,10 @@ int moldyn_integrate(t_moldyn *moldyn) {
                /* integration step */
                moldyn->integrate(moldyn);
 
+               /* p/t scaling */
+               if(moldyn->pt_scale&(T_SCALE_BERENDSEN|T_SCALE_DIRECT))
+                       scale_velocity(moldyn,FALSE);
+
                /* increase absolute time */
                moldyn->time+=moldyn->tau;
 
@@ -700,7 +764,7 @@ int velocity_verlet(t_moldyn *moldyn) {
                v3_add(&(atom[i].r),&(atom[i].r),&delta);
                v3_scale(&delta,&(atom[i].f),0.5*tau_square/atom[i].mass);
                v3_add(&(atom[i].r),&(atom[i].r),&delta);
-               v3_per_bound(&(atom[i].r),&(moldyn->dim));
+               check_per_bound(moldyn,&(atom[i].r));
 
                /* velocities */
                v3_scale(&delta,&(atom[i].f),0.5*tau/atom[i].mass);
@@ -735,106 +799,152 @@ int velocity_verlet(t_moldyn *moldyn) {
 int potential_force_calc(t_moldyn *moldyn) {
 
        int i,j,k,count;
-       t_atom *atom,*btom,*ktom;
+       t_atom *itom,*jtom,*ktom;
        t_linkcell *lc;
-       t_list neighbour[27];
-       t_list *this,*thisk,*neighbourk;
-       u8 bc,bck;
+       t_list neighbour_i[27];
+       t_list neighbour_i2[27];
+       //t_list neighbour_j[27];
+       t_list *this,*that;
+       u8 bc_ij,bc_ijk;
        int countn,dnlc;
 
        count=moldyn->count;
-       atom=moldyn->atom;
+       itom=moldyn->atom;
        lc=&(moldyn->lc);
 
        /* reset energy */
        moldyn->energy=0.0;
 
-printf("DEBUG: count = %d\n",count);
        for(i=0;i<count;i++) {
-       
+
                /* reset force */
-               v3_zero(&(atom[i].f));
+               v3_zero(&(itom[i].f));
 
                /* single particle potential/force */
-               if(atom[i].attr&ATOM_ATTR_1BP)
-                       moldyn->func1b(moldyn,&(atom[i]));
+               if(itom[i].attr&ATOM_ATTR_1BP)
+                       moldyn->func1b(moldyn,&(itom[i]));
 
                /* 2 body pair potential/force */
-               if(atom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)) {
+               if(itom[i].attr&(ATOM_ATTR_2BP|ATOM_ATTR_3BP)) {
        
-printf("DEBUG: processing atom %d\n",i);
                        link_cell_neighbour_index(moldyn,
-                               (atom[i].r.x+moldyn->dim.x/2)/lc->x,
-                               (atom[i].r.y+moldyn->dim.y/2)/lc->y,
-                               (atom[i].r.z+moldyn->dim.z/2)/lc->z,
-                               neighbour);
+                               (itom[i].r.x+moldyn->dim.x/2)/lc->x,
+                               (itom[i].r.y+moldyn->dim.y/2)/lc->y,
+                               (itom[i].r.z+moldyn->dim.z/2)/lc->z,
+                               neighbour_i);
 
                        countn=lc->countn;
                        dnlc=lc->dnlc;
 
-printf("DEBUG: countn = %d - dnslc = %d\n",countn,dnlc);
                        for(j=0;j<countn;j++) {
 
-                               this=&(neighbour[j]);
+                               this=&(neighbour_i[j]);
                                list_reset(this);
 
                                if(this->start==NULL)
                                        continue;
 
-                               bc=(j<dnlc)?0:1;
+                               bc_ij=(j<dnlc)?0:1;
 
                                do {
-                                       btom=this->current->data;
+                                       jtom=this->current->data;
 
-                                       if(btom==&(atom[i]))
+                                       if(jtom==&(itom[i]))
                                                continue;
 
-                                       if((btom->attr&ATOM_ATTR_2BP)&
-                                          (atom[i].attr&ATOM_ATTR_2BP))
-printf("DEBUG: calling func2b\n");
+                                       if((jtom->attr&ATOM_ATTR_2BP)&
+                                          (itom[i].attr&ATOM_ATTR_2BP))
                                                moldyn->func2b(moldyn,
-                                                              &(atom[i]),
-                                                              btom,
-                                                              bc);
+                                                              &(itom[i]),
+                                                              jtom,
+                                                              bc_ij);
 
                                        /* 3 body potential/force */
 
-                                       if(!(atom[i].attr&ATOM_ATTR_3BP)||
-                                          !(btom->attr&ATOM_ATTR_3BP))
+                                       if(!(itom[i].attr&ATOM_ATTR_3BP)||
+                                          !(jtom->attr&ATOM_ATTR_3BP))
                                                continue;
 
-                                       link_cell_neighbour_index(moldyn,
-                                          (btom->r.x+moldyn->dim.x/2)/lc->x,
-                                          (btom->r.y+moldyn->dim.y/2)/lc->y,
-                                          (btom->r.z+moldyn->dim.z/2)/lc->z,
-                                          neighbourk);
-
-                                       for(k=0;k<lc->countn;k++) {
-
-                                               thisk=&(neighbourk[k]);
-                                               list_reset(thisk);
+                       /*
+                        * according to mr. nordlund, we dont need to take the 
+                        * sum over all atoms now, as 'this is centered' around
+                        * atom i ...
+                        * i am not quite sure though! there is a not vanishing
+                        * part even if f_c_ik is zero ...
+                        * this analytical potentials suck!
+                        * switching from mc to md to dft soon!
+                        */
+
+                       //              link_cell_neighbour_index(moldyn,
+                       //                 (jtom->r.x+moldyn->dim.x/2)/lc->x,
+                       //                 (jtom->r.y+moldyn->dim.y/2)/lc->y,
+                       //                 (jtom->r.z+moldyn->dim.z/2)/lc->z,
+                       //                 neighbour_j);
+
+//                                     /* neighbours of j */
+//                                     for(k=0;k<lc->countn;k++) {
+//
+//                                             that=&(neighbour_j[k]);
+//                                             list_reset(that);
+//                                     
+//                                             if(that->start==NULL)
+//                                                     continue;
+//
+//                                             bc_ijk=(k<lc->dnlc)?0:1;
+//
+//                                             do {
+//
+//                     ktom=that->current->data;
+//
+//                     if(!(ktom->attr&ATOM_ATTR_3BP))
+//                             continue;
+//
+//                     if(ktom==jtom)
+//                             continue;
+//
+//                     if(ktom==&(itom[i]))
+//                             continue;
+//
+//                     moldyn->func3b(moldyn,&(itom[i]),jtom,ktom,bc_ijk);
+//
+/*                                             } while(list_next(that)!=\ */
+//                                                     L_NO_NEXT_ELEMENT);
+//
+//                                     }
+                       
+                                       /* copy the neighbour lists */
+                                       memcpy(neighbour_i2,neighbour_i,
+                                              27*sizeof(t_list));
+
+                                       /* get neighbours of i */
+                                       for(k=0;k<countn;k++) {
+
+                                               that=&(neighbour_i2[k]);
+                                               list_reset(that);
                                        
-                                               if(thisk->start==NULL)
+                                               if(that->start==NULL)
                                                        continue;
 
-                                               bck=(k<lc->dnlc)?0:1;
+                                               bc_ijk=(k<dnlc)?0:1;
 
                                                do {
 
-                       ktom=thisk->current->data;
+                       ktom=that->current->data;
 
                        if(!(ktom->attr&ATOM_ATTR_3BP))
                                continue;
 
-                       if(ktom==btom)
+                       if(ktom==jtom)
                                continue;
 
-                       if(ktom==&(atom[i]))
+                       if(ktom==&(itom[i]))
                                continue;
 
-                       moldyn->func3b(moldyn,&(atom[i]),btom,ktom,bck);
+printf("Debug: atom %d before 3bp: %08x %08x %08x | %.15f %.15f %.15f\n",i,&itom[i],jtom,ktom,itom[i].r.x,itom[i].f.x,itom[i].v.x);
+                       moldyn->func3b(moldyn,&(itom[i]),jtom,ktom,bc_ijk);
+printf("Debug: atom %d after 3bp: %08x %08x %08x | %.15f %.15f %.15f\n",i,&itom[i],jtom,ktom,itom[i].r.x,itom[i].f.x,itom[i].v.x);
 
-                                               } while(list_next(thisk)!=\
+                                               } while(list_next(that)!=\
                                                        L_NO_NEXT_ELEMENT);
 
                                        }
@@ -898,7 +1008,6 @@ int harmonic_oscillator(t_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc) {
 
        v3_sub(&distance,&(ai->r),&(aj->r));
        
-       v3_per_bound(&distance,&(moldyn->dim));
        if(bc) check_per_bound(moldyn,&distance);
        d=v3_norm(&distance);
        if(d<=moldyn->cutoff) {
@@ -942,7 +1051,7 @@ int lennard_jones(t_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc) {
                d=+h1-h2;
                d*=eps;
                v3_scale(&force,&distance,d);
-               v3_add(&(ai->f),&(aj->f),&force);
+               v3_add(&(ai->f),&(ai->f),&force);
        }
 
        return 0;
@@ -952,6 +1061,20 @@ int lennard_jones(t_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc) {
  * tersoff potential & force for 2 sorts of atoms
  */
 
+/* create mixed terms from parameters and set them */
+int tersoff_mult_complete_params(t_tersoff_mult_params *p) {
+
+       printf("[moldyn] tersoff parameter completion\n");
+       p->Smixed=sqrt(p->S[0]*p->S[1]);
+       p->Rmixed=sqrt(p->R[0]*p->R[1]);
+       p->Amixed=sqrt(p->A[0]*p->A[1]);
+       p->Bmixed=sqrt(p->B[0]*p->B[1]);
+       p->lambda_m=0.5*(p->lambda[0]+p->lambda[1]);
+       p->mu_m=0.5*(p->mu[0]+p->mu[1]);
+
+       return 0;
+}
+
 /* tersoff 1 body part */
 int tersoff_mult_1bp(t_moldyn *moldyn,t_atom *ai) {
 
@@ -1013,8 +1136,11 @@ int tersoff_mult_2bp(t_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc) {
 
        if(bc) check_per_bound(moldyn,&dist_ij);
 
-       /* save for use in 3bp */ /* REALLY ?!?!?! */
-       exchange->dist_ij=dist_ij;
+       d_ij=v3_norm(&dist_ij);
+
+       /* save for use in 3bp */
+       exchange->dist_ij=dist_ij; /* <- needed ? */
+       exchange->d_ij=d_ij;
 
        /* constants */
        if(num==aj->bnum) {
@@ -1040,11 +1166,6 @@ int tersoff_mult_2bp(t_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc) {
                params->exchange.chi=params->chi;
        }
 
-       d_ij=v3_norm(&dist_ij);
-
-       /* save for use in 3bp */
-       exchange->d_ij=d_ij;
-
        if(d_ij>S)
                return 0;
 
@@ -1072,8 +1193,8 @@ int tersoff_mult_2bp(t_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc) {
 
        /* add forces */
        v3_add(&(ai->f),&(ai->f),&force);
-       /* energy is 0.5 f_r f_c, but we will sum it up twice ... */
-       moldyn->energy+=(0.25*f_r*f_c);
+       /* energy is 0.5 f_r f_c ... */
+       moldyn->energy+=(0.5*f_r*f_c);
 
        /* save for use in 3bp */
        exchange->f_c=f_c;