Merge branch 'master' into write-bonus-data

This commit is contained in:
Axel Kohlmeyer 2020-07-16 12:15:59 -04:00
commit 35b030015d
No known key found for this signature in database
GPG Key ID: D9B44E93BF0C375A
8 changed files with 249 additions and 10 deletions

View File

@ -42,6 +42,7 @@ OPT.
* :doc:`bocs <fix_bocs>`
* :doc:`bond/break <fix_bond_break>`
* :doc:`bond/create <fix_bond_create>`
* :doc:`bond/create/angle <fix_bond_create>`
* :doc:`bond/react <fix_bond_react>`
* :doc:`bond/swap <fix_bond_swap>`
* :doc:`box/relax <fix_box_relax>`

View File

@ -185,6 +185,7 @@ accelerated styles exist.
* :doc:`bocs <fix_bocs>` - NPT style time integration with pressure correction
* :doc:`bond/break <fix_bond_break>` - break bonds on the fly
* :doc:`bond/create <fix_bond_create>` - create bonds on the fly
* :doc:`bond/create/angle <fix_bond_create>` - create bonds on the fly with angle constraints
* :doc:`bond/react <fix_bond_react>` - apply topology changes to model reactions
* :doc:`bond/swap <fix_bond_swap>` - Monte Carlo bond swapping
* :doc:`box/relax <fix_box_relax>` - relax box size during energy minimization

View File

