Restarting additional BPM settings, adding virial contributions from tangential bonds

This commit is contained in:
jtclemm
2022-08-30 13:25:36 -06:00
parent 80257099de
commit 8fafd4d8fb
14 changed files with 5424 additions and 4519 deletions

View File

@ -14,6 +14,7 @@
#include "bond_bpm.h"
#include "atom.h"
#include "comm.h"
#include "domain.h"
#include "error.h"
#include "fix_bond_history.h"
@ -286,6 +287,33 @@ double BondBPM::equilibrium_distance(int /*i*/)
return max_stretch * r0_max_estimate / 1.5;
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void BondBPM::write_restart(FILE *fp)
{
fwrite(&max_stretch, sizeof(double), 1, fp);
fwrite(&r0_max_estimate, sizeof(double), 1, fp);
fwrite(&overlay_flag, sizeof(int), 1, fp);
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void BondBPM::read_restart(FILE *fp)
{
if (comm->me == 0) {
utils::sfread(FLERR, &max_stretch, sizeof(double), 1, fp, nullptr, error);
utils::sfread(FLERR, &r0_max_estimate, sizeof(double), 1, fp, nullptr, error);
utils::sfread(FLERR, &overlay_flag, sizeof(int), 1, fp, nullptr, error);
}
MPI_Bcast(&max_stretch, 1, MPI_DOUBLE, 0, world);
MPI_Bcast(&r0_max_estimate, 1, MPI_DOUBLE, 0, world);
MPI_Bcast(&overlay_flag, 1, MPI_INT, 0, world);
}
/* ---------------------------------------------------------------------- */
void BondBPM::process_broken(int i, int j)

View File

@ -29,9 +29,8 @@ class BondBPM : public Bond {
void init_style() override;
void settings(int, char **) override;
double equilibrium_distance(int) override;
void write_restart(FILE *) override{};
void read_restart(FILE *) override{};
void write_data(FILE *) override{};
void write_restart(FILE *) override;
void read_restart(FILE *) override;
double single(int, double, int, int, double &) override = 0;
protected:

View File

@ -35,18 +35,6 @@ using namespace LAMMPS_NS;
BondBPMRotational::BondBPMRotational(LAMMPS *_lmp) : BondBPM(_lmp)
{
Kr = nullptr;
Ks = nullptr;
Kt = nullptr;
Kb = nullptr;
Fcr = nullptr;
Fcs = nullptr;
Tct = nullptr;
Tcb = nullptr;
gnorm = nullptr;
gslide = nullptr;
groll = nullptr;
gtwist = nullptr;
partial_flag = 1;
smooth_flag = 1;
@ -195,7 +183,7 @@ void BondBPMRotational::store_data()
2) P. Mora & Y. Wang Advances in Geomcomputing 2009
---------------------------------------------------------------------- */
double BondBPMRotational::elastic_forces(int i1, int i2, int type, double &Fr, double r_mag,
double BondBPMRotational::elastic_forces(int i1, int i2, int type, double r_mag,
double r0_mag, double r_mag_inv, double * /*rhat*/,
double *r, double *r0, double *force1on2,
double *torque1on2, double *torque2on1)
@ -208,7 +196,7 @@ double BondBPMRotational::elastic_forces(int i1, int i2, int type, double &Fr, d
double q1[4], q2[4];
double q2inv[4], mq[4], mqinv[4], qp21[4], q21[4], qtmp[4];
double rb[3], rb_x_r0[3], s[3], t[3];
double Fs[3], Fsp[3], F_rot[3], Ftmp[3];
double Fr, Fs[3], Fsp[3], F_rot[3], Ftmp[3];
double Ts[3], Tb[3], Tt[3], Tbp[3], Ttp[3], Tsp[3], T_rot[3], Ttmp[3];
double **quat = atom->quat;
@ -377,7 +365,7 @@ double BondBPMRotational::elastic_forces(int i1, int i2, int type, double &Fr, d
Note: n points towards 1 vs pointing towards 2
---------------------------------------------------------------------- */
void BondBPMRotational::damping_forces(int i1, int i2, int type, double &Fr, double *rhat,
void BondBPMRotational::damping_forces(int i1, int i2, int type, double *rhat,
double *r, double *force1on2, double *torque1on2,
double *torque2on1)
{
@ -398,7 +386,6 @@ void BondBPMRotational::damping_forces(int i1, int i2, int type, double &Fr, dou
MathExtra::sub3(vn1, vn2, tmp);
MathExtra::scale3(gnorm[type], tmp);
Fr = MathExtra::lensq3(tmp);
MathExtra::add3(force1on2, tmp, force1on2);
// Damp tangential objective velocities
@ -464,7 +451,7 @@ void BondBPMRotational::compute(int eflag, int vflag)
int i1, i2, itmp, n, type;
double r[3], r0[3], rhat[3];
double rsq, r0_mag, r_mag, r_mag_inv;
double Fr, breaking, smooth;
double breaking, smooth;
double force1on2[3], torque1on2[3], torque2on1[3];
ev_init(eflag, vflag);
@ -520,7 +507,7 @@ void BondBPMRotational::compute(int eflag, int vflag)
// Calculate forces, check if bond breaks
// ------------------------------------------------------//
breaking = elastic_forces(i1, i2, type, Fr, r_mag, r0_mag, r_mag_inv, rhat, r, r0, force1on2,
breaking = elastic_forces(i1, i2, type, r_mag, r0_mag, r_mag_inv, rhat, r, r0, force1on2,
torque1on2, torque2on1);
if (breaking >= 1.0) {
@ -529,7 +516,7 @@ void BondBPMRotational::compute(int eflag, int vflag)
continue;
}
damping_forces(i1, i2, type, Fr, rhat, r, force1on2, torque1on2, torque2on1);
damping_forces(i1, i2, type, rhat, r, force1on2, torque1on2, torque2on1);
if (smooth_flag) {
smooth = breaking * breaking;
@ -562,7 +549,7 @@ void BondBPMRotational::compute(int eflag, int vflag)
torque[i2][2] += torque1on2[2] * smooth;
}
if (evflag) ev_tally(i1, i2, nlocal, newton_bond, 0.0, Fr * smooth, r[0], r[1], r[2]);
if (evflag) ev_tally_xyz(i1, i2, nlocal, newton_bond, 0.0, -force1on2[0] * smooth, -force1on2[1] * smooth, -force1on2[2] * smooth, r[0], r[1], r[2]);
}
}
@ -688,6 +675,9 @@ void BondBPMRotational::settings(int narg, char **arg)
void BondBPMRotational::write_restart(FILE *fp)
{
BondBPM::write_restart(fp);
write_restart_settings(fp);
fwrite(&Kr[1], sizeof(double), atom->nbondtypes, fp);
fwrite(&Ks[1], sizeof(double), atom->nbondtypes, fp);
fwrite(&Kt[1], sizeof(double), atom->nbondtypes, fp);
@ -708,6 +698,8 @@ void BondBPMRotational::write_restart(FILE *fp)
void BondBPMRotational::read_restart(FILE *fp)
{
BondBPM::read_restart(fp);
read_restart_settings(fp);
allocate();
if (comm->me == 0) {
@ -740,6 +732,26 @@ void BondBPMRotational::read_restart(FILE *fp)
for (int i = 1; i <= atom->nbondtypes; i++) setflag[i] = 1;
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void BondBPMRotational::write_restart_settings(FILE *fp)
{
fwrite(&smooth_flag, sizeof(int), 1, fp);
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void BondBPMRotational::read_restart_settings(FILE *fp)
{
if (comm->me == 0)
utils::sfread(FLERR, &smooth_flag, sizeof(int), 1, fp, nullptr, error);
MPI_Bcast(&smooth_flag, 1, MPI_INT, 0, world);
}
/* ----------------------------------------------------------------------
proc 0 writes to data file
------------------------------------------------------------------------- */
@ -785,14 +797,13 @@ double BondBPMRotational::single(int type, double rsq, int i, int j, double &ffo
r_mag_inv = 1.0 / r_mag;
MathExtra::scale3(r_mag_inv, r, rhat);
double breaking, smooth, Fr;
double force1on2[3], torque1on2[3], torque2on1[3];
breaking = elastic_forces(i, j, type, Fr, r_mag, r0_mag, r_mag_inv, rhat, r, r0, force1on2,
double breaking = elastic_forces(i, j, type, r_mag, r0_mag, r_mag_inv, rhat, r, r0, force1on2,
torque1on2, torque2on1);
fforce = Fr;
damping_forces(i, j, type, Fr, rhat, r, force1on2, torque1on2, torque2on1);
fforce += Fr;
damping_forces(i, j, type, rhat, r, force1on2, torque1on2, torque2on1);
fforce = MathExtra::len3(force1on2);
double smooth = 1.0;
if (smooth_flag) {
smooth = breaking * breaking;
smooth = 1.0 - smooth * smooth;
@ -806,16 +817,16 @@ double BondBPMRotational::single(int type, double rsq, int i, int j, double &ffo
svector[1] = -r0[0];
svector[2] = -r0[1];
svector[3] = -r0[2];
svector[4] = force1on2[0];
svector[5] = force1on2[1];
svector[6] = force1on2[2];
svector[4] = force1on2[0] * smooth;
svector[5] = force1on2[1] * smooth;
svector[6] = force1on2[2] * smooth;
} else {
svector[1] = r0[0];
svector[2] = r0[1];
svector[3] = r0[2];
svector[4] = -force1on2[0];
svector[5] = -force1on2[1];
svector[6] = -force1on2[2];
svector[4] = -force1on2[0] * smooth;
svector[5] = -force1on2[1] * smooth;
svector[6] = -force1on2[2] * smooth;
}

View File

@ -34,6 +34,8 @@ class BondBPMRotational : public BondBPM {
void settings(int, char **) override;
void write_restart(FILE *) override;
void read_restart(FILE *) override;
void write_restart_settings(FILE *) override;
void read_restart_settings(FILE *) override;
void write_data(FILE *) override;
double single(int, double, int, int, double &) override;
@ -44,9 +46,9 @@ class BondBPMRotational : public BondBPM {
double acos_limit(double);
double elastic_forces(int, int, int, double &, double, double, double, double *, double *,
double elastic_forces(int, int, int, double, double, double, double *, double *,
double *, double *, double *, double *);
void damping_forces(int, int, int, double &, double *, double *, double *, double *, double *);
void damping_forces(int, int, int, double *, double *, double *, double *, double *);
void allocate();
void store_data();

View File

@ -311,6 +311,8 @@ void BondBPMSpring::settings(int narg, char **arg)
void BondBPMSpring::write_restart(FILE *fp)
{
write_restart_settings(fp);
fwrite(&k[1], sizeof(double), atom->nbondtypes, fp);
fwrite(&ecrit[1], sizeof(double), atom->nbondtypes, fp);
fwrite(&gamma[1], sizeof(double), atom->nbondtypes, fp);
@ -322,6 +324,7 @@ void BondBPMSpring::write_restart(FILE *fp)
void BondBPMSpring::read_restart(FILE *fp)
{
read_restart_settings(fp);
allocate();
if (comm->me == 0) {
@ -336,6 +339,28 @@ void BondBPMSpring::read_restart(FILE *fp)
for (int i = 1; i <= atom->nbondtypes; i++) setflag[i] = 1;
}
/* ----------------------------------------------------------------------
proc 0 writes to restart file
------------------------------------------------------------------------- */
void BondBPMSpring::write_restart_settings(FILE *fp)
{
BondBPM::write_restart_settings(fp);
fwrite(&smooth_flag, sizeof(int), 1, fp);
}
/* ----------------------------------------------------------------------
proc 0 reads from restart file, bcasts
------------------------------------------------------------------------- */
void BondBPMSpring::read_restart_settings(FILE *fp)
{
BondBPM::read_restart_settings(fp);
if (comm->me == 0)
utils::sfread(FLERR, &smooth_flag, sizeof(int), 1, fp, nullptr, error);
MPI_Bcast(&smooth_flag, 1, MPI_INT, 0, world);
}
/* ----------------------------------------------------------------------
proc 0 writes to data file
------------------------------------------------------------------------- */

View File

@ -34,6 +34,8 @@ class BondBPMSpring : public BondBPM {
void settings(int, char **) override;
void write_restart(FILE *) override;
void read_restart(FILE *) override;
void write_restart_settings(FILE *) override;
void read_restart_settings(FILE *) override;
void write_data(FILE *) override;
double single(int, double, int, int, double &) override;

View File

@ -250,6 +250,91 @@ void Bond::ev_tally(int i, int j, int nlocal, int newton_bond, double ebond, dou
}
}
/* ----------------------------------------------------------------------
tally energy and virial into global or per-atom accumulators
for virial, have delx,dely,delz and fx,fy,fz
------------------------------------------------------------------------- */
void Bond::ev_tally_xyz(int i, int j, int nlocal, int newton_bond,
double ebond, double fx, double fy, double fz,
double delx, double dely, double delz)
{
double ebondhalf,v[6];
if (eflag_either) {
if (eflag_global) {
if (newton_bond) {
energy += ebond;
} else {
ebondhalf = 0.5 * ebond;
if (i < nlocal) energy += ebondhalf;
if (j < nlocal) energy += ebondhalf;
}
}
if (eflag_atom) {
ebondhalf = 0.5 * ebond;
if (newton_bond || i < nlocal) eatom[i] += ebondhalf;
if (newton_bond || j < nlocal) eatom[j] += ebondhalf;
}
}
if (vflag_either) {
v[0] = delx * fx;
v[1] = dely * fy;
v[2] = delz * fz;
v[3] = delx * fy;
v[4] = delx * fz;
v[5] = dely * fz;
if (vflag_global) {
if (newton_bond) {
virial[0] += v[0];
virial[1] += v[1];
virial[2] += v[2];
virial[3] += v[3];
virial[4] += v[4];
virial[5] += v[5];
} else {
if (i < nlocal) {
virial[0] += 0.5 * v[0];
virial[1] += 0.5 * v[1];
virial[2] += 0.5 * v[2];
virial[3] += 0.5 * v[3];
virial[4] += 0.5 * v[4];
virial[5] += 0.5 * v[5];
}
if (j < nlocal) {
virial[0] += 0.5 * v[0];
virial[1] += 0.5 * v[1];
virial[2] += 0.5 * v[2];
virial[3] += 0.5 * v[3];
virial[4] += 0.5 * v[4];
virial[5] += 0.5 * v[5];
}
}
}
if (vflag_atom) {
if (newton_bond || i < nlocal) {
vatom[i][0] += 0.5 * v[0];
vatom[i][1] += 0.5 * v[1];
vatom[i][2] += 0.5 * v[2];
vatom[i][3] += 0.5 * v[3];
vatom[i][4] += 0.5 * v[4];
vatom[i][5] += 0.5 * v[5];
}
if (newton_bond || j < nlocal) {
vatom[j][0] += 0.5 * v[0];
vatom[j][1] += 0.5 * v[1];
vatom[j][2] += 0.5 * v[2];
vatom[j][3] += 0.5 * v[3];
vatom[j][4] += 0.5 * v[4];
vatom[j][5] += 0.5 * v[5];
}
}
}
}
/* ----------------------------------------------------------------------
write a table of bond potential energy/force vs distance to a file
------------------------------------------------------------------------- */

View File

@ -102,6 +102,7 @@ class Bond : protected Pointers {
}
void ev_setup(int, int, int alloc = 1);
void ev_tally(int, int, int, int, double, double, double, double, double);
void ev_tally_xyz(int, int, int, int, double, double, double, double, double, double, double);
};
} // namespace LAMMPS_NS

View File

@ -17,6 +17,7 @@
#include "atom_vec.h"
#include "error.h"
#include "force.h"
#include "modify.h"
#include "neigh_list.h"
#include "neighbor.h"
#include "pair.h"
@ -50,6 +51,11 @@ int FixUpdateSpecialBonds::setmask()
void FixUpdateSpecialBonds::setup(int /*vflag*/)
{
// error if more than one fix update/special/bonds
if (modify->get_fix_by_style("UPDATE_SPECIAL_BONDS").size() > 1)
error->all(FLERR,"More than one fix update/special/bonds");
// Require atoms know about all of their bonds and if they break
if (force->newton_bond) error->all(FLERR, "Fix update/special/bonds requires Newton bond off");