add support for determining damping scale factor from atom-style variable

This commit is contained in:
Axel Kohlmeyer
2022-03-23 10:43:26 -04:00
parent 1fd699d279
commit ae41996967
6 changed files with 199 additions and 46 deletions

View File

@ -19,9 +19,10 @@ Syntax
.. parsed-literal::
keyword = *scale*
*scale* values = type ratio
*scale* values = *type ratio* or *v_name*
type = atom type (1-N)
ratio = factor to scale the damping coefficients by
v_name = reference to atom style variable *name*
Examples
""""""""
@ -30,6 +31,7 @@ Examples
fix 1 all damping/cundall 0.8 0.8
fix 1 all damping/cundall 0.8 0.5 scale 3 2.5
fix a all damping/cundall 0.8 0.5 scale v_radscale
Description
"""""""""""
@ -58,9 +60,13 @@ particle. Damping is applied component-by-component in each direction
{T_d}_k = - \gamma_a \, T_k \, \mathrm{sign}(T_k \omega_k)
The larger the coefficients, the faster the kinetic energy is reduced.
If the optional keyword *scale* is used, :math:`\gamma_l` and :math:`\gamma_a`
can be scaled up or down by the specified factor for atoms of that type.
It can be used multiple times to adjust :math:`\gamma` for several atom types.
can be scaled up or down by the specified factor for atoms. This factor can be
set for different atom types and thus the *scale* keyword used multiple times
followed by the atom type and the associated scale factor. Alternately the
scaling factor can be computed for each atom (e.g. based on its radius) by
using an :doc:`atom-style variable <variable>`.
.. Note::

View File

@ -18,9 +18,10 @@ Syntax
.. parsed-literal::
keyword = *scale*
*scale* values = type ratio
*scale* values = *type ratio* or *v_name*
type = atom type (1-N)
ratio = factor to scale the damping coefficient by
ratio = factor to scale the damping coefficients by
v_name = reference to atom style variable *name*
Examples
""""""""
@ -29,6 +30,7 @@ Examples
fix 1 flow viscous/sphere 0.1
fix 1 damp viscous/sphere 0.5 scale 3 2.5
fix 1 damp viscous/sphere 0.5 scale v_radscale
Description
"""""""""""
@ -43,9 +45,13 @@ technique.
The damping torque :math:`T_i` is given by :math:`T_i = - \gamma \omega_i`.
The larger the coefficient, the faster the rotational kinetic energy is reduced.
If the optional keyword *scale* is used, :math:`\gamma` can be scaled up or
down by the specified factor for atoms of that type. It can be used
multiple times to adjust :math:`\gamma` for several atom types.
If the optional keyword *scale* is used, :math:`\gamma` can be scaled up
or down by the specified factor for atoms. This factor can be set for
different atom types and thus the *scale* keyword used multiple times
followed by the atom type and the associated scale factor. Alternately
the scaling factor can be computed for each atom (e.g. based on its
radius) by using an :doc:`atom-style variable <variable>`.
.. note::

View File

