661 lines
20 KiB
C++
661 lines
20 KiB
C++
/* ----------------------------------------------------------------------
|
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
|
https://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: K. Michael Salerno (NRL)
|
|
Based on tabulated dihedral (dihedral_table.cpp) by Andrew Jewett
|
|
------------------------------------------------------------------------- */
|
|
|
|
#include "dihedral_table_cut.h"
|
|
|
|
#include "atom.h"
|
|
#include "citeme.h"
|
|
#include "comm.h"
|
|
#include "error.h"
|
|
#include "force.h"
|
|
#include "math_const.h"
|
|
#include "memory.h"
|
|
#include "neighbor.h"
|
|
#include "update.h"
|
|
|
|
#include <cmath>
|
|
#include <fstream> // IWYU pragma: keep
|
|
#include <sstream> // IWYU pragma: keep
|
|
|
|
using namespace LAMMPS_NS;
|
|
using namespace MathConst;
|
|
|
|
static const char cite_dihedral_tablecut[] =
|
|
"dihedral_style table/cut command:\n\n"
|
|
"@Article{Salerno17,\n"
|
|
" author = {K. M. Salerno and N. Bernstein},\n"
|
|
" title = {Persistence Length, End-to-End Distance, and Structure of Coarse-Grained Polymers},\n"
|
|
" journal = {J.~Chem.~Theory Comput.},\n"
|
|
" year = 2018,\n"
|
|
" DOI = 10.1021/acs.jctc.7b01229"
|
|
"}\n\n";
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
#define TOLERANCE 0.05
|
|
#define SMALL 0.0000001
|
|
|
|
// ------------------------------------------------------------------------
|
|
// The following auxiliary functions were left out of the
|
|
// DihedralTable class either because they require template parameters,
|
|
// or because they have nothing to do with dihedral angles.
|
|
// ------------------------------------------------------------------------
|
|
|
|
// -------------------------------------------------------------------
|
|
// --------- The function was taken verbatim from the ---------
|
|
// --------- GNU Scientific Library (GSL, version 1.15) ---------
|
|
// -------------------------------------------------------------------
|
|
|
|
/* Author: Gerard Jungman */
|
|
/* for description of method see [Engeln-Mullges + Uhlig, p. 96]
|
|
*
|
|
* diag[0] offdiag[0] 0 ..... offdiag[N-1]
|
|
* offdiag[0] diag[1] offdiag[1] .....
|
|
* 0 offdiag[1] diag[2]
|
|
* 0 0 offdiag[2] .....
|
|
* ... ...
|
|
* offdiag[N-1] ...
|
|
*
|
|
*/
|
|
// -- (A non-symmetric version of this function is also available.) --
|
|
|
|
enum { //GSL status return codes.
|
|
GSL_FAILURE = -1,
|
|
GSL_SUCCESS = 0,
|
|
GSL_ENOMEM = 8,
|
|
GSL_EZERODIV = 12,
|
|
GSL_EBADLEN = 19
|
|
};
|
|
|
|
|
|
|
|
// cyc_splintD(): Evaluate the deriviative of a cyclic spline at position x,
|
|
// with n control points at xa[], ya[], with parameters y2a[].
|
|
// The xa[] must be monotonically increasing and their
|
|
// range should not exceed period (ie xa[n-1] < xa[0] + period).
|
|
// x must lie in the range: [(xa[n-1]-period), (xa[0]+period)]
|
|
// "period" is typically 2*PI.
|
|
static double cyc_splintD(double const *xa,
|
|
double const *ya,
|
|
double const *y2a,
|
|
int n,
|
|
double period,
|
|
double x)
|
|
{
|
|
int klo = -1;
|
|
int khi = n; // (not n-1)
|
|
int k;
|
|
double xlo = xa[n-1] - period;
|
|
double xhi = xa[0] + period;
|
|
while (khi-klo > 1) {
|
|
k = (khi+klo) >> 1; //(k=(khi+klo)/2)
|
|
if (xa[k] > x) {
|
|
khi = k;
|
|
xhi = xa[k];
|
|
}
|
|
else {
|
|
klo = k;
|
|
xlo = xa[k];
|
|
}
|
|
}
|
|
|
|
if (khi == n) khi = 0;
|
|
if (klo ==-1) klo = n-1;
|
|
|
|
double yhi = ya[khi];
|
|
double ylo = ya[klo];
|
|
double h = xhi-xlo;
|
|
double g = yhi-ylo;
|
|
double a = (xhi-x) / h;
|
|
double b = (x-xlo) / h;
|
|
// Formula below taken from equation 3.3.5 of "numerical recipes in c"
|
|
// "yD" = the derivative of y
|
|
double yD = g/h - ( (3.0*a*a-1.0)*y2a[klo] - (3.0*b*b-1.0)*y2a[khi] ) * h/6.0;
|
|
// For rerefence: y = a*ylo + b*yhi +
|
|
// ((a*a*a-a)*y2a[klo] + (b*b*b-b)*y2a[khi]) * (h*h)/6.0;
|
|
|
|
return yD;
|
|
|
|
} // cyc_splintD()
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
DihedralTableCut::DihedralTableCut(LAMMPS *lmp) : DihedralTable(lmp)
|
|
{
|
|
if (lmp->citeme) lmp->citeme->add(cite_dihedral_tablecut);
|
|
aat_k = aat_theta0_1 = aat_theta0_2 = nullptr;
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
DihedralTableCut::~DihedralTableCut()
|
|
{
|
|
memory->destroy(aat_k);
|
|
memory->destroy(aat_theta0_1);
|
|
memory->destroy(aat_theta0_2);
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
void DihedralTableCut::compute(int eflag, int vflag)
|
|
{
|
|
|
|
int i1,i2,i3,i4,i,j,k,n,type;
|
|
double edihedral;
|
|
double vb1x,vb1y,vb1z,vb2x,vb2y,vb2z,vb3x,vb3y,vb3z,vb2xm,vb2ym,vb2zm;
|
|
double fphi,fpphi;
|
|
double r1mag2,r1,r2mag2,r2,r3mag2,r3;
|
|
double sb1,rb1,sb2,rb2,sb3,rb3,c0,r12c1;
|
|
double r12c2,costh12,costh13,costh23,sc1,sc2,s1,s2,c;
|
|
double phi,sinphi,a11,a22,a33,a12,a13,a23,sx1,sx2;
|
|
double sx12,sy1,sy2,sy12,sz1,sz2,sz12;
|
|
double t1,t2,t3,t4;
|
|
double da1,da2;
|
|
double s12,sin2;
|
|
double dcosphidr[4][3],dphidr[4][3],dthetadr[2][4][3];
|
|
double fabcd[4][3];
|
|
|
|
edihedral = 0.0;
|
|
ev_init(eflag,vflag);
|
|
|
|
double **x = atom->x;
|
|
double **f = atom->f;
|
|
int **dihedrallist = neighbor->dihedrallist;
|
|
int ndihedrallist = neighbor->ndihedrallist;
|
|
int nlocal = atom->nlocal;
|
|
int newton_bond = force->newton_bond;
|
|
|
|
for (n = 0; n < ndihedrallist; n++) {
|
|
i1 = dihedrallist[n][0];
|
|
i2 = dihedrallist[n][1];
|
|
i3 = dihedrallist[n][2];
|
|
i4 = dihedrallist[n][3];
|
|
type = dihedrallist[n][4];
|
|
|
|
// 1st bond
|
|
|
|
vb1x = x[i1][0] - x[i2][0];
|
|
vb1y = x[i1][1] - x[i2][1];
|
|
vb1z = x[i1][2] - x[i2][2];
|
|
|
|
// 2nd bond
|
|
|
|
vb2x = x[i3][0] - x[i2][0];
|
|
vb2y = x[i3][1] - x[i2][1];
|
|
vb2z = x[i3][2] - x[i2][2];
|
|
|
|
vb2xm = -vb2x;
|
|
vb2ym = -vb2y;
|
|
vb2zm = -vb2z;
|
|
|
|
// 3rd bond
|
|
|
|
vb3x = x[i4][0] - x[i3][0];
|
|
vb3y = x[i4][1] - x[i3][1];
|
|
vb3z = x[i4][2] - x[i3][2];
|
|
|
|
// distances
|
|
|
|
r1mag2 = vb1x*vb1x + vb1y*vb1y + vb1z*vb1z;
|
|
r1 = sqrt(r1mag2);
|
|
r2mag2 = vb2x*vb2x + vb2y*vb2y + vb2z*vb2z;
|
|
r2 = sqrt(r2mag2);
|
|
r3mag2 = vb3x*vb3x + vb3y*vb3y + vb3z*vb3z;
|
|
r3 = sqrt(r3mag2);
|
|
|
|
sb1 = 1.0/r1mag2;
|
|
rb1 = 1.0/r1;
|
|
sb2 = 1.0/r2mag2;
|
|
rb2 = 1.0/r2;
|
|
sb3 = 1.0/r3mag2;
|
|
rb3 = 1.0/r3;
|
|
|
|
c0 = (vb1x*vb3x + vb1y*vb3y + vb1z*vb3z) * rb1*rb3;
|
|
|
|
// angles
|
|
|
|
r12c1 = rb1*rb2;
|
|
r12c2 = rb2*rb3;
|
|
costh12 = (vb1x*vb2x + vb1y*vb2y + vb1z*vb2z) * r12c1;
|
|
costh13 = c0;
|
|
costh23 = (vb2xm*vb3x + vb2ym*vb3y + vb2zm*vb3z) * r12c2;
|
|
|
|
// cos and sin of 2 angles and final c
|
|
|
|
sin2 = MAX(1.0 - costh12*costh12,0.0);
|
|
sc1 = sqrt(sin2);
|
|
if (sc1 < SMALL) sc1 = SMALL;
|
|
sc1 = 1.0/sc1;
|
|
|
|
sin2 = MAX(1.0 - costh23*costh23,0.0);
|
|
sc2 = sqrt(sin2);
|
|
if (sc2 < SMALL) sc2 = SMALL;
|
|
sc2 = 1.0/sc2;
|
|
|
|
s1 = sc1 * sc1;
|
|
s2 = sc2 * sc2;
|
|
s12 = sc1 * sc2;
|
|
c = (c0 + costh12*costh23) * s12;
|
|
|
|
// error check
|
|
|
|
if (c > 1.0 + TOLERANCE || c < (-1.0 - TOLERANCE)) {
|
|
int me;
|
|
MPI_Comm_rank(world,&me);
|
|
if (screen) {
|
|
char str[128];
|
|
sprintf(str,"Dihedral problem: %d " BIGINT_FORMAT " "
|
|
TAGINT_FORMAT " " TAGINT_FORMAT " "
|
|
TAGINT_FORMAT " " TAGINT_FORMAT,
|
|
me,update->ntimestep,
|
|
atom->tag[i1],atom->tag[i2],atom->tag[i3],atom->tag[i4]);
|
|
error->warning(FLERR,str,0);
|
|
fprintf(screen," 1st atom: %d %g %g %g\n",
|
|
me,x[i1][0],x[i1][1],x[i1][2]);
|
|
fprintf(screen," 2nd atom: %d %g %g %g\n",
|
|
me,x[i2][0],x[i2][1],x[i2][2]);
|
|
fprintf(screen," 3rd atom: %d %g %g %g\n",
|
|
me,x[i3][0],x[i3][1],x[i3][2]);
|
|
fprintf(screen," 4th atom: %d %g %g %g\n",
|
|
me,x[i4][0],x[i4][1],x[i4][2]);
|
|
}
|
|
}
|
|
|
|
if (c > 1.0) c = 1.0;
|
|
if (c < -1.0) c = -1.0;
|
|
double phil = acos(c);
|
|
phi = acos(c);
|
|
|
|
sinphi = sqrt(1.0 - c*c);
|
|
sinphi = MAX(sinphi,SMALL);
|
|
|
|
// n123 = vb1 x vb2
|
|
|
|
double n123x = vb1y*vb2z - vb1z*vb2y;
|
|
double n123y = vb1z*vb2x - vb1x*vb2z;
|
|
double n123z = vb1x*vb2y - vb1y*vb2x;
|
|
double n123_dot_vb3 = n123x*vb3x + n123y*vb3y + n123z*vb3z;
|
|
if (n123_dot_vb3 > 0.0) {
|
|
phil = -phil;
|
|
phi = -phi;
|
|
sinphi = -sinphi;
|
|
}
|
|
|
|
a11 = -c*sb1*s1;
|
|
a22 = sb2 * (2.0*costh13*s12 - c*(s1+s2));
|
|
a33 = -c*sb3*s2;
|
|
a12 = r12c1 * (costh12*c*s1 + costh23*s12);
|
|
a13 = rb1*rb3*s12;
|
|
a23 = r12c2 * (-costh23*c*s2 - costh12*s12);
|
|
|
|
sx1 = a11*vb1x + a12*vb2x + a13*vb3x;
|
|
sx2 = a12*vb1x + a22*vb2x + a23*vb3x;
|
|
sx12 = a13*vb1x + a23*vb2x + a33*vb3x;
|
|
sy1 = a11*vb1y + a12*vb2y + a13*vb3y;
|
|
sy2 = a12*vb1y + a22*vb2y + a23*vb3y;
|
|
sy12 = a13*vb1y + a23*vb2y + a33*vb3y;
|
|
sz1 = a11*vb1z + a12*vb2z + a13*vb3z;
|
|
sz2 = a12*vb1z + a22*vb2z + a23*vb3z;
|
|
sz12 = a13*vb1z + a23*vb2z + a33*vb3z;
|
|
|
|
// set up d(cos(phi))/d(r) and dphi/dr arrays
|
|
|
|
dcosphidr[0][0] = -sx1;
|
|
dcosphidr[0][1] = -sy1;
|
|
dcosphidr[0][2] = -sz1;
|
|
dcosphidr[1][0] = sx2 + sx1;
|
|
dcosphidr[1][1] = sy2 + sy1;
|
|
dcosphidr[1][2] = sz2 + sz1;
|
|
dcosphidr[2][0] = sx12 - sx2;
|
|
dcosphidr[2][1] = sy12 - sy2;
|
|
dcosphidr[2][2] = sz12 - sz2;
|
|
dcosphidr[3][0] = -sx12;
|
|
dcosphidr[3][1] = -sy12;
|
|
dcosphidr[3][2] = -sz12;
|
|
|
|
for (i = 0; i < 4; i++)
|
|
for (j = 0; j < 3; j++)
|
|
dphidr[i][j] = -dcosphidr[i][j] / sinphi;
|
|
|
|
|
|
for (i = 0; i < 4; i++)
|
|
for (j = 0; j < 3; j++)
|
|
fabcd[i][j] = 0;
|
|
edihedral = 0;
|
|
|
|
|
|
// set up d(theta)/d(r) array
|
|
// dthetadr(i,j,k) = angle i, atom j, coordinate k
|
|
|
|
for (i = 0; i < 2; i++)
|
|
for (j = 0; j < 4; j++)
|
|
for (k = 0; k < 3; k++)
|
|
dthetadr[i][j][k] = 0.0;
|
|
|
|
t1 = costh12 / r1mag2;
|
|
t2 = costh23 / r2mag2;
|
|
t3 = costh12 / r2mag2;
|
|
t4 = costh23 / r3mag2;
|
|
|
|
// angle12
|
|
|
|
dthetadr[0][0][0] = sc1 * ((t1 * vb1x) - (vb2x * r12c1));
|
|
dthetadr[0][0][1] = sc1 * ((t1 * vb1y) - (vb2y * r12c1));
|
|
dthetadr[0][0][2] = sc1 * ((t1 * vb1z) - (vb2z * r12c1));
|
|
|
|
dthetadr[0][1][0] = sc1 * ((-t1 * vb1x) + (vb2x * r12c1) +
|
|
(-t3 * vb2x) + (vb1x * r12c1));
|
|
dthetadr[0][1][1] = sc1 * ((-t1 * vb1y) + (vb2y * r12c1) +
|
|
(-t3 * vb2y) + (vb1y * r12c1));
|
|
dthetadr[0][1][2] = sc1 * ((-t1 * vb1z) + (vb2z * r12c1) +
|
|
(-t3 * vb2z) + (vb1z * r12c1));
|
|
|
|
dthetadr[0][2][0] = sc1 * ((t3 * vb2x) - (vb1x * r12c1));
|
|
dthetadr[0][2][1] = sc1 * ((t3 * vb2y) - (vb1y * r12c1));
|
|
dthetadr[0][2][2] = sc1 * ((t3 * vb2z) - (vb1z * r12c1));
|
|
|
|
// angle23
|
|
|
|
dthetadr[1][1][0] = sc2 * ((t2 * vb2x) + (vb3x * r12c2));
|
|
dthetadr[1][1][1] = sc2 * ((t2 * vb2y) + (vb3y * r12c2));
|
|
dthetadr[1][1][2] = sc2 * ((t2 * vb2z) + (vb3z * r12c2));
|
|
|
|
dthetadr[1][2][0] = sc2 * ((-t2 * vb2x) - (vb3x * r12c2) +
|
|
(t4 * vb3x) + (vb2x * r12c2));
|
|
dthetadr[1][2][1] = sc2 * ((-t2 * vb2y) - (vb3y * r12c2) +
|
|
(t4 * vb3y) + (vb2y * r12c2));
|
|
dthetadr[1][2][2] = sc2 * ((-t2 * vb2z) - (vb3z * r12c2) +
|
|
(t4 * vb3z) + (vb2z * r12c2));
|
|
|
|
dthetadr[1][3][0] = -sc2 * ((t4 * vb3x) + (vb2x * r12c2));
|
|
dthetadr[1][3][1] = -sc2 * ((t4 * vb3y) + (vb2y * r12c2));
|
|
dthetadr[1][3][2] = -sc2 * ((t4 * vb3z) + (vb2z * r12c2));
|
|
|
|
// angle/angle/torsion cutoff
|
|
|
|
da1 = acos(costh12) - aat_theta0_1[type] ;
|
|
da2 = acos(costh23) - aat_theta0_1[type] ;
|
|
double dtheta = aat_theta0_2[type]-aat_theta0_1[type];
|
|
|
|
fphi = 0.0;
|
|
fpphi = 0.0;
|
|
if (phil < 0) phil +=MY_2PI;
|
|
uf_lookup(type, phil, fphi, fpphi);
|
|
|
|
double gt = aat_k[type];
|
|
double gtt = aat_k[type];
|
|
double gpt = 0;
|
|
double gptt = 0;
|
|
|
|
if ( acos(costh12) > aat_theta0_1[type]) {
|
|
gt *= 1-da1*da1/dtheta/dtheta;
|
|
gpt = -aat_k[type]*2*da1/dtheta/dtheta;
|
|
}
|
|
|
|
if ( acos(costh23) > aat_theta0_1[type]) {
|
|
gtt *= 1-da2*da2/dtheta/dtheta;
|
|
gptt = -aat_k[type]*2*da2/dtheta/dtheta;
|
|
}
|
|
|
|
if (eflag) edihedral = gt*gtt*fphi;
|
|
|
|
for (i = 0; i < 4; i++)
|
|
for (j = 0; j < 3; j++)
|
|
fabcd[i][j] -= gt*gtt*fpphi*dphidr[i][j]
|
|
- gt*gptt*fphi*dthetadr[1][i][j] + gpt*gtt*fphi*dthetadr[0][i][j];
|
|
|
|
// apply force to each of 4 atoms
|
|
|
|
if (newton_bond || i1 < nlocal) {
|
|
f[i1][0] += fabcd[0][0];
|
|
f[i1][1] += fabcd[0][1];
|
|
f[i1][2] += fabcd[0][2];
|
|
}
|
|
|
|
if (newton_bond || i2 < nlocal) {
|
|
f[i2][0] += fabcd[1][0];
|
|
f[i2][1] += fabcd[1][1];
|
|
f[i2][2] += fabcd[1][2];
|
|
}
|
|
|
|
if (newton_bond || i3 < nlocal) {
|
|
f[i3][0] += fabcd[2][0];
|
|
f[i3][1] += fabcd[2][1];
|
|
f[i3][2] += fabcd[2][2];
|
|
}
|
|
|
|
if (newton_bond || i4 < nlocal) {
|
|
f[i4][0] += fabcd[3][0];
|
|
f[i4][1] += fabcd[3][1];
|
|
f[i4][2] += fabcd[3][2];
|
|
}
|
|
|
|
if (evflag)
|
|
ev_tally(i1,i2,i3,i4,nlocal,newton_bond,edihedral,
|
|
fabcd[0],fabcd[2],fabcd[3],
|
|
vb1x,vb1y,vb1z,vb2x,vb2y,vb2z,vb3x,vb3y,vb3z);
|
|
}
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
void DihedralTableCut::allocate()
|
|
{
|
|
allocated = 1;
|
|
int n = atom->ndihedraltypes;
|
|
|
|
memory->create(aat_k,n+1,"dihedral:aat_k");
|
|
memory->create(aat_theta0_1,n+1,"dihedral:aat_theta0_1");
|
|
memory->create(aat_theta0_2,n+1,"dihedral:aat_theta0_2");
|
|
|
|
memory->create(tabindex,n+1,"dihedral:tabindex");
|
|
memory->create(setflag,n+1,"dihedral:setflag");
|
|
for (int i = 1; i <= n; i++) setflag[i] = 0;
|
|
}
|
|
|
|
/* ----------------------------------------------------------------------
|
|
set coeffs for one or more types
|
|
arg1 = "aat" -> AngleAngleTorsion coeffs
|
|
arg1 -> Dihedral coeffs
|
|
------------------------------------------------------------------------- */
|
|
|
|
void DihedralTableCut::coeff(int narg, char **arg)
|
|
{
|
|
if (narg != 7) error->all(FLERR,"Incorrect args for dihedral coefficients");
|
|
if (!allocated) allocate();
|
|
|
|
int ilo,ihi;
|
|
utils::bounds(FLERR,arg[0],1,atom->ndihedraltypes,ilo,ihi,error);
|
|
|
|
double k_one = utils::numeric(FLERR,arg[2],false,lmp);
|
|
double theta0_1_one = utils::numeric(FLERR,arg[3],false,lmp);
|
|
double theta0_2_one = utils::numeric(FLERR,arg[4],false,lmp);
|
|
|
|
// convert theta0's from degrees to radians
|
|
|
|
for (int i = ilo; i <= ihi; i++) {
|
|
aat_k[i] = k_one;
|
|
aat_theta0_1[i] = theta0_1_one/180.0 * MY_PI;
|
|
aat_theta0_2[i] = theta0_2_one/180.0 * MY_PI;
|
|
}
|
|
|
|
int me;
|
|
MPI_Comm_rank(world,&me);
|
|
tables = (Table *)
|
|
memory->srealloc(tables,(ntables+1)*sizeof(Table), "dihedral:tables");
|
|
Table *tb = &tables[ntables];
|
|
null_table(tb);
|
|
if (me == 0) read_table(tb,arg[5],arg[6]);
|
|
bcast_table(tb);
|
|
|
|
// --- check the angle data for range errors ---
|
|
// --- and resolve issues with periodicity ---
|
|
|
|
if (tb->ninput < 2)
|
|
error->all(FLERR,fmt::format("Invalid dihedral table length: {}",arg[5]));
|
|
else if ((tb->ninput == 2) && (tabstyle == SPLINE))
|
|
error->all(FLERR,fmt::format("Invalid dihedral spline table length: {} "
|
|
"(Try linear)",arg[5]));
|
|
|
|
// check for monotonicity
|
|
for (int i=0; i < tb->ninput-1; i++) {
|
|
if (tb->phifile[i] >= tb->phifile[i+1]) {
|
|
auto err_msg =fmt::format("Dihedral table values are not increasing "
|
|
"({}, entry {})",arg[5],i+1);
|
|
if (i==0)
|
|
err_msg += "\n(This is probably a mistake with your table format.)\n";
|
|
error->all(FLERR,err_msg);
|
|
}
|
|
}
|
|
|
|
// check the range of angles
|
|
double philo = tb->phifile[0];
|
|
double phihi = tb->phifile[tb->ninput-1];
|
|
if (tb->use_degrees) {
|
|
if ((phihi - philo) >= 360)
|
|
error->all(FLERR,fmt::format("Dihedral table angle range must be < 360 "
|
|
"degrees ({})",arg[5]));
|
|
} else {
|
|
if ((phihi - philo) >= MY_2PI)
|
|
error->all(FLERR,fmt::format("Dihedral table angle range must be < 2*PI "
|
|
"radians ({})",arg[5]));
|
|
}
|
|
|
|
// convert phi from degrees to radians
|
|
if (tb->use_degrees) {
|
|
for (int i=0; i < tb->ninput; i++) {
|
|
tb->phifile[i] *= MY_PI/180.0;
|
|
// I assume that if angles are in degrees, then the forces (f=dU/dphi)
|
|
// are specified with "phi" in degrees as well.
|
|
tb->ffile[i] *= 180.0/MY_PI;
|
|
}
|
|
}
|
|
|
|
// We want all the phi dihedral angles to lie in the range from 0 to 2*PI.
|
|
// But I don't want to restrict users to input their data in this range.
|
|
// We also want the angles to be sorted in increasing order.
|
|
// This messy code fixes these problems with the user's data:
|
|
{
|
|
double *phifile_tmp = new double [tb->ninput]; //temporary arrays
|
|
double *ffile_tmp = new double [tb->ninput]; //used for sorting
|
|
double *efile_tmp = new double [tb->ninput];
|
|
|
|
// After re-imaging, does the range of angles cross the 0 or 2*PI boundary?
|
|
// If so, find the discontinuity:
|
|
int i_discontinuity = tb->ninput;
|
|
for (int i=0; i < tb->ninput; i++) {
|
|
double phi = tb->phifile[i];
|
|
// Add a multiple of 2*PI to phi until it lies in the range [0, 2*PI).
|
|
phi -= MY_2PI * floor(phi/MY_2PI);
|
|
phifile_tmp[i] = phi;
|
|
efile_tmp[i] = tb->efile[i];
|
|
ffile_tmp[i] = tb->ffile[i];
|
|
if ((i>0) && (phifile_tmp[i] < phifile_tmp[i-1])) {
|
|
//There should only be at most one discontinuity, because we have
|
|
//insured that the data was sorted before imaging, and because the
|
|
//range of angle values does not exceed 2*PI.
|
|
i_discontinuity = i;
|
|
}
|
|
}
|
|
|
|
int I = 0;
|
|
for (int i = i_discontinuity; i < tb->ninput; i++) {
|
|
tb->phifile[I] = phifile_tmp[i];
|
|
tb->efile[I] = efile_tmp[i];
|
|
tb->ffile[I] = ffile_tmp[i];
|
|
I++;
|
|
}
|
|
for (int i = 0; i < i_discontinuity; i++) {
|
|
tb->phifile[I] = phifile_tmp[i];
|
|
tb->efile[I] = efile_tmp[i];
|
|
tb->ffile[I] = ffile_tmp[i];
|
|
I++;
|
|
}
|
|
|
|
// clean up temporary storage
|
|
delete[] phifile_tmp;
|
|
delete[] ffile_tmp;
|
|
delete[] efile_tmp;
|
|
}
|
|
|
|
// spline read-in and compute r,e,f vectors within table
|
|
|
|
spline_table(tb);
|
|
compute_table(tb);
|
|
|
|
// Optional: allow the user to print out the interpolated spline tables
|
|
|
|
if (me == 0) {
|
|
if (!checkU_fname.empty()) {
|
|
std::ofstream checkU_file;
|
|
checkU_file.open(checkU_fname, std::ios::out);
|
|
for (int i=0; i < tablength; i++) {
|
|
double phi = i*MY_2PI/tablength;
|
|
double u = tb->e[i];
|
|
if (tb->use_degrees)
|
|
phi *= 180.0/MY_PI;
|
|
checkU_file << phi << " " << u << "\n";
|
|
}
|
|
checkU_file.close();
|
|
}
|
|
if (!checkF_fname.empty()) {
|
|
std::ofstream checkF_file;
|
|
checkF_file.open(checkF_fname, std::ios::out);
|
|
for (int i=0; i < tablength; i++) {
|
|
double phi = i*MY_2PI/tablength;
|
|
double f;
|
|
if ((tabstyle == SPLINE) && (tb->f_unspecified)) {
|
|
double dU_dphi =
|
|
// (If the user did not specify the forces now, AND the user
|
|
// selected the "spline" option, (as opposed to "linear")
|
|
// THEN the tb->f array is uninitialized, so there's
|
|
// no point to print out the contents of the tb->f[] array.
|
|
// Instead, later on, we will calculate the force using the
|
|
// -cyc_splintD() routine to calculate the derivative of the
|
|
// energy spline, using the energy data (tb->e[]).
|
|
// To be nice and report something, I do the same thing here.)
|
|
cyc_splintD(tb->phi, tb->e, tb->e2, tablength, MY_2PI,phi);
|
|
f = -dU_dphi;
|
|
} else
|
|
// Otherwise we calculated the tb->f[] array. Report its contents.
|
|
f = tb->f[i];
|
|
if (tb->use_degrees) {
|
|
phi *= 180.0/MY_PI;
|
|
// If the user wants degree angle units, we should convert our
|
|
// internal force tables (in energy/radians) to (energy/degrees)
|
|
f *= MY_PI/180.0;
|
|
}
|
|
checkF_file << phi << " " << f << "\n";
|
|
}
|
|
checkF_file.close();
|
|
} // if (checkF_fname && (strlen(checkF_fname) != 0))
|
|
} // if (me == 0)
|
|
|
|
// store ptr to table in tabindex
|
|
int count = 0;
|
|
for (int i = ilo; i <= ihi; i++) {
|
|
tabindex[i] = ntables;
|
|
//phi0[i] = tb->phi0; <- equilibrium dihedral angles not supported
|
|
setflag[i] = 1;
|
|
count++;
|
|
}
|
|
ntables++;
|
|
|
|
if (count == 0) error->all(FLERR,"Incorrect args for dihedral coefficients");
|
|
}
|