Commit 29cd4eb5 authored by Andrew Schultz's avatar Andrew Schultz
Browse files

Add HMA compute

parent 3367a408
Loading
Loading
Loading
Loading
+4 −67
Original line number Diff line number Diff line
@@ -173,65 +173,43 @@ There are also additional accelerated compute styles included in the
LAMMPS distribution for faster performance on CPUs, GPUs, and KNLs.
The individual style names on the "Commands
compute"_Commands_compute.html doc page are followed by one or more of
(g,i,k,o,t) to indicate which accelerated styles exist.
(g,i,k,o,t) to indicate which accerlerated styles exist.

"ackland/atom"_compute_ackland_atom.html - 
"aggregate/atom"_compute_cluster_atom.html - aggregate ID for each atom
"angle"_compute_angle.html - 
"angle/local"_compute_angle_local.html - 
"angle/local"_compute_bond_local.html - theta and energy of each angle
"angmom/chunk"_compute_angmom_chunk.html - angular momentum for each chunk
"basal/atom"_compute_basal_atom.html - 
"body/local"_compute_body_local.html - attributes of body sub-particles
"bond"_compute_bond.html - values computed by a bond style
"bond/local"_compute_bond_local.html - distance and energy of each bond
"centro/atom"_compute_centro_atom.html - centro-symmetry parameter for each atom
"chunk/atom"_compute_chunk_atom.html - assign chunk IDs to each atom
"chunk/spread/atom"_compute_chunk_spread_atom.html - spreads chunk values to each atom in chunk
"cluster/atom"_compute_cluster_atom.html - cluster ID for each atom
"cna/atom"_compute_cna_atom.html - common neighbor analysis (CNA) for each atom
"cnp/atom"_compute_cnp_atom.html - 
"com"_compute_com.html - center-of-mass of group of atoms
"com/chunk"_compute_com_chunk.html - center-of-mass for each chunk
"contact/atom"_compute_contact_atom.html - contact count for each spherical particle
"coord/atom"_compute_coord_atom.html - coordination number for each atom
"damage/atom"_compute_damage_atom.html - Peridynamic damage for each atom
"dihedral"_compute_dihedral.html - 
"dihedral/local"_compute_dihedral_local.html - angle of each dihedral
"dilatation/atom"_compute_dilatation_atom.html - Peridynamic dilatation for each atom
"dipole/chunk"_compute_dipole_chunk.html - 
"displace/atom"_compute_displace_atom.html - displacement of each atom
"dpd"_compute_dpd.html - 
"dpd/atom"_compute_dpd_atom.html - 
"edpd/temp/atom"_compute_edpd_temp_atom.html - 
"entropy/atom"_compute_entropy_atom.html - 
"erotate/asphere"_compute_erotate_asphere.html - rotational energy of aspherical particles
"erotate/rigid"_compute_erotate_rigid.html - rotational energy of rigid bodies
"erotate/sphere"_compute_erotate_sphere.html - rotational energy of spherical particles
"erotate/sphere/atom"_compute_erotate_sphere.html - rotational energy for each spherical particle
"erotate/sphere/atom"_compute_erotate_sphere_atom.html - 
"event/displace"_compute_event_displace.html - detect event on atom displacement
"fep"_compute_fep.html - 
"force/tally"_compute_tally.html - 
"fragment/atom"_compute_cluster_atom.html - fragment ID for each atom
"global/atom"_compute_global_atom.html - 
"group/group"_compute_group_group.html - energy/force between two groups of atoms
"gyration"_compute_gyration.html - radius of gyration of group of atoms
"gyration/chunk"_compute_gyration_chunk.html - radius of gyration for each chunk
"heat/flux"_compute_heat_flux.html - heat flux through a group of atoms
"heat/flux/tally"_compute_tally.html - 
"hexorder/atom"_compute_hexorder_atom.html - bond orientational order parameter q6
"improper"_compute_improper.html - 
"hma"_compute_hma.html - harmonically mapped averaging for atomic crystals
"improper/local"_compute_improper_local.html - angle of each improper
"inertia/chunk"_compute_inertia_chunk.html - inertia tensor for each chunk
"ke"_compute_ke.html - translational kinetic energy
"ke/atom"_compute_ke_atom.html - kinetic energy for each atom
"ke/atom/eff"_compute_ke_atom_eff.html - 
"ke/eff"_compute_ke_eff.html - 
"ke/rigid"_compute_ke_rigid.html - translational kinetic energy of rigid bodies
"meso/e/atom"_compute_meso_e_atom.html - 
"meso/rho/atom"_compute_meso_rho_atom.html - 
"meso/t/atom"_compute_meso_t_atom.html - 
"msd"_compute_msd.html - mean-squared displacement of group of atoms
"msd/chunk"_compute_msd_chunk.html - mean-squared displacement for each chunk
"msd/nongauss"_compute_msd_nongauss.html - MSD and non-Gaussian parameter of group of atoms
@@ -241,77 +219,36 @@ compute"_Commands_compute.html doc page are followed by one or more of
"pair/local"_compute_pair_local.html - distance/energy/force of each pairwise interaction
"pe"_compute_pe.html - potential energy
"pe/atom"_compute_pe_atom.html - potential energy for each atom
"pe/mol/tally"_compute_tally.html - 
"pe/tally"_compute_tally.html - 
"plasticity/atom"_compute_plasticity_atom.html - Peridynamic plasticity for each atom
"pressure"_compute_pressure.html - total pressure and pressure tensor
"pressure/cylinder"_compute_pressure_cylinder.html - 
"pressure/uef"_compute_pressure_uef.html - 
"property/atom"_compute_property_atom.html - convert atom attributes to per-atom vectors/arrays
"property/chunk"_compute_property_chunk.html - extract various per-chunk attributes
"property/local"_compute_property_local.html - convert local attributes to localvectors/arrays
"ptm/atom"_compute_ptm_atom.html - 
"property/chunk"_compute_property_chunk.html - extract various per-chunk attributes
"rdf"_compute_rdf.html - radial distribution function g(r) histogram of group of atoms
"reduce"_compute_reduce.html - combine per-atom quantities into a single global value
"reduce/chunk"_compute_reduce_chunk.html - reduce per-atom quantities within each chunk
"reduce/region"_compute_reduce.html - same as compute reduce, within a region
"rigid/local"_compute_rigid_local.html - extract rigid body attributes
"saed"_compute_saed.html - 
"slice"_compute_slice.html - extract values from global vector or array
"smd/contact/radius"_compute_smd_contact_radius.html - 
"smd/damage"_compute_smd_damage.html - 
"smd/hourglass/error"_compute_smd_hourglass_error.html - 
"smd/internal/energy"_compute_smd_internal_energy.html - 
"smd/plastic/strain"_compute_smd_plastic_strain.html - 
"smd/plastic/strain/rate"_compute_smd_plastic_strain_rate.html - 
"smd/rho"_compute_smd_rho.html - 
"smd/tlsph/defgrad"_compute_smd_tlsph_defgrad.html - 
"smd/tlsph/dt"_compute_smd_tlsph_dt.html - 
"smd/tlsph/num/neighs"_compute_smd_tlsph_num_neighs.html - 
"smd/tlsph/shape"_compute_smd_tlsph_shape.html - 
"smd/tlsph/strain"_compute_smd_tlsph_strain.html - 
"smd/tlsph/strain/rate"_compute_smd_tlsph_strain_rate.html - 
"smd/tlsph/stress"_compute_smd_tlsph_stress.html - 
"smd/triangle/vertices"_compute_smd_triangle_vertices.html -
"smd/triangle/vertices"_compute_smd_triangle_vertices.html - 
"smd/ulsph/num/neighs"_compute_smd_ulsph_num_neighs.html - 
"smd/ulsph/strain"_compute_smd_ulsph_strain.html - 
"smd/ulsph/strain/rate"_compute_smd_ulsph_strain_rate.html - 
"smd/ulsph/stress"_compute_smd_ulsph_stress.html - 
"smd/vol"_compute_smd_vol.html - 
"sna/atom"_compute_sna_atom.html - calculate bispectrum coefficients for each atom
"snad/atom"_compute_sna_atom.html - derivative of bispectrum coefficients for each atom
"snav/atom"_compute_sna_atom.html - virial contribution from bispectrum coefficients for each atom
"spin"_compute_spin.html - 
"stress/atom"_compute_stress_atom.html - stress tensor for each atom
"stress/mop"_compute_stress_mop.html - 
"stress/mop/profile"_compute_stress_mop.html - 
"stress/tally"_compute_tally.html - 
"tdpd/cc/atom"_compute_tdpd_cc_atom.html - 
"temp"_compute_temp.html - temperature of group of atoms
"temp/asphere"_compute_temp_asphere.html - temperature of aspherical particles
"temp/body"_compute_temp_body.html - temperature of body particles
"temp/chunk"_compute_temp_chunk.html - temperature of each chunk
"temp/com"_compute_temp_com.html - temperature after subtracting center-of-mass velocity
"temp/cs"_compute_temp_cs.html - 
"temp/deform"_compute_temp_deform.html - temperature excluding box deformation velocity
"temp/deform/eff"_compute_temp_deform_eff.html - 
"temp/drude"_compute_temp_drude.html - 
"temp/eff"_compute_temp_eff.html - 
"temp/partial"_compute_temp_partial.html - temperature excluding one or more dimensions of velocity
"temp/profile"_compute_temp_profile.html - temperature excluding a binned velocity profile
"temp/ramp"_compute_temp_ramp.html - temperature excluding ramped velocity component
"temp/region"_compute_temp_region.html - temperature of a region of atoms
"temp/region/eff"_compute_temp_region_eff.html - 
"temp/rotate"_compute_temp_rotate.html - 
"temp/sphere"_compute_temp_sphere.html - temperature of spherical particles
"temp/uef"_compute_temp_uef.html - 
"ti"_compute_ti.html - thermodynamic integration free energy values
"torque/chunk"_compute_torque_chunk.html - torque applied on each chunk
"vacf"_compute_vacf.html - velocity-autocorrelation function of group of atoms
"vcm/chunk"_compute_vcm_chunk.html - velocity of center-of-mass for each chunk
"voronoi/atom"_compute_voronoi_atom.html - Voronoi volume and neighbors for each atom
"xrd"_compute_xrd.html - :ul
"voronoi/atom"_compute_voronoi_atom.html - Voronoi volume and neighbors for each atom :ul

