Initial commit for pair style srp-react files
This commit is contained in:
@ -1,18 +1,23 @@
|
|||||||
.. index:: pair_style srp
|
.. index:: pair_style srp
|
||||||
|
.. index:: pair_style srp/react
|
||||||
|
|
||||||
pair_style srp command
|
pair_style srp command
|
||||||
======================
|
======================
|
||||||
|
|
||||||
|
pair_style srp/react command
|
||||||
|
============================
|
||||||
Syntax
|
Syntax
|
||||||
""""""
|
""""""
|
||||||
|
|
||||||
.. code-block:: LAMMPS
|
.. code-block:: LAMMPS
|
||||||
|
|
||||||
pair_style srp cutoff btype dist keyword value ...
|
pair_style srp cutoff btype dist keyword value ...
|
||||||
|
pair_style srp/react cutoff btype dist react-id keyword value ...
|
||||||
|
|
||||||
* cutoff = global cutoff for SRP interactions (distance units)
|
* cutoff = global cutoff for SRP interactions (distance units)
|
||||||
* btype = bond type to apply SRP interactions to (can be wildcard, see below)
|
* btype = bond type to apply SRP interactions to (can be wildcard, see below)
|
||||||
* distance = *min* or *mid*
|
* distance = *min* or *mid*
|
||||||
|
* react-id = id of either fix bond/break or fix bond/create
|
||||||
* zero or more keyword/value pairs may be appended
|
* zero or more keyword/value pairs may be appended
|
||||||
* keyword = *exclude*
|
* keyword = *exclude*
|
||||||
|
|
||||||
@ -36,13 +41,19 @@ Examples
|
|||||||
pair_coeff 1 2 none
|
pair_coeff 1 2 none
|
||||||
pair_coeff 2 2 srp 40.0
|
pair_coeff 2 2 srp 40.0
|
||||||
|
|
||||||
|
fix create all bond/create 100 1 2 1.0 1 prob 0.2 19852
|
||||||
|
pair_style hybrid dpd 1.0 1.0 12345 srp/react 0.8 * min create exclude yes
|
||||||
|
pair_coeff 1 1 dpd 60.0 50 1.0
|
||||||
|
pair_coeff 1 2 none
|
||||||
|
pair_coeff 2 2 srp/react 40.0
|
||||||
|
|
||||||
pair_style hybrid srp 0.8 2 mid
|
pair_style hybrid srp 0.8 2 mid
|
||||||
pair_coeff 1 1 none
|
pair_coeff 1 1 none
|
||||||
pair_coeff 1 2 none
|
pair_coeff 1 2 none
|
||||||
pair_coeff 2 2 srp 100.0 0.8
|
pair_coeff 2 2 srp 100.0 0.8
|
||||||
|
|
||||||
Description
|
Description
|
||||||
"""""""""""
|
|
||||||
|
|
||||||
Style *srp* computes a soft segmental repulsive potential (SRP) that
|
Style *srp* computes a soft segmental repulsive potential (SRP) that
|
||||||
acts between pairs of bonds. This potential is useful for preventing
|
acts between pairs of bonds. This potential is useful for preventing
|
||||||
@ -121,6 +132,18 @@ at the cutoff distance :math:`r_c`.
|
|||||||
|
|
||||||
----------
|
----------
|
||||||
|
|
||||||
|
Pair style *srp/react* interfaces the pair style *srp* with the
|
||||||
|
bond breaking and formation mechanisms provided by fix *bond/break*
|
||||||
|
and fix *bond/create*, respectively. When using this pair style, whenever a
|
||||||
|
bond breaking (or formation) reaction occurs, the corresponding fictitious
|
||||||
|
particle is deleted (or inserted) during the same simulation time step as
|
||||||
|
the reaction. This is useful in the simulation of reactive systems involving
|
||||||
|
large polymeric molecules :ref:`(Palkar) <Palkar>` where the segmental repulsive
|
||||||
|
potential is necessary to minimize topological violations, and also needs to be
|
||||||
|
turned on and off according to the progress of the reaction.
|
||||||
|
|
||||||
|
----------
|
||||||
|
|
||||||
Mixing, shift, table, tail correction, restart, rRESPA info
|
Mixing, shift, table, tail correction, restart, rRESPA info
|
||||||
"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
|
"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
|
||||||
|
|
||||||
@ -178,3 +201,8 @@ The default keyword value is exclude = yes.
|
|||||||
|
|
||||||
**(Sirk)** Sirk TW, Sliozberg YR, Brennan JK, Lisal M, Andzelm JW, J
|
**(Sirk)** Sirk TW, Sliozberg YR, Brennan JK, Lisal M, Andzelm JW, J
|
||||||
Chem Phys, 136 (13) 134903, 2012.
|
Chem Phys, 136 (13) 134903, 2012.
|
||||||
|
|
||||||
|
.. _Palkar:
|
||||||
|
|
||||||
|
**(Palkar)** Palkar V, Kuksenok O, J. Phys. Chem. B, 126 (1), 336-346, 2022
|
||||||
|
|
||||||
|
|||||||
@ -347,6 +347,7 @@ accelerated styles exist.
|
|||||||
* :doc:`spin/magelec <pair_spin_magelec>` -
|
* :doc:`spin/magelec <pair_spin_magelec>` -
|
||||||
* :doc:`spin/neel <pair_spin_neel>` -
|
* :doc:`spin/neel <pair_spin_neel>` -
|
||||||
* :doc:`srp <pair_srp>` -
|
* :doc:`srp <pair_srp>` -
|
||||||
|
* :doc:`srp/react <pair_srp>` -
|
||||||
* :doc:`sw <pair_sw>` - Stillinger-Weber 3-body potential
|
* :doc:`sw <pair_sw>` - Stillinger-Weber 3-body potential
|
||||||
* :doc:`sw/angle/table <pair_sw_angle_table>` - Stillinger-Weber potential with tabulated angular term
|
* :doc:`sw/angle/table <pair_sw_angle_table>` - Stillinger-Weber potential with tabulated angular term
|
||||||
* :doc:`sw/mod <pair_sw>` - modified Stillinger-Weber 3-body potential
|
* :doc:`sw/mod <pair_sw>` - modified Stillinger-Weber 3-body potential
|
||||||
|
|||||||
16
examples/PACKAGES/srp_react/README
Normal file
16
examples/PACKAGES/srp_react/README
Normal file
@ -0,0 +1,16 @@
|
|||||||
|
This directory contains an input script for performing
|
||||||
|
simulations with the srp/react pair style. The pair style
|
||||||
|
srp/react interfaces fix bond/break and fix bond/create commands
|
||||||
|
with the segmental repulsive potential. This is useful in simulating
|
||||||
|
reactions with soft potentials such as DPD where minimizing
|
||||||
|
topological violations is important.
|
||||||
|
|
||||||
|
The input script in.srp_react is an example of a simulation of
|
||||||
|
a degrading nanogel particle. An initial equilibrated structure
|
||||||
|
of a nanogel particle (prior to degradation) is read from the
|
||||||
|
restart file. The degradation reaction is simulated via
|
||||||
|
the fix bond/break command. The simulation will generate the
|
||||||
|
file bonds_broken.txt containing the number of bonds
|
||||||
|
broken and fraction of bonds intact over the simulation time.
|
||||||
|
For more details see the LAMMPS online documentation and the
|
||||||
|
paper: Palkar, V., & Kuksenok, O. (2022). JPC B, 126, 336.
|
||||||
22959
examples/PACKAGES/srp_react/gel_equil.dat
Normal file
22959
examples/PACKAGES/srp_react/gel_equil.dat
Normal file
File diff suppressed because it is too large
Load Diff
95
examples/PACKAGES/srp_react/in.srp_react
Normal file
95
examples/PACKAGES/srp_react/in.srp_react
Normal file
@ -0,0 +1,95 @@
|
|||||||
|
## srp_react example script
|
||||||
|
## Author: Vaibhav Palkar
|
||||||
|
##
|
||||||
|
## Simulates controlled degradation of a nanogel particle
|
||||||
|
## in a simulation box and prints statistics regarding
|
||||||
|
## the fraction of bonds broken over time.
|
||||||
|
|
||||||
|
variable rseeddpd equal 26817
|
||||||
|
variable rseedvel equal 5991
|
||||||
|
variable breakstep equal 10
|
||||||
|
variable probbreak equal 0.0009
|
||||||
|
variable rseedbreak equal 6777
|
||||||
|
|
||||||
|
# simulation time
|
||||||
|
#***********************************************
|
||||||
|
variable mainsteps equal 10000
|
||||||
|
|
||||||
|
# simulation setup
|
||||||
|
#***********************************************
|
||||||
|
units lj
|
||||||
|
atom_style molecular
|
||||||
|
boundary p p p
|
||||||
|
bond_style harmonic
|
||||||
|
#lattice fcc 3.0
|
||||||
|
comm_modify cutoff 4.0 vel yes
|
||||||
|
|
||||||
|
# initial nanogel
|
||||||
|
#***********************************************
|
||||||
|
read_data gel_equil.dat
|
||||||
|
|
||||||
|
# define groups, create solvent atoms
|
||||||
|
#**************************************************
|
||||||
|
group Npoly type 1 2 3 4
|
||||||
|
group water type 5
|
||||||
|
group N_all type 1 2 3 4 5
|
||||||
|
|
||||||
|
# density check
|
||||||
|
#***********************************************
|
||||||
|
variable N_atoms equal count(all)
|
||||||
|
variable tdens equal count(all)/vol
|
||||||
|
print "The system density is now ${tdens}"
|
||||||
|
|
||||||
|
# bond break settings
|
||||||
|
# #***********************************************
|
||||||
|
fix break N_all bond/break ${breakstep} 2 0 prob ${probbreak} ${rseedbreak}
|
||||||
|
|
||||||
|
# interaction parameter setting
|
||||||
|
#***********************************************
|
||||||
|
mass * 1.0
|
||||||
|
bond_coeff * 500.0 0.70
|
||||||
|
special_bonds lj 1 1 1
|
||||||
|
newton on
|
||||||
|
|
||||||
|
pair_style hybrid dpd 1.0 1.0 ${rseeddpd} srp/react 0.8 * mid break
|
||||||
|
#***********************************************
|
||||||
|
pair_coeff *5 *5 dpd 78.000 4.5 1.0
|
||||||
|
pair_coeff *4 5 dpd 79.500 4.5 1.0
|
||||||
|
pair_coeff *5 6 none
|
||||||
|
pair_coeff 6 6 srp/react 80.00 0.800
|
||||||
|
|
||||||
|
# initial velocity
|
||||||
|
#***********************************************
|
||||||
|
velocity all create 1.0 ${rseedvel} dist gaussian mom yes
|
||||||
|
|
||||||
|
# integrator control
|
||||||
|
#***********************************************
|
||||||
|
neighbor 0.3 bin
|
||||||
|
neigh_modify every 1 delay 5 check no
|
||||||
|
timestep 0.02
|
||||||
|
|
||||||
|
# Access variables of fix bond/break
|
||||||
|
#**********************************************
|
||||||
|
variable Nbreak equal f_break[2] # Number of bonds broken
|
||||||
|
variable TIME equal time
|
||||||
|
|
||||||
|
# ensemble setting
|
||||||
|
#***********************************************
|
||||||
|
fix 1 all nve
|
||||||
|
|
||||||
|
# print bonds breaking stats
|
||||||
|
# ***********************************************
|
||||||
|
variable TotBreak equal 100 # total breakable bonds in current system
|
||||||
|
fix print_bonds_broken all print 100 "${TIME} ${Nbreak} $((v_TotBreak-v_Nbreak)/v_TotBreak)" file bonds_broken.txt screen no title "time bonds_broken fraction_bonds_intact"
|
||||||
|
|
||||||
|
# thermo output
|
||||||
|
#***********************************************
|
||||||
|
thermo 100
|
||||||
|
thermo_style custom step temp pe ke etotal epair
|
||||||
|
thermo_modify flush yes norm no
|
||||||
|
|
||||||
|
reset_timestep 0
|
||||||
|
#dump 1 Npoly custom 5000 traj.lammpstrj id type x y z
|
||||||
|
#*********************************************
|
||||||
|
run ${mainsteps}
|
||||||
|
#***********************************************
|
||||||
4
src/.gitignore
vendored
4
src/.gitignore
vendored
@ -1522,6 +1522,8 @@
|
|||||||
/fix_smd_wall_surface.h
|
/fix_smd_wall_surface.h
|
||||||
/fix_srp.cpp
|
/fix_srp.cpp
|
||||||
/fix_srp.h
|
/fix_srp.h
|
||||||
|
/fix_srp_react.cpp
|
||||||
|
/fix_srp_react.h
|
||||||
/fix_tfmc.cpp
|
/fix_tfmc.cpp
|
||||||
/fix_tfmc.h
|
/fix_tfmc.h
|
||||||
/fix_ttm.cpp
|
/fix_ttm.cpp
|
||||||
@ -1562,6 +1564,8 @@
|
|||||||
/pair_smd_ulsph.h
|
/pair_smd_ulsph.h
|
||||||
/pair_srp.cpp
|
/pair_srp.cpp
|
||||||
/pair_srp.h
|
/pair_srp.h
|
||||||
|
/pair_srp_react.cpp
|
||||||
|
/pair_srp_react.h
|
||||||
/pair_thole.cpp
|
/pair_thole.cpp
|
||||||
/pair_thole.h
|
/pair_thole.h
|
||||||
/pair_buck_mdf.cpp
|
/pair_buck_mdf.cpp
|
||||||
|
|||||||
@ -25,6 +25,7 @@ FixStyle(bond/break,FixBondBreak);
|
|||||||
namespace LAMMPS_NS {
|
namespace LAMMPS_NS {
|
||||||
|
|
||||||
class FixBondBreak : public Fix {
|
class FixBondBreak : public Fix {
|
||||||
|
friend class FixSRPREACT;
|
||||||
public:
|
public:
|
||||||
FixBondBreak(class LAMMPS *, int, char **);
|
FixBondBreak(class LAMMPS *, int, char **);
|
||||||
~FixBondBreak() override;
|
~FixBondBreak() override;
|
||||||
|
|||||||
@ -25,6 +25,7 @@ FixStyle(bond/create,FixBondCreate);
|
|||||||
namespace LAMMPS_NS {
|
namespace LAMMPS_NS {
|
||||||
|
|
||||||
class FixBondCreate : public Fix {
|
class FixBondCreate : public Fix {
|
||||||
|
friend class FixSRPREACT;
|
||||||
public:
|
public:
|
||||||
FixBondCreate(class LAMMPS *, int, char **);
|
FixBondCreate(class LAMMPS *, int, char **);
|
||||||
~FixBondCreate() override;
|
~FixBondCreate() override;
|
||||||
|
|||||||
@ -54,7 +54,7 @@ class FixSRP : public Fix {
|
|||||||
|
|
||||||
double **array;
|
double **array;
|
||||||
|
|
||||||
private:
|
protected:
|
||||||
int btype;
|
int btype;
|
||||||
int bptype;
|
int bptype;
|
||||||
};
|
};
|
||||||
|
|||||||
186
src/MISC/fix_srp_react.cpp
Normal file
186
src/MISC/fix_srp_react.cpp
Normal file
@ -0,0 +1,186 @@
|
|||||||
|
// clang-format off
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
||||||
|
https://www.lammps.org/, 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: Vaibhav Palkar (Kuksenok Lab, Clemson University)
|
||||||
|
based on the pair srp code by Timothy Sirk (ARL)
|
||||||
|
|
||||||
|
This pair style interfaces the pair style srp with the stochastic reaction
|
||||||
|
fixes bond/break and bond/create by updating pseudo beads corresponding to
|
||||||
|
bonds as bond breaking and formation takes place. This is useful in
|
||||||
|
simulation of reactions in polymers with soft potentials such as DPD.
|
||||||
|
|
||||||
|
See the doc page for pair_style srp/react command for usage instructions.
|
||||||
|
|
||||||
|
There is an example script for this package in examples/PACKAGES/srp_react/.
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
#include "fix_srp_react.h"
|
||||||
|
|
||||||
|
#include "atom.h"
|
||||||
|
#include "atom_vec.h"
|
||||||
|
#include "comm.h"
|
||||||
|
#include "domain.h"
|
||||||
|
#include "error.h"
|
||||||
|
#include "force.h"
|
||||||
|
#include "memory.h"
|
||||||
|
#include "modify.h"
|
||||||
|
#include "neighbor.h"
|
||||||
|
|
||||||
|
#include "fix_bond_break.h"
|
||||||
|
#include "fix_bond_create.h"
|
||||||
|
|
||||||
|
#include <cmath>
|
||||||
|
#include <cstring>
|
||||||
|
|
||||||
|
using namespace LAMMPS_NS;
|
||||||
|
using namespace FixConst;
|
||||||
|
|
||||||
|
/* ---------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
FixSRPREACT::FixSRPREACT(LAMMPS *lmp, int narg, char **arg) : FixSRP(lmp, narg, arg)
|
||||||
|
{
|
||||||
|
// default idbreak and idcreate = NULL
|
||||||
|
idbreak = nullptr;
|
||||||
|
idcreate = nullptr;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ---------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
int FixSRPREACT::setmask()
|
||||||
|
{
|
||||||
|
int mask = 0;
|
||||||
|
mask |= PRE_FORCE;
|
||||||
|
mask |= PRE_EXCHANGE;
|
||||||
|
mask |= POST_RUN;
|
||||||
|
mask |= POST_NEIGHBOR;
|
||||||
|
|
||||||
|
return mask;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ---------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void FixSRPREACT::init()
|
||||||
|
{
|
||||||
|
if (force->pair_match("hybrid",1) == nullptr && force->pair_match("hybrid/overlay",1) == nullptr)
|
||||||
|
error->all(FLERR,"Cannot use pair srp without pair_style hybrid");
|
||||||
|
|
||||||
|
int has_rigid = 0;
|
||||||
|
for (int i = 0; i < modify->nfix; i++)
|
||||||
|
if (utils::strmatch(modify->fix[i]->style,"^rigid")) ++has_rigid;
|
||||||
|
|
||||||
|
if (has_rigid > 0)
|
||||||
|
error->all(FLERR,"Pair srp is not compatible with rigid fixes.");
|
||||||
|
|
||||||
|
if ((bptype < 1) || (bptype > atom->ntypes))
|
||||||
|
error->all(FLERR,"Illegal bond particle type");
|
||||||
|
|
||||||
|
// this fix must come before any fix which migrates atoms in its pre_exchange()
|
||||||
|
// because this fix's pre_exchange() creates per-atom data structure
|
||||||
|
// that data must be current for atom migration to carry it along
|
||||||
|
|
||||||
|
for (int i = 0; i < modify->nfix; i++) {
|
||||||
|
if (modify->fix[i] == this) break;
|
||||||
|
if (modify->fix[i]->pre_exchange_migrate)
|
||||||
|
error->all(FLERR,"Fix SRP comes after a fix which "
|
||||||
|
"migrates atoms in pre_exchange");
|
||||||
|
}
|
||||||
|
|
||||||
|
// setup neigh exclusions for diff atom types
|
||||||
|
// bond particles do not interact with other types
|
||||||
|
// type bptype only interacts with itself
|
||||||
|
|
||||||
|
for (int z = 1; z < atom->ntypes; z++) {
|
||||||
|
if (z == bptype)
|
||||||
|
continue;
|
||||||
|
neighbor->modify_params(fmt::format("exclude type {} {}",z,bptype));
|
||||||
|
}
|
||||||
|
|
||||||
|
// find fix bond break
|
||||||
|
if( idbreak != nullptr )
|
||||||
|
f_bb = (FixBondBreak *) modify->get_fix_by_id(idbreak);
|
||||||
|
|
||||||
|
// find fix bond create
|
||||||
|
if( idcreate != nullptr )
|
||||||
|
f_bc = (FixBondCreate *) modify->get_fix_by_id(idcreate);
|
||||||
|
|
||||||
|
// free memory
|
||||||
|
delete [] idbreak;
|
||||||
|
delete [] idcreate;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
rebuild bond particle array
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
void FixSRPREACT::post_neighbor()
|
||||||
|
{
|
||||||
|
// store ncalls as it is reset in fix srp setup pre force
|
||||||
|
int ncalls = neighbor->ncalls;
|
||||||
|
|
||||||
|
if( idbreak != nullptr)
|
||||||
|
if(f_bb->breakcount)
|
||||||
|
{
|
||||||
|
setup_pre_force(0);
|
||||||
|
|
||||||
|
//reset break count before exiting
|
||||||
|
// not reseting breakcount would lead to redundant rebuilds
|
||||||
|
f_bb->breakcount=0;
|
||||||
|
|
||||||
|
// count additional call during setup_pre_force
|
||||||
|
neighbor->ncalls = ncalls+1;
|
||||||
|
}
|
||||||
|
if( idcreate != nullptr)
|
||||||
|
if(f_bc->createcount)
|
||||||
|
{
|
||||||
|
setup_pre_force(0);
|
||||||
|
|
||||||
|
//reset create count before exiting
|
||||||
|
// not reseting createcount would lead to redundant rebuilds
|
||||||
|
f_bc->createcount=0;
|
||||||
|
|
||||||
|
// count additional call during setup_pre_force
|
||||||
|
neighbor->ncalls = ncalls+1;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
interface with pair class
|
||||||
|
in addition to bond type and bond particle type,
|
||||||
|
pair srp react sets id of either fix bond break or bond create
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
int FixSRPREACT::modify_param(int /*narg*/, char **arg)
|
||||||
|
{
|
||||||
|
if (strcmp(arg[0],"btype") == 0) {
|
||||||
|
btype = utils::inumeric(FLERR,arg[1],false,lmp);
|
||||||
|
return 2;
|
||||||
|
}
|
||||||
|
if (strcmp(arg[0],"bptype") == 0) {
|
||||||
|
bptype = utils::inumeric(FLERR,arg[1],false,lmp);
|
||||||
|
return 2;
|
||||||
|
}
|
||||||
|
if (strcmp(arg[0],"bond/break") == 0) {
|
||||||
|
int n = strlen(arg[1]) + 1;
|
||||||
|
idbreak = new char[n];
|
||||||
|
strcpy(idbreak,arg[1]);
|
||||||
|
return 2;
|
||||||
|
}
|
||||||
|
if (strcmp(arg[0],"bond/create") == 0) {
|
||||||
|
int n = strlen(arg[1]) + 1;
|
||||||
|
idcreate = new char[n];
|
||||||
|
strcpy(idcreate,arg[1]);
|
||||||
|
return 2;
|
||||||
|
}
|
||||||
|
return 0;
|
||||||
|
}
|
||||||
50
src/MISC/fix_srp_react.h
Normal file
50
src/MISC/fix_srp_react.h
Normal file
@ -0,0 +1,50 @@
|
|||||||
|
/* -*- c++ -*- ----------------------------------------------------------
|
||||||
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
||||||
|
https://www.lammps.org/, 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.
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
#ifdef FIX_CLASS
|
||||||
|
// clang-format off
|
||||||
|
FixStyle(SRPREACT,FixSRPREACT);
|
||||||
|
// clang-format on
|
||||||
|
#else
|
||||||
|
|
||||||
|
#ifndef LMP_FIX_SRP_REACT_H
|
||||||
|
#define LMP_FIX_SRP_REACT_H
|
||||||
|
|
||||||
|
#include "fix_srp.h"
|
||||||
|
|
||||||
|
namespace LAMMPS_NS {
|
||||||
|
|
||||||
|
class FixSRPREACT : public FixSRP {
|
||||||
|
public:
|
||||||
|
FixSRPREACT(class LAMMPS *, int, char **);
|
||||||
|
int setmask() override;
|
||||||
|
void init() override;
|
||||||
|
void post_neighbor() override;
|
||||||
|
int modify_param(int, char **) override;
|
||||||
|
|
||||||
|
private:
|
||||||
|
|
||||||
|
class FixBondBreak* f_bb;
|
||||||
|
char *idbreak;
|
||||||
|
class FixBondCreate* f_bc;
|
||||||
|
char *idcreate;
|
||||||
|
};
|
||||||
|
|
||||||
|
} // namespace LAMMPS_NS
|
||||||
|
|
||||||
|
#endif
|
||||||
|
#endif
|
||||||
|
|
||||||
|
/* ERROR/WARNING messages:
|
||||||
|
|
||||||
|
*/
|
||||||
@ -126,7 +126,7 @@ PairSRP::~PairSRP()
|
|||||||
}
|
}
|
||||||
|
|
||||||
// check nfix in case all fixes have already been deleted
|
// check nfix in case all fixes have already been deleted
|
||||||
if (modify->nfix) modify->delete_fix(f_srp->id);
|
if (modify->nfix && modify->get_fix_by_id(f_srp->id)!=nullptr) modify->delete_fix(f_srp->id);
|
||||||
}
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
|
|||||||
@ -40,7 +40,7 @@ class PairSRP : public Pair {
|
|||||||
void write_restart_settings(FILE *) override;
|
void write_restart_settings(FILE *) override;
|
||||||
void read_restart_settings(FILE *) override;
|
void read_restart_settings(FILE *) override;
|
||||||
|
|
||||||
private:
|
protected:
|
||||||
inline void onetwoexclude(int *&, int &, int *&, int *&, int **&);
|
inline void onetwoexclude(int *&, int &, int *&, int *&, int **&);
|
||||||
inline void remapBonds(int &);
|
inline void remapBonds(int &);
|
||||||
void allocate();
|
void allocate();
|
||||||
|
|||||||
239
src/MISC/pair_srp_react.cpp
Normal file
239
src/MISC/pair_srp_react.cpp
Normal file
@ -0,0 +1,239 @@
|
|||||||
|
// clang-format off
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
||||||
|
https://www.lammps.org/, 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: Vaibhav Palkar (Kuksenok Lab, Clemson University)
|
||||||
|
based on the pair srp code by Timothy Sirk (ARL)
|
||||||
|
|
||||||
|
This pair style interfaces the pair style srp with the stochastic reaction
|
||||||
|
fixes bond/break and bond/create by updating pseudo beads corresponding to
|
||||||
|
bonds as bond breaking and formation takes place. This is useful in
|
||||||
|
simulation of reactions in polymers with soft potentials such as DPD.
|
||||||
|
|
||||||
|
See the doc page for pair_style srp/react command for usage instructions.
|
||||||
|
|
||||||
|
There is an example script for this package in examples/PACKAGES/srp_react/.
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
#include "pair_srp_react.h"
|
||||||
|
|
||||||
|
#include "atom.h"
|
||||||
|
#include "citeme.h"
|
||||||
|
#include "comm.h"
|
||||||
|
#include "domain.h"
|
||||||
|
#include "error.h"
|
||||||
|
#include "fix.h"
|
||||||
|
#include "fix_srp_react.h"
|
||||||
|
#include "force.h"
|
||||||
|
#include "memory.h"
|
||||||
|
#include "modify.h"
|
||||||
|
#include "neigh_list.h"
|
||||||
|
#include "neighbor.h"
|
||||||
|
#include "output.h"
|
||||||
|
#include "thermo.h"
|
||||||
|
|
||||||
|
#include <cmath>
|
||||||
|
#include <cstring>
|
||||||
|
|
||||||
|
using namespace LAMMPS_NS;
|
||||||
|
|
||||||
|
static const char cite_srpreact[] =
|
||||||
|
"@Article{palkar2022\n"
|
||||||
|
" author = {Palkar, Vaibhav and Kuksenok, Olga},\n"
|
||||||
|
" title = {Controlling Degradation and Erosion of Polymer Networks: Insights from Mesoscale Modeling},\n"
|
||||||
|
" journal = {J. Phys. Chem. B},\n"
|
||||||
|
" year = {2022},\n"
|
||||||
|
" volume = {126},\n"
|
||||||
|
" pages = {336-346}\n"
|
||||||
|
"}\n\n";
|
||||||
|
|
||||||
|
static int srp_instance = 0;
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
constructor
|
||||||
|
---------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
PairSRPREACT::PairSRPREACT(LAMMPS *lmp) : PairSRP(lmp)
|
||||||
|
{
|
||||||
|
|
||||||
|
if (lmp->citeme) lmp->citeme->add(cite_srpreact);
|
||||||
|
|
||||||
|
// pair srp/react has its own fix, hence delete fix srp instance
|
||||||
|
// created in the constructor of pair srp
|
||||||
|
for( int ifix = 0; ifix<modify->nfix; ifix++)
|
||||||
|
if( strcmp(modify->get_fix_by_index(ifix)->style, "SRP") == 0)
|
||||||
|
modify->delete_fix(ifix);
|
||||||
|
|
||||||
|
// similar to fix SRP, create fix SRP REACT instance here with unique fix id
|
||||||
|
f_srp = (FixSRPREACT *) modify->add_fix(fmt::format("{:02d}_FIX_SRP_REACT all SRPREACT",srp_instance));
|
||||||
|
++srp_instance;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
global settings
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void PairSRPREACT::settings(int narg, char **arg)
|
||||||
|
{
|
||||||
|
if (narg < 3 || narg > 8)
|
||||||
|
error->all(FLERR,"Illegal pair_style command");
|
||||||
|
|
||||||
|
if (atom->tag_enable == 0)
|
||||||
|
error->all(FLERR,"Pair_style srp requires atom IDs");
|
||||||
|
|
||||||
|
cut_global = utils::numeric(FLERR,arg[0],false,lmp);
|
||||||
|
// wildcard
|
||||||
|
if (strcmp(arg[1],"*") == 0) {
|
||||||
|
btype = 0;
|
||||||
|
} else {
|
||||||
|
btype = utils::inumeric(FLERR,arg[1],false,lmp);
|
||||||
|
if ((btype > atom->nbondtypes) || (btype <= 0))
|
||||||
|
error->all(FLERR,"Illegal pair_style command");
|
||||||
|
}
|
||||||
|
|
||||||
|
// settings
|
||||||
|
midpoint = false;
|
||||||
|
min = false;
|
||||||
|
|
||||||
|
if (strcmp(arg[2],"min") == 0) min = true;
|
||||||
|
else if (strcmp(arg[2],"mid") == 0) midpoint = true;
|
||||||
|
else
|
||||||
|
error->all(FLERR,"Illegal pair_style command");
|
||||||
|
|
||||||
|
// default for bond/break and bond/create settings
|
||||||
|
bond_create=false;
|
||||||
|
bond_break=false;
|
||||||
|
idbreak = nullptr;
|
||||||
|
idcreate= nullptr;
|
||||||
|
|
||||||
|
// find whether id is of bond/break or bond/create
|
||||||
|
const char* reactid = arg[3];
|
||||||
|
if(strcmp(modify->get_fix_by_id(reactid)->style,"bond/break") == 0)
|
||||||
|
{
|
||||||
|
bond_break = true;
|
||||||
|
int n = strlen(reactid) + 1;
|
||||||
|
idbreak = new char[n];
|
||||||
|
strcpy(idbreak,reactid);
|
||||||
|
}
|
||||||
|
else if(strcmp(modify->get_fix_by_id(reactid)->style,"bond/create") == 0)
|
||||||
|
{
|
||||||
|
bond_create = true;
|
||||||
|
int n = strlen(reactid) + 1;
|
||||||
|
idcreate = new char[n];
|
||||||
|
strcpy(idcreate,reactid);
|
||||||
|
}
|
||||||
|
else
|
||||||
|
error->all(FLERR,"Illegal pair_style command");
|
||||||
|
|
||||||
|
int iarg = 4;
|
||||||
|
// default exclude 1-2
|
||||||
|
// scaling for 1-2, etc not supported
|
||||||
|
exclude = 1;
|
||||||
|
|
||||||
|
// use last atom type by default for bond particles
|
||||||
|
bptype = atom->ntypes;
|
||||||
|
|
||||||
|
while (iarg < narg) {
|
||||||
|
if (strcmp(arg[iarg],"exclude") == 0) {
|
||||||
|
if (iarg+2 > narg) error->all(FLERR,"Illegal pair srp command");
|
||||||
|
exclude = utils::logical(FLERR, arg[iarg+1], false, lmp);
|
||||||
|
if (min && !exclude) error->all(FLERR,"Illegal exclude option in pair srp command");
|
||||||
|
iarg += 2;
|
||||||
|
} else if (strcmp(arg[iarg],"bptype") == 0) {
|
||||||
|
if (iarg+2 > narg) error->all(FLERR,"Illegal pair srp command");
|
||||||
|
bptype = utils::inumeric(FLERR,arg[iarg+1],false,lmp);
|
||||||
|
if ((bptype < 1) || (bptype > atom->ntypes))
|
||||||
|
error->all(FLERR,"Illegal bond particle type for srp");
|
||||||
|
iarg += 2;
|
||||||
|
} else error->all(FLERR,"Illegal pair srp command");
|
||||||
|
}
|
||||||
|
|
||||||
|
// reset cutoffs if explicitly set
|
||||||
|
if (allocated) {
|
||||||
|
int i,j;
|
||||||
|
for (i = 1; i <= bptype; i++)
|
||||||
|
for (j = i; j <= bptype; j++)
|
||||||
|
if (setflag[i][j]) cut[i][j] = cut_global;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
init specific to this pair style
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void PairSRPREACT::init_style()
|
||||||
|
{
|
||||||
|
if (!force->newton_pair)
|
||||||
|
error->all(FLERR,"PairSRPREACT: Pair srp requires newton pair on");
|
||||||
|
|
||||||
|
// verify that fix SRP is still defined and has not been changed.
|
||||||
|
|
||||||
|
if (strcmp(f_srp->style,"SRPREACT") != 0)
|
||||||
|
error->all(FLERR,"Fix SRP has been changed unexpectedly");
|
||||||
|
|
||||||
|
if (comm->me == 0)
|
||||||
|
utils::logmesg(lmp,"Using type {} for bond particles\n",bptype);
|
||||||
|
|
||||||
|
// set bond and bond particle types in fix srp
|
||||||
|
// bonds of this type will be represented by bond particles
|
||||||
|
// if bond type is 0, then all bonds have bond particles
|
||||||
|
// btype = bond type
|
||||||
|
|
||||||
|
char c0[20];
|
||||||
|
char* arg0[2];
|
||||||
|
sprintf(c0, "%d", btype);
|
||||||
|
arg0[0] = (char *) "btype";
|
||||||
|
arg0[1] = c0;
|
||||||
|
f_srp->modify_params(2, arg0);
|
||||||
|
|
||||||
|
// bptype = bond particle type
|
||||||
|
sprintf(c0, "%d", bptype);
|
||||||
|
arg0[0] = (char *) "bptype";
|
||||||
|
arg0[1] = c0;
|
||||||
|
f_srp->modify_params(2, arg0);
|
||||||
|
|
||||||
|
// if using fix bond/break, set id of fix bond/break in fix srp
|
||||||
|
// idbreak = id of fix bond break
|
||||||
|
if( bond_break )
|
||||||
|
{
|
||||||
|
sprintf(c0, "%s", idbreak);
|
||||||
|
arg0[0] = (char *) "bond/break";
|
||||||
|
arg0[1] = c0;
|
||||||
|
f_srp->modify_params(2, arg0);
|
||||||
|
delete [] idbreak;
|
||||||
|
}
|
||||||
|
|
||||||
|
// if using fix bond/create, set id of fix bond/create in fix srp
|
||||||
|
// idcreate = id of fix bond break
|
||||||
|
if( bond_create )
|
||||||
|
{
|
||||||
|
sprintf(c0, "%s", idcreate);
|
||||||
|
arg0[0] = (char *) "bond/create";
|
||||||
|
arg0[1] = c0;
|
||||||
|
f_srp->modify_params(2, arg0);
|
||||||
|
delete [] idcreate;
|
||||||
|
}
|
||||||
|
|
||||||
|
// bond particles do not contribute to energy or virial
|
||||||
|
// bond particles do not belong to group all
|
||||||
|
// but thermo normalization is by nall
|
||||||
|
// therefore should turn off normalization
|
||||||
|
char *arg1[2];
|
||||||
|
arg1[0] = (char *) "norm";
|
||||||
|
arg1[1] = (char *) "no";
|
||||||
|
output->thermo->modify_params(2, arg1);
|
||||||
|
if (comm->me == 0) error->message(FLERR,"Thermo normalization turned off by pair srp");
|
||||||
|
|
||||||
|
neighbor->request(this,instance_me);
|
||||||
|
}
|
||||||
46
src/MISC/pair_srp_react.h
Normal file
46
src/MISC/pair_srp_react.h
Normal file
@ -0,0 +1,46 @@
|
|||||||
|
/* -*- c++ -*- ----------------------------------------------------------
|
||||||
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
||||||
|
https://www.lammps.org/, 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.
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
#ifdef PAIR_CLASS
|
||||||
|
// clang-format off
|
||||||
|
PairStyle(srp/react,PairSRPREACT);
|
||||||
|
// clang-format on
|
||||||
|
#else
|
||||||
|
|
||||||
|
#ifndef LMP_PAIR_SRP_REACT_H
|
||||||
|
#define LMP_PAIR_SRP_REACT_H
|
||||||
|
|
||||||
|
#include "pair_srp.h"
|
||||||
|
|
||||||
|
namespace LAMMPS_NS {
|
||||||
|
|
||||||
|
class PairSRPREACT : public PairSRP {
|
||||||
|
public:
|
||||||
|
PairSRPREACT(class LAMMPS *);
|
||||||
|
void settings(int, char **) override;
|
||||||
|
void init_style() override;
|
||||||
|
|
||||||
|
private:
|
||||||
|
char *idbreak;
|
||||||
|
char *idcreate;
|
||||||
|
bool bond_create, bond_break;
|
||||||
|
};
|
||||||
|
|
||||||
|
} // namespace LAMMPS_NS
|
||||||
|
|
||||||
|
#endif
|
||||||
|
#endif
|
||||||
|
|
||||||
|
/* ERROR/WARNING messages:
|
||||||
|
|
||||||
|
*/
|
||||||
Reference in New Issue
Block a user