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