#include "gmxpre.h"
#include "ewald.h"
#include <cmath>
#include <cstdio>
#include <cstdlib>
#include <algorithm>
#include <array>
#include <filesystem>
#include "gromacs/ewald/ewald_utils.h"
#include "gromacs/math/functions.h"
#include "gromacs/math/gmxcomplex.h"
#include "gromacs/math/units.h"
#include "gromacs/math/utilities.h"
#include "gromacs/math/vec.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/mdtypes/commrec.h"
#include "gromacs/mdtypes/forcerec.h"
#include "gromacs/mdtypes/inputrec.h"
#include "gromacs/mdtypes/interaction_const.h"
#include "gromacs/mdtypes/md_enums.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/smalloc.h"
This file contains function definitions necessary for computing energies and forces for the plain-Ewald long-ranged part, and the correction for overall system charge for all Ewald-family methods.
- Author
- David van der Spoel david.nosp@m..van.nosp@m.dersp.nosp@m.oel@.nosp@m.icm.u.nosp@m.u.se
-
Mark Abraham mark..nosp@m.j.ab.nosp@m.raham.nosp@m.@gma.nosp@m.il.co.nosp@m.m
|
using | cvec = std::array< t_complex, DIM > |
|
|
static void | calc_lll (const rvec box, rvec lll) |
| Calculates wave vectors.
|
|
static void | tabulateStructureFactors (int natom, gmx::ArrayRef< const gmx::RVec > x, int kmax, cvec **eir, const rvec lll) |
| Make tables for the structure factor parts.
|
|
real | do_ewald (bool havePbcXY2Walls, real wallEwaldZfac, real epsilonR, FreeEnergyPerturbationType freeEnergyPerturbationType, gmx::ArrayRef< const gmx::RVec > coords, gmx::ArrayRef< gmx::RVec > forces, gmx::ArrayRef< const real > chargeA, gmx::ArrayRef< const real > chargeB, const matrix box, const t_commrec *commrec, int natoms, matrix lrvir, real ewaldcoeff, real lambda, real *dvdlambda, gmx_ewald_tab_t *et) |
| Do the long-ranged part of an Ewald calculation.
|
|
real | ewald_charge_correction (const t_commrec *commrec, const real epsilonR, const real ewaldcoeffQ, gmx::ArrayRef< const double > qsum, const real lambda, const matrix box, real *dvdlambda, tensor vir) |
| Calculate the correction to the Ewald sum, due to a net system charge. More...
|
|
real ewald_charge_correction |
( |
const t_commrec * |
commrec, |
|
|
real |
epsilonR, |
|
|
real |
ewaldcoeffQ, |
|
|
gmx::ArrayRef< const double > |
qsum, |
|
|
real |
lambda, |
|
|
const matrix |
box, |
|
|
real * |
dvdlambda, |
|
|
tensor |
vir |
|
) |
| |
Calculate the correction to the Ewald sum, due to a net system charge.
Should only be called on one thread.