Gromacs  2019-beta1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Macros | Typedefs | Enumerations | Functions
forcetable.h File Reference
#include <cstdio>
#include "gromacs/mdtypes/fcdata.h"
#include "gromacs/mdtypes/forcerec.h"
#include "gromacs/mdtypes/interaction_const.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for forcetable.h:

Description

Old routines for table generation (will eventually be replaced)

Author
Erik Lindahl erik..nosp@m.lind.nosp@m.ahl@g.nosp@m.mail.nosp@m..com

Macros

#define GMX_MAKETABLES_FORCEUSER   (1<<0)
 Flag to select user tables for make_tables.
 
#define GMX_MAKETABLES_14ONLY   (1<<1)
 Flag to only make 1,4 pair tables for make_tables.
 

Typedefs

typedef double(* real_space_grid_contribution_computer )(double, double)
 Function pointer to calculate the grid contribution for coulomb/LJ. More...
 

Enumerations

enum  { etiCOUL, etiLJ6, etiLJ12, etiNR }
 Enumerated type to describe the interaction types in a table. More...
 

Functions

void table_spline3_fill_ewald_lr (real *table_F, real *table_V, real *table_FDV0, int ntab, double dx, real beta, real_space_grid_contribution_computer v_lr)
 Fill tables with the Ewald long-range force interaction. More...
 
real ewald_spline3_table_scale (const interaction_const_t *ic)
 Compute scaling for the Ewald quadratic spline tables. More...
 
double v_q_ewald_lr (double beta, double r)
 Return the real space grid contribution for Ewald. More...
 
double v_lj_ewald_lr (double beta, double r)
 Return the real space grid contribution for LJ-Ewald. More...
 
t_forcetable * make_tables (FILE *fp, const interaction_const_t *ic, const char *fn, real rtab, int flags)
 Return tables for inner loops. More...
 
bondedtable_t make_bonded_table (FILE *fplog, const char *fn, int angle)
 Return a table for bonded interactions,. More...
 
t_forcetable * makeDispersionCorrectionTable (FILE *fp, const interaction_const_t *ic, real rtab, const char *tabfn)
 Construct and return tabulated dispersion and repulsion interactions. More...
 

Typedef Documentation

typedef double(* real_space_grid_contribution_computer)(double, double)

Function pointer to calculate the grid contribution for coulomb/LJ.

Used to tell table_spline3_fill_ewald_lr whether it should calculate the grid contribution for electrostatics or LJ.

Enumeration Type Documentation

anonymous enum

Enumerated type to describe the interaction types in a table.

Enumerator
etiCOUL 

Coulomb.

etiLJ6 

Dispersion.

etiLJ12 

Repulsion.

etiNR 

Total number of interaction types.

Function Documentation

real ewald_spline3_table_scale ( const interaction_const_t *  ic)

Compute scaling for the Ewald quadratic spline tables.

Parameters
icPointer to interaction constant structure
Returns
The scaling factor
bondedtable_t make_bonded_table ( FILE *  fplog,
const char *  fn,
int  angle 
)

Return a table for bonded interactions,.

Parameters
fplogPointer to log file
fnFile name
angleType of angle: bonds 0, angles 1, dihedrals 2
Returns
New bonded table datatype
t_forcetable* make_tables ( FILE *  fp,
const interaction_const_t *  ic,
const char *  fn,
real  rtab,
int  flags 
)

Return tables for inner loops.

Parameters
fpLog file pointer
icNon-bonded interaction constants
fnFile name from which to read user tables
rtabLargest interaction distance to tabulate
flagsFlags to select table settings
Returns
Pointer to inner loop table structure
t_forcetable* makeDispersionCorrectionTable ( FILE *  fp,
const interaction_const_t *  ic,
real  rtab,
const char *  tabfn 
)

Construct and return tabulated dispersion and repulsion interactions.

This table can be used to compute long-range dispersion corrections

void table_spline3_fill_ewald_lr ( real table_F,
real table_V,
real table_FDV0,
int  ntab,
double  dx,
real  beta,
real_space_grid_contribution_computer  v_lr 
)

Fill tables with the Ewald long-range force interaction.

Fill tables of ntab points with spacing dr with the ewald long-range (mesh) force. There are three separate tables with format FDV0, F, and V. This function interpolates the Ewald mesh potential contribution with coefficient beta using a quadratic spline. The force can then be interpolated linearly.

Parameters
table_FForce table
table_VPotential table
table_FDV0Combined table optimized for SIMD loads
ntabNumber of points in tables
dxSpacing
betaEwald splitting paramter
v_lrPointer to function calculating real-space grid contribution
double v_lj_ewald_lr ( double  beta,
double  r 
)

Return the real space grid contribution for LJ-Ewald.

Parameters
betaEwald splitting parameter
rDistance for which to calculate the real-space contrib
Returns
Real space grid contribution for Ewald Lennard-Jones interaction
double v_q_ewald_lr ( double  beta,
double  r 
)

Return the real space grid contribution for Ewald.

Parameters
betaEwald splitting parameter
rDistance for which to calculate the real-space contrib
Returns
Real space grid contribution for Ewald electrostatics