Gromacs  2024.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Enumerations | Functions | Variables
#include <string>
#include "gromacs/libgromacs_export.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/topology/ifunc.h"
#include "gromacs/utility/basedefinitions.h"
+ Include dependency graph for bonded.h:

Description

This file contains declarations necessary for low-level functions for computing energies and forces for bonded interactions.

Author
Mark Abraham mark..nosp@m.j.ab.nosp@m.raham.nosp@m.@gma.nosp@m.il.co.nosp@m.m

Classes

struct  gmx::EnumerationArray< EnumType, DataType, ArraySize >
 Wrapper for a C-style array with size and indexing defined by an enum. Useful for declaring arrays of enum names for debug or other printing. An ArrayRef<DataType> may be constructed from an object of this type. More...
 
class  gmx::ArrayRef< typename >
 STL-like interface to a C array of T (or part of a std container of T). More...
 

Enumerations

enum  BondedKernelFlavor {
  BondedKernelFlavor::ForcesSimdWhenAvailable, BondedKernelFlavor::ForcesNoSimd, BondedKernelFlavor::ForcesAndVirialAndEnergy, BondedKernelFlavor::ForcesAndEnergy,
  BondedKernelFlavor::Count
}
 For selecting which flavor of bonded kernel is used for simple bonded types. More...
 

Functions

real bond_angle (const rvec xi, const rvec xj, const rvec xk, const struct t_pbc *pbc, rvec r_ij, rvec r_kj, real *costh, int *t1, int *t2)
 Calculate bond-angle. No PBC is taken into account (use mol-shift)
 
real dih_angle (const rvec xi, const rvec xj, const rvec xk, const rvec xl, const struct t_pbc *pbc, rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n, int *t1, int *t2, int *t3)
 Calculate dihedral-angle. No PBC is taken into account (use mol-shift)
 
void do_dih_fup (int i, int j, int k, int l, real ddphi, rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n, rvec4 f[], rvec fshift[], const struct t_pbc *pbc, const rvec *x, int t1, int t2, int t3)
 Do an update of the forces for dihedral potentials.
 
void make_dp_periodic (real *dp)
 Make a dihedral fall in the range (-pi,pi)
 
real cmap_dihs (int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], const gmx_cmap_t *cmap_grid, const rvec x[], rvec4 f[], rvec fshift[], const struct t_pbc *pbc, real gmx_unused lambda, real gmx_unused *dvdlambda, gmx::ArrayRef< const real >, t_fcdata gmx_unused *fcd, t_disresdata gmx_unused *disresdata, t_oriresdata gmx_unused *oriresdata, int gmx_unused *global_atom_index)
 Compute CMAP dihedral energies and forces.
 
static constexpr bool computeEnergy (const BondedKernelFlavor flavor)
 Returns whether the energy should be computed.
 
static constexpr bool computeVirial (const BondedKernelFlavor flavor)
 Returns whether the virial should be computed.
 
static constexpr bool computeEnergyOrVirial (const BondedKernelFlavor flavor)
 Returns whether the energy and/or virial should be computed.
 
real calculateSimpleBond (int ftype, int numForceatoms, const t_iatom forceatoms[], const t_iparams forceparams[], const rvec x[], rvec4 f[], rvec fshift[], const struct t_pbc *pbc, real lambda, real *dvdlambda, gmx::ArrayRef< const real > charge, t_fcdata *fcd, t_disresdata *disresdata, t_oriresdata *oriresdata, int gmx_unused *global_atom_index, BondedKernelFlavor bondedKernelFlavor)
 Calculates bonded interactions for simple bonded types. More...
 
int nrnbIndex (int ftype)
 Getter for finding the flop count for an ftype interaction.
 

Variables

LIBGROMACS_EXPORT const
gmx::EnumerationArray
< BondedKernelFlavor,
std::string,
BondedKernelFlavor::Count
c_bondedKernelFlavorStrings
 Helper strings for human-readable messages.
 

Enumeration Type Documentation

enum BondedKernelFlavor
strong

For selecting which flavor of bonded kernel is used for simple bonded types.

Enumerator
ForcesSimdWhenAvailable 

Compute only forces, use SIMD when available; should not be used with perturbed parameters.

ForcesNoSimd 

Compute only forces, do not use SIMD.

ForcesAndVirialAndEnergy 

Compute forces, virial and energy (no SIMD)

ForcesAndEnergy 

Compute forces and energy (no SIMD)

Count 

The number of flavors.

Function Documentation

real calculateSimpleBond ( int  ftype,
int  numForceatoms,
const t_iatom  forceatoms[],
const t_iparams  forceparams[],
const rvec  x[],
rvec4  f[],
rvec  fshift[],
const struct t_pbc pbc,
real  lambda,
real dvdlambda,
gmx::ArrayRef< const real charge,
t_fcdata *  fcd,
t_disresdata *  disresdata,
t_oriresdata oriresdata,
int gmx_unused global_atom_index,
BondedKernelFlavor  bondedKernelFlavor 
)

Calculates bonded interactions for simple bonded types.

Exits with an error when the bonded type is not simple All pointers should be non-null, except for pbc and g which can be nullptr.

Returns
the energy or 0 when bondedKernelFlavor did not request the energy.