\(\renewcommand{\AA}{\text{Å}}\)
1.1.3. System properties
This section documents the following functions:
The library interface allows the extraction of different kinds of information about the active simulation instance and also - in some cases - to apply modifications to it. This enables combining of a LAMMPS simulation with other processing and simulation methods computed by the calling code, or by another code that is coupled to LAMMPS via the library interface. In some cases the data returned is direct reference to the original data inside LAMMPS, cast to a void pointer. In that case the data needs to be cast to a suitable pointer for the calling program to access it, and you may need to know the correct dimensions and lengths. This also means you can directly change those value(s) from the calling program (e.g., to modify atom positions). Of course, changing values should be done with care. When accessing per-atom data, please note that these data are the per-processor local data and are indexed accordingly. Per-atom data can change sizes and ordering at every neighbor list rebuild or atom sort event as atoms migrate between subdomains and processors.
#include "library.h"
#include <stdio.h>
int main(int argc, char **argv)
{
  void *handle;
  int i;
  handle = lammps_open_no_mpi(0, NULL, NULL);
  lammps_file(handle,"in.sysinit");
  printf("Running a simulation with %g atoms.\n",
         lammps_get_natoms(handle));
  printf(" %d local and %d ghost atoms. %d atom types\n",
         lammps_extract_setting(handle,"nlocal"),
         lammps_extract_setting(handle,"nghost"),
         lammps_extract_setting(handle,"ntypes"));
  double  *dt = (double *)lammps_extract_global(handle,"dt");
  printf("Changing timestep from %g to 0.5\n", *dt);
  *dt = 0.5;
  lammps_command(handle,"run 1000 post no");
  for (i=0; i < 10; ++i) {
    lammps_command(handle,"run 100 pre no post no");
    printf("PE = %g\nKE = %g\n",
           lammps_get_thermo(handle,"pe"),
           lammps_get_thermo(handle,"ke"));
  }
  lammps_close(handle);
  return 0;
}
- 
double lammps_get_natoms(void *handle)
- Return the total number of atoms in the system. - This number may be very large when running large simulations across multiple processes. Depending on compile time choices, LAMMPS may be using either 32-bit or a 64-bit integer to store this number. For portability this function returns thus a double precision floating point number, which can represent up to a 53-bit signed integer exactly (\(\approx 10^{16}\)). - As an alternative, you can use - lammps_extract_global()and cast the resulting pointer to an integer pointer of the correct size and dereference it. The size of that integer (in bytes) can be queried by calling- lammps_extract_setting()to return the size of a- bigintinteger.- Changed in version 18Sep2020: The type of the return value was changed from - intto- doubleto accommodate reporting atom counts for larger systems that would overflow a 32-bit int without having to depend on a 64-bit bit integer type definition.- Parameters:
- handle – pointer to a previously created LAMMPS instance 
- Returns:
- total number of atoms or 0 if value is too large 
 
- 
double lammps_get_thermo(void *handle, const char *keyword)
- Evaluate a thermo keyword. - This function returns the current value of a thermo keyword. Unlike - lammps_extract_global()it does not give access to the storage of the desired data but returns its value as a- double, so it can also return information that is computed on-the-fly. Use- lammps_last_thermo()to get access to the cached data from the last thermo output.- Parameters:
- handle – pointer to a previously created LAMMPS instance 
- keyword – string with the name of the thermo keyword 
 
- Returns:
- value of the requested thermo property or 0.0 
 
