// clang-format off /* ---------------------------------------------------------------------- LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator https://www.lammps.org/, Sandia National Laboratories LAMMPS development team: developers@lammps.org 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. ------------------------------------------------------------------------- */ #include "pair_lj_class2_coul_long_soft.h" #include "atom.h" #include "comm.h" #include "error.h" #include "force.h" #include "kspace.h" #include "math_const.h" #include "memory.h" #include "neigh_list.h" #include "neighbor.h" #include #include using namespace LAMMPS_NS; using namespace MathConst; #define EWALD_F 1.12837917 #define EWALD_P 0.3275911 #define A1 0.254829592 #define A2 -0.284496736 #define A3 1.421413741 #define A4 -1.453152027 #define A5 1.061405429 /* ---------------------------------------------------------------------- */ PairLJClass2CoulLongSoft::PairLJClass2CoulLongSoft(LAMMPS *lmp) : Pair(lmp) { ewaldflag = pppmflag = 1; writedata = 1; } /* ---------------------------------------------------------------------- */ PairLJClass2CoulLongSoft::~PairLJClass2CoulLongSoft() { if (copymode) return; if (allocated) { memory->destroy(setflag); memory->destroy(cutsq); memory->destroy(cut_lj); memory->destroy(cut_ljsq); memory->destroy(epsilon); memory->destroy(sigma); memory->destroy(lambda); memory->destroy(lj1); memory->destroy(lj2); memory->destroy(lj3); memory->destroy(lj4); memory->destroy(offset); } } /* ---------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::compute(int eflag, int vflag) { int i,j,ii,jj,inum,jnum,itype,jtype; double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,evdwl,ecoul,fpair; double rsq,r,forcecoul,forcelj; double grij,expm2,prefactor,t,erfc; double factor_coul,factor_lj; double denc, denlj, r4sig6; int *ilist,*jlist,*numneigh,**firstneigh; evdwl = ecoul = 0.0; ev_init(eflag,vflag); double **x = atom->x; double **f = atom->f; double *q = atom->q; int *type = atom->type; int nlocal = atom->nlocal; double *special_coul = force->special_coul; double *special_lj = force->special_lj; int newton_pair = force->newton_pair; double qqrd2e = force->qqrd2e; 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]; qtmp = q[i]; 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]; factor_lj = special_lj[sbmask(j)]; factor_coul = special_coul[sbmask(j)]; j &= NEIGHMASK; 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_coulsq) { r = sqrt(rsq); grij = g_ewald * r; expm2 = exp(-grij*grij); t = 1.0 / (1.0 + EWALD_P*grij); erfc = t * (A1+t*(A2+t*(A3+t*(A4+t*A5)))) * expm2; denc = sqrt(lj4[itype][jtype] + rsq); prefactor = qqrd2e * lj1[itype][jtype] * qtmp*q[j] / (denc*denc*denc); forcecoul = prefactor * (erfc + EWALD_F*grij*expm2); if (factor_coul < 1.0) forcecoul -= (1.0-factor_coul)*prefactor; } else forcecoul = 0.0; if (rsq < cut_ljsq[itype][jtype]) { r4sig6 = rsq*rsq / lj2[itype][jtype]; denlj = lj3[itype][jtype] + rsq*r4sig6; forcelj = lj1[itype][jtype] * epsilon[itype][jtype] * (18.0*r4sig6*pow(denlj, -2.5) - 18.0*r4sig6*pow(denlj, -2)); } else forcelj = 0.0; fpair = forcecoul + factor_lj*forcelj; 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) { if (rsq < cut_coulsq) { prefactor = qqrd2e * lj1[itype][jtype] * qtmp*q[j] / denc; ecoul = prefactor*erfc; if (factor_coul < 1.0) ecoul -= (1.0-factor_coul)*prefactor; } else ecoul = 0.0; if (rsq < cut_ljsq[itype][jtype]) { r4sig6 = rsq*rsq / lj2[itype][jtype]; denlj = lj3[itype][jtype] + rsq*r4sig6; evdwl = lj1[itype][jtype] * epsilon[itype][jtype] * (2.0/(denlj*sqrt(denlj)) - 3.0/denlj) - offset[itype][jtype]; evdwl *= factor_lj; } else evdwl = 0.0; } if (evflag) ev_tally(i,j,nlocal,newton_pair, evdwl,ecoul,fpair,delx,dely,delz); } } } if (vflag_fdotr) virial_fdotr_compute(); } /* ---------------------------------------------------------------------- allocate all arrays ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::allocate() { allocated = 1; int n = atom->ntypes; memory->create(setflag,n+1,n+1,"pair:setflag"); for (int i = 1; i <= n; i++) for (int j = i; j <= n; j++) setflag[i][j] = 0; memory->create(cutsq,n+1,n+1,"pair:cutsq"); memory->create(cut_lj,n+1,n+1,"pair:cut_lj"); memory->create(cut_ljsq,n+1,n+1,"pair:cut_ljsq"); memory->create(epsilon,n+1,n+1,"pair:epsilon"); memory->create(sigma,n+1,n+1,"pair:sigma"); memory->create(lambda,n+1,n+1,"pair:lambda"); memory->create(lj1,n+1,n+1,"pair:lj1"); memory->create(lj2,n+1,n+1,"pair:lj2"); memory->create(lj3,n+1,n+1,"pair:lj3"); memory->create(lj4,n+1,n+1,"pair:lj4"); memory->create(offset,n+1,n+1,"pair:offset"); } /* ---------------------------------------------------------------------- global settings ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::settings(int narg, char **arg) { if (narg < 4 || narg > 5) error->all(FLERR,"Illegal pair_style command"); nlambda = utils::numeric(FLERR,arg[0],false,lmp); alphalj = utils::numeric(FLERR,arg[1],false,lmp); alphac = utils::numeric(FLERR,arg[2],false,lmp); cut_lj_global = utils::numeric(FLERR,arg[3],false,lmp); if (narg == 4) cut_coul = cut_lj_global; else cut_coul = utils::numeric(FLERR,arg[4],false,lmp); // reset cutoffs that have been explicitly set if (allocated) { int i,j; for (i = 1; i <= atom->ntypes; i++) for (j = i; j <= atom->ntypes; j++) if (setflag[i][j]) cut_lj[i][j] = cut_lj_global; } } /* ---------------------------------------------------------------------- set coeffs for one or more type pairs ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::coeff(int narg, char **arg) { if (narg < 5 || narg > 6) error->all(FLERR,"Incorrect args for pair coefficients"); if (!allocated) allocate(); int ilo,ihi,jlo,jhi; utils::bounds(FLERR,arg[0],1,atom->ntypes,ilo,ihi,error); utils::bounds(FLERR,arg[1],1,atom->ntypes,jlo,jhi,error); double epsilon_one = utils::numeric(FLERR,arg[2],false,lmp); double sigma_one = utils::numeric(FLERR,arg[3],false,lmp); double lambda_one = utils::numeric(FLERR,arg[4],false,lmp); if (sigma_one <= 0.0) error->all(FLERR,"Incorrect args for pair coefficients"); double cut_lj_one = cut_lj_global; if (narg == 6) cut_lj_one = utils::numeric(FLERR,arg[5],false,lmp); 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; lambda[i][j] = lambda_one; cut_lj[i][j] = cut_lj_one; setflag[i][j] = 1; count++; } } if (count == 0) error->all(FLERR,"Incorrect args for pair coefficients"); } /* ---------------------------------------------------------------------- init specific to this pair style ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::init_style() { if (!atom->q_flag) error->all(FLERR, "Pair style lj/class2/coul/long/soft requires atom attribute q"); neighbor->add_request(this); cut_coulsq = cut_coul * cut_coul; // ensure use of KSpace long-range solver, set g_ewald if (force->kspace == nullptr) error->all(FLERR,"Pair style requires a KSpace style"); g_ewald = force->kspace->g_ewald; } /* ---------------------------------------------------------------------- init for one type pair i,j and corresponding j,i ------------------------------------------------------------------------- */ double PairLJClass2CoulLongSoft::init_one(int i, int j) { // always mix epsilon,sigma via sixthpower rules // mix distance via user-defined rule if (setflag[i][j] == 0) { epsilon[i][j] = 2.0 * sqrt(epsilon[i][i]*epsilon[j][j]) * pow(sigma[i][i],3.0) * pow(sigma[j][j],3.0) / (pow(sigma[i][i],6.0) + pow(sigma[j][j],6.0)); sigma[i][j] = pow((0.5 * (pow(sigma[i][i],6.0) + pow(sigma[j][j],6.0))),1.0/6.0); if (lambda[i][i] != lambda[j][j]) error->all(FLERR,"Pair lj/class2/coul/long/soft different lambda values in mix"); lambda[i][j] = lambda[i][i]; cut_lj[i][j] = mix_distance(cut_lj[i][i],cut_lj[j][j]); } double cut = MAX(cut_lj[i][j],cut_coul); cut_ljsq[i][j] = cut_lj[i][j] * cut_lj[i][j]; lj1[i][j] = pow(lambda[i][j], nlambda); lj2[i][j] = pow(sigma[i][j], 6.0); lj3[i][j] = alphalj * (1.0 - lambda[i][j])*(1.0 - lambda[i][j]); lj4[i][j] = alphac * (1.0 - lambda[i][j])*(1.0 - lambda[i][j]); if (offset_flag && (cut_lj[i][j] > 0.0)) { double denlj = lj3[i][j] + pow(cut_lj[i][j] / sigma[i][j], 6.0); offset[i][j] = epsilon[i][j] * (2.0*pow(denlj,-1.5) - 3.0*pow(denlj,-1.0)); } else offset[i][j] = 0.0; epsilon[j][i] = epsilon[i][j]; sigma[j][i] = sigma[i][j]; lambda[j][i] = lambda[i][j]; cut_ljsq[j][i] = cut_ljsq[i][j]; 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]; // 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 sig3 = sigma[i][j]*sigma[i][j]*sigma[i][j]; double sig6 = sig3*sig3; double rc3 = cut_lj[i][j]*cut_lj[i][j]*cut_lj[i][j]; double rc6 = rc3*rc3; etail_ij = 2.0*MY_PI*all[0]*all[1]*lj1[i][j] *epsilon[i][j] * sig6 * (sig3 - 3.0*rc3) / (3.0*rc6); ptail_ij = 2.0*MY_PI*all[0]*all[1]*lj1[i][j] *epsilon[i][j] * sig6 * (sig3 - 2.0*rc3) / rc6; } return cut; } /* ---------------------------------------------------------------------- proc 0 writes to restart file ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::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(&lambda[i][j],sizeof(double),1,fp); fwrite(&cut_lj[i][j],sizeof(double),1,fp); } } } /* ---------------------------------------------------------------------- proc 0 reads from restart file, bcasts ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::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) utils::sfread(FLERR,&setflag[i][j],sizeof(int),1,fp,nullptr,error); MPI_Bcast(&setflag[i][j],1,MPI_INT,0,world); if (setflag[i][j]) { if (me == 0) { utils::sfread(FLERR,&epsilon[i][j],sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&sigma[i][j],sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&lambda[i][j],sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&cut_lj[i][j],sizeof(double),1,fp,nullptr,error); } MPI_Bcast(&epsilon[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&sigma[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&lambda[i][j],1,MPI_DOUBLE,0,world); MPI_Bcast(&cut_lj[i][j],1,MPI_DOUBLE,0,world); } } } /* ---------------------------------------------------------------------- proc 0 writes to restart file ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::write_restart_settings(FILE *fp) { fwrite(&nlambda,sizeof(double),1,fp); fwrite(&alphalj,sizeof(double),1,fp); fwrite(&alphac,sizeof(double),1,fp); fwrite(&cut_lj_global,sizeof(double),1,fp); fwrite(&cut_coul,sizeof(double),1,fp); fwrite(&offset_flag,sizeof(int),1,fp); fwrite(&mix_flag,sizeof(int),1,fp); fwrite(&tail_flag,sizeof(int),1,fp); } /* ---------------------------------------------------------------------- proc 0 reads from restart file, bcasts ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::read_restart_settings(FILE *fp) { if (comm->me == 0) { utils::sfread(FLERR,&nlambda,sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&alphalj,sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&alphac,sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&cut_lj_global,sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&cut_coul,sizeof(double),1,fp,nullptr,error); utils::sfread(FLERR,&offset_flag,sizeof(int),1,fp,nullptr,error); utils::sfread(FLERR,&mix_flag,sizeof(int),1,fp,nullptr,error); utils::sfread(FLERR,&tail_flag,sizeof(int),1,fp,nullptr,error); } MPI_Bcast(&nlambda,1,MPI_DOUBLE,0,world); MPI_Bcast(&alphalj,1,MPI_DOUBLE,0,world); MPI_Bcast(&alphac,1,MPI_DOUBLE,0,world); MPI_Bcast(&cut_lj_global,1,MPI_DOUBLE,0,world); MPI_Bcast(&cut_coul,1,MPI_DOUBLE,0,world); MPI_Bcast(&offset_flag,1,MPI_INT,0,world); MPI_Bcast(&mix_flag,1,MPI_INT,0,world); MPI_Bcast(&tail_flag,1,MPI_INT,0,world); } /* ---------------------------------------------------------------------- proc 0 writes to data file ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::write_data(FILE *fp) { for (int i = 1; i <= atom->ntypes; i++) fprintf(fp,"%d %g %g %g\n",i,epsilon[i][i],sigma[i][i],lambda[i][i]); } /* ---------------------------------------------------------------------- proc 0 writes all pairs to data file ------------------------------------------------------------------------- */ void PairLJClass2CoulLongSoft::write_data_all(FILE *fp) { for (int i = 1; i <= atom->ntypes; i++) for (int j = i; j <= atom->ntypes; j++) fprintf(fp,"%d %d %g %g %g %g\n",i,j,epsilon[i][j],sigma[i][j], lambda[i][j],cut_lj[i][j]); } /* ---------------------------------------------------------------------- */ double PairLJClass2CoulLongSoft::single(int i, int j, int itype, int jtype, double rsq, double factor_coul, double factor_lj, double &fforce) { double denc,r,denlj,r4sig6,grij,expm2,t,erfc,prefactor; double forcecoul,forcelj,phicoul,philj; if (rsq < cut_coulsq) { r = sqrt(rsq); grij = g_ewald * r; expm2 = exp(-grij*grij); t = 1.0 / (1.0 + EWALD_P*grij); erfc = t * (A1+t*(A2+t*(A3+t*(A4+t*A5)))) * expm2; denc = sqrt(lj4[itype][jtype] + rsq); prefactor = force->qqrd2e * lj1[itype][jtype] * atom->q[i]*atom->q[j] / (denc*denc*denc); forcecoul = prefactor * (erfc + EWALD_F*grij*expm2); if (factor_coul < 1.0) forcecoul -= (1.0-factor_coul)*prefactor; } else forcecoul = 0.0; if (rsq < cut_ljsq[itype][jtype]) { r4sig6 = rsq*rsq / lj2[itype][jtype]; denlj = lj3[itype][jtype] + rsq*r4sig6; forcelj = lj1[itype][jtype] * epsilon[itype][jtype] * (18.0*r4sig6/(denlj*denlj*sqrt(denlj)) - 18.0*r4sig6/(denlj*denlj)); } else forcelj = 0.0; fforce = forcecoul + factor_lj*forcelj; double eng = 0.0; if (rsq < cut_coulsq) { prefactor = force->qqrd2e * lj1[itype][jtype] * atom->q[i]*atom->q[j] / denc; phicoul = prefactor*erfc; if (factor_coul < 1.0) phicoul -= (1.0-factor_coul)*prefactor; eng += phicoul; } if (rsq < cut_ljsq[itype][jtype]) { philj = lj1[itype][jtype] * epsilon[itype][jtype] * (2.0/(denlj*sqrt(denlj)) - 3.0/denlj) - offset[itype][jtype]; eng += factor_lj*philj; } return eng; } /* ---------------------------------------------------------------------- */ void *PairLJClass2CoulLongSoft::extract(const char *str, int &dim) { dim = 0; if (strcmp(str,"cut_coul") == 0) return (void *) &cut_coul; dim = 2; if (strcmp(str,"epsilon") == 0) return (void *) epsilon; if (strcmp(str,"sigma") == 0) return (void *) sigma; if (strcmp(str,"lambda") == 0) return (void *) lambda; return nullptr; }