Page MenuHomec4science

pair_lj_cut_tip4p_cut.cpp
No OneTemporary

File Metadata

Created
Sat, Jul 6, 04:17

pair_lj_cut_tip4p_cut.cpp

/* ----------------------------------------------------------------------
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: Pavel Elkind (Gothenburg University)
------------------------------------------------------------------------- */
#include "math.h"
#include "stdlib.h"
#include "string.h"
#include "pair_lj_cut_tip4p_cut.h"
#include "atom.h"
#include "force.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "domain.h"
#include "angle.h"
#include "bond.h"
#include "comm.h"
#include "math_const.h"
#include "memory.h"
#include "error.h"
using namespace LAMMPS_NS;
using namespace MathConst;
/* ---------------------------------------------------------------------- */
PairLJCutTIP4PCut::PairLJCutTIP4PCut(LAMMPS *lmp) : Pair(lmp)
{
single_enable = 0;
writedata = 1;
nmax = 0;
hneigh = NULL;
newsite = NULL;
// TIP4P cannot compute virial as F dot r
// due to finding bonded H atoms which are not near O atom
no_virial_fdotr_compute = 1;
}
/* ---------------------------------------------------------------------- */
PairLJCutTIP4PCut::~PairLJCutTIP4PCut()
{
if (allocated) {
memory->destroy(setflag);
memory->destroy(cutsq);
memory->destroy(cut_lj);
memory->destroy(cut_ljsq);
memory->destroy(epsilon);
memory->destroy(sigma);
memory->destroy(lj1);
memory->destroy(lj2);
memory->destroy(lj3);
memory->destroy(lj4);
memory->destroy(offset);
}
memory->destroy(hneigh);
memory->destroy(newsite);
}
/* ---------------------------------------------------------------------- */
void PairLJCutTIP4PCut::compute(int eflag, int vflag)
{
int i,j,ii,jj,inum,jnum,itype,jtype;
double qtmp,xtmp,ytmp,ztmp,delx,dely,delz,evdwl,ecoul;
double rsq,r2inv,r6inv,forcecoul,forcelj,factor_lj,factor_coul;
int *ilist,*jlist,*numneigh,**firstneigh;
int key;
int n,vlist[6];
int iH1,iH2,jH1,jH2;
double cforce;
double fO[3],fH[3],fd[3],v[6],xH1[3],xH2[3];
double *x1,*x2;
evdwl = ecoul = 0.0;
if (eflag || vflag) ev_setup(eflag,vflag);
else evflag = vflag_fdotr = 0;
// reallocate hneigh & newsite if necessary
// initialize hneigh[0] to -1 on steps when reneighboring occurred
// initialize hneigh[2] to 0 every step
int nlocal = atom->nlocal;
int nall = nlocal + atom->nghost;
if (atom->nmax > nmax) {
nmax = atom->nmax;
memory->destroy(hneigh);
memory->create(hneigh,nmax,3,"pair:hneigh");
memory->destroy(newsite);
memory->create(newsite,nmax,3,"pair:newsite");
}
if (neighbor->ago == 0)
for (i = 0; i < nall; i++) hneigh[i][0] = -1;
for (i = 0; i < nall; i++) hneigh[i][2] = 0;
double **f = atom->f;
double **x = atom->x;
double *q = atom->q;
tagint *tag = atom->tag;
int *type = atom->type;
double *special_lj = force->special_lj;
double *special_coul = force->special_coul;
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];
if (itype == typeO) {
if (hneigh[i][0] < 0) {
hneigh[i][0] = iH1 = atom->map(tag[i] + 1);
hneigh[i][1] = iH2 = atom->map(tag[i] + 2);
hneigh[i][2] = 1;
if (iH1 == -1 || iH2 == -1)
error->one(FLERR,"TIP4P hydrogen is missing");
if (atom->type[iH1] != typeH || atom->type[iH2] != typeH)
error->one(FLERR,"TIP4P hydrogen has incorrect atom type");
compute_newsite(x[i],x[iH1],x[iH2],newsite[i]);
} else {
iH1 = hneigh[i][0];
iH2 = hneigh[i][1];
if (hneigh[i][2] == 0) {
hneigh[i][2] = 1;
compute_newsite(x[i],x[iH1],x[iH2],newsite[i]);
}
}
x1 = newsite[i];
} else x1 = x[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];
// LJ interaction based on true rsq
if (rsq < cut_ljsq[itype][jtype]) {
r2inv = 1.0/rsq;
r6inv = r2inv*r2inv*r2inv;
forcelj = r6inv * (lj1[itype][jtype]*r6inv - lj2[itype][jtype]);
forcelj *= factor_lj * r2inv;
f[i][0] += delx*forcelj;
f[i][1] += dely*forcelj;
f[i][2] += delz*forcelj;
f[j][0] -= delx*forcelj;
f[j][1] -= dely*forcelj;
f[j][2] -= delz*forcelj;
if (eflag) {
evdwl = r6inv*(lj3[itype][jtype]*r6inv-lj4[itype][jtype]) -
offset[itype][jtype];
evdwl *= factor_lj;
} else evdwl = 0.0;
if (evflag) ev_tally(i,j,nlocal,newton_pair,
evdwl,0.0,forcelj,delx,dely,delz);
}
// adjust rsq and delxyz for off-site O charge(s) if necessary
// but only if they are within reach
if (rsq < cut_coulsqplus) {
if (itype == typeO || jtype == typeO) {
// if atom J = water O, set x2 = offset charge site
// else x2 = x of atom J
if (jtype == typeO) {
if (hneigh[j][0] < 0) {
hneigh[j][0] = jH1 = atom->map(tag[j] + 1);
hneigh[j][1] = jH2 = atom->map(tag[j] + 2);
hneigh[j][2] = 1;
if (jH1 == -1 || jH2 == -1)
error->one(FLERR,"TIP4P hydrogen is missing");
if (atom->type[jH1] != typeH || atom->type[jH2] != typeH)
error->one(FLERR,"TIP4P hydrogen has incorrect atom type");
compute_newsite(x[j],x[jH1],x[jH2],newsite[j]);
} else {
jH1 = hneigh[j][0];
jH2 = hneigh[j][1];
if (hneigh[j][2] == 0) {
hneigh[j][2] = 1;
compute_newsite(x[j],x[jH1],x[jH2],newsite[j]);
}
}
x2 = newsite[j];
} else x2 = x[j];
delx = x1[0] - x2[0];
dely = x1[1] - x2[1];
delz = x1[2] - x2[2];
rsq = delx*delx + dely*dely + delz*delz;
}
// Coulombic interaction based on modified rsq
if (rsq < cut_coulsq) {
r2inv = 1.0 / rsq;
forcecoul = qqrd2e * qtmp * q[j] * sqrt(r2inv);
cforce = factor_coul * forcecoul * r2inv;
// if i,j are not O atoms, force is applied directly;
// if i or j are O atoms, force is on fictitious atom & partitioned
// force partitioning due to Feenstra, J Comp Chem, 20, 786 (1999)
// f_f = fictitious force, fO = f_f (1 - 2 alpha), fH = alpha f_f
// preserves total force and torque on water molecule
// virial = sum(r x F) where each water's atoms are near xi and xj
// vlist stores 2,4,6 atoms whose forces contribute to virial
n = 0;
key = 0;
if (itype != typeO) {
f[i][0] += delx * cforce;
f[i][1] += dely * cforce;
f[i][2] += delz * cforce;
if (vflag) {
v[0] = x[i][0] * delx * cforce;
v[1] = x[i][1] * dely * cforce;
v[2] = x[i][2] * delz * cforce;
v[3] = x[i][0] * dely * cforce;
v[4] = x[i][0] * delz * cforce;
v[5] = x[i][1] * delz * cforce;
}
vlist[n++] = i;
} else {
key++;
fd[0] = delx*cforce;
fd[1] = dely*cforce;
fd[2] = delz*cforce;
fO[0] = fd[0]*(1.0 - alpha);
fO[1] = fd[1]*(1.0 - alpha);
fO[2] = fd[2]*(1.0 - alpha);
fH[0] = 0.5 * alpha * fd[0];
fH[1] = 0.5 * alpha * fd[1];
fH[2] = 0.5 * alpha * fd[2];
f[i][0] += fO[0];
f[i][1] += fO[1];
f[i][2] += fO[2];
f[iH1][0] += fH[0];
f[iH1][1] += fH[1];
f[iH1][2] += fH[2];
f[iH2][0] += fH[0];
f[iH2][1] += fH[1];
f[iH2][2] += fH[2];
if(vflag) {
domain->closest_image(x[i],x[iH1],xH1);
domain->closest_image(x[i],x[iH2],xH2);
v[0] = x[i][0]*fO[0] + xH1[0]*fH[0] + xH2[0]*fH[0];
v[1] = x[i][1]*fO[1] + xH1[1]*fH[1] + xH2[1]*fH[1];
v[2] = x[i][2]*fO[2] + xH1[2]*fH[2] + xH2[2]*fH[2];
v[3] = x[i][0]*fO[1] + xH1[0]*fH[1] + xH2[0]*fH[1];
v[4] = x[i][0]*fO[2] + xH1[0]*fH[2] + xH2[0]*fH[2];
v[5] = x[i][1]*fO[2] + xH1[1]*fH[2] + xH2[1]*fH[2];
}
vlist[n++] = i;
vlist[n++] = iH1;
vlist[n++] = iH2;
}
if (jtype != typeO) {
f[j][0] -= delx * cforce;
f[j][1] -= dely * cforce;
f[j][2] -= delz * cforce;
if (vflag) {
v[0] -= x[j][0] * delx * cforce;
v[1] -= x[j][1] * dely * cforce;
v[2] -= x[j][2] * delz * cforce;
v[3] -= x[j][0] * dely * cforce;
v[4] -= x[j][0] * delz * cforce;
v[5] -= x[j][1] * delz * cforce;
}
vlist[n++] = j;
} else {
key += 2;
fd[0] = -delx*cforce;
fd[1] = -dely*cforce;
fd[2] = -delz*cforce;
fO[0] = fd[0]*(1 - alpha);
fO[1] = fd[1]*(1 - alpha);
fO[2] = fd[2]*(1 - alpha);
fH[0] = 0.5 * alpha * fd[0];
fH[1] = 0.5 * alpha * fd[1];
fH[2] = 0.5 * alpha * fd[2];
f[j][0] += fO[0];
f[j][1] += fO[1];
f[j][2] += fO[2];
f[jH1][0] += fH[0];
f[jH1][1] += fH[1];
f[jH1][2] += fH[2];
f[jH2][0] += fH[0];
f[jH2][1] += fH[1];
f[jH2][2] += fH[2];
if (vflag) {
domain->closest_image(x[j],x[jH1],xH1);
domain->closest_image(x[j],x[jH2],xH2);
v[0] += x[j][0]*fO[0] + xH1[0]*fH[0] + xH2[0]*fH[0];
v[1] += x[j][1]*fO[1] + xH1[1]*fH[1] + xH2[1]*fH[1];
v[2] += x[j][2]*fO[2] + xH1[2]*fH[2] + xH2[2]*fH[2];
v[3] += x[j][0]*fO[1] + xH1[0]*fH[1] + xH2[0]*fH[1];
v[4] += x[j][0]*fO[2] + xH1[0]*fH[2] + xH2[0]*fH[2];
v[5] += x[j][1]*fO[2] + xH1[1]*fH[2] + xH2[1]*fH[2];
}
vlist[n++] = j;
vlist[n++] = jH1;
vlist[n++] = jH2;
}
if (eflag) {
ecoul = qqrd2e * qtmp * q[j] * sqrt(r2inv);
ecoul *= factor_coul;
} else ecoul = 0.0;
if (evflag) ev_tally_tip4p(key,vlist,v,ecoul,alpha);
}
}
}
}
}
/* ----------------------------------------------------------------------
allocate all arrays
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::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(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 PairLJCutTIP4PCut::settings(int narg, char **arg)
{
if (narg < 6 || narg > 7) error->all(FLERR,"Illegal pair_style command");
typeO = force->inumeric(FLERR,arg[0]);
typeH = force->inumeric(FLERR,arg[1]);
typeB = force->inumeric(FLERR,arg[2]);
typeA = force->inumeric(FLERR,arg[3]);
qdist = force->numeric(FLERR,arg[4]);
cut_lj_global = force->numeric(FLERR,arg[5]);
if (narg == 6) cut_coul = cut_lj_global;
else cut_coul = force->numeric(FLERR,arg[6]);
cut_coulsq = cut_coul * cut_coul;
cut_coulsqplus = (cut_coul + 2.0*qdist) * (cut_coul + 2.0*qdist);
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_lj[i][j] = cut_lj_global;
}
}
/* ----------------------------------------------------------------------
set coeffs for one or more type pairs
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::coeff(int narg, char **arg)
{
if (narg < 4 || narg > 5)
error->all(FLERR,"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 = force->numeric(FLERR,arg[2]);
double sigma_one = force->numeric(FLERR,arg[3]);
double cut_lj_one = cut_lj_global;
if (narg == 5) cut_lj_one = force->numeric(FLERR,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_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 PairLJCutTIP4PCut::init_style()
{
if (atom->tag_enable == 0)
error->all(FLERR,"Pair style lj/cut/tip4p/cut requires atom IDs");
if (!force->newton_pair)
error->all(FLERR,
"Pair style lj/cut/tip4p/cut requires newton pair on");
if (!atom->q_flag)
error->all(FLERR,
"Pair style lj/cut/tip4p/cut requires atom attribute q");
if (force->bond == NULL)
error->all(FLERR,"Must use a bond style with TIP4P potential");
if (force->angle == NULL)
error->all(FLERR,"Must use an angle style with TIP4P potential");
neighbor->request(this,instance_me);
// set alpha parameter
double theta = force->angle->equilibrium_angle(typeA);
double blen = force->bond->equilibrium_distance(typeB);
alpha = qdist / (cos(0.5*theta) * blen);
}
/* ----------------------------------------------------------------------
init for one type pair i,j and corresponding j,i
------------------------------------------------------------------------- */
double PairLJCutTIP4PCut::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_lj[i][j] = mix_distance(cut_lj[i][i],cut_lj[j][j]);
}
// include TIP4P qdist in full cutoff, qdist = 0.0 if not TIP4P
double cut = MAX(cut_lj[i][j],cut_coul+2.0*qdist);
cut_ljsq[i][j] = cut_lj[i][j] * cut_lj[i][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_lj[i][j];
offset[i][j] = 4.0 * epsilon[i][j] * (pow(ratio,12.0) - pow(ratio,6.0));
} else offset[i][j] = 0.0;
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 sig2 = sigma[i][j]*sigma[i][j];
double sig6 = sig2*sig2*sig2;
double rc3 = cut_lj[i][j]*cut_lj[i][j]*cut_lj[i][j];
double rc6 = rc3*rc3;
double rc9 = rc3*rc6;
etail_ij = 8.0*MY_PI*all[0]*all[1]*epsilon[i][j] *
sig6 * (sig6 - 3.0*rc6) / (9.0*rc9);
ptail_ij = 16.0*MY_PI*all[0]*all[1]*epsilon[i][j] *
sig6 * (2.0*sig6 - 3.0*rc6) / (9.0*rc9);
}
// check that LJ epsilon = 0.0 for water H
// set LJ cutoff to 0.0 for any interaction involving water H
// so LJ term isn't calculated in compute()
if ((i == typeH && epsilon[i][i] != 0.0) ||
(j == typeH && epsilon[j][j] != 0.0))
error->all(FLERR,"Water H epsilon must be 0.0 for "
"pair style lj/cut/tip4p/cut");
if (i == typeH || j == typeH)
cut_ljsq[j][i] = cut_ljsq[i][j] = 0.0;
return cut;
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::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_lj[i][j],sizeof(double),1,fp);
}
}
}
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::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_lj[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_lj[i][j],1,MPI_DOUBLE,0,world);
}
}
}
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::write_restart_settings(FILE *fp)
{
fwrite(&typeO,sizeof(int),1,fp);
fwrite(&typeH,sizeof(int),1,fp);
fwrite(&typeB,sizeof(int),1,fp);
fwrite(&typeA,sizeof(int),1,fp);
fwrite(&qdist,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 PairLJCutTIP4PCut::read_restart_settings(FILE *fp)
{
if (comm->me == 0) {
fread(&typeO,sizeof(int),1,fp);
fread(&typeH,sizeof(int),1,fp);
fread(&typeB,sizeof(int),1,fp);
fread(&typeA,sizeof(int),1,fp);
fread(&qdist,sizeof(double),1,fp);
fread(&cut_lj_global,sizeof(double),1,fp);
fread(&cut_coul,sizeof(double),1,fp);
fread(&offset_flag,sizeof(int),1,fp);
fread(&mix_flag,sizeof(int),1,fp);
fread(&tail_flag,sizeof(int),1,fp);
}
MPI_Bcast(&typeO,1,MPI_INT,0,world);
MPI_Bcast(&typeH,1,MPI_INT,0,world);
MPI_Bcast(&typeB,1,MPI_INT,0,world);
MPI_Bcast(&typeA,1,MPI_INT,0,world);
MPI_Bcast(&qdist,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);
cut_coulsq = cut_coul * cut_coul;
cut_coulsqplus = (cut_coul + 2.0*qdist) * (cut_coul + 2.0*qdist);
}
/* ----------------------------------------------------------------------
proc 0 writes to data file
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::write_data(FILE *fp)
{
for (int i = 1; i <= atom->ntypes; i++)
fprintf(fp,"%d %g %g\n",i,epsilon[i][i],sigma[i][i]);
}
/* ----------------------------------------------------------------------
proc 0 writes all pairs to data file
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::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\n",i,j,epsilon[i][j],sigma[i][j],cut_lj[i][j]);
}
/* ----------------------------------------------------------------------
compute position xM of fictitious charge site for O atom and 2 H atoms
return it as xM
------------------------------------------------------------------------- */
void PairLJCutTIP4PCut::compute_newsite(double *xO, double *xH1,
double *xH2, double *xM)
{
double delx1 = xH1[0] - xO[0];
double dely1 = xH1[1] - xO[1];
double delz1 = xH1[2] - xO[2];
domain->minimum_image(delx1,dely1,delz1);
double delx2 = xH2[0] - xO[0];
double dely2 = xH2[1] - xO[1];
double delz2 = xH2[2] - xO[2];
domain->minimum_image(delx2,dely2,delz2);
xM[0] = xO[0] + alpha * 0.5 * (delx1 + delx2);
xM[1] = xO[1] + alpha * 0.5 * (dely1 + dely2);
xM[2] = xO[2] + alpha * 0.5 * (delz1 + delz2);
}
/* ---------------------------------------------------------------------- */
void *PairLJCutTIP4PCut::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;
return NULL;
}
/* ----------------------------------------------------------------------
memory usage of hneigh
------------------------------------------------------------------------- */
double PairLJCutTIP4PCut::memory_usage()
{
double bytes = maxeatom * sizeof(double);
bytes += maxvatom*6 * sizeof(double);
bytes += 2 * nmax * sizeof(double);
return bytes;
}

Event Timeline