- 
void *lammps_last_thermo(void *handle, const char *what, int index)
- Access cached data from last thermo output - Added in version 15Jun2023. - This function provides access to cached data from the last thermo output. This differs from - lammps_get_thermo()in that it does not trigger an evaluation. Instead it provides direct access to a read-only location of the last thermo output data and the corresponding keyword strings. How to handle the return value depends on the value of the what argument string. When accessing the data from a concurrent thread while LAMMPS is running, the cache needs to be locked first and then unlocked after the data is obtained, so that the data is not corrupted while reading in case LAMMPS wants to update it at the same time. Outside of a run, the lock/unlock calls have no effect.- Value of what - Description of return value - Data type - Uses index - setup - 1 if setup is not completed and thus thermo data invalid, 0 otherwise - pointer to int - no - line - line number (0-based) of current line in current file or buffer - pointer to int - no - imagename - file name of the last dump image file written - pointer to 0-terminated const char array - no - step - timestep when the last thermo output was generated or -1 - pointer to bigint - no - num - number of fields in thermo output - pointer to int - no - keyword - column keyword for thermo output - pointer to 0-terminated const char array - yes - type - data type of thermo output column; see - _LMP_DATATYPE_CONST- pointer to int - yes - data - actual field data for column - pointer to int, int64_t or double - yes - lock - acquires lock to thermo data cache - NULL pointer - no - unlock - releases lock to thermo data cache - NULL pointer - no - Note - The type property points to a static location that is reassigned with every call, so the returned pointer should be recast, dereferenced, and assigned immediately. Otherwise, its value may be changed with the next invocation of the function. - Parameters:
- handle – pointer to a previously created LAMMPS instance 
- what – string with the kind of data requested 
- index – integer with index into data arrays, ignored for scalar data 
 
- Returns:
- pointer to location of requested data cast to void or NULL 
 
- 
void lammps_extract_box(void *handle, double *boxlo, double *boxhi, double *xy, double *yz, double *xz, int *pflags, int *boxflag)
- Extract simulation box parameters. - This function (re-)initializes the simulation box and boundary information and then assign the designated data to the locations in the pointers passed as arguments. Any argument (except the first) may be a NULL pointer and then will not be assigned. - Parameters:
- handle – pointer to a previously created LAMMPS instance 
- boxlo – pointer to 3 doubles where the lower box boundary is stored 
- boxhi – pointer to 3 doubles where the upper box boundary is stored 
- xy – pointer to a double where the xy tilt factor is stored 
- yz – pointer to a double where the yz tilt factor is stored 
- xz – pointer to a double where the xz tilt factor is stored 
- pflags – pointer to 3 ints, set to 1 for periodic boundaries and 0 for non-periodic 
- boxflag – pointer to an int, which is set to 1 if the box will be changed during a simulation by a fix and 0 if not. 
 
 
- 
void lammps_reset_box(void *handle, double *boxlo, double *boxhi, double xy, double yz, double xz)
- Reset simulation box parameters. - This function sets the simulation box dimensions (upper and lower bounds and tilt factors) from the provided data and then re-initializes the box information and all derived settings. It may only be called before atoms are created. - Parameters:
- handle – pointer to a previously created LAMMPS instance 
- boxlo – pointer to 3 doubles containing the lower box boundary 
- boxhi – pointer to 3 doubles containing the upper box boundary 
- xy – xy tilt factor 
- yz – yz tilt factor 
- xz – xz tilt factor 
 
 
- 
void lammps_memory_usage(void *handle, double *meminfo)
- Get memory usage information - Added in version 18Sep2020. - This function will retrieve memory usage information for the current LAMMPS instance or process. The meminfo buffer will be filled with 3 different numbers (if supported by the operating system). The first is the tally (in MBytes) of all large memory allocations made by LAMMPS. This is a lower boundary of how much memory is requested and does not account for memory allocated on the stack or allocations via - new. The second number is the current memory allocation of the current process as returned by a memory allocation reporting in the system library. The third number is the maximum amount of RAM (not swap) used by the process so far. If any of the two latter parameters is not supported by the operating system it will be set to zero.- Parameters:
- handle – pointer to a previously created LAMMPS instance 
- meminfo – buffer with space for at least 3 double to store data in. 
 
 
- 
int lammps_get_mpi_comm(void *handle)
- Return current LAMMPS world communicator as integer - Added in version 18Sep2020. - This will take the LAMMPS “world” communicator and convert it to an integer using - MPI_Comm_c2f(), so it is equivalent to the corresponding MPI communicator in Fortran. This way it can be safely passed around between different programming languages. To convert it to the C language representation use- MPI_Comm_f2c().- If LAMMPS was compiled with MPI_STUBS, this function returns -1. - See also
 - Parameters:
