print single warning when some rho[i] exceeded rhomax of the current EAM potential
This commit is contained in:
@ -140,6 +140,20 @@ The OpenKIM Project at
|
|||||||
provides EAM potentials that can be used directly in LAMMPS with the
|
provides EAM potentials that can be used directly in LAMMPS with the
|
||||||
:doc:`kim command <kim_commands>` interface.
|
:doc:`kim command <kim_commands>` interface.
|
||||||
|
|
||||||
|
.. warning::
|
||||||
|
|
||||||
|
The EAM potential files tabulate the embedding energy as a function
|
||||||
|
of the local electron density :math:`\rho`. When atoms get too
|
||||||
|
close, this electron density may exceed the range for which the
|
||||||
|
embedding energy was tabulated for. For simplicity and to avoid
|
||||||
|
errors during equilibration of randomized geometries, LAMMPS will
|
||||||
|
assume a linearly increasing embedding energy for electron densities
|
||||||
|
beyond the maximum tabulated value. This usually means that the EAM
|
||||||
|
model is not a good model for the kind of system under investigation.
|
||||||
|
LAMMPS will print a single warning when this happens. It may be
|
||||||
|
harmless at the beginning of an equilibration but would be a big
|
||||||
|
concern for accuracy if it happens during production runs.
|
||||||
|
|
||||||
----------
|
----------
|
||||||
|
|
||||||
For style *eam*, potential values are read from a file that is in the
|
For style *eam*, potential values are read from a file that is in the
|
||||||
|
|||||||
@ -234,7 +234,7 @@ void PairEAMIntel::eval(const int offload, const int vflag,
|
|||||||
const int istride = fc.rhor_istride();
|
const int istride = fc.rhor_istride();
|
||||||
const int jstride = fc.rhor_jstride();
|
const int jstride = fc.rhor_jstride();
|
||||||
const int fstride = fc.frho_stride();
|
const int fstride = fc.frho_stride();
|
||||||
|
int beyond_rhomax = 0;
|
||||||
{
|
{
|
||||||
#if defined(__MIC__) && defined(_LMP_INTEL_OFFLOAD)
|
#if defined(__MIC__) && defined(_LMP_INTEL_OFFLOAD)
|
||||||
*timer_compute = MIC_Wtime();
|
*timer_compute = MIC_Wtime();
|
||||||
@ -453,7 +453,10 @@ void PairEAMIntel::eval(const int offload, const int vflag,
|
|||||||
if (EFLAG) {
|
if (EFLAG) {
|
||||||
flt_t phi = ((frho_spline_e[ioff].a*p + frho_spline_e[ioff].b)*p +
|
flt_t phi = ((frho_spline_e[ioff].a*p + frho_spline_e[ioff].b)*p +
|
||||||
frho_spline_e[ioff].c)*p + frho_spline_e[ioff].d;
|
frho_spline_e[ioff].c)*p + frho_spline_e[ioff].d;
|
||||||
if (rho[i] > frhomax) phi += fp_f[i] * (rho[i]-frhomax);
|
if (rho[i] > frhomax) {
|
||||||
|
phi += fp_f[i] * (rho[i]-frhomax);
|
||||||
|
beyond_rhomax = 1;
|
||||||
|
}
|
||||||
if (!ONETYPE) {
|
if (!ONETYPE) {
|
||||||
const int ptr_off=itype*ntypes + itype;
|
const int ptr_off=itype*ntypes + itype;
|
||||||
oscale = scale_f[ptr_off];
|
oscale = scale_f[ptr_off];
|
||||||
@ -653,6 +656,16 @@ void PairEAMIntel::eval(const int offload, const int vflag,
|
|||||||
else
|
else
|
||||||
fix->stop_watch(TIME_HOST_PAIR);
|
fix->stop_watch(TIME_HOST_PAIR);
|
||||||
|
|
||||||
|
if (EFLAG && (exceeded_rhomax >= 0)) {
|
||||||
|
MPI_Allreduce(&beyond_rhomax, &exceeded_rhomax, 1, MPI_INT, MPI_MAX, world);
|
||||||
|
if (exceeded_rhomax > 0) {
|
||||||
|
if (comm->me == 0)
|
||||||
|
error->warning(FLERR, "Local rho[i] exceeded rhomax of EAM potential table. "
|
||||||
|
"Computed embedding term is unreliable.");
|
||||||
|
exceeded_rhomax = -1;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
if (EFLAG || vflag)
|
if (EFLAG || vflag)
|
||||||
fix->add_result_array(f_start, ev_global, offload, eatom, 0, vflag);
|
fix->add_result_array(f_start, ev_global, offload, eatom, 0, vflag);
|
||||||
else
|
else
|
||||||
@ -849,4 +862,3 @@ void PairEAMIntel::unpack_forward_comm(int n, int first, double *buf,
|
|||||||
last = first + n;
|
last = first + n;
|
||||||
for (i = first; i < last; i++) fp_f[i] = buf[m++];
|
for (i = first; i < last; i++) fp_f[i] = buf[m++];
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@ -43,6 +43,7 @@ PairEAM::PairEAM(LAMMPS *lmp) : Pair(lmp)
|
|||||||
unit_convert_flag = utils::get_supported_conversions(utils::ENERGY);
|
unit_convert_flag = utils::get_supported_conversions(utils::ENERGY);
|
||||||
|
|
||||||
nmax = 0;
|
nmax = 0;
|
||||||
|
exceeded_rhomax = 0;
|
||||||
rho = nullptr;
|
rho = nullptr;
|
||||||
fp = nullptr;
|
fp = nullptr;
|
||||||
numforce = nullptr;
|
numforce = nullptr;
|
||||||
@ -145,6 +146,7 @@ void PairEAM::compute(int eflag, int vflag)
|
|||||||
double *coeff;
|
double *coeff;
|
||||||
int *ilist,*jlist,*numneigh,**firstneigh;
|
int *ilist,*jlist,*numneigh,**firstneigh;
|
||||||
|
|
||||||
|
int beyond_rhomax = 0;
|
||||||
evdwl = 0.0;
|
evdwl = 0.0;
|
||||||
ev_init(eflag,vflag);
|
ev_init(eflag,vflag);
|
||||||
|
|
||||||
@ -237,7 +239,9 @@ void PairEAM::compute(int eflag, int vflag)
|
|||||||
fp[i] = (coeff[0]*p + coeff[1])*p + coeff[2];
|
fp[i] = (coeff[0]*p + coeff[1])*p + coeff[2];
|
||||||
if (eflag) {
|
if (eflag) {
|
||||||
phi = ((coeff[3]*p + coeff[4])*p + coeff[5])*p + coeff[6];
|
phi = ((coeff[3]*p + coeff[4])*p + coeff[5])*p + coeff[6];
|
||||||
if (rho[i] > rhomax) phi += fp[i] * (rho[i]-rhomax);
|
if (rho[i] > rhomax) {
|
||||||
|
phi += fp[i] * (rho[i]-rhomax);
|
||||||
|
}
|
||||||
phi *= scale[type[i]][type[i]];
|
phi *= scale[type[i]][type[i]];
|
||||||
if (eflag_global) eng_vdwl += phi;
|
if (eflag_global) eng_vdwl += phi;
|
||||||
if (eflag_atom) eatom[i] += phi;
|
if (eflag_atom) eatom[i] += phi;
|
||||||
@ -322,6 +326,16 @@ void PairEAM::compute(int eflag, int vflag)
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
if (eflag && (exceeded_rhomax >= 0)) {
|
||||||
|
MPI_Allreduce(&beyond_rhomax, &exceeded_rhomax, 1, MPI_INT, MPI_MAX, world);
|
||||||
|
if (exceeded_rhomax > 0) {
|
||||||
|
if (comm->me == 0)
|
||||||
|
error->warning(FLERR, "Local rho[i] exceeded rhomax of EAM potential table. "
|
||||||
|
"Computed embedding term is unreliable.");
|
||||||
|
exceeded_rhomax = -1;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
if (vflag_fdotr) virial_fdotr_compute();
|
if (vflag_fdotr) virial_fdotr_compute();
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -424,6 +438,8 @@ void PairEAM::init_style()
|
|||||||
|
|
||||||
neighbor->add_request(this);
|
neighbor->add_request(this);
|
||||||
embedstep = -1;
|
embedstep = -1;
|
||||||
|
|
||||||
|
exceeded_rhomax = 0;
|
||||||
}
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
|
|||||||
@ -63,7 +63,8 @@ class PairEAM : public Pair {
|
|||||||
void swap_eam(double *, double **) override;
|
void swap_eam(double *, double **) override;
|
||||||
|
|
||||||
protected:
|
protected:
|
||||||
int nmax; // allocated size of per-atom arrays
|
int nmax; // allocated size of per-atom arrays
|
||||||
|
int exceeded_rhomax; // 0 or 1 when rho[i] exceeded rhomax, -1 when not to check
|
||||||
double cutforcesq;
|
double cutforcesq;
|
||||||
double **scale;
|
double **scale;
|
||||||
bigint embedstep; // timestep, the embedding term was computed
|
bigint embedstep; // timestep, the embedding term was computed
|
||||||
|
|||||||
@ -17,6 +17,7 @@
|
|||||||
|
|
||||||
#include "atom.h"
|
#include "atom.h"
|
||||||
#include "comm.h"
|
#include "comm.h"
|
||||||
|
#include "error.h"
|
||||||
#include "force.h"
|
#include "force.h"
|
||||||
#include "memory.h"
|
#include "memory.h"
|
||||||
#include "neigh_list.h"
|
#include "neigh_list.h"
|
||||||
@ -102,6 +103,7 @@ void PairEAMOMP::eval(int iifrom, int iito, ThrData * const thr)
|
|||||||
double *coeff;
|
double *coeff;
|
||||||
int *ilist,*jlist,*numneigh,**firstneigh;
|
int *ilist,*jlist,*numneigh,**firstneigh;
|
||||||
|
|
||||||
|
int beyond_rhomax = 0;
|
||||||
evdwl = 0.0;
|
evdwl = 0.0;
|
||||||
|
|
||||||
const auto * _noalias const x = (dbl3_t *) atom->x[0];
|
const auto * _noalias const x = (dbl3_t *) atom->x[0];
|
||||||
@ -203,7 +205,10 @@ void PairEAMOMP::eval(int iifrom, int iito, ThrData * const thr)
|
|||||||
fp[i] = (coeff[0]*p + coeff[1])*p + coeff[2];
|
fp[i] = (coeff[0]*p + coeff[1])*p + coeff[2];
|
||||||
if (EFLAG) {
|
if (EFLAG) {
|
||||||
phi = ((coeff[3]*p + coeff[4])*p + coeff[5])*p + coeff[6];
|
phi = ((coeff[3]*p + coeff[4])*p + coeff[5])*p + coeff[6];
|
||||||
if (rho[i] > rhomax) phi += fp[i] * (rho[i]-rhomax);
|
if (rho[i] > rhomax) {
|
||||||
|
phi += fp[i] * (rho[i]-rhomax);
|
||||||
|
beyond_rhomax = 1;
|
||||||
|
}
|
||||||
e_tally_thr(this, i, i, nlocal, NEWTON_PAIR, scale[type[i]][type[i]]*phi, 0.0, thr);
|
e_tally_thr(this, i, i, nlocal, NEWTON_PAIR, scale[type[i]][type[i]]*phi, 0.0, thr);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -298,6 +303,16 @@ void PairEAMOMP::eval(int iifrom, int iito, ThrData * const thr)
|
|||||||
f[i].y += fytmp;
|
f[i].y += fytmp;
|
||||||
f[i].z += fztmp;
|
f[i].z += fztmp;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
if (EFLAG && (exceeded_rhomax >= 0)) {
|
||||||
|
MPI_Allreduce(&beyond_rhomax, &exceeded_rhomax, 1, MPI_INT, MPI_MAX, world);
|
||||||
|
if (exceeded_rhomax > 0) {
|
||||||
|
if (comm->me == 0)
|
||||||
|
error->warning(FLERR, "Local rho[i] exceeded rhomax of EAM potential table. "
|
||||||
|
"Computed embedding term is unreliable.");
|
||||||
|
exceeded_rhomax = -1;
|
||||||
|
}
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
|
|||||||
@ -23,6 +23,7 @@
|
|||||||
|
|
||||||
#include "atom.h"
|
#include "atom.h"
|
||||||
#include "comm.h"
|
#include "comm.h"
|
||||||
|
#include "error.h"
|
||||||
#include "force.h"
|
#include "force.h"
|
||||||
#include "memory.h"
|
#include "memory.h"
|
||||||
#include "neigh_list.h"
|
#include "neigh_list.h"
|
||||||
@ -118,6 +119,7 @@ template <int EVFLAG, int EFLAG, int NEWTON_PAIR> void PairEAMOpt::eval()
|
|||||||
|
|
||||||
int ntypes = atom->ntypes;
|
int ntypes = atom->ntypes;
|
||||||
int ntypes2 = ntypes * ntypes;
|
int ntypes2 = ntypes * ntypes;
|
||||||
|
int beyond_rhomax = 0;
|
||||||
|
|
||||||
auto *_noalias fast_alpha =
|
auto *_noalias fast_alpha =
|
||||||
(fast_alpha_t *) malloc((size_t) ntypes2 * (nr + 1) * sizeof(fast_alpha_t));
|
(fast_alpha_t *) malloc((size_t) ntypes2 * (nr + 1) * sizeof(fast_alpha_t));
|
||||||
@ -251,7 +253,10 @@ template <int EVFLAG, int EFLAG, int NEWTON_PAIR> void PairEAMOpt::eval()
|
|||||||
fp[i] = (coeff[0] * p + coeff[1]) * p + coeff[2];
|
fp[i] = (coeff[0] * p + coeff[1]) * p + coeff[2];
|
||||||
if (EFLAG) {
|
if (EFLAG) {
|
||||||
double phi = ((coeff[3] * p + coeff[4]) * p + coeff[5]) * p + coeff[6];
|
double phi = ((coeff[3] * p + coeff[4]) * p + coeff[5]) * p + coeff[6];
|
||||||
if (rho[i] > rhomax) phi += fp[i] * (rho[i] - rhomax);
|
if (rho[i] > rhomax) {
|
||||||
|
phi += fp[i] * (rho[i] - rhomax);
|
||||||
|
beyond_rhomax = 1;
|
||||||
|
}
|
||||||
phi *= scale[type[i]][type[i]];
|
phi *= scale[type[i]][type[i]];
|
||||||
if (eflag_global) eng_vdwl += phi;
|
if (eflag_global) eng_vdwl += phi;
|
||||||
if (eflag_atom) eatom[i] += phi;
|
if (eflag_atom) eatom[i] += phi;
|
||||||
@ -361,5 +366,15 @@ template <int EVFLAG, int EFLAG, int NEWTON_PAIR> void PairEAMOpt::eval()
|
|||||||
free(fast_gamma);
|
free(fast_gamma);
|
||||||
fast_gamma = nullptr;
|
fast_gamma = nullptr;
|
||||||
|
|
||||||
|
if (EFLAG && (exceeded_rhomax >= 0)) {
|
||||||
|
MPI_Allreduce(&beyond_rhomax, &exceeded_rhomax, 1, MPI_INT, MPI_MAX, world);
|
||||||
|
if (exceeded_rhomax > 0) {
|
||||||
|
if (comm->me == 0)
|
||||||
|
error->warning(FLERR, "Local rho[i] exceeded rhomax of EAM potential table. "
|
||||||
|
"Computed embedding term is unreliable.");
|
||||||
|
exceeded_rhomax = -1;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
if (vflag_fdotr) virial_fdotr_compute();
|
if (vflag_fdotr) virial_fdotr_compute();
|
||||||
}
|
}
|
||||||
|
|||||||
Reference in New Issue
Block a user