523 lines
17 KiB
C++
523 lines
17 KiB
C++
/* ----------------------------------------------------------------------
|
|
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.
|
|
------------------------------------------------------------------------- */
|
|
|
|
/* ----------------------------------------------------------------------
|
|
Contributing author: Eduardo Bringa (LLNL)
|
|
------------------------------------------------------------------------- */
|
|
|
|
#include "pair_buck_coul_cut.h"
|
|
|
|
#include "atom.h"
|
|
#include "comm.h"
|
|
#include "error.h"
|
|
#include "force.h"
|
|
#include "math_const.h"
|
|
#include "memory.h"
|
|
#include "neigh_list.h"
|
|
#include "neighbor.h"
|
|
|
|
#include <cmath>
|
|
#include <cstring>
|
|
|
|
using namespace LAMMPS_NS;
|
|
using namespace MathConst;
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
PairBuckCoulCut::PairBuckCoulCut(LAMMPS *lmp) : Pair(lmp)
|
|
{
|
|
born_matrix_enable = 1;
|
|
writedata = 1;
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
PairBuckCoulCut::~PairBuckCoulCut()
|
|
{
|
|
if (copymode) return;
|
|
|
|
if (allocated) {
|
|
memory->destroy(setflag);
|
|
memory->destroy(cutsq);
|
|
|
|
memory->destroy(cut_lj);
|
|
memory->destroy(cut_ljsq);
|
|
memory->destroy(cut_coul);
|
|
memory->destroy(cut_coulsq);
|
|
memory->destroy(a);
|
|
memory->destroy(rho);
|
|
memory->destroy(c);
|
|
memory->destroy(rhoinv);
|
|
memory->destroy(buck1);
|
|
memory->destroy(buck2);
|
|
memory->destroy(offset);
|
|
}
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::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, r2inv, r6inv, r, rexp, forcecoul, forcebuck, factor_coul, factor_lj;
|
|
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]) {
|
|
r2inv = 1.0 / rsq;
|
|
r = sqrt(rsq);
|
|
|
|
if (rsq < cut_coulsq[itype][jtype])
|
|
forcecoul = qqrd2e * qtmp * q[j] / r;
|
|
else
|
|
forcecoul = 0.0;
|
|
|
|
if (rsq < cut_ljsq[itype][jtype]) {
|
|
r6inv = r2inv * r2inv * r2inv;
|
|
rexp = exp(-r * rhoinv[itype][jtype]);
|
|
forcebuck = buck1[itype][jtype] * r * rexp - buck2[itype][jtype] * r6inv;
|
|
} else
|
|
forcebuck = 0.0;
|
|
|
|
fpair = (factor_coul * forcecoul + factor_lj * forcebuck) * 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) {
|
|
if (rsq < cut_coulsq[itype][jtype])
|
|
ecoul = factor_coul * qqrd2e * qtmp * q[j] / r;
|
|
else
|
|
ecoul = 0.0;
|
|
if (rsq < cut_ljsq[itype][jtype]) {
|
|
evdwl = a[itype][jtype] * rexp - c[itype][jtype] * r6inv - 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 PairBuckCoulCut::allocate()
|
|
{
|
|
allocated = 1;
|
|
int np1 = atom->ntypes + 1;
|
|
|
|
memory->create(setflag, np1, np1, "pair:setflag");
|
|
|
|
for (int i = 1; i < np1; i++)
|
|
for (int j = i; j < np1; j++) setflag[i][j] = 0;
|
|
|
|
memory->create(cutsq, np1, np1, "pair:cutsq");
|
|
|
|
memory->create(cut_lj, np1, np1, "pair:cut_lj");
|
|
memory->create(cut_ljsq, np1, np1, "pair:cut_ljsq");
|
|
memory->create(cut_coul, np1, np1, "pair:cut_coul");
|
|
memory->create(cut_coulsq, np1, np1, "pair:cut_coulsq");
|
|
memory->create(a, np1, np1, "pair:a");
|
|
memory->create(rho, np1, np1, "pair:rho");
|
|
memory->create(c, np1, np1, "pair:c");
|
|
memory->create(rhoinv, np1, np1, "pair:rhoinv");
|
|
memory->create(buck1, np1, np1, "pair:buck1");
|
|
memory->create(buck2, np1, np1, "pair:buck2");
|
|
memory->create(offset, np1, np1, "pair:offset");
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
global settings
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::settings(int narg, char **arg)
|
|
{
|
|
if (narg < 1 || narg > 2) error->all(FLERR, "Illegal pair_style command");
|
|
|
|
cut_lj_global = utils::numeric(FLERR, arg[0], false, lmp);
|
|
if (narg == 1)
|
|
cut_coul_global = cut_lj_global;
|
|
else
|
|
cut_coul_global = utils::numeric(FLERR, arg[1], 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;
|
|
cut_coul[i][j] = cut_coul_global;
|
|
}
|
|
}
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
set coeffs for one or more type pairs
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::coeff(int narg, char **arg)
|
|
{
|
|
if (narg < 5 || narg > 7) 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 a_one = utils::numeric(FLERR, arg[2], false, lmp);
|
|
double rho_one = utils::numeric(FLERR, arg[3], false, lmp);
|
|
if (rho_one <= 0) error->all(FLERR, "Incorrect args for pair coefficients");
|
|
double c_one = utils::numeric(FLERR, arg[4], false, lmp);
|
|
|
|
double cut_lj_one = cut_lj_global;
|
|
double cut_coul_one = cut_coul_global;
|
|
if (narg >= 6) cut_coul_one = cut_lj_one = utils::numeric(FLERR, arg[5], false, lmp);
|
|
if (narg == 7) cut_coul_one = utils::numeric(FLERR, arg[6], false, lmp);
|
|
|
|
int count = 0;
|
|
for (int i = ilo; i <= ihi; i++) {
|
|
for (int j = MAX(jlo, i); j <= jhi; j++) {
|
|
a[i][j] = a_one;
|
|
rho[i][j] = rho_one;
|
|
c[i][j] = c_one;
|
|
cut_lj[i][j] = cut_lj_one;
|
|
cut_coul[i][j] = cut_coul_one;
|
|
setflag[i][j] = 1;
|
|
count++;
|
|
}
|
|
}
|
|
|
|
if (count == 0) error->all(FLERR, "Incorrect args for pair coefficients");
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
init specific to this pair style
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::init_style()
|
|
{
|
|
if (!atom->q_flag) error->all(FLERR, "Pair style buck/coul/cut requires atom attribute q");
|
|
|
|
neighbor->add_request(this);
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
init for one type pair i,j and corresponding j,i
|
|
------------------------------------------------------------------------- */
|
|
|
|
double PairBuckCoulCut::init_one(int i, int j)
|
|
{
|
|
if (setflag[i][j] == 0) error->all(FLERR, "All pair coeffs are not set");
|
|
|
|
double cut = MAX(cut_lj[i][j], cut_coul[i][j]);
|
|
cut_ljsq[i][j] = cut_lj[i][j] * cut_lj[i][j];
|
|
cut_coulsq[i][j] = cut_coul[i][j] * cut_coul[i][j];
|
|
|
|
rhoinv[i][j] = 1.0 / rho[i][j];
|
|
buck1[i][j] = a[i][j] / rho[i][j];
|
|
buck2[i][j] = 6.0 * c[i][j];
|
|
|
|
if (offset_flag && (cut_lj[i][j] > 0.0)) {
|
|
double rexp = exp(-cut_lj[i][j] / rho[i][j]);
|
|
offset[i][j] = a[i][j] * rexp - c[i][j] / pow(cut_lj[i][j], 6.0);
|
|
} else
|
|
offset[i][j] = 0.0;
|
|
|
|
cut_ljsq[j][i] = cut_ljsq[i][j];
|
|
cut_coulsq[j][i] = cut_coulsq[i][j];
|
|
a[j][i] = a[i][j];
|
|
c[j][i] = c[i][j];
|
|
rhoinv[j][i] = rhoinv[i][j];
|
|
buck1[j][i] = buck1[i][j];
|
|
buck2[j][i] = buck2[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 rho1 = rho[i][j];
|
|
double rho2 = rho1 * rho1;
|
|
double rho3 = rho2 * rho1;
|
|
double rc = cut_lj[i][j];
|
|
double rc2 = rc * rc;
|
|
double rc3 = rc2 * rc;
|
|
etail_ij = 2.0 * MY_PI * all[0] * all[1] *
|
|
(a[i][j] * exp(-rc / rho1) * rho1 * (rc2 + 2.0 * rho1 * rc + 2.0 * rho2) -
|
|
c[i][j] / (3.0 * rc3));
|
|
ptail_ij = (-1 / 3.0) * 2.0 * MY_PI * all[0] * all[1] *
|
|
(-a[i][j] * exp(-rc / rho1) * (rc3 + 3.0 * rho1 * rc2 + 6.0 * rho2 * rc + 6.0 * rho3) +
|
|
2.0 * c[i][j] / rc3);
|
|
}
|
|
|
|
return cut;
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
proc 0 writes to restart file
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::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(&rho[i][j], sizeof(double), 1, fp);
|
|
fwrite(&c[i][j], sizeof(double), 1, fp);
|
|
fwrite(&cut_lj[i][j], sizeof(double), 1, fp);
|
|
fwrite(&cut_coul[i][j], sizeof(double), 1, fp);
|
|
}
|
|
}
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
proc 0 reads from restart file, bcasts
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::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, &a[i][j], sizeof(double), 1, fp, nullptr, error);
|
|
utils::sfread(FLERR, &rho[i][j], sizeof(double), 1, fp, nullptr, error);
|
|
utils::sfread(FLERR, &c[i][j], sizeof(double), 1, fp, nullptr, error);
|
|
utils::sfread(FLERR, &cut_lj[i][j], sizeof(double), 1, fp, nullptr, error);
|
|
utils::sfread(FLERR, &cut_coul[i][j], sizeof(double), 1, fp, nullptr, error);
|
|
}
|
|
MPI_Bcast(&a[i][j], 1, MPI_DOUBLE, 0, world);
|
|
MPI_Bcast(&rho[i][j], 1, MPI_DOUBLE, 0, world);
|
|
MPI_Bcast(&c[i][j], 1, MPI_DOUBLE, 0, world);
|
|
MPI_Bcast(&cut_lj[i][j], 1, MPI_DOUBLE, 0, world);
|
|
MPI_Bcast(&cut_coul[i][j], 1, MPI_DOUBLE, 0, world);
|
|
}
|
|
}
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
proc 0 writes to restart file
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::write_restart_settings(FILE *fp)
|
|
{
|
|
fwrite(&cut_lj_global, sizeof(double), 1, fp);
|
|
fwrite(&cut_coul_global, 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 PairBuckCoulCut::read_restart_settings(FILE *fp)
|
|
{
|
|
if (comm->me == 0) {
|
|
utils::sfread(FLERR, &cut_lj_global, sizeof(double), 1, fp, nullptr, error);
|
|
utils::sfread(FLERR, &cut_coul_global, 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(&cut_lj_global, 1, MPI_DOUBLE, 0, world);
|
|
MPI_Bcast(&cut_coul_global, 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 PairBuckCoulCut::write_data(FILE *fp)
|
|
{
|
|
for (int i = 1; i <= atom->ntypes; i++)
|
|
fprintf(fp, "%d %g %g %g\n", i, a[i][i], rho[i][i], c[i][i]);
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
proc 0 writes all pairs to data file
|
|
------------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::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 %g\n", i, j, a[i][j], rho[i][j], c[i][j], cut_lj[i][j],
|
|
cut_coul[i][j]);
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
double PairBuckCoulCut::single(int i, int j, int itype, int jtype, double rsq, double factor_coul,
|
|
double factor_lj, double &fforce)
|
|
{
|
|
double r2inv, r6inv, r, rexp, forcecoul, forcebuck, phicoul, phibuck;
|
|
|
|
r2inv = 1.0 / rsq;
|
|
if (rsq < cut_coulsq[itype][jtype])
|
|
forcecoul = force->qqrd2e * atom->q[i] * atom->q[j] * sqrt(r2inv);
|
|
else
|
|
forcecoul = 0.0;
|
|
if (rsq < cut_ljsq[itype][jtype]) {
|
|
r6inv = r2inv * r2inv * r2inv;
|
|
r = sqrt(rsq);
|
|
rexp = exp(-r * rhoinv[itype][jtype]);
|
|
forcebuck = buck1[itype][jtype] * r * rexp - buck2[itype][jtype] * r6inv;
|
|
} else
|
|
forcebuck = 0.0;
|
|
fforce = (factor_coul * forcecoul + factor_lj * forcebuck) * r2inv;
|
|
|
|
double eng = 0.0;
|
|
if (rsq < cut_coulsq[itype][jtype]) {
|
|
phicoul = force->qqrd2e * atom->q[i] * atom->q[j] * sqrt(r2inv);
|
|
eng += factor_coul * phicoul;
|
|
}
|
|
if (rsq < cut_ljsq[itype][jtype]) {
|
|
phibuck = a[itype][jtype] * rexp - c[itype][jtype] * r6inv - offset[itype][jtype];
|
|
eng += factor_lj * phibuck;
|
|
}
|
|
return eng;
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
void PairBuckCoulCut::born_matrix(int i, int j, int itype, int jtype, double rsq,
|
|
double factor_coul, double factor_lj, double &dupair,
|
|
double &du2pair)
|
|
{
|
|
double rinv, r2inv, r3inv, r6inv, r7inv, r8inv, r, rexp;
|
|
double du_lj, du2_lj, du_coul, du2_coul;
|
|
|
|
double *q = atom->q;
|
|
double qqrd2e = force->qqrd2e;
|
|
|
|
r = sqrt(rsq);
|
|
rexp = exp(-r * rhoinv[itype][jtype]);
|
|
|
|
r2inv = 1.0 / rsq;
|
|
rinv = sqrt(r2inv);
|
|
r3inv = r2inv * rinv;
|
|
r6inv = r2inv * r2inv * r2inv;
|
|
r7inv = r6inv * rinv;
|
|
r8inv = r6inv * r2inv;
|
|
|
|
// Reminder: buck1[itype][jtype] = a[itype][jtype]/rho[itype][jtype];
|
|
// Reminder: buck2[itype][jtype] = 6.0*c[itype][jtype];
|
|
|
|
du_lj = buck2[itype][jtype] * r7inv - buck1[itype][jtype] * rexp;
|
|
du2_lj = (buck1[itype][jtype] / rho[itype][jtype]) * rexp - 7 * buck2[itype][jtype] * r8inv;
|
|
|
|
// Reminder: qqrd2e converts q^2/r to energy w/ dielectric constant
|
|
|
|
du_coul = -qqrd2e * q[i] * q[j] * r2inv;
|
|
du2_coul = 2.0 * qqrd2e * q[i] * q[j] * r3inv;
|
|
|
|
dupair = factor_lj * du_lj + factor_coul * du_coul;
|
|
du2pair = factor_lj * du2_lj + factor_coul * du2_coul;
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
void *PairBuckCoulCut::extract(const char *str, int &dim)
|
|
{
|
|
dim = 2;
|
|
if (strcmp(str, "a") == 0) return (void *) a;
|
|
if (strcmp(str, "c") == 0) return (void *) c;
|
|
return nullptr;
|
|
}
|