[Restrictions:] none

+147 −0
Original line number Diff line number Diff line
"LAMMPS WWW Site"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c

:link(lws,http://lammps.sandia.gov)
:link(ld,Manual.html)
:link(lc,Section_commands.html#comm)

:line

compute hma command :h3

[Syntax:]

compute ID group-ID hma temp-ID keyword ... :pre

ID, group-ID are documented in "compute"_compute.html command :l
hma = style name of this compute command :l
temp-ID = ID of fix that specifies the set temperature during canonical simulation :l
keyword = {anharmonic} {u} {p Pharm} {cv} :l
  {anharmonic} = compute will return anharmonic property values
  {u} = compute will return potential energy
  {p} = compute will return pressure.  the following keyword must be the difference between the harmonic pressure and lattice pressure as described below
  {cv} = compute will return the heat capacity :pre
:ule

[Examples:]

compute 2 all hma 1 u
compute 2 all hma 1 anharmonic u p 0.9
compute 2 all hma 1 u cv :pre



[Description:]

Define a computation that calculates the properties of a solid (potential
energy, pressure or heat capacity), using the harmonically-mapped averaging
(HMA) method. 
This command yields much higher precision than the equivalent compute commands
("compute pe"_compute_pe.html, "compute pressure"_compute_pressure.html, etc.)
commands during a canonical simulation of an atomic crystal.


In this method, the analytically known harmonic behavior of a crystal is removed from the traditional ensemble
averages, which leads to an accurate and precise measurement of the anharmonic contributions without contamination 
by noise produced by the already-known harmonic behavior. 
A detailed description of this method can be found in ("Moustafa"_#hma-Moustafa). The potential energy is computed by the formula:

\begin\{equation\}
\left< U\right>_\{HMA\} = \frac\{d(N-1)\}\{2\beta\}  + \left< U + \frac\{1\}\{2\} F\bullet\Delta r \right>
\end\{equation\}

where \(N\) is the number of atoms in the system, \(\beta\) is the reciprocal of the thermodynamic temperature, \(d\) is the
dimensionality of the system (2 or 3 for 2d/3d), \(F\bullet\Delta r\) is the sum of dot products of the 
atomic force vectors and displacement (from lattice sites) vectors, and \(U\) is the sum of 
pair, bond, angle, dihedral, improper, kspace (long-range), and fix energies. 

The pressure is computed by the formula:

\begin\{equation\}
\left< P\right>_\{HMA\} = \Delta \hat P + \left< P_\{vir\} + \frac\{\beta \Delta \hat P - \rho\}\{d(N-1)\} F\bullet\Delta r \right>
\end\{equation\}

where \(\rho\) is the number density of the system, \(\Delta \hat P\) is the
difference between theh harmonic and lattice pressure, and \(P_\{vir\}\) is
the virial pressure computed as the sum of pair, bond, angle, dihedral,
improper, kspace (long-range), and fix contributions to the force on each
atom.  Although the method will work for any value of \(\Delta \hat P\)
specified (use pressure "units"_units.html), the precision of the resultant
pressure is sensitive to \(\Delta \hat P\); the precision tends to be
best when \(\Delta \hat P\) is the actual the difference between the lattice
pressure and harmonic pressure.

\begin\{equation\}
\left<C_V \right>_\{HMA\} = \frac\{d k_B (N-1)\}\{2\} + \beta \left( \left<
U_\{HMA\}^2 \right> - \left<U_\{HMA\}\right>^2 \right)/T + \frac\{1\}\{4 T\}
\left< F\bullet\Delta r + \Delta r \bullet \Phi \bullet \Delta r \right>
\end\{equation\}

where \(\Phi\) is the Hessian of second derivatives. The compute hma command
computes the full expression for \(C_V\) except for the
\(\left<U_\{HMA\}^2\right>^2\) in the variance term, which can be obtained by
passing the {u} keyword; you must add this extra contribution to the \(C_V\)
value reported by this compute.  The variance term can cause significant
roundoff error when computing \(C_V\).  To address this, the {anharmonic}
keyword can be passed and/or the output format can be speicified with more
digits.

thermo_modify format float '%22.15e' :pre

The {anharmonic} keyword will instruct the compute to return anharmonic
properties rather than the full properties (lattice, harmonic and anharmonic).
When using this keyword, the compute must be first active (it must be included
via a "thermo_style custom"_thermo_style.html command) while the atoms are
still at their lattice sites (before equilibration).

The temp-ID specified with compute hma command should be same as the fix-ID of Nose-Hoover ("fix nvt"_fix_nh.html) or 
Berendsen ("fix temp/berendsen"_fix_temp_berendsen.html) thermostat used for the simulation. While using this command, Langevin thermostat 
("fix langevin"_fix_langevin.html) 
should be avoided as its extra forces interfere with the HMA implementation. 

 

NOTE: Compute hma command should be used right after the energy minimization, when the atoms are at their lattice sites. 
The simulation should not be started before this command has been used in the input script.


The following example illustrates the placement of this command in the input script:


min_style cg 
minimize 1e-35 1e-15 50000 500000 
compute 1 all hma thermostatid u
fix thermostatid all nvt temp 600.0 600.0 100.0 :pre  



NOTE: Compute hma should be used when the atoms of the solid do not diffuse. Diffusion will reduce the precision in the potential energy computation.

 
NOTE: The "fix_modify energy yes"_fix_modify.html command must also be specified if a fix is to contribute potential energy to this command.

:line

[Output info:]

This compute calculates a global vector that includes the n properties
requested as arguments to the command (the potential energy, pressure or heat
capacity).  The elements of the vector can be accessed by indices 1-n by any
command that uses global vector values as input.  See the "Howto
output"_Howto_output.html doc page for an overview of LAMMPS output options.

The vector values calculated by this compute are "extensive".  The
scalar value will be in energy "units"_units.html.

[Restrictions:] Usage restricted to canonical (NVT) ensemble simulation only.

[Related commands:]

"compute pe"_compute_pe.html, "compute pressure"_compute_pressure.html

[Default:] none

:line

:link(hma-Moustafa)
[(Moustafa)] Sabry G. Moustafa, Andrew J. Schultz, and David A. Kofke, {Very fast averaging of thermal properties of crystals by molecular simulation}, 
"Phys. Rev. E \[92\], 043303 (2015)"_https://link.aps.org/doi/10.1103/PhysRevE.92.043303
+2 −7
Original line number Diff line number Diff line
@@ -6,7 +6,6 @@ Computes :h1
   :maxdepth: 1

   compute_ackland_atom
   compute_adf
   compute_angle
   compute_angle_local
   compute_angmom_chunk
@@ -16,7 +15,6 @@ Computes :h1
   compute_bond_local
   compute_centro_atom
   compute_chunk_atom
   compute_chunk_spread_atom
   compute_cluster_atom
   compute_cna_atom
   compute_cnp_atom
@@ -46,6 +44,7 @@ Computes :h1
   compute_gyration_chunk
   compute_heat_flux
   compute_hexorder_atom
   compute_hma
   compute_improper
   compute_improper_local
   compute_inertia_chunk
@@ -68,15 +67,12 @@ Computes :h1
   compute_pe_atom
   compute_plasticity_atom
   compute_pressure
   compute_pressure_cylinder
   compute_pressure_uef
   compute_property_atom
   compute_property_chunk
   compute_property_local
   compute_ptm_atom
   compute_rdf
   compute_reduce
   compute_reduce_chunk
   compute_rigid_local
   compute_saed
   compute_slice
@@ -94,7 +90,7 @@ Computes :h1
   compute_smd_tlsph_strain
   compute_smd_tlsph_strain_rate
   compute_smd_tlsph_stress
   compute_smd_triangle_vertices
   compute_smd_triangle_mesh_vertices
   compute_smd_ulsph_num_neighs
   compute_smd_ulsph_strain
   compute_smd_ulsph_strain_rate
@@ -103,7 +99,6 @@ Computes :h1
   compute_sna_atom
   compute_spin
   compute_stress_atom
   compute_stress_mop
   compute_tally
   compute_tdpd_cc_atom
   compute_temp

src/compute_hma.cpp

0 → 100644
+470 −0

File added.

Preview size limit exceeded, changes collapsed.

src/compute_hma.h

0 → 100644
+66 −0
Original line number Diff line number Diff line
/* -*- 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 COMPUTE_CLASS

ComputeStyle(HMA,ComputeHMA)

#else

#ifndef LMP_COMPUTE_HMA_H
#define LMP_COMPUTE_HMA_H

#include "compute.h"

namespace LAMMPS_NS {

class ComputeHMA : public Compute {
 public:
  ComputeHMA(class LAMMPS *, int, char **);
  ~ComputeHMA();
  void setup();
  void init();
  void init_list(int, class NeighList *);
  void compute_vector();
  void set_arrays(int);
  int pack_forward_comm(int, int *, double *, int, int *);
  void unpack_forward_comm(int, int, double *);

 private:
  int nmax;
  int atomsingroup;
  char *id_fix;
  char *id_temp;
  double finaltemp;
  class FixStore *fix;
  double boltz, nktv2p, inv_volume;
  double deltaPcap;
  double virial_compute(int);
  static double sumVirial(int n, double* v) {
    double x = 0;
    for (int i=0; i<n; i++) x += v[i];
    return x;
  }
  int computeU, computeP, computeCv;
  class NeighList *list; // half neighbor list
  double **deltaR;
  int returnAnharmonic;
  double uLat, pLat;
};

}

#endif
#endif