/* ---------------------------------------------------------------------- LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator http://lammps.sandia.gov, Sandia National Laboratories Steve Plimpton, sjplimp@sandia.gov Copyright (2003) Sandia Corporation. Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains certain rights in this software. This software is distributed under the GNU General Public License. See the README file in the top-level LAMMPS directory. ------------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- Contributing author: Sai Jayaraman (Sandia) ------------------------------------------------------------------------- */ #include "math.h" #include "stdio.h" #include "stdlib.h" #include "string.h" #include "pair_gauss.h" #include "atom.h" #include "comm.h" #include "force.h" #include "neighbor.h" #include "neigh_list.h" #include "memory.h" #include "error.h" using namespace LAMMPS_NS; #define MIN(a,b) ((a) < (b) ? (a) : (b)) #define MAX(a,b) ((a) > (b) ? (a) : (b)) #define EPSILON 1.0e-10 /* ---------------------------------------------------------------------- */ PairGauss::PairGauss(LAMMPS *lmp) :Pair(lmp) { nextra = 1; pvector = new double[1]; } /* ---------------------------------------------------------------------- */ PairGauss::~PairGauss() { delete [] pvector; if (allocated) { memory->destroy_2d_int_array(setflag); memory->destroy_2d_double_array(cutsq); memory->destroy_2d_double_array(cut); memory->destroy_2d_double_array(a); memory->destroy_2d_double_array(b); memory->destroy_2d_double_array(offset); } } /* ---------------------------------------------------------------------- */ void PairGauss::compute(int eflag, int vflag) { int i,j,ii,jj,inum,jnum,itype,jtype; double xtmp,ytmp,ztmp,delx,dely,delz,evdwl,fpair; double r,rsq,r2inv,r6inv,forcelj,rexp; int *ilist,*jlist,*numneigh,**firstneigh; evdwl = 0.0; if (eflag || vflag) ev_setup(eflag,vflag); else evflag = vflag_fdotr = 0; int occ = 0; double **x = atom->x; double **f = atom->f; int *type = atom->type; int nlocal = atom->nlocal; int nall = nlocal + atom->nghost; int newton_pair = force->newton_pair; inum = list->inum; ilist = list->ilist; numneigh = list->numneigh; firstneigh = list->firstneigh; // loop over neighbors of my atoms for (ii = 0; ii < inum; ii++) { i = ilist[ii]; xtmp = x[i][0]; ytmp = x[i][1]; ztmp = x[i][2]; itype = type[i]; jlist = firstneigh[i]; jnum = numneigh[i]; for (jj = 0; jj < jnum; jj++) { j = jlist[jj]; delx = xtmp - x[j][0]; dely = ytmp - x[j][1]; delz = ztmp - x[j][2]; rsq = delx*delx + dely*dely + delz*delz; jtype = type[j]; // define a Gaussian well to be occupied if // the site it interacts with is within the force maximum if (eflag_global && rsq < 0.5/b[itype][jtype]) occ++; if (rsq < cutsq[itype][jtype]) { r2inv = 1.0/rsq; r = sqrt(rsq); forcelj = - 2.0*a[itype][jtype]*b[itype][jtype] * rsq * exp(-b[itype][jtype]*rsq); fpair = forcelj*r2inv; f[i][0] += delx*fpair; f[i][1] += dely*fpair; f[i][2] += delz*fpair; if (newton_pair || j < nlocal) { f[j][0] -= delx*fpair; f[j][1] -= dely*fpair; f[j][2] -= delz*fpair; } if (eflag) evdwl = -(a[itype][jtype]*exp(-b[itype][jtype]*rsq) - offset[itype][jtype]); if (evflag) ev_tally(i,j,nlocal,newton_pair, evdwl,0.0,fpair,delx,dely,delz); } } } if (eflag_global) pvector[0] = occ; if (vflag_fdotr) virial_compute(); } /* ---------------------------------------------------------------------- allocate all arrays ------------------------------------------------------------------------- */ void PairGauss::allocate() { allocated = 1; int n = atom->ntypes; setflag = memory->create_2d_int_array(n+1,n+1,"pair:setflag"); for (int i = 1; i <= n; i++) for (int j = 1; j <= n; j++) setflag[i][j] = 0; cutsq = memory->create_2d_double_array(n+1,n+1,"pair:cutsq"); cut = memory->create_2d_double_array(n+1,n+1,"pair:cut_gauss"); a = memory->create_2d_double_array(n+1,n+1,"pair:a"); b = memory->create_2d_double_array(n+1,n+1,"pair:b"); offset = memory->create_2d_double_array(n+1,n+1,"pair:offset"); } /* ---------------------------------------------------------------------- global settings ------------------------------------------------------------------------- */ void PairGauss::settings(int narg, char **arg) { if (narg != 1) error->all("Illegal pair_style command"); cut_global = atof(arg[0]); // reset cutoffs that have been explicity set if (allocated) { int i,j; for (i = 1; i <= atom->ntypes; i++) for (j = i+1; j <= atom->ntypes; j++) if (setflag[i][j]) cut[i][j] = cut_global; } } /* ---------------------------------------------------------------------- set coeffs for one or more type pairs ------------------------------------------------------------------------- */ void PairGauss::coeff(int narg, char **arg) { if (narg < 4 || narg > 5) error->all("Incorrect args for pair coefficients"); if (!allocated) allocate(); int ilo, ihi, jlo, jhi; force->bounds(arg[0],atom->ntypes,ilo,ihi); force->bounds(arg[1],atom->ntypes,jlo,jhi); double a_one = atof(arg[2]); double b_one = atof(arg[3]); double cut_one = cut_global; if (narg == 5) cut_one = atof(arg[4]); int count = 0; for (int i = ilo; i <= ihi; i++) { for (int j = MAX(jlo,i); j<=jhi; j++) { a[i][j] = a_one; b[i][j] = b_one; cut[i][j] = cut_one; setflag[i][j] = 1; count++ ; } } if (count == 0) error->all("Incorrect args for pair coefficients"); } /* ---------------------------------------------------------------------- init for one type pair i,j and corresponding j,i ------------------------------------------------------------------------- */ double PairGauss::init_one(int i, int j) { if (setflag[i][j] == 0) error->all("All pair coeffs are not set"); if (offset_flag) offset[i][j] = a[i][j]*exp(-b[i][j]*cut[i][j]*cut[i][j]); else offset[i][j] = 0.0; a[j][i] = a[i][j]; b[j][i] = b[i][j]; offset[j][i] = offset[i][j]; return cut[i][j]; } /* ---------------------------------------------------------------------- proc 0 writes to restart file ------------------------------------------------------------------------- */ void PairGauss::write_restart(FILE *fp) { write_restart_settings(fp); int i,j; for (i = 1; i <= atom->ntypes; i++) for (j = i; j <= atom->ntypes; j++) { fwrite(&setflag[i][j],sizeof(int),1,fp); if (setflag[i][j]) { fwrite(&a[i][j],sizeof(double),1,fp); fwrite(&b[i][j],sizeof(double),1,fp); fwrite(&cut[i][j],sizeof(double),1,fp); } } } /* ---------------------------------------------------------------------- proc 0 reads from restart file, bcasts ------------------------------------------------------------------------- */ void PairGauss::read_restart(FILE *fp) { read_restart_settings(fp); allocate(); int i,j; int me = comm->me; for (i = 1; i <= atom->ntypes; i++) for (j = i; j <= atom->ntypes; j++) { if (me == 0) fread(&setflag[i][j],sizeof(int),1,fp); MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world); if (setflag[i][j]) { if (me == 0) { fread(&a[i][j],sizeof(double),1,fp); fread(&b[i][j],sizeof(double),1,fp); fread(&cut[i][j],sizeof(double),1,fp); } MPI_Bcast(&a[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&b[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&cut[i][j],1,MPI_DOUBLE,0,world); } } } /* ---------------------------------------------------------------------- proc 0 writes to restart file ------------------------------------------------------------------------- */ void PairGauss::write_restart_settings(FILE *fp) { fwrite(&cut_global,sizeof(double),1,fp); fwrite(&offset_flag,sizeof(int),1,fp); fwrite(&mix_flag,sizeof(int),1,fp); } /* ---------------------------------------------------------------------- proc 0 reads from restart file, bcasts ------------------------------------------------------------------------- */ void PairGauss::read_restart_settings(FILE *fp) { if (comm->me == 0) { fread(&cut_global,sizeof(double),1,fp); fread(&offset_flag,sizeof(int),1,fp); fread(&mix_flag,sizeof(int),1,fp); } MPI_Bcast(&cut_global,1,MPI_DOUBLE,0,world); MPI_Bcast(&offset_flag,1,MPI_INT,0,world); MPI_Bcast(&mix_flag,1,MPI_INT,0,world); } /* ---------------------------------------------------------------------- */ double PairGauss::single(int i, int j, int itype, int jtype, double rsq, double factor_coul, double factor_lj, double &fforce) { double r2inv,r6inv,forcelj,philj,r; r = sqrt(rsq); r2inv = 1.0/rsq; philj = -(a[itype][jtype]*exp(-b[itype][jtype]*rsq) - offset[itype][jtype]); forcelj = -2.0*a[itype][jtype]*b[itype][jtype]*rsq*exp(-b[itype][jtype]*rsq); fforce = forcelj*r2inv; return philj; } /* ---------------------------------------------------------------------- */ void *PairGauss::extract(char *str, int &dim) { dim = 2; if (strcmp(str,"a") == 0) return (void *) a; return NULL; }