Gromacs
2020.4
|
#include "gmxpre.h"
#include "pairs.h"
#include <cmath>
#include "gromacs/listed_forces/bonded.h"
#include "gromacs/math/functions.h"
#include "gromacs/math/vec.h"
#include "gromacs/mdtypes/group.h"
#include "gromacs/mdtypes/md_enums.h"
#include "gromacs/mdtypes/nblist.h"
#include "gromacs/mdtypes/simulation_workload.h"
#include "gromacs/pbcutil/ishift.h"
#include "gromacs/pbcutil/mshift.h"
#include "gromacs/pbcutil/pbc.h"
#include "gromacs/pbcutil/pbc_simd.h"
#include "gromacs/simd/simd.h"
#include "gromacs/simd/simd_math.h"
#include "gromacs/simd/vector_operations.h"
#include "gromacs/tables/forcetable.h"
#include "gromacs/utility/basedefinitions.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxassert.h"
#include "listed_internal.h"
This file defines functions for "pair" interactions (i.e. listed non-bonded interactions, e.g. 1-4 interactions)
Functions | |
static void | warning_rlimit (const rvec *x, int ai, int aj, int *global_atom_index, real r, real rlimit) |
Issue a warning if a listed interaction is beyond a table limit. | |
static real | evaluate_single (real r2, real tabscale, const real *vftab, real tableStride, real qq, real c6, real c12, real *velec, real *vvdw) |
Compute the energy and force for a single pair interaction. | |
static real | free_energy_evaluate_single (real r2, real sc_r_power, real alpha_coul, real alpha_vdw, real tabscale, const real *vftab, real tableStride, real qqA, real c6A, real c12A, real qqB, real c6B, real c12B, const real LFC[2], const real LFV[2], const real DLF[2], const real lfac_coul[2], const real lfac_vdw[2], const real dlfac_coul[2], const real dlfac_vdw[2], real sigma6_def, real sigma6_min, real sigma2_def, real sigma2_min, real *velectot, real *vvdwtot, real *dvdl) |
Compute the energy and force for a single pair interaction under FEP. | |
template<BondedKernelFlavor flavor> | |
static real | do_pairs_general (int ftype, int nbonds, const t_iatom iatoms[], const t_iparams iparams[], const rvec x[], rvec4 f[], rvec fshift[], const struct t_pbc *pbc, const struct t_graph *g, const real *lambda, real *dvdl, const t_mdatoms *md, const t_forcerec *fr, gmx_grppairener_t *grppener, int *global_atom_index) |
Calculate pair interactions, supports all types and conditions. | |
template<typename T , int pack_size, typename pbc_type > | |
static void | do_pairs_simple (int nbonds, const t_iatom iatoms[], const t_iparams iparams[], const rvec x[], rvec4 f[], const pbc_type pbc, const t_mdatoms *md, const real scale_factor) |
Calculate pairs, only for plain-LJ + plain Coulomb normal type. More... | |
void | do_pairs (int ftype, int nbonds, const t_iatom iatoms[], const t_iparams iparams[], const rvec x[], rvec4 f[], rvec fshift[], const struct t_pbc *pbc, const struct t_graph *g, const real *lambda, real *dvdl, const t_mdatoms *md, const t_forcerec *fr, const bool havePerturbedInteractions, const gmx::StepWorkload &stepWork, gmx_grppairener_t *grppener, int *global_atom_index) |
Calculate all listed pair interactions. More... | |
void do_pairs | ( | int | ftype, |
int | nbonds, | ||
const t_iatom | iatoms[], | ||
const t_iparams | iparams[], | ||
const rvec | x[], | ||
rvec4 | f[], | ||
rvec | fshift[], | ||
const struct t_pbc * | pbc, | ||
const struct t_graph * | g, | ||
const real * | lambda, | ||
real * | dvdl, | ||
const t_mdatoms * | md, | ||
const t_forcerec * | fr, | ||
const bool | havePerturbedInteractions, | ||
const gmx::StepWorkload & | stepWork, | ||
gmx_grppairener_t * | grppener, | ||
int * | global_atom_index | ||
) |
Calculate all listed pair interactions.
Calculate VdW/charge listed pair interactions (usually 1-4 interactions).
|
static |
Calculate pairs, only for plain-LJ + plain Coulomb normal type.
This function is templated for real/SimdReal and for optimization.