- handle – pointer to a previously created LAMMPS instance 
- Returns:
- Fortran representation of the LAMMPS world communicator 
 
- 
int lammps_extract_setting(void *handle, const char *keyword)
- Query LAMMPS about global settings. - This function will retrieve or compute global properties. In contrast to - lammps_get_thermo()this function returns an- int. The following tables list the currently supported keyword. If a keyword is not recognized, the function returns -1. The integer sizes functions may be called without a valid LAMMPS object handle (it is ignored).- Integer sizes - Keyword - Description / Return value - bigint - size of the - bigintinteger type, 4 or 8 bytes. Set at compile time.- tagint - size of the - tagintinteger type, 4 or 8 bytes. Set at compile time.- imageint - size of the - imageintinteger type, 4 or 8 bytes. Set at compile time.- Image masks - These settings are related to how LAMMPS stores and interprets periodic images. The values are used internally by the Fortran interface and are not likely to be useful to users. - Keyword - Description / Return value - IMGMASK - Bit-mask used to convert image flags to a single integer - IMGMAX - Maximum allowed image number for a particular atom - IMGBITS - Bits used in image counts - IMG2BITS - Second bitmask used in image counts - System status - Keyword - Description / Return value - dimension - Number of dimensions: 2 or 3. See dimension command. - box_exist - 1 if the simulation box is defined, 0 if not. See create_box command. - kokkos_active - 1 if the KOKKOS package is compiled in and activated, 0 if not. See KOKKOS package. - kokkos_nthreads - Number of Kokkos threads per MPI process, 0 if Kokkos is not active. See KOKKOS package. - kokkos_ngpus - Number of Kokkos gpus per physical node, 0 if Kokkos is not active or no GPU support. See KOKKOS package. - nthreads - Number of requested OpenMP threads per MPI process for LAMMPS’ execution - newton_bond - 1 if Newton’s 3rd law is applied to bonded interactions, 0 if not. - newton_pair - 1 if Newton’s 3rd law is applied to non-bonded interactions, 0 if not. - triclinic - 1 if the the simulation box is triclinic, 0 if orthogonal. See change_box command. - Communication status - Keyword - Description / Return value - universe_rank - MPI rank on LAMMPS’ universe communicator (0 <= universe_rank < universe_size) - universe_size - Number of ranks on LAMMPS’ universe communicator (world_size <= universe_size) - world_rank - MPI rank on LAMMPS’ world communicator (0 <= world_rank < world_size, = comm->me) - world_size - Number of ranks on LAMMPS’ world communicator (aka comm->nprocs) - comm_style - communication style (0 = BRICK, 1 = TILED) - comm_layout - communication layout (0 = LAYOUT_UNIFORM, 1 = LAYOUT_NONUNIFORM, 2 = LAYOUT_TILED) - comm_mode - communication mode (0 = SINGLE, 1 = MULTI, 2 = MULTIOLD) - ghost_velocity - whether velocities are communicated for ghost atoms (0 = no, 1 = yes) - System sizes - Keyword - Description / Return value - nlocal - number of “owned” atoms of the current MPI rank. - nghost - number of “ghost” atoms of the current MPI rank. - nall - number of all “owned” and “ghost” atoms of the current MPI rank. - nmax - maximum of nlocal+nghost across all MPI ranks (for per-atom data array size). - ntypes - number of atom types - nbondtypes - number of bond types - nangletypes - number of angle types - ndihedraltypes - number of dihedral types - nimpropertypes - number of improper types - bond_per_atom - size of per-atom bond data arrays - angle_per_atom - size of per-atom angle data arrays - dihedral_per_atom - size of per-atom dihedral data arrays - improper_per_atom - size of per-atom improper data arrays - maxspecial - size of per-atom special data array - nellipsoids - number of atoms that have ellipsoid data - nlines - number of atoms that have line data (see pair style line/lj) - ntris - number of atoms that have triangle data (see pair style tri/lj) - nbodies - number of atoms that have body data (see the Body particle HowTo) - Neighbor list settings - neigh_every - neighbor lists are rebuild every this many steps - neigh_delay - neighbor lists are rebuild delayed this many steps - neigh_dist_check - 0 if always rebuild, 1 rebuild after 1/2 skin - neigh_ago - neighbor lists were rebuilt this many steps ago - nbondlist - number of entries in bondlist (get list with lammps_extract_global()) - nanglelist - number of entries in anglelist (get list with lammps_extract_global()) - ndihedrallist - number of entries in dihedrallist (get list with lammps_extract_global()) - nimproperlist - number of entries in improperlist (get list with lammps_extract_global()) - Atom style flags - Keyword - Description / Return value - molecule_flag - 1 if the atom style includes molecular topology data. See atom_style command. - q_flag - 1 if the atom style includes point charges. See atom_style command. - mu_flag - 1 if the atom style includes point dipoles. See atom_style command. - rmass_flag - 1 if the atom style includes per-atom masses, 0 if there are per-type masses. See atom_style command. - radius_flag - 1 if the atom style includes a per-atom radius. See atom_style command. - ellipsoid_flag - 1 if the atom style describes extended particles that may be ellipsoidal. See atom_style command. - omega_flag - 1 if the atom style can store per-atom rotational velocities. See atom_style command. - torque_flag - 1 if the atom style can store per-atom torques. See atom_style command. - angmom_flag - 1 if the atom style can store per-atom angular momentum. See atom_style command. - Thermo settings - Keyword - Description / Return value - thermo_every - The current interval of thermo output. See thermo command. - thermo_norm - 1 if the thermo output is normalized. See thermo_modify command. - See also
 - Parameters:
