61669276407f61e036a5b7cc39e0b54b46479dda
[alexxy/gromacs.git] / src / gromacs / ewald / ewald.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \libinternal \defgroup module_ewald Ewald-family treatments of long-ranged forces
38  * \ingroup group_mdrun
39  *
40  * \brief Computes energies and forces for long-ranged interactions
41  * using the Ewald decomposition. Includes plain Ewald, PME, P3M for
42  * Coulomb, PME for Lennard-Jones, load-balancing for PME, and
43  * supporting code.
44  *
45  * \author Berk Hess <hess@kth.se>
46  * \author Erik Lindahl <erik@kth.se>
47  * \author Roland Schulz <roland@rschulz.eu>
48  * \author Mark Abraham <mark.j.abraham@gmail.com>
49  * \author Christan Wennberg <cwennberg@kth.se>
50  */
51
52 /*! \libinternal \file
53  *
54  * \brief This file contains function declarations necessary for
55  * computing energies and forces for the plain-Ewald long-ranged part,
56  * and the correction for overall system charge for all Ewald-family
57  * methods.
58  *
59  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
60  * \author Mark Abraham <mark.j.abraham@gmail.com>
61  * \inlibraryapi
62  * \ingroup module_ewald
63  */
64
65 #ifndef GMX_EWALD_EWALD_H
66 #define GMX_EWALD_EWALD_H
67
68 #include <stdio.h>
69
70 #include "gromacs/math/vectypes.h"
71 #include "gromacs/utility/real.h"
72
73 struct t_commrec;
74 struct t_forcerec;
75 struct t_inputrec;
76
77 /* Forward declaration of type for managing Ewald tables */
78 struct gmx_ewald_tab_t;
79
80 /*! \brief Initialize the tables used in the Ewald long-ranged part */
81 void init_ewald_tab(struct gmx_ewald_tab_t** et, const t_inputrec* ir, FILE* fp);
82
83 /*! \brief Do the long-ranged part of an Ewald calculation */
84 real do_ewald(const t_inputrec* ir,
85               const rvec        x[],
86               rvec              f[],
87               const real        chargeA[],
88               const real        chargeB[],
89               const matrix      box,
90               const t_commrec*  cr,
91               int               natoms,
92               matrix            lrvir,
93               real              ewaldcoeff,
94               real              lambda,
95               real*             dvdlambda,
96               gmx_ewald_tab_t*  et);
97
98 /*! \brief Calculate the correction to the Ewald sum, due to a net system
99  * charge.
100  *
101  * Should only be called on one thread. */
102 real ewald_charge_correction(const t_commrec*  cr,
103                              const t_forcerec* fr,
104                              real              lambda,
105                              const matrix      box,
106                              real*             dvdlambda,
107                              tensor            vir);
108
109 #endif