@ -14,19 +14,29 @@
#include "fix_viscous_sphere.h"
#include "atom.h"
#include "comm.h"
#include "error.h"
#include "input.h"
#include "memory.h"
#include "modify.h"
#include "respa.h"
#include "update.h"
#include "variable.h"
#include <cstring>
using namespace LAMMPS_NS;
using namespace FixConst;
// type of scaling
enum { NONE, TYPE, VARIABLE };
/* ---------------------------------------------------------------------- */
FixViscousSphere::FixViscousSphere(LAMMPS *lmp, int narg, char **arg) :
Fix(lmp, narg, arg), gamma(nullptr)
FixViscousSphere::FixViscousSphere(LAMMPS *_lmp, int narg, char **arg) :
Fix(_lmp, narg, arg), scalegamma(nullptr), scaleval(nullptr), scalevarid(nullptr),
scalestyle(NONE)
{
dynamic_group_allow = 1;
@ -34,22 +44,41 @@ FixViscousSphere::FixViscousSphere(LAMMPS *lmp, int narg, char **arg) :
if (narg < 4) error->all(FLERR, "Illegal fix viscous/sphere command");
double gamma_one = utils::numeric(FLERR, arg[3], false, lmp);
gamma = new double[atom->ntypes + 1];
for (int i = 1; i <= atom->ntypes; i++) gamma[i] = gamma_one;
gamma = utils::numeric(FLERR, arg[3], false, lmp);
// optional args
int iarg = 4;
while (iarg < narg) {
if (strcmp(arg[iarg], "scale") == 0) {
if (iarg + 3 > narg) error->all(FLERR, "Illegal fix viscous/sphere command");
int itype = utils::inumeric(FLERR, arg[iarg + 1], false, lmp);
double scale = utils::numeric(FLERR, arg[iarg + 2], false, lmp);
if (itype <= 0 || itype > atom->ntypes)
error->all(FLERR, "Illegal fix viscous/sphere command");
gamma[itype] = gamma_one * scale;
iarg += 3;
if (iarg + 2 > narg) error->all(FLERR, "Illegal fix viscous/sphere command");
if (utils::strmatch(arg[iarg + 1], "^v_")) {
if (scalestyle != NONE) error->all(FLERR, "Must use only one style of scaling");
scalevarid = utils::strdup(arg[iarg + 1] + 2);
int ivariable = input->variable->find(scalevarid);
if (ivariable < 0)
error->all(FLERR, "Variable name {} for fix viscous/sphere does not exist", scalevarid);
if (input->variable->atomstyle(ivariable) == 0)
error->all(FLERR, "Fix viscous/scale variable {} is not atom-style variable", scalevarid);
scalestyle = VARIABLE;
memory->destroy(scaleval);
memory->create(scaleval, atom->nmax, "fix_viscous/sphere:scaleval");
iarg += 2;
} else {
if (scalestyle == VARIABLE) error->all(FLERR, "Must use only one style of scaling");
if (iarg + 3 > narg) error->all(FLERR, "Illegal fix viscous/sphere command");
if (!scalegamma) {
scalegamma = new double[atom->ntypes + 1];
for (int i = 1; i <= atom->ntypes; i++) scalegamma[i] = 1.0;
}
int itype = utils::inumeric(FLERR, arg[iarg + 1], false, lmp);
double scale = utils::numeric(FLERR, arg[iarg + 2], false, lmp);
if ((itype < 1) || (itype > atom->ntypes))
error->all(FLERR, "Atom type {} out of range for fix viscous/sphere command:", itype);
scalegamma[itype] = scale;
scalestyle = TYPE;
iarg += 3;
}
} else
error->all(FLERR, "Illegal fix viscous/sphere command");
}
@ -62,7 +91,9 @@ FixViscousSphere::FixViscousSphere(LAMMPS *lmp, int narg, char **arg) :
FixViscousSphere::~FixViscousSphere()
{
delete[] gamma;
memory->destroy(scaleval);
delete[] scalegamma;
delete[] scalevarid;
}
/* ---------------------------------------------------------------------- */
@ -86,6 +117,22 @@ void FixViscousSphere::init()
ilevel_respa = max_respa = ((Respa *) update->integrate)->nlevels - 1;
if (respa_level >= 0) ilevel_respa = MIN(respa_level, max_respa);
}
bool fflag = false;
for (auto ifix : modify->get_fix_list()) {
if (fflag && (comm->me == 0) && (ifix->setmask() & POST_FORCE))
error->warning(FLERR, "Fix {} alters forces after fix viscous/sphere", ifix->id);
if (ifix == this) fflag = true;
}
if (scalestyle == VARIABLE) {
int ivariable = input->variable->find(scalevarid);
if (ivariable < 0)
error->all(FLERR, "Variable name {} for fix viscous/sphere does not exist", scalevarid);
if (input->variable->atomstyle(ivariable) == 0)
error->all(FLERR, "Fix viscous/sphere variable {} is not atom-style variable", scalevarid);
scalevar = ivariable;
}
}
/* ---------------------------------------------------------------------- */
@ -123,10 +170,20 @@ void FixViscousSphere::post_force(int /*vflag*/)
int nlocal = atom->nlocal;
double drag;
if (scalestyle == VARIABLE) {
memory->grow(scaleval, atom->nmax, "fix_viscous/sphere:scaleval");
input->variable->compute_atom(scalevar, igroup, scaleval, 1, 0);
}
for (int i = 0; i < nlocal; i++)
if (mask[i] & groupbit) {
drag = gamma[type[i]];
if (scalestyle == TYPE) {
drag = gamma * scalegamma[type[i]];
} else if (scalestyle == VARIABLE) {
drag = gamma * scaleval[i];
} else { // scalestyle == NONE
drag = gamma;
}
torque[i][0] -= drag * omega[i][0];
torque[i][1] -= drag * omega[i][1];
torque[i][2] -= drag * omega[i][2];
@ -146,3 +203,15 @@ void FixViscousSphere::min_post_force(int vflag)
{
post_force(vflag);
}
/* ---------------------------------------------------------------------- */
double FixViscousSphere::memory_usage()
{
if (scalestyle == VARIABLE)
return (double) sizeof(double) * atom->nmax;
else if (scalestyle == TYPE)
return (double) sizeof(double) * atom->ntypes;
else
return 0.0;
}

View File

@ -35,9 +35,12 @@ class FixViscousSphere : public Fix {
void post_force(int);
void post_force_respa(int, int, int);
void min_post_force(int);
double memory_usage();
protected:
double *gamma;
double gamma, *scalegamma, *scaleval;
const char *scalevarid;
int scalestyle, scalevar;
int ilevel_respa;
};

View File

@ -14,19 +14,29 @@
#include "fix_damping_cundall.h"
#include "atom.h"
#include "comm.h"
#include "error.h"
#include "input.h"
#include "memory.h"
#include "modify.h"
#include "respa.h"
#include "update.h"
#include "variable.h"
#include <cstring>
using namespace LAMMPS_NS;
using namespace FixConst;
// type of scaling
enum { NONE, TYPE, VARIABLE };
/* ---------------------------------------------------------------------- */
FixDampingCundall::FixDampingCundall(LAMMPS *lmp, int narg, char **arg) :
Fix(lmp, narg, arg), gamma_lin(nullptr), gamma_ang(nullptr)
FixDampingCundall::FixDampingCundall(LAMMPS *_lmp, int narg, char **arg) :
Fix(_lmp, narg, arg), scalegamma(nullptr), scaleval(nullptr), scalevarid(nullptr),
scalestyle(NONE)
{
dynamic_group_allow = 1;
@ -34,28 +44,42 @@ FixDampingCundall::FixDampingCundall(LAMMPS *lmp, int narg, char **arg) :
if (narg < 5) error->all(FLERR, "Illegal fix damping/cundall command");
double gamma_lin_one = utils::numeric(FLERR, arg[3], false, lmp);
double gamma_ang_one = utils::numeric(FLERR, arg[4], false, lmp);
gamma_lin = new double[atom->ntypes + 1];
gamma_ang = new double[atom->ntypes + 1];
for (int i = 1; i <= atom->ntypes; i++) {
gamma_lin[i] = gamma_lin_one;
gamma_ang[i] = gamma_ang_one;
}
gamma_lin = utils::numeric(FLERR, arg[3], false, lmp);
gamma_ang = utils::numeric(FLERR, arg[4], false, lmp);
// optional args
int iarg = 5;
while (iarg < narg) {
if (strcmp(arg[iarg], "scale") == 0) {
if (iarg + 3 > narg) error->all(FLERR, "Illegal fix damping/cundall command");
int itype = utils::inumeric(FLERR, arg[iarg + 1], false, lmp);
double scale = utils::numeric(FLERR, arg[iarg + 2], false, lmp);
if (itype <= 0 || itype > atom->ntypes)
error->all(FLERR, "Illegal fix damping/cundall command");
gamma_lin[itype] = gamma_lin_one * scale;
gamma_ang[itype] = gamma_ang_one * scale;
iarg += 3;
if (iarg + 2 > narg) error->all(FLERR, "Illegal fix damping/cundall command");
if (utils::strmatch(arg[iarg + 1], "^v_")) {
if (scalestyle != NONE) error->all(FLERR, "Must use only one style of scaling");
scalevarid = utils::strdup(arg[iarg + 1] + 2);
int ivariable = input->variable->find(scalevarid);
if (ivariable < 0)
error->all(FLERR, "Variable name {} for fix damping/cundall does not exist", scalevarid);
if (input->variable->atomstyle(ivariable) == 0)
error->all(FLERR, "Fix viscous/scale variable {} is not atom-style variable", scalevarid);
scalestyle = VARIABLE;
memory->destroy(scaleval);
memory->create(scaleval, atom->nmax, "fix_damping/cundall:scaleval");
iarg += 2;
} else {
if (scalestyle == VARIABLE) error->all(FLERR, "Must use only one style of scaling");
if (iarg + 3 > narg) error->all(FLERR, "Illegal fix damping/cundall command");
if (!scalegamma) {
scalegamma = new double[atom->ntypes + 1];
for (int i = 1; i <= atom->ntypes; i++) scalegamma[i] = 1.0;
}
int itype = utils::inumeric(FLERR, arg[iarg + 1], false, lmp);
double scale = utils::numeric(FLERR, arg[iarg + 2], false, lmp);
if ((itype < 1) || (itype > atom->ntypes))
error->all(FLERR, "Atom type {} out of range for fix damping/cundall command:", itype);
scalegamma[itype] = scale;
scalestyle = TYPE;
iarg += 3;
}
} else
error->all(FLERR, "Illegal fix damping/cundall command");
}
@ -68,8 +92,9 @@ FixDampingCundall::FixDampingCundall(LAMMPS *lmp, int narg, char **arg) :
FixDampingCundall::~FixDampingCundall()
{
delete[] gamma_lin;
delete[] gamma_ang;
memory->destroy(scaleval);
delete[] scalegamma;
delete[] scalevarid;
}
/* ---------------------------------------------------------------------- */
@ -93,6 +118,22 @@ void FixDampingCundall::init()
ilevel_respa = max_respa = ((Respa *) update->integrate)->nlevels - 1;
if (respa_level >= 0) ilevel_respa = MIN(respa_level, max_respa);
}
bool fflag = false;
for (auto ifix : modify->get_fix_list()) {
if (fflag && (comm->me == 0) && (ifix->setmask() & POST_FORCE))
error->warning(FLERR, "Fix {} alters forces after fix damping/cundall", ifix->id);
if (ifix == this) fflag = true;
}
if (scalestyle == VARIABLE) {
int ivariable = input->variable->find(scalevarid);
if (ivariable < 0)
error->all(FLERR, "Variable name {} for fix damping/cundall does not exist", scalevarid);
if (input->variable->atomstyle(ivariable) == 0)
error->all(FLERR, "Fix damping/cundall variable {} is not atom-style variable", scalevarid);
scalevar = ivariable;
}
}
/* ---------------------------------------------------------------------- */
@ -140,10 +181,23 @@ void FixDampingCundall::post_force(int /*vflag*/)
int signf0, signf1, signf2;
int signt0, signt1, signt2;
if (scalestyle == VARIABLE) {
memory->grow(scaleval, atom->nmax, "fix_damping/cundall:scaleval");
input->variable->compute_atom(scalevar, igroup, scaleval, 1, 0);
}
for (int i = 0; i < nlocal; i++)
if (mask[i] & groupbit) {
gamma_l = gamma_lin[type[i]];
gamma_a = gamma_ang[type[i]];
if (scalestyle == TYPE) {
gamma_l = gamma_lin * scalegamma[type[i]];
gamma_a = gamma_ang * scalegamma[type[i]];
} else if (scalestyle == VARIABLE) {
gamma_l = gamma_lin * scaleval[i];
gamma_a = gamma_ang * scaleval[i];
} else { // scalestyle NONE
gamma_l = gamma_lin;
gamma_a = gamma_ang;
}
signf0 = (f[i][0] * v[i][0] >= 0.0) ? 1.0 : -1.0;
signf1 = (f[i][1] * v[i][1] >= 0.0) ? 1.0 : -1.0;
@ -174,3 +228,15 @@ void FixDampingCundall::min_post_force(int vflag)
{
post_force(vflag);
}
/* ---------------------------------------------------------------------- */
double FixDampingCundall::memory_usage()
{
if (scalestyle == VARIABLE)
return (double) sizeof(double) * atom->nmax;
else if (scalestyle == TYPE)
return 2.0 * sizeof(double) * atom->ntypes;
else
return 0.0;
}

View File

@ -35,9 +35,12 @@ class FixDampingCundall : public Fix {
void post_force(int);
void post_force_respa(int, int, int);
void min_post_force(int);
double memory_usage();
protected:
double *gamma_lin, *gamma_ang;
double gamma_lin, gamma_ang, *scalegamma, *scaleval;
const char *scalevarid;
int scalestyle, scalevar;
int ilevel_respa;
};