- handle – pointer to a previously created LAMMPS instance 
- keyword – string with the name of the thermo keyword 
 
- Returns:
- value of the queried setting or -1 if unknown 
 
- 
int lammps_extract_global_datatype(void *handle, const char *name)
- Get data type of internal global LAMMPS variables or arrays. - Added in version 18Sep2020. - This function returns an integer that encodes the data type of the global property with the specified name. See - _LMP_DATATYPE_CONSTfor valid values. Callers of- lammps_extract_global()can use this information to then decide how to cast the- void *pointer and access the data.- Parameters:
- handle – pointer to a previously created LAMMPS instance (unused) 
- name – string with the name of the extracted property 
 
- Returns:
- integer constant encoding the data type of the property or -1 if not found. 
 
- 
void *lammps_extract_global(void *handle, const char *name)
- Get pointer to internal global LAMMPS variables or arrays. - This function returns a pointer to the location of some global property stored in one of the constituent classes of a LAMMPS instance. The returned pointer is cast to - void *and needs to be cast to a pointer of the type that the entity represents. The pointers returned by this function are generally persistent; therefore it is not necessary to call the function again, unless a clear command command is issued which wipes out and recreates the contents of the- LAMMPSclass.- Please also see - lammps_extract_setting(),- lammps_get_thermo(), and- lammps_extract_box().- The following tables list the supported names, their data types, length of the data area, and a short description. The data type can also be queried through calling - lammps_extract_global_datatype(). The- biginttype may be defined to be either an- intor an- int64_t. This is set at compile time of the LAMMPS library and can be queried through calling- lammps_extract_setting(). The function- lammps_extract_global_datatype()will directly report the “native” data type. The following tables are provided:- Timestep settings - Name - Type - Length - Description - dt - double - 1 - length of the time step. See timestep command. - ntimestep - bigint - 1 - current time step number. See reset_timestep command. - atime - double - 1 - accumulated simulation time in time units. - atimestep - bigint - 1 - the number of the timestep when “atime” was last updated. - respa_levels - int - 1 - \(N_{respa}\) = number of r-RESPA levels. See run_style command. - respa_dt - double - \(N_{respa}\) - length of the time steps with r-RESPA. See run_style command. - Simulation box settings - Name - Type - Length - Description - boxlo - double - 3 - lower box boundaries in x-, y-, and z-direction; see create_box command. - boxhi - double - 3 - lower box boundaries in x-, y-, and z-direction; see create_box command. - boxxlo - double - 1 - lower box boundary in x-direction; see create_box command. - boxxhi - double - 1 - upper box boundary in x-direction; see create_box command. - boxylo - double - 1 - lower box boundary in y-direction; see create_box command. - boxyhi - double - 1 - upper box boundary in y-direction; see create_box command. - boxzlo - double - 1 - lower box boundary in z-direction; see create_box command. - boxzhi - double - 1 - upper box boundary in z-direction; see create_box command. - sublo - double - 3 - subbox lower boundaries in x-, y-, and z-direction - subhi - double - 3 - subbox upper boundaries in x-, y-, and z-direction - sublo_lambda - double - 3 - subbox lower boundaries in fractional coordinates (for triclinic cells only) - subhi_lambda - double - 3 - subbox upper boundaries in fractional coordinates (for triclinic cells only) - periodicity - int - 3 - 0 if non-periodic, 1 if periodic for x, y, and z; see boundary command. - triclinic - int - 1 - 1 if box is triclinic, 0 if orthogonal; see change_box command. - xy - double - 1 - triclinic tilt factor; see Triclinic (non-orthogonal) simulation boxes. - yz - double - 1 - triclinic tilt factor; see Triclinic (non-orthogonal) simulation boxes. - xz - double - 1 - triclinic tilt factor; see Triclinic (non-orthogonal) simulation boxes. - xlattice - double - 1 - lattice spacing in x-direction; see lattice command. - ylattice - double - 1 - lattice spacing in y-direction; see lattice command. - zlattice - double - 1 - lattice spacing in z-direction; see lattice command. - procgrid - int - 3 - processor count in x-, y-, and z- direction; see processors command. - System property settings - Name - Type - Length - Description - natoms - bigint - 1 - total number of atoms in the simulation. - nbonds - bigint - 1 - total number of bonds in the simulation. - nangles - bigint - 1 - total number of angles in the simulation. - ndihedrals - bigint - 1 - total number of dihedrals in the simulation. - nimpropers - bigint - 1 - total number of impropers in the simulation. - nlocal - int - 1 - number of “owned” atoms of the current MPI rank. - nghost - int - 1 - number of “ghost” atoms of the current MPI rank. - nmax - int - 1 - maximum of nlocal+nghost across all MPI ranks (for per-atom data array size). - ntypes - int - 1 - number of atom types - special_lj - double - 4 - special pair weighting factors for LJ interactions (first element is always 1.0) - special_coul - double - 4 - special pair weighting factors for Coulomb interactions (first element is always 1.0) - map_style - int - 1 - atom map setting: 0 = none, 1 = array, 2 = hash, 3 = yes - map_tag_max - int/bigint - 1 - largest atom ID that can be mapped to a local index (bigint with -DLAMMPS_BIGBIG) - sametag - int - variable - index of next local atom with the same ID in ascending order. -1 signals end. - sortfreq - int - 1 - frequency of atom sorting. 0 means sorting is off. - nextsort - bigint - 1 - timestep when atoms are sorted next - q_flag - int - 1 - deprecated. Use - lammps_extract_setting()instead.- atom_style - char * - 1 - string with the current atom style. - pair_style - char * - 1 - string with the current pair style. - bond_style - char * - 1 - string with the current bond style. - angle_style - char * - 1 - string with the current angle style. - dihedral_style - char * - 1 - string with the current dihedral style. - improper_style - char * - 1 - string with the current improper style. - kspace_style - char * - 1 - string with the current KSpace style. - Neighbor topology data - Get length of lists with lammps_extract_setting(). - Name - Type - Length - Description - neigh_bondlist - 2d int - nbondlist - list of bonds (atom1, atom2, type) - neigh_anglelist - 2d int - nanglelist - list of angles (atom1, atom2, atom3, type) - neigh_dihedrallist - 2d int - ndihedrallist - list of dihedrals (atom1, atom2, atom3, atom4, type) - neigh_improperlist - 2d int - nimproperlist - list of impropers (atom1, atom2, atom3, atom4, type) - Energy and virial tally settings - Name - Type - Length - Description - eflag_global - bigint - 1 - timestep global energy is tallied on - eflag_atom - bigint - 1 - timestep per-atom energy is tallied on - vflag_global - bigint - 1 - timestep global virial is tallied on - vflag_atom - bigint - 1 - timestep per-atom virial is tallied on - Git revision and version settings - Name - Type - Length - Description - git_commit - const char * - 1 - Git commit hash for the LAMMPS version. - git_branch - const char * - 1 - Git branch for the LAMMPS version. - git_descriptor - const char * - 1 - Combined descriptor for the git revision - lammps_version - const char * - 1 - LAMMPS version string. - Unit settings - Name - Type - Length - Description - units - char * - 1 - string with the current unit style. See units command. - boltz - double - 1 - value of the “boltz” constant. See units command. - hplanck - double - 1 - value of the “hplanck” constant. See units command. - mvv2e - double - 1 - factor to convert \(\frac{1}{2}mv^2\) for a particle to the current energy unit; See units command. - ftm2v - double - 1 - (description missing) See units command. - mv2d - double - 1 - (description missing) See units command. - nktv2p - double - 1 - (description missing) See units command. - qqr2e - double - 1 - factor to convert \(\frac{q_i q_j}{r}\) to energy units; See units command. - qe2f - double - 1 - (description missing) See units command. - vxmu2f - double - 1 - (description missing) See units command. - xxt2kmu - double - 1 - (description missing) See units command. - dielectric - double - 1 - value of the dielectric constant. See dielectric command. - qqrd2e - double - 1 - (description missing) See units command. - e_mass - double - 1 - (description missing) See units command. - hhmrr2e - double - 1 - (description missing) See units command. - mvh2r - double - 1 - (description missing) See units command. - angstrom - double - 1 - constant to convert current length unit to angstroms; 1.0 for reduced (aka “lj”) units. See units command. - femtosecond - double - 1 - constant to convert current time unit to femtoseconds; 1.0 for reduced (aka “lj”) units - qelectron - double - 1 - (description missing) See units command. - Warning - Modifying the data in the location pointed to by the returned pointer may lead to inconsistent internal data and thus may cause failures or crashes or bogus simulations. In general it is thus usually better to use a LAMMPS input command that sets or changes these parameters. Those will take care of all side effects and necessary updates of settings derived from such settings. Where possible, a reference to such a command or a relevant section of the manual is given below. - Parameters:
- handle – pointer to a previously created LAMMPS instance 
- name – string with the name of the extracted property 
 
