Gromacs  2025.0-dev-20241011-013a99c
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Typedefs | Functions | Variables
anonymous_namespace{bonded.cpp} Namespace Reference

Typedefs

using BondedFunction = real(*)(int nbonds, const t_iatom iatoms[], const t_iparams iparams[], const rvec x[], rvec4 f[], rvec fshift[], const t_pbc *pbc, real lambda, real *dvdlambda, gmx::ArrayRef< const real > charge, t_fcdata *fcd, t_disresdata *disresdata, t_oriresdata *oriresdata, int *ddgatindex)
 Type of CPU function to compute a bonded interaction.
 
using CmapForceStructure = std::array< real, 4 >
 

Functions

int pbc_rvec_sub (const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
 Compute dx = xi - xj, modulo PBC if non-NULL. More...
 
template<BondedKernelFlavor flavor>
void spreadBondForces (const real bondForce, const rvec dx, const int ai, const int aj, rvec4 *f, int shiftIndex, rvec *fshift)
 Spreads and accumulates the bonded forces to the two atoms and adds the virial contribution when needed. More...
 
template<BondedKernelFlavor flavor>
real morse_bonds (int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], const rvec x[], rvec4 f[], rvec fshift[], const t_pbc *pbc, real lambda, real *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)
 Morse potential bond. More...
 
int cmap_setup_grid_index (int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
 Mysterious undocumented function.
 
gmx::RVec processCmapForceComponent (const gmx::RVec a, const gmx::RVec b, const real df, const real gaa, const real fga, const real gbb, const real hgb, const int dim)
 
CmapForceStructure applyCmapForceComponent (const gmx::RVec forceComponent)
 
void accumulateCmapForces (const rvec x[], rvec4 f[], rvec fshift[], const t_pbc *pbc, const gmx::RVec r_ij, const gmx::RVec r_kj, const gmx::RVec r_kl, const gmx::RVec a, const gmx::RVec b, gmx::RVec h, const real ra2r, const real rb2r, const real rgr, const real rg, const int ai, const int aj, const int ak, const int al, const real df, const int t1, const int t2)
 

Variables

const int cmap_coeff_matrix []
 Mysterious CMAP coefficient matrix. More...
 

Function Documentation

template<BondedKernelFlavor flavor>
real anonymous_namespace{bonded.cpp}::morse_bonds ( int  nbonds,
const t_iatom  forceatoms[],
const t_iparams  forceparams[],
const rvec  x[],
rvec4  f[],
rvec  fshift[],
const t_pbc pbc,
real  lambda,
real 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 
)

Morse potential bond.

By Frank Everdij. Three parameters needed:

b0 = equilibrium distance in nm be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e)) cb = well depth in kJ/mol

Note: the potential is referenced to be +cb at infinite separation and zero at the equilibrium distance!

int anonymous_namespace{bonded.cpp}::pbc_rvec_sub ( const t_pbc pbc,
const rvec  xi,
const rvec  xj,
rvec  dx 
)

Compute dx = xi - xj, modulo PBC if non-NULL.

Todo:
This kind of code appears in many places. Consolidate it
template<BondedKernelFlavor flavor>
void anonymous_namespace{bonded.cpp}::spreadBondForces ( const real  bondForce,
const rvec  dx,
const int  ai,
const int  aj,
rvec4 *  f,
int  shiftIndex,
rvec fshift 
)
inline

Spreads and accumulates the bonded forces to the two atoms and adds the virial contribution when needed.

shiftIndex is used as the periodic shift.

Variable Documentation

const int anonymous_namespace{bonded.cpp}::cmap_coeff_matrix[]
Initial value:
= {
1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4, 0, 0, 0, 0, 0, 0, 0, 0,
3, 0, -9, 6, -2, 0, 6, -4, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4, 0, 0, 0, 0, 1, 0, -3, 2,
-2, 0, 6, -4, 1, 0, -3, 2, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, 3, -2,
0, 0, -6, 4, 0, 0, 3, -2, 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2, 0, 0, 0, 0, 0, 0, 0, 0,
0, 0, -3, 3, 0, 0, 2, -2, 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1, 0, 0, 0, 0, 0, 0, 0, 0,
0, -1, 2, -1, 0, 1, -2, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
}

Mysterious CMAP coefficient matrix.