@ -3,6 +3,9 @@
fix bond/create command
=======================
fix bond/create/angle command
=============================
Syntax
""""""
@ -17,7 +20,7 @@ Syntax
* Rmin = 2 atoms separated by less than Rmin can bond (distance units)
* bondtype = type of created bonds
* zero or more keyword/value pairs may be appended to args
* keyword = *iparam* or *jparam* or *prob* or *atype* or *dtype* or *itype*
* keyword = *iparam* or *jparam* or *prob* or *atype* or *dtype* or *itype* or *aconstrain*
.. parsed-literal::
@ -36,6 +39,9 @@ Syntax
dihedraltype = type of created dihedrals
*itype* value = impropertype
impropertype = type of created impropers
*aconstrain* value = amin amax
amin = minimal angle at which new bonds can be created
amax = maximal angle at which new bonds can be created
Examples
""""""""
@ -45,6 +51,7 @@ Examples
fix 5 all bond/create 10 1 2 0.8 1
fix 5 all bond/create 1 3 3 0.8 1 prob 0.5 85784 iparam 2 3
fix 5 all bond/create 1 3 3 0.8 1 prob 0.5 85784 iparam 2 3 atype 1 dtype 2
fix 5 all bond/create/angle 10 1 2 1.122 1 aconstrain 120 180 prob 1 4928459 iparam 2 1 jparam 2 2
Description
"""""""""""
@ -110,7 +117,16 @@ actually created. The *fraction* setting must be a value between 0.0
and 1.0. A uniform random number between 0.0 and 1.0 is generated and
the eligible bond is only created if the random number < fraction.
Any bond that is created is assigned a bond type of *bondtype*
The *aconstrain* keyword is only available with the fix
bond/create/angle command. It allows to specify a minimal and maximal
angle *amin* and *amax* between the two prospective bonding partners and
a third particle that is already bonded to one of the two partners.
Such a criterion can be important when new angles are defined together
with the formation of a new bond. Without a restriction on the
permissible angle, and for stiffer angle potentials, very large energies
can arise and lead to uncontrolled behavior.
Any bond that is created is assigned a bond type of *bondtype*.
When a bond is created, data structures within LAMMPS that store bond
topology are updated to reflect the creation. If the bond is part of
@ -218,12 +234,14 @@ You can dump out snapshots of the current bond topology via the :doc:`dump local
**Restart, fix_modify, output, run start/stop, minimize info:**
No information about this fix is written to :doc:`binary restart files <restart>`. None of the :doc:`fix_modify <fix_modify>` options
are relevant to this fix.
No information about this fix is written to :doc:`binary restart files
<restart>`. None of the :doc:`fix_modify <fix_modify>` options are
relevant to this fix.
This fix computes two statistics which it stores in a global vector of
length 2, which can be accessed by various :doc:`output commands <Howto_output>`. The vector values calculated by this fix
are "intensive".
length 2, which can be accessed by various :doc:`output commands
<Howto_output>`. The vector values calculated by this fix are
"intensive".
These are the 2 quantities:

View File

@ -27,9 +27,11 @@
#include "random_mars.h"
#include "memory.h"
#include "error.h"
#include "math_const.h"
using namespace LAMMPS_NS;
using namespace FixConst;
using namespace MathConst;
#define BIG 1.0e20
#define DELTA 16
@ -79,6 +81,11 @@ FixBondCreate::FixBondCreate(LAMMPS *lmp, int narg, char **arg) :
int seed = 12345;
atype = dtype = itype = 0;
constrainflag = 0;
constrainpass = 0;
amin = 0;
amax = 180;
int iarg = 8;
while (iarg < narg) {
if (strcmp(arg[iarg],"iparam") == 0) {
@ -120,6 +127,22 @@ FixBondCreate::FixBondCreate(LAMMPS *lmp, int narg, char **arg) :
itype = force->inumeric(FLERR,arg[iarg+1]);
if (itype < 0) error->all(FLERR,"Illegal fix bond/create command");
iarg += 2;
} else if (strcmp(arg[iarg],"aconstrain") == 0 &&
strcmp(style,"bond/create/angle") == 0) {
if (iarg+3 > narg)
error->all(FLERR,"Illegal fix bond/create/angle command");
amin = force->numeric(FLERR,arg[iarg+1]);
amax = force->inumeric(FLERR,arg[iarg+2]);
if (amin >= amax)
error->all(FLERR,"Illegal fix bond/create/angle command");
if (amin < 0 || amin > 180)
error->all(FLERR,"Illegal fix bond/create/angle command");
if (amax < 0 || amax > 180)
error->all(FLERR,"Illegal fix bond/create/angle command");
amin = (MY_PI/180.0) * amin;
amax = (MY_PI/180.0) * amax;
constrainflag = 1;
iarg += 3;
} else error->all(FLERR,"Illegal fix bond/create command");
}
@ -434,6 +457,11 @@ void FixBondCreate::post_integrate()
rsq = delx*delx + dely*dely + delz*delz;
if (rsq >= cutsq) continue;
if (constrainflag) {
constrainpass = constrain(i,j,amin,amax);
if (!constrainpass) continue;
}
if (rsq < distsq[i]) {
partner[i] = tag[j];
distsq[i] = rsq;
@ -442,6 +470,7 @@ void FixBondCreate::post_integrate()
partner[j] = tag[i];
distsq[j] = rsq;
}
}
}

View File

@ -27,7 +27,7 @@ namespace LAMMPS_NS {
class FixBondCreate : public Fix {
public:
FixBondCreate(class LAMMPS *, int, char **);
~FixBondCreate();
virtual ~FixBondCreate();
int setmask();
void init();
void init_list(int, class NeighList *);
@ -46,15 +46,18 @@ class FixBondCreate : public Fix {
double compute_vector(int);
double memory_usage();
private:
protected:
int me;
int iatomtype,jatomtype;
int btype,seed;
int imaxbond,jmaxbond;
int inewtype,jnewtype;
int constrainflag,constrainpass;
double amin,amax;
double cutsq,fraction;
int atype,dtype,itype;
int angleflag,dihedralflag,improperflag;
int overflow;
tagint lastcheck;
@ -84,6 +87,8 @@ class FixBondCreate : public Fix {
void create_impropers(int);
int dedup(int, int, tagint *);
virtual int constrain(int, int, double, double) {return 1;}
// DEBUG
void print_bb();

View File

@ -0,0 +1,147 @@
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://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 authors:
Joao Gregorio, Oliver Henrich (University of Strathclyde, Glasgow, UK)
------------------------------------------------------------------------- */
#include "fix_bond_create_angle.h"
#include "atom.h"
#include <cmath>
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
FixBondCreateAngle::FixBondCreateAngle(LAMMPS *lmp, int narg, char **arg) :
FixBondCreate(lmp, narg, arg)
{
}
/* ---------------------------------------------------------------------- */
int FixBondCreateAngle::constrain(int i, int j, double amin, double amax)
{
double **x = atom->x;
tagint *tag = atom->tag;
tagint **bond_atom = atom->bond_atom;
int *num_bond = atom->num_bond;
int **nspecial = atom->nspecial;
tagint **special = atom->special;
int *mask = atom->mask;
int *type = atom->type;
double v1x = 0.0;
double v1y = 0.0;
double v1z = 0.0;
double v2x = 0.0;
double v2y = 0.0;
double v2z = 0.0;
double angle1,angle2;
int flag = 0;
// pass if both atoms have no neighbors: bond is always formed
if (nspecial[i][0] == 0 && nspecial[j][0] == 0) flag = 1;
// pass if i has at least one neighbor and angle constraint is met
if (nspecial[i][0] != 0 && nspecial[j][0] == 0) {
// calculate first vector between i and j
// calculate second vector between i and its first neighbor
v1x = x[i][0] - x[j][0];
v1y = x[i][1] - x[j][1];
v1z = x[i][2] - x[j][2];
v2x = x[i][0] - x[atom->map(special[i][0])][0];
v2y = x[i][1] - x[atom->map(special[i][0])][1];
v2z = x[i][2] - x[atom->map(special[i][0])][2];
// calculate angle between both vectors
// set flag if the angle constraint is met
angle1 = acos((v1x*v2x + v1y*v2y + v1z*v2z)/
(sqrt(v1x*v1x + v1y*v1y + v1z*v1z)*sqrt(v2x*v2x + v2y*v2y + v2z*v2z)));
if (amin <= angle1 && angle1 <= amax) flag = 1;
}
// pass if j has at least one neighbor and angle constraint is met
if (nspecial[j][0] != 0 && nspecial[i][0] == 0) {
// calculate first vector between j and i
// calculate second vector between j and its first neighbor
v1x = x[j][0] - x[i][0];
v1y = x[j][1] - x[i][1];
v1z = x[j][2] - x[i][2];
v2x = x[j][0] - x[atom->map(special[j][0])][0];
v2y = x[j][1] - x[atom->map(special[j][0])][1];
v2z = x[j][2] - x[atom->map(special[j][0])][2];
// calculate angle between both vectors
// set flag if angle constraint is met
angle1 = acos((v1x*v2x + v1y*v2y + v1z*v2z) /
(sqrt(v1x*v1x + v1y*v1y + v1z*v1z)*sqrt(v2x*v2x + v2y*v2y + v2z*v2z)));
if (amin <= angle1 && angle1 <= amax) flag = 1;
}
// pass if both i and j have at least one neighbor and angle constraint is met
if (nspecial[i][0] != 0 && nspecial[j][0] != 0) {
// Calculate 1st angle
// calculate first vector between i and j
// calculate second vector between i and its first neighbor
v1x = x[i][0] - x[j][0];
v1y = x[i][1] - x[j][1];
v1z = x[i][2] - x[j][2];
v2x = x[i][0] - x[atom->map(special[i][0])][0];
v2y = x[i][1] - x[atom->map(special[i][0])][1];
v2z = x[i][2] - x[atom->map(special[i][0])][2];
// calculate angle between both vectors
angle1 = acos((v1x*v2x + v1y*v2y + v1z*v2z) /
(sqrt(v1x*v1x + v1y*v1y + v1z*v1z)*sqrt(v2x*v2x + v2y*v2y + v2z*v2z)));
// Calculate 2nd angle
// calculate first vector between i and j
// calculate second vector between i and its first neighbor
v1x = x[j][0] - x[i][0];
v1y = x[j][1] - x[i][1];
v1z = x[j][2] - x[i][2];
v2x = x[j][0] - x[atom->map(special[j][0])][0];
v2y = x[j][1] - x[atom->map(special[j][0])][1];
v2z = x[j][2] - x[atom->map(special[j][0])][2];
// calculate angle between both vectors
angle2 = acos((v1x*v2x + v1y*v2y + v1z*v2z) /
(sqrt(v1x*v1x + v1y*v1y + v1z*v1z)*sqrt(v2x*v2x + v2y*v2y + v2z*v2z)));
// set flag if the angle constraint is met
if ((amin <= angle1 && angle1 <= amax) && (amin <= angle2 && angle2 <= amax))
flag = 1;
}
return flag;
}

View File

@ -0,0 +1,38 @@
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://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.
------------------------------------------------------------------------- */
#ifdef FIX_CLASS
FixStyle(bond/create/angle,FixBondCreateAngle)
#else
#ifndef LMP_FIX_BOND_CREATE_ANGLE_H
#define LMP_FIX_BOND_CREATE_ANGLE_H
#include "fix_bond_create.h"
namespace LAMMPS_NS {
class FixBondCreateAngle : public FixBondCreate {
public:
FixBondCreateAngle(class LAMMPS *, int, char **);
private:
int constrain(int, int, double, double);
};
}
#endif
#endif

View File

@ -7,9 +7,9 @@ SHELL = /bin/sh
ROOT = lmp
EXE = $(ROOT)_$@
SRC = angle.cpp angle_charmm.cpp angle_cosine.cpp angle_cosine_delta.cpp angle_cosine_squared.cpp angle_harmonic.cpp angle_hybrid.cpp angle_table.cpp atom.cpp atom_vec.cpp atom_vec_angle.cpp atom_vec_atomic.cpp atom_vec_bond.cpp atom_vec_charge.cpp atom_vec_full.cpp atom_vec_hybrid.cpp atom_vec_molecular.cpp bond.cpp bond_fene.cpp bond_fene_expand.cpp bond_harmonic.cpp bond_hybrid.cpp bond_morse.cpp bond_nonlinear.cpp bond_quartic.cpp bond_table.cpp change_box.cpp comm.cpp compute.cpp compute_angle_local.cpp compute_bond_local.cpp compute_centro_atom.cpp compute_cna_atom.cpp compute_com.cpp compute_com_molecule.cpp compute_coord_atom.cpp compute_dihedral_local.cpp compute_displace_atom.cpp compute_erotate_sphere.cpp compute_group_group.cpp compute_gyration.cpp compute_gyration_molecule.cpp compute_heat_flux.cpp compute_improper_local.cpp compute_ke.cpp compute_ke_atom.cpp compute_msd.cpp compute_msd_molecule.cpp compute_pair_local.cpp compute_pe.cpp compute_pe_atom.cpp compute_pressure.cpp compute_property_atom.cpp compute_property_local.cpp compute_property_molecule.cpp compute_rdf.cpp compute_reduce.cpp compute_reduce_region.cpp compute_stress_atom.cpp compute_temp.cpp compute_temp_com.cpp compute_temp_deform.cpp compute_temp_partial.cpp compute_temp_profile.cpp compute_temp_ramp.cpp compute_temp_region.cpp compute_temp_sphere.cpp create_atoms.cpp create_box.cpp delete_atoms.cpp delete_bonds.cpp dihedral.cpp dihedral_charmm.cpp dihedral_harmonic.cpp dihedral_helix.cpp dihedral_hybrid.cpp dihedral_multi_harmonic.cpp dihedral_opls.cpp displace_atoms.cpp displace_box.cpp domain.cpp dump.cpp dump_atom.cpp dump_cfg.cpp dump_custom.cpp dump_dcd.cpp dump_local.cpp dump_xyz.cpp error.cpp ewald.cpp fft3d.cpp fft3d_wrap.cpp finish.cpp fix.cpp fix_adapt.cpp fix_addforce.cpp fix_ave_atom.cpp fix_ave_correlate.cpp fix_ave_histo.cpp fix_ave_spatial.cpp fix_ave_time.cpp fix_aveforce.cpp fix_bond_break.cpp fix_bond_create.cpp fix_bond_swap.cpp fix_box_relax.cpp fix_deform.cpp fix_deposit.cpp fix_drag.cpp fix_dt_reset.cpp fix_efield.cpp fix_enforce2d.cpp fix_evaporate.cpp fix_gravity.cpp fix_heat.cpp fix_indent.cpp fix_langevin.cpp fix_lineforce.cpp fix_minimize.cpp fix_momentum.cpp fix_move.cpp fix_nh.cpp fix_nh_sphere.cpp fix_nph.cpp fix_nph_sphere.cpp fix_npt.cpp fix_npt_sphere.cpp fix_nve.cpp fix_nve_limit.cpp fix_nve_noforce.cpp fix_nve_sphere.cpp fix_nvt.cpp fix_nvt_sllod.cpp fix_nvt_sphere.cpp fix_orient_fcc.cpp fix_planeforce.cpp fix_press_berendsen.cpp fix_print.cpp fix_qeq_comb.cpp fix_recenter.cpp fix_respa.cpp fix_rigid.cpp fix_rigid_nve.cpp fix_rigid_nvt.cpp fix_setforce.cpp fix_shake.cpp fix_shear_history.cpp fix_spring.cpp fix_spring_rg.cpp fix_spring_self.cpp fix_store_force.cpp fix_store_state.cpp fix_temp_berendsen.cpp fix_temp_rescale.cpp fix_thermal_conductivity.cpp fix_tmd.cpp fix_ttm.cpp fix_viscosity.cpp fix_viscous.cpp fix_wall.cpp fix_wall_harmonic.cpp fix_wall_lj126.cpp fix_wall_lj93.cpp fix_wall_reflect.cpp fix_wall_region.cpp force.cpp group.cpp improper.cpp improper_cvff.cpp improper_harmonic.cpp improper_hybrid.cpp input.cpp integrate.cpp kspace.cpp lammps.cpp lattice.cpp library.cpp main.cpp memory.cpp min.cpp min_cg.cpp min_hftn.cpp min_linesearch.cpp min_sd.cpp minimize.cpp modify.cpp neigh_bond.cpp neigh_derive.cpp neigh_full.cpp neigh_gran.cpp neigh_half_bin.cpp neigh_half_multi.cpp neigh_half_nsq.cpp neigh_list.cpp neigh_request.cpp neigh_respa.cpp neigh_stencil.cpp neighbor.cpp output.cpp pack.cpp pair.cpp pair_airebo.cpp pair_born_coul_long.cpp pair_buck.cpp pair_buck_coul_cut.cpp pair_buck_coul_long.cpp pair_comb.cpp pair_coul_cut.cpp pair_coul_debye.cpp pair_coul_long.cpp pair_dpd.cpp pair_dpd_tstat.cpp pair_eam.cpp pair_eam_alloy.cpp pair_eam_fs.cpp pair_eim.cpp pair_hybrid.cpp pair_hybrid_overlay.cpp pair_lj96_cut.cpp pair_lj_charmm_coul_charmm.cpp pair_lj_charmm_coul_charmm_implicit.cpp pair_lj_charmm_coul_long.cpp pair_lj_cut.cpp pair_lj_cut_coul_cut.cpp pair_lj_cut_coul_debye.cpp pair_lj_cut_coul_long.cpp pair_lj_cut_coul_long_tip4p.cpp pair_lj_expand.cpp pair_lj_gromacs.cpp pair_lj_gromacs_coul_gromacs.cpp pair_lj_smooth.cpp pair_morse.cpp pair_soft.cpp pair_sw.cpp pair_table.cpp pair_tersoff.cpp pair_tersoff_zbl.cpp pair_yukawa.cpp pppm.cpp pppm_tip4p.cpp random_mars.cpp random_park.cpp read_data.cpp read_restart.cpp region.cpp region_block.cpp region_cone.cpp region_cylinder.cpp region_intersect.cpp region_plane.cpp region_prism.cpp region_sphere.cpp region_union.cpp remap.cpp remap_wrap.cpp replicate.cpp respa.cpp run.cpp set.cpp shell.cpp special.cpp temper.cpp thermo.cpp timer.cpp universe.cpp update.cpp variable.cpp velocity.cpp verlet.cpp write_restart.cpp
SRC = angle.cpp angle_charmm.cpp angle_cosine.cpp angle_cosine_delta.cpp angle_cosine_squared.cpp angle_harmonic.cpp angle_hybrid.cpp angle_table.cpp atom.cpp atom_vec.cpp atom_vec_angle.cpp atom_vec_atomic.cpp atom_vec_bond.cpp atom_vec_charge.cpp atom_vec_full.cpp atom_vec_hybrid.cpp atom_vec_molecular.cpp bond.cpp bond_fene.cpp bond_fene_expand.cpp bond_harmonic.cpp bond_hybrid.cpp bond_morse.cpp bond_nonlinear.cpp bond_quartic.cpp bond_table.cpp change_box.cpp comm.cpp compute.cpp compute_angle_local.cpp compute_bond_local.cpp compute_centro_atom.cpp compute_cna_atom.cpp compute_com.cpp compute_com_molecule.cpp compute_coord_atom.cpp compute_dihedral_local.cpp compute_displace_atom.cpp compute_erotate_sphere.cpp compute_group_group.cpp compute_gyration.cpp compute_gyration_molecule.cpp compute_heat_flux.cpp compute_improper_local.cpp compute_ke.cpp compute_ke_atom.cpp compute_msd.cpp compute_msd_molecule.cpp compute_pair_local.cpp compute_pe.cpp compute_pe_atom.cpp compute_pressure.cpp compute_property_atom.cpp compute_property_local.cpp compute_property_molecule.cpp compute_rdf.cpp compute_reduce.cpp compute_reduce_region.cpp compute_stress_atom.cpp compute_temp.cpp compute_temp_com.cpp compute_temp_deform.cpp compute_temp_partial.cpp compute_temp_profile.cpp compute_temp_ramp.cpp compute_temp_region.cpp compute_temp_sphere.cpp create_atoms.cpp create_box.cpp delete_atoms.cpp delete_bonds.cpp dihedral.cpp dihedral_charmm.cpp dihedral_harmonic.cpp dihedral_helix.cpp dihedral_hybrid.cpp dihedral_multi_harmonic.cpp dihedral_opls.cpp displace_atoms.cpp displace_box.cpp domain.cpp dump.cpp dump_atom.cpp dump_cfg.cpp dump_custom.cpp dump_dcd.cpp dump_local.cpp dump_xyz.cpp error.cpp ewald.cpp fft3d.cpp fft3d_wrap.cpp finish.cpp fix.cpp fix_adapt.cpp fix_addforce.cpp fix_ave_atom.cpp fix_ave_correlate.cpp fix_ave_histo.cpp fix_ave_spatial.cpp fix_ave_time.cpp fix_aveforce.cpp fix_bond_break.cpp fix_bond_create.cpp fix_bond_create_angle.cpp fix_bond_swap.cpp fix_box_relax.cpp fix_deform.cpp fix_deposit.cpp fix_drag.cpp fix_dt_reset.cpp fix_efield.cpp fix_enforce2d.cpp fix_evaporate.cpp fix_gravity.cpp fix_heat.cpp fix_indent.cpp fix_langevin.cpp fix_lineforce.cpp fix_minimize.cpp fix_momentum.cpp fix_move.cpp fix_nh.cpp fix_nh_sphere.cpp fix_nph.cpp fix_nph_sphere.cpp fix_npt.cpp fix_npt_sphere.cpp fix_nve.cpp fix_nve_limit.cpp fix_nve_noforce.cpp fix_nve_sphere.cpp fix_nvt.cpp fix_nvt_sllod.cpp fix_nvt_sphere.cpp fix_orient_fcc.cpp fix_planeforce.cpp fix_press_berendsen.cpp fix_print.cpp fix_qeq_comb.cpp fix_recenter.cpp fix_respa.cpp fix_rigid.cpp fix_rigid_nve.cpp fix_rigid_nvt.cpp fix_setforce.cpp fix_shake.cpp fix_shear_history.cpp fix_spring.cpp fix_spring_rg.cpp fix_spring_self.cpp fix_store_force.cpp fix_store_state.cpp fix_temp_berendsen.cpp fix_temp_rescale.cpp fix_thermal_conductivity.cpp fix_tmd.cpp fix_ttm.cpp fix_viscosity.cpp fix_viscous.cpp fix_wall.cpp fix_wall_harmonic.cpp fix_wall_lj126.cpp fix_wall_lj93.cpp fix_wall_reflect.cpp fix_wall_region.cpp force.cpp group.cpp improper.cpp improper_cvff.cpp improper_harmonic.cpp improper_hybrid.cpp input.cpp integrate.cpp kspace.cpp lammps.cpp lattice.cpp library.cpp main.cpp memory.cpp min.cpp min_cg.cpp min_hftn.cpp min_linesearch.cpp min_sd.cpp minimize.cpp modify.cpp neigh_bond.cpp neigh_derive.cpp neigh_full.cpp neigh_gran.cpp neigh_half_bin.cpp neigh_half_multi.cpp neigh_half_nsq.cpp neigh_list.cpp neigh_request.cpp neigh_respa.cpp neigh_stencil.cpp neighbor.cpp output.cpp pack.cpp pair.cpp pair_airebo.cpp pair_born_coul_long.cpp pair_buck.cpp pair_buck_coul_cut.cpp pair_buck_coul_long.cpp pair_comb.cpp pair_coul_cut.cpp pair_coul_debye.cpp pair_coul_long.cpp pair_dpd.cpp pair_dpd_tstat.cpp pair_eam.cpp pair_eam_alloy.cpp pair_eam_fs.cpp pair_eim.cpp pair_hybrid.cpp pair_hybrid_overlay.cpp pair_lj96_cut.cpp pair_lj_charmm_coul_charmm.cpp pair_lj_charmm_coul_charmm_implicit.cpp pair_lj_charmm_coul_long.cpp pair_lj_cut.cpp pair_lj_cut_coul_cut.cpp pair_lj_cut_coul_debye.cpp pair_lj_cut_coul_long.cpp pair_lj_cut_coul_long_tip4p.cpp pair_lj_expand.cpp pair_lj_gromacs.cpp pair_lj_gromacs_coul_gromacs.cpp pair_lj_smooth.cpp pair_morse.cpp pair_soft.cpp pair_sw.cpp pair_table.cpp pair_tersoff.cpp pair_tersoff_zbl.cpp pair_yukawa.cpp pppm.cpp pppm_tip4p.cpp random_mars.cpp random_park.cpp read_data.cpp read_restart.cpp region.cpp region_block.cpp region_cone.cpp region_cylinder.cpp region_intersect.cpp region_plane.cpp region_prism.cpp region_sphere.cpp region_union.cpp remap.cpp remap_wrap.cpp replicate.cpp respa.cpp run.cpp set.cpp shell.cpp special.cpp temper.cpp thermo.cpp timer.cpp universe.cpp update.cpp variable.cpp velocity.cpp verlet.cpp write_restart.cpp
INC = angle.h angle_charmm.h angle_cosine.h angle_cosine_delta.h angle_cosine_squared.h angle_harmonic.h angle_hybrid.h angle_table.h atom.h atom_vec.h atom_vec_angle.h atom_vec_atomic.h atom_vec_bond.h atom_vec_charge.h atom_vec_full.h atom_vec_hybrid.h atom_vec_molecular.h bond.h bond_fene.h bond_fene_expand.h bond_harmonic.h bond_hybrid.h bond_morse.h bond_nonlinear.h bond_quartic.h bond_table.h change_box.h comm.h compute.h compute_angle_local.h compute_bond_local.h compute_centro_atom.h compute_cna_atom.h compute_com.h compute_com_molecule.h compute_coord_atom.h compute_dihedral_local.h compute_displace_atom.h compute_erotate_sphere.h compute_group_group.h compute_gyration.h compute_gyration_molecule.h compute_heat_flux.h compute_improper_local.h compute_ke.h compute_ke_atom.h compute_msd.h compute_msd_molecule.h compute_pair_local.h compute_pe.h compute_pe_atom.h compute_pressure.h compute_property_atom.h compute_property_local.h compute_property_molecule.h compute_rdf.h compute_reduce.h compute_reduce_region.h compute_stress_atom.h compute_temp.h compute_temp_com.h compute_temp_deform.h compute_temp_partial.h compute_temp_profile.h compute_temp_ramp.h compute_temp_region.h compute_temp_sphere.h create_atoms.h create_box.h delete_atoms.h delete_bonds.h dihedral.h dihedral_charmm.h dihedral_harmonic.h dihedral_helix.h dihedral_hybrid.h dihedral_multi_harmonic.h dihedral_opls.h displace_atoms.h displace_box.h domain.h dump.h dump_atom.h dump_cfg.h dump_custom.h dump_dcd.h dump_local.h dump_xyz.h error.h ewald.h fft3d.h fft3d_wrap.h finish.h fix.h fix_adapt.h fix_addforce.h fix_ave_atom.h fix_ave_correlate.h fix_ave_histo.h fix_ave_spatial.h fix_ave_time.h fix_aveforce.h fix_bond_break.h fix_bond_create.h fix_bond_swap.h fix_box_relax.h fix_deform.h fix_deposit.h fix_drag.h fix_dt_reset.h fix_efield.h fix_enforce2d.h fix_evaporate.h fix_gravity.h fix_heat.h fix_indent.h fix_langevin.h fix_lineforce.h fix_minimize.h fix_momentum.h fix_move.h fix_nh.h fix_nh_sphere.h fix_nph.h fix_nph_sphere.h fix_npt.h fix_npt_sphere.h fix_nve.h fix_nve_limit.h fix_nve_noforce.h fix_nve_sphere.h fix_nvt.h fix_nvt_sllod.h fix_nvt_sphere.h fix_orient_fcc.h fix_planeforce.h fix_press_berendsen.h fix_print.h fix_qeq_comb.h fix_recenter.h fix_respa.h fix_rigid.h fix_rigid_nve.h fix_rigid_nvt.h fix_setforce.h fix_shake.h fix_shear_history.h fix_spring.h fix_spring_rg.h fix_spring_self.h fix_store_force.h fix_store_state.h fix_temp_berendsen.h fix_temp_rescale.h fix_thermal_conductivity.h fix_tmd.h fix_ttm.h fix_viscosity.h fix_viscous.h fix_wall.h fix_wall_harmonic.h fix_wall_lj126.h fix_wall_lj93.h fix_wall_reflect.h fix_wall_region.h force.h group.h improper.h improper_cvff.h improper_harmonic.h improper_hybrid.h input.h integrate.h kspace.h lammps.h lattice.h library.h math_extra.h memory.h min.h min_cg.h min_hftn.h min_linesearch.h min_sd.h minimize.h modify.h neigh_list.h neigh_request.h neighbor.h output.h pack.h pair.h pair_airebo.h pair_born_coul_long.h pair_buck.h pair_buck_coul_cut.h pair_buck_coul_long.h pair_comb.h pair_coul_cut.h pair_coul_debye.h pair_coul_long.h pair_dpd.h pair_dpd_tstat.h pair_eam.h pair_eam_alloy.h pair_eam_fs.h pair_eim.h pair_hybrid.h pair_hybrid_overlay.h pair_lj96_cut.h pair_lj_charmm_coul_charmm.h pair_lj_charmm_coul_charmm_implicit.h pair_lj_charmm_coul_long.h pair_lj_cut.h pair_lj_cut_coul_cut.h pair_lj_cut_coul_debye.h pair_lj_cut_coul_long.h pair_lj_cut_coul_long_tip4p.h pair_lj_expand.h pair_lj_gromacs.h pair_lj_gromacs_coul_gromacs.h pair_lj_smooth.h pair_morse.h pair_soft.h pair_sw.h pair_table.h pair_tersoff.h pair_tersoff_zbl.h pair_yukawa.h pointers.h pppm.h pppm_tip4p.h random_mars.h random_park.h read_data.h read_restart.h region.h region_block.h region_cone.h region_cylinder.h region_intersect.h region_plane.h region_prism.h region_sphere.h region_union.h remap.h remap_wrap.h replicate.h respa.h run.h set.h shell.h special.h style_angle.h style_atom.h style_bond.h style_command.h style_compute.h style_dihedral.h style_dump.h style_fix.h style_improper.h style_integrate.h style_kspace.h style_minimize.h style_pair.h style_region.h temper.h thermo.h timer.h universe.h update.h variable.h velocity.h verlet.h version.h write_restart.h
INC = angle.h angle_charmm.h angle_cosine.h angle_cosine_delta.h angle_cosine_squared.h angle_harmonic.h angle_hybrid.h angle_table.h atom.h atom_vec.h atom_vec_angle.h atom_vec_atomic.h atom_vec_bond.h atom_vec_charge.h atom_vec_full.h atom_vec_hybrid.h atom_vec_molecular.h bond.h bond_fene.h bond_fene_expand.h bond_harmonic.h bond_hybrid.h bond_morse.h bond_nonlinear.h bond_quartic.h bond_table.h change_box.h comm.h compute.h compute_angle_local.h compute_bond_local.h compute_centro_atom.h compute_cna_atom.h compute_com.h compute_com_molecule.h compute_coord_atom.h compute_dihedral_local.h compute_displace_atom.h compute_erotate_sphere.h compute_group_group.h compute_gyration.h compute_gyration_molecule.h compute_heat_flux.h compute_improper_local.h compute_ke.h compute_ke_atom.h compute_msd.h compute_msd_molecule.h compute_pair_local.h compute_pe.h compute_pe_atom.h compute_pressure.h compute_property_atom.h compute_property_local.h compute_property_molecule.h compute_rdf.h compute_reduce.h compute_reduce_region.h compute_stress_atom.h compute_temp.h compute_temp_com.h compute_temp_deform.h compute_temp_partial.h compute_temp_profile.h compute_temp_ramp.h compute_temp_region.h compute_temp_sphere.h create_atoms.h create_box.h delete_atoms.h delete_bonds.h dihedral.h dihedral_charmm.h dihedral_harmonic.h dihedral_helix.h dihedral_hybrid.h dihedral_multi_harmonic.h dihedral_opls.h displace_atoms.h displace_box.h domain.h dump.h dump_atom.h dump_cfg.h dump_custom.h dump_dcd.h dump_local.h dump_xyz.h error.h ewald.h fft3d.h fft3d_wrap.h finish.h fix.h fix_adapt.h fix_addforce.h fix_ave_atom.h fix_ave_correlate.h fix_ave_histo.h fix_ave_spatial.h fix_ave_time.h fix_aveforce.h fix_bond_break.h fix_bond_create.h fix_bond_create_angle.h fix_bond_swap.h fix_box_relax.h fix_deform.h fix_deposit.h fix_drag.h fix_dt_reset.h fix_efield.h fix_enforce2d.h fix_evaporate.h fix_gravity.h fix_heat.h fix_indent.h fix_langevin.h fix_lineforce.h fix_minimize.h fix_momentum.h fix_move.h fix_nh.h fix_nh_sphere.h fix_nph.h fix_nph_sphere.h fix_npt.h fix_npt_sphere.h fix_nve.h fix_nve_limit.h fix_nve_noforce.h fix_nve_sphere.h fix_nvt.h fix_nvt_sllod.h fix_nvt_sphere.h fix_orient_fcc.h fix_planeforce.h fix_press_berendsen.h fix_print.h fix_qeq_comb.h fix_recenter.h fix_respa.h fix_rigid.h fix_rigid_nve.h fix_rigid_nvt.h fix_setforce.h fix_shake.h fix_shear_history.h fix_spring.h fix_spring_rg.h fix_spring_self.h fix_store_force.h fix_store_state.h fix_temp_berendsen.h fix_temp_rescale.h fix_thermal_conductivity.h fix_tmd.h fix_ttm.h fix_viscosity.h fix_viscous.h fix_wall.h fix_wall_harmonic.h fix_wall_lj126.h fix_wall_lj93.h fix_wall_reflect.h fix_wall_region.h force.h group.h improper.h improper_cvff.h improper_harmonic.h improper_hybrid.h input.h integrate.h kspace.h lammps.h lattice.h library.h math_extra.h memory.h min.h min_cg.h min_hftn.h min_linesearch.h min_sd.h minimize.h modify.h neigh_list.h neigh_request.h neighbor.h output.h pack.h pair.h pair_airebo.h pair_born_coul_long.h pair_buck.h pair_buck_coul_cut.h pair_buck_coul_long.h pair_comb.h pair_coul_cut.h pair_coul_debye.h pair_coul_long.h pair_dpd.h pair_dpd_tstat.h pair_eam.h pair_eam_alloy.h pair_eam_fs.h pair_eim.h pair_hybrid.h pair_hybrid_overlay.h pair_lj96_cut.h pair_lj_charmm_coul_charmm.h pair_lj_charmm_coul_charmm_implicit.h pair_lj_charmm_coul_long.h pair_lj_cut.h pair_lj_cut_coul_cut.h pair_lj_cut_coul_debye.h pair_lj_cut_coul_long.h pair_lj_cut_coul_long_tip4p.h pair_lj_expand.h pair_lj_gromacs.h pair_lj_gromacs_coul_gromacs.h pair_lj_smooth.h pair_morse.h pair_soft.h pair_sw.h pair_table.h pair_tersoff.h pair_tersoff_zbl.h pair_yukawa.h pointers.h pppm.h pppm_tip4p.h random_mars.h random_park.h read_data.h read_restart.h region.h region_block.h region_cone.h region_cylinder.h region_intersect.h region_plane.h region_prism.h region_sphere.h region_union.h remap.h remap_wrap.h replicate.h respa.h run.h set.h shell.h special.h style_angle.h style_atom.h style_bond.h style_command.h style_compute.h style_dihedral.h style_dump.h style_fix.h style_improper.h style_integrate.h style_kspace.h style_minimize.h style_pair.h style_region.h temper.h thermo.h timer.h universe.h update.h variable.h velocity.h verlet.h version.h write_restart.h
OBJ = $(SRC:.cpp=.o)