/* ---------------------------------------------------------------------- LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator www.cs.sandia.gov/~sjplimp/lammps.html Steve Plimpton, sjplimp@sandia.gov, Sandia National Laboratories 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: Paul Crozier (SNL) ------------------------------------------------------------------------- */ #include "math.h" #include "stdio.h" #include "stdlib.h" #include "string.h" #include "pair_lj_cut.h" #include "atom.h" #include "comm.h" #include "force.h" #include "neighbor.h" #include "update.h" #include "integrate.h" #include "respa.h" #include "memory.h" #include "error.h" #define MIN(a,b) ((a) < (b) ? (a) : (b)) #define MAX(a,b) ((a) > (b) ? (a) : (b)) /* ---------------------------------------------------------------------- */ PairLJCut::PairLJCut() { respa_enable = 1; } /* ---------------------------------------------------------------------- free all arrays ------------------------------------------------------------------------- */ PairLJCut::~PairLJCut() { 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(epsilon); memory->destroy_2d_double_array(sigma); memory->destroy_2d_double_array(lj1); memory->destroy_2d_double_array(lj2); memory->destroy_2d_double_array(lj3); memory->destroy_2d_double_array(lj4); memory->destroy_2d_double_array(offset); } } /* ---------------------------------------------------------------------- */ void PairLJCut::compute(int eflag, int vflag) { int i,j,k,numneigh,itype,jtype; double xtmp,ytmp,ztmp,delx,dely,delz; double rsq,r2inv,r6inv,forcelj,fforce,factor_lj,philj; int *neighs; double **f; eng_vdwl = 0.0; if (vflag) for (i = 0; i < 6; i++) virial[i] = 0.0; if (vflag == 2) f = update->f_pair; else f = atom->f; double **x = atom->x; int *type = atom->type; int nlocal = atom->nlocal; int nall = atom->nlocal + atom->nghost; double *special_lj = force->special_lj; int newton_pair = force->newton_pair; // loop over neighbors of my atoms for (i = 0; i < nlocal; i++) { xtmp = x[i][0]; ytmp = x[i][1]; ztmp = x[i][2]; itype = type[i]; neighs = neighbor->firstneigh[i]; numneigh = neighbor->numneigh[i]; for (k = 0; k < numneigh; k++) { j = neighs[k]; if (j < nall) factor_lj = 1.0; else { factor_lj = special_lj[j/nall]; j %= nall; } 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]; if (rsq < cutsq[itype][jtype]) { r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]); fforce = factor_lj*forcelj*r2inv; f[i][0] += delx*fforce; f[i][1] += dely*fforce; f[i][2] += delz*fforce; if (newton_pair || j < nlocal) { f[j][0] -= delx*fforce; f[j][1] -= dely*fforce; f[j][2] -= delz*fforce; } if (eflag) { philj = r6inv*(lj3[itype][jtype]*r6inv-lj4[itype][jtype]) - offset[itype][jtype]; if (newton_pair || j < nlocal) eng_vdwl += factor_lj*philj; else eng_vdwl += 0.5*factor_lj*philj; } if (vflag == 1) { if (newton_pair || j < nlocal) { virial[0] += delx*delx*fforce; virial[1] += dely*dely*fforce; virial[2] += delz*delz*fforce; virial[3] += delx*dely*fforce; virial[4] += delx*delz*fforce; virial[5] += dely*delz*fforce; } else { virial[0] += 0.5*delx*delx*fforce; virial[1] += 0.5*dely*dely*fforce; virial[2] += 0.5*delz*delz*fforce; virial[3] += 0.5*delx*dely*fforce; virial[4] += 0.5*delx*delz*fforce; virial[5] += 0.5*dely*delz*fforce; } } } } } if (vflag == 2) virial_compute(); } /* ---------------------------------------------------------------------- */ void PairLJCut::compute_inner() { int i,j,k,numneigh,itype,jtype; double xtmp,ytmp,ztmp,delx,dely,delz; double rsq,r2inv,r6inv,forcelj,fforce,factor_lj; double rsw; int *neighs; double **f; f = atom->f; double **x = atom->x; int *type = atom->type; int nlocal = atom->nlocal; int nall = atom->nlocal + atom->nghost; double *special_lj = force->special_lj; int newton_pair = force->newton_pair; double cut_out_on = cut_respa[0]; double cut_out_off = cut_respa[1]; double cut_out_diff = cut_out_off - cut_out_on; double cut_out_on_sq = cut_out_on*cut_out_on; double cut_out_off_sq = cut_out_off*cut_out_off; // loop over neighbors of my atoms for (i = 0; i < nlocal; i++) { xtmp = x[i][0]; ytmp = x[i][1]; ztmp = x[i][2]; itype = type[i]; neighs = neighbor->firstneigh_inner[i]; numneigh = neighbor->numneigh_inner[i]; for (k = 0; k < numneigh; k++) { j = neighs[k]; if (j < nall) factor_lj = 1.0; else { factor_lj = special_lj[j/nall]; j %= nall; } delx = xtmp - x[j][0]; dely = ytmp - x[j][1]; delz = ztmp - x[j][2]; rsq = delx*delx + dely*dely + delz*delz; if (rsq < cut_out_off_sq) { r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; jtype = type[j]; forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]); fforce = factor_lj*forcelj*r2inv; if (rsq > cut_out_on_sq) { rsw = (sqrt(rsq) - cut_out_on)/cut_out_diff; fforce *= 1.0 - rsw*rsw*(3.0 - 2.0*rsw); } f[i][0] += delx*fforce; f[i][1] += dely*fforce; f[i][2] += delz*fforce; if (newton_pair || j < nlocal) { f[j][0] -= delx*fforce; f[j][1] -= dely*fforce; f[j][2] -= delz*fforce; } } } } } /* ---------------------------------------------------------------------- */ void PairLJCut::compute_middle() { int i,j,k,numneigh,itype,jtype; double xtmp,ytmp,ztmp,delx,dely,delz; double rsq,r2inv,r6inv,forcelj,fforce,factor_lj; double rsw; int *neighs; double **f; f = atom->f; double **x = atom->x; int *type = atom->type; int nlocal = atom->nlocal; int nall = atom->nlocal + atom->nghost; double *special_lj = force->special_lj; int newton_pair = force->newton_pair; double cut_in_off = cut_respa[0]; double cut_in_on = cut_respa[1]; double cut_out_on = cut_respa[2]; double cut_out_off = cut_respa[3]; double cut_in_diff = cut_in_on - cut_in_off; double cut_out_diff = cut_out_off - cut_out_on; double cut_in_off_sq = cut_in_off*cut_in_off; double cut_in_on_sq = cut_in_on*cut_in_on; double cut_out_on_sq = cut_out_on*cut_out_on; double cut_out_off_sq = cut_out_off*cut_out_off; // loop over neighbors of my atoms for (i = 0; i < nlocal; i++) { xtmp = x[i][0]; ytmp = x[i][1]; ztmp = x[i][2]; itype = type[i]; neighs = neighbor->firstneigh_middle[i]; numneigh = neighbor->numneigh_middle[i]; for (k = 0; k < numneigh; k++) { j = neighs[k]; if (j < nall) factor_lj = 1.0; else { factor_lj = special_lj[j/nall]; j %= nall; } delx = xtmp - x[j][0]; dely = ytmp - x[j][1]; delz = ztmp - x[j][2]; rsq = delx*delx + dely*dely + delz*delz; if (rsq < cut_out_off_sq && rsq > cut_in_off_sq) { r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; jtype = type[j]; forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]); fforce = factor_lj*forcelj*r2inv; if (rsq < cut_in_on_sq) { rsw = (sqrt(rsq) - cut_in_off)/cut_in_diff; fforce *= rsw*rsw*(3.0 - 2.0*rsw); } if (rsq > cut_out_on_sq) { rsw = (sqrt(rsq) - cut_out_on)/cut_out_diff; fforce *= 1.0 + rsw*rsw*(2.0*rsw - 3.0); } f[i][0] += delx*fforce; f[i][1] += dely*fforce; f[i][2] += delz*fforce; if (newton_pair || j < nlocal) { f[j][0] -= delx*fforce; f[j][1] -= dely*fforce; f[j][2] -= delz*fforce; } } } } } /* ---------------------------------------------------------------------- */ void PairLJCut::compute_outer(int eflag, int vflag) { int i,j,k,numneigh,itype,jtype; double xtmp,ytmp,ztmp,delx,dely,delz; double rsq,r2inv,r6inv,forcelj,fforce,factor_lj,philj; double rsw; int *neighs; double **f; eng_vdwl = 0.0; if (vflag) for (i = 0; i < 6; i++) virial[i] = 0.0; f = atom->f; double **x = atom->x; int *type = atom->type; int nlocal = atom->nlocal; int nall = atom->nlocal + atom->nghost; double *special_lj = force->special_lj; int newton_pair = force->newton_pair; double cut_in_off = cut_respa[2]; double cut_in_on = cut_respa[3]; double cut_in_diff = cut_in_on - cut_in_off; double cut_in_off_sq = cut_in_off*cut_in_off; double cut_in_on_sq = cut_in_on*cut_in_on; // loop over neighbors of my atoms for (i = 0; i < nlocal; i++) { xtmp = x[i][0]; ytmp = x[i][1]; ztmp = x[i][2]; itype = type[i]; neighs = neighbor->firstneigh[i]; numneigh = neighbor->numneigh[i]; for (k = 0; k < numneigh; k++) { j = neighs[k]; if (j < nall) factor_lj = 1.0; else { factor_lj = special_lj[j/nall]; j %= nall; } 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]; if (rsq < cutsq[itype][jtype]) { if (rsq > cut_in_off_sq) { r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]); fforce = factor_lj*forcelj*r2inv; if (rsq < cut_in_on_sq) { rsw = (sqrt(rsq) - cut_in_off)/cut_in_diff; fforce *= rsw*rsw*(3.0 - 2.0*rsw); } f[i][0] += delx*fforce; f[i][1] += dely*fforce; f[i][2] += delz*fforce; if (newton_pair || j < nlocal) { f[j][0] -= delx*fforce; f[j][1] -= dely*fforce; f[j][2] -= delz*fforce; } } if (eflag) { r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; philj = r6inv*(lj3[itype][jtype]*r6inv-lj4[itype][jtype]) - offset[itype][jtype]; if (newton_pair || j < nlocal) eng_vdwl += factor_lj*philj; else eng_vdwl += 0.5*factor_lj*philj; } if (vflag) { if (rsq <= cut_in_off_sq) { r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]); fforce = factor_lj*forcelj*r2inv; } else if (rsq < cut_in_on_sq) fforce = factor_lj*forcelj*r2inv; if (newton_pair || j < nlocal) { virial[0] += delx*delx*fforce; virial[1] += dely*dely*fforce; virial[2] += delz*delz*fforce; virial[3] += delx*dely*fforce; virial[4] += delx*delz*fforce; virial[5] += dely*delz*fforce; } else { virial[0] += 0.5*delx*delx*fforce; virial[1] += 0.5*dely*dely*fforce; virial[2] += 0.5*delz*delz*fforce; virial[3] += 0.5*delx*dely*fforce; virial[4] += 0.5*delx*delz*fforce; virial[5] += 0.5*dely*delz*fforce; } } } } } } /* ---------------------------------------------------------------------- allocate all arrays ------------------------------------------------------------------------- */ void PairLJCut::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 = i; 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"); epsilon = memory->create_2d_double_array(n+1,n+1,"pair:epsilon"); sigma = memory->create_2d_double_array(n+1,n+1,"pair:sigma"); lj1 = memory->create_2d_double_array(n+1,n+1,"pair:lj1"); lj2 = memory->create_2d_double_array(n+1,n+1,"pair:lj2"); lj3 = memory->create_2d_double_array(n+1,n+1,"pair:lj3"); lj4 = memory->create_2d_double_array(n+1,n+1,"pair:lj4"); offset = memory->create_2d_double_array(n+1,n+1,"pair:offset"); } /* ---------------------------------------------------------------------- global settings ------------------------------------------------------------------------- */ void PairLJCut::settings(int narg, char **arg) { if (narg != 1) error->all("Illegal pair_style command"); cut_global = atof(arg[0]); // reset cutoffs that have been explicitly 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 PairLJCut::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 epsilon_one = atof(arg[2]); double sigma_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++) { epsilon[i][j] = epsilon_one; sigma[i][j] = sigma_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 PairLJCut::init_one(int i, int j) { if (setflag[i][j] == 0) { epsilon[i][j] = mix_energy(epsilon[i][i],epsilon[j][j], sigma[i][i],sigma[j][j]); sigma[i][j] = mix_distance(sigma[i][i],sigma[j][j]); cut[i][j] = mix_distance(cut[i][i],cut[j][j]); } lj1[i][j] = 48.0 * epsilon[i][j] * pow(sigma[i][j],12.0); lj2[i][j] = 24.0 * epsilon[i][j] * pow(sigma[i][j],6.0); lj3[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],12.0); lj4[i][j] = 4.0 * epsilon[i][j] * pow(sigma[i][j],6.0); if (offset_flag) { double ratio = sigma[i][j] / cut[i][j]; offset[i][j] = 4.0 * epsilon[i][j] * (pow(ratio,12.0) - pow(ratio,6.0)); } else offset[i][j] = 0.0; lj1[j][i] = lj1[i][j]; lj2[j][i] = lj2[i][j]; lj3[j][i] = lj3[i][j]; lj4[j][i] = lj4[i][j]; offset[j][i] = offset[i][j]; // set & error check interior rRESPA cutoff if (strcmp(update->integrate_style,"respa") == 0) { if (((Respa *) update->integrate)->level_inner >= 0) { cut_respa = ((Respa *) update->integrate)->cutoff; if (cut[i][j] < cut_respa[3]) error->all("Pair cutoff < Respa interior cutoff"); } } else cut_respa = NULL; // compute I,J contribution to long-range tail correction // count total # of atoms of type I and J via Allreduce if (tail_flag) { int *type = atom->type; int nlocal = atom->nlocal; double count[2],all[2]; count[0] = count[1] = 0.0; for (int k = 0; k < nlocal; k++) { if (type[k] == i) count[0] += 1.0; if (type[k] == j) count[1] += 1.0; } MPI_Allreduce(count,all,2,MPI_DOUBLE,MPI_SUM,world); double PI = 4.0*atan(1.0); double sig2 = sigma[i][j]*sigma[i][j]; double sig6 = sig2*sig2*sig2; double rc3 = cut[i][j]*cut[i][j]*cut[i][j]; double rc6 = rc3*rc3; double rc9 = rc3*rc6; etail_ij = 8.0*PI*all[0]*all[1]*epsilon[i][j] * sig6 * (sig6 - 3.0*rc6) / (9.0*rc9); ptail_ij = 16.0*PI*all[0]*all[1]*epsilon[i][j] * sig6 * (2.0*sig6 - 3.0*rc6) / (9.0*rc9); } return cut[i][j]; } /* ---------------------------------------------------------------------- proc 0 writes to restart file ------------------------------------------------------------------------- */ void PairLJCut::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(&epsilon[i][j],sizeof(double),1,fp); fwrite(&sigma[i][j],sizeof(double),1,fp); fwrite(&cut[i][j],sizeof(double),1,fp); } } } /* ---------------------------------------------------------------------- proc 0 reads from restart file, bcasts ------------------------------------------------------------------------- */ void PairLJCut::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(&epsilon[i][j],sizeof(double),1,fp); fread(&sigma[i][j],sizeof(double),1,fp); fread(&cut[i][j],sizeof(double),1,fp); } MPI_Bcast(&epsilon[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&sigma[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&cut[i][j],1,MPI_DOUBLE,0,world); } } } /* ---------------------------------------------------------------------- proc 0 writes to restart file ------------------------------------------------------------------------- */ void PairLJCut::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 PairLJCut::read_restart_settings(FILE *fp) { int me = comm->me; if (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); } /* ---------------------------------------------------------------------- */ void PairLJCut::single(int i, int j, int itype, int jtype, double rsq, double factor_coul, double factor_lj, int eflag, One &one) { double r2inv,r6inv,forcelj,philj; r2inv = 1.0/rsq; r6inv = r2inv*r2inv*r2inv; forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]); one.fforce = factor_lj*forcelj*r2inv; if (eflag) { philj = r6inv*(lj3[itype][jtype]*r6inv-lj4[itype][jtype]) - offset[itype][jtype]; one.eng_vdwl = factor_lj*philj; one.eng_coul = 0.0; } }