7144dd0f125d71e60fd6f5501d1d97f64150c698
[alexxy/gromacs.git] / src / gromacs / listed_forces / disre.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,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,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 \file
39  * \brief
40  * Declares functions for handling distance restraints.
41  *
42  * \inlibraryapi
43  * \ingroup module_listed_forces
44  */
45 #ifndef GMX_LISTED_FORCES_DISRE_H
46 #define GMX_LISTED_FORCES_DISRE_H
47
48 #include <cstdio>
49
50 #include "gromacs/topology/ifunc.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/gmxmpi.h"
53
54 struct gmx_mtop_t;
55 struct gmx_multisim_t;
56 class history_t;
57 struct t_commrec;
58 struct t_disresdata;
59 struct t_fcdata;
60 struct t_inputrec;
61 struct t_pbc;
62 class t_state;
63 enum class DDRole;
64 enum class NumRanks;
65
66 //! Whether distance restraints are called from mdrun or from an analysis tool
67 enum class DisResRunMode
68 {
69     MDRun,
70     AnalysisTool
71 };
72
73 /*! \brief
74  * Initiates *disresdata.
75  *
76  * Must be called once, nbonds is the number
77  * of iatoms in the ilist of the idef struct.
78  * When time averaging is used, the history is initialized in state,
79  * unless it was read before from a checkpoint file.
80  * The implementation of distance restraints with -multidir
81  * must differ according to whether REMD is active.
82  */
83 void init_disres(FILE*                 fplog,
84                  const gmx_mtop_t*     mtop,
85                  t_inputrec*           ir,
86                  DisResRunMode         disResRunMode,
87                  DDRole                ddRole,
88                  NumRanks              numRanks,
89                  MPI_Comm              communicator,
90                  const gmx_multisim_t* ms,
91                  t_disresdata*         disresdata,
92                  t_state*              state,
93                  gmx_bool              bIsREMD);
94
95 /*! \brief
96  * Calculates r and r^-3 (inst. and time averaged) for all pairs
97  * and the ensemble averaged r^-6 (inst. and time averaged) for all restraints
98  */
99 void calc_disres_R_6(const t_commrec*      cr,
100                      const gmx_multisim_t* ms,
101                      int                   nfa,
102                      const t_iatom*        fa,
103                      const rvec*           x,
104                      const t_pbc*          pbc,
105                      t_disresdata*         disresdata,
106                      history_t*            hist);
107
108 //! Calculates the distance restraint forces, return the potential.
109 real ta_disres(int              nfa,
110                const t_iatom    forceatoms[],
111                const t_iparams  ip[],
112                const rvec       x[],
113                rvec4            f[],
114                rvec             fshift[],
115                const t_pbc*     pbc,
116                real             lambda,
117                real*            dvdlambda,
118                const t_mdatoms* md,
119                t_fcdata*        fcdata,
120                int*             global_atom_index);
121
122 //! Copies the new time averages that have been calculated in calc_disres_R_6.
123 void update_disres_history(const t_disresdata& disresdata, history_t* hist);
124
125 #endif