- Returns:
- pointer (cast to - void *) to the location of the requested property. NULL if name is not known.
 
- 
int lammps_extract_pair_dimension(void *handle, const char *name)
- Get data dimension of pair style data accessible via Pair::extract(). - Added in version 29Aug2024. - This function returns an integer that specified the dimensionality of the data that can be extracted from the current pair style with - Pair::extract(). Callers of- lammps_extract_pair()can use this information to then decide how to cast the- void *pointer and access the data.- Parameters:
- handle – pointer to a previously created LAMMPS instance 
- name – string with the name of the extracted property 
 
- Returns:
- integer constant encoding the dimensionality of the extractable pair style property or -1 if not found. 
 
- 
void *lammps_extract_pair(void *handle, const char *name)
- Get extract pair style data accessible via Pair::extract(). - Added in version 29Aug2024. - This function returns a pointer to data available from the current pair style with - Pair::extract(). The dimensionality of the returned pointer can be determined with- lammps_extract_pair_dimension().- Parameters:
- handle – pointer to a previously created LAMMPS instance 
- name – string with the name of the extracted property 
 
- Returns:
- pointer (cast to - void *) to the location of the requested property. NULL if name is not known.
 
- 
int lammps_map_atom(void *handle, const void *id)
- Map global atom ID to local atom index - Added in version 27June2024. - This function returns an integer that corresponds to the local atom index for an atom with the global atom ID id. The atom ID is passed as a void pointer so that it can use the same interface for either a 32-bit or 64-bit tagint. The size of the tagint can be determined using - lammps_extract_setting().- Parameters:
- handle – pointer to a previously created LAMMPS instance 
- id – void pointer to the atom ID (of data type tagint, i.e. 32-bit or 64-bit integer) 
 
- Returns:
- local atom index or -1 if the atom is not found or no map exists 
 
