Gromacs  5.1.5
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions | Variables
#include "gmxpre.h"
#include "bonded.h"
#include "config.h"
#include <assert.h>
#include <cmath>
#include <algorithm>
#include "gromacs/math/utilities.h"
#include "gromacs/math/vec.h"
#include "gromacs/pbcutil/ishift.h"
#include "gromacs/pbcutil/mshift.h"
#include "gromacs/pbcutil/pbc-simd.h"
#include "gromacs/pbcutil/pbc.h"
#include "gromacs/simd/simd.h"
#include "gromacs/simd/simd_math.h"
#include "gromacs/simd/vector_operations.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/smalloc.h"
#include "listed-internal.h"
#include "pairs.h"
#include "restcbt.h"
+ Include dependency graph for bonded.cpp:

Description

This file defines low-level functions necessary for computing energies and forces for listed interactions.

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

Functions

static void gmx_simdcall gmx_hack_simd_gather_rvec_dist_two_index (const rvec *v, const int *index0, const int *index1, real *buf, gmx_simd_float_t *dx, gmx_simd_float_t *dy, gmx_simd_float_t *dz)
 Store differences between indexed rvecs in SIMD registers. More...
 
static 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...
 
real morse_bonds (int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], const rvec x[], rvec f[], rvec fshift[], const t_pbc *pbc, const t_graph *g, real lambda, real *dvdlambda, const t_mdatoms *md, t_fcdata *fcd, int *global_atom_index)
 Morse potential bond. More...
 
static int cmap_setup_grid_index (int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
 Mysterious undocumented function.
 
real cmap_dihs (int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], const gmx_cmap_t *cmap_grid, const rvec x[], rvec f[], rvec fshift[], const struct t_pbc *pbc, const struct t_graph *g, real lambda, real *dvdlambda, const t_mdatoms *md, t_fcdata *fcd, int *global_atom_index)
 Compute CMAP dihedral energies and forces.
 

Variables

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

Function Documentation

static void gmx_simdcall gmx_hack_simd_gather_rvec_dist_two_index ( const rvec *  v,
const int *  index0,
const int *  index1,
real buf,
gmx_simd_float_t dx,
gmx_simd_float_t dy,
gmx_simd_float_t dz 
)
inlinestatic

Store differences between indexed rvecs in SIMD registers.

Returns SIMD register with the difference vectors: v[index0[i]] - v[index1[i]]

Parameters
[in]vArray of rvecs
[in]index0Index into the vector array
[in]index1Index into the vector array
[in,out]bufAligned tmp buffer of size 3*GMX_SIMD_REAL_WIDTH
[out]dxSIMD register with x difference
[out]dySIMD register with y difference
[out]dzSIMD register with z difference
real morse_bonds ( int  nbonds,
const t_iatom  forceatoms[],
const t_iparams  forceparams[],
const rvec  x[],
rvec  f[],
rvec  fshift[],
const t_pbc *  pbc,
const t_graph *  g,
real  lambda,
real dvdlambda,
const t_mdatoms *  md,
t_fcdata *  fcd,
int *  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!

static int pbc_rvec_sub ( const t_pbc *  pbc,
const rvec  xi,
const rvec  xj,
rvec  dx 
)
static

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

Todo:
This kind of code appears in many places. Consolidate it

Variable Documentation

const int 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.