#include <cstdio>
#include <memory>
#include "gromacs/mdtypes/fcdata.h"
#include "gromacs/mdtypes/forcerec.h"
#include "gromacs/mdtypes/interaction_const.h"
#include "gromacs/utility/real.h"
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
|
EwaldCorrectionTables | generateEwaldCorrectionTables (int numPoints, double tableScaling, real beta, real_space_grid_contribution_computer v_lr) |
| Construct tables with the Ewald long-range force interaction. More...
|
|
real | ewald_spline3_table_scale (const interaction_const_t &ic, bool generateCoulombTables, bool generateVdwTables) |
| 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...
|
|
std::unique_ptr< 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 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.
Enumerated type to describe the interaction types in a table.
Enumerator |
---|
etiCOUL |
Coulomb.
|
etiLJ6 |
Dispersion.
|
etiLJ12 |
Repulsion.
|
etiNR |
Total number of interaction types.
|
real ewald_spline3_table_scale |
( |
const interaction_const_t & |
ic, |
|
|
bool |
generateCoulombTables, |
|
|
bool |
generateVdwTables |
|
) |
| |
Compute scaling for the Ewald quadratic spline tables.
- Parameters
-
ic | Pointer to interaction constant structure |
generateCoulombTables | Take the spacing for Coulomb Ewald corrections into account |
generateVdwTables | Take the spacing for Van der Waals Ewald corrections into account |
- Returns
- The scaling factor in units 1/nm
Construct tables with the Ewald long-range force interaction.
Creates and fills tables of numPoints points with the spacing set to 1/tableScaling with the Ewald long-range (mesh) force. There are three separate tables with format F, V, FDV0. This function interpolates the Ewald mesh potential contribution with coefficient beta using a quadratic spline. The force is then be interpolated linearly.
- Parameters
-
numPoints | Number of points in the tables |
tableScaling | 1/spacing, units 1/nm |
beta | Ewald splitting parameter, units 1/nm |
v_lr | Pointer to function calculating real-space grid contribution |
- Returns
- a set of Ewald correction tables
bondedtable_t make_bonded_table |
( |
FILE * |
fplog, |
|
|
const char * |
fn, |
|
|
int |
angle |
|
) |
| |
Return a table for bonded interactions,.
- Parameters
-
fplog | Pointer to log file |
fn | File name |
angle | Type 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
-
fp | Log file pointer |
ic | Non-bonded interaction constants |
fn | File name from which to read user tables |
rtab | Largest interaction distance to tabulate |
flags | Flags to select table settings |
- Returns
- Pointer to inner loop table structure
std::unique_ptr<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. Returns pointer owning nothing when tabfn=nullptr.
double v_lj_ewald_lr |
( |
double |
beta, |
|
|
double |
r |
|
) |
| |
Return the real space grid contribution for LJ-Ewald.
- Parameters
-
beta | Ewald splitting parameter |
r | Distance 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
-
beta | Ewald splitting parameter |
r | Distance for which to calculate the real-space contrib |
- Returns
- Real space grid contribution for Ewald electrostatics