37065e8be559164e7f12ac6204bc5e3322edde0e
[alexxy/gromacs.git] / src / gromacs / listed_forces / orires.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) 2010,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 \file
39  * \brief
40  * Declares functions for handling orientation restraints.
41  *
42  * \inlibraryapi
43  * \ingroup module_listed_forces
44  */
45 #ifndef GMX_LISTED_FORCES_ORIRES_H
46 #define GMX_LISTED_FORCES_ORIRES_H
47
48 #include <cstdio>
49
50 #include "gromacs/topology/ifunc.h"
51
52 struct gmx_mtop_t;
53 struct gmx_multisim_t;
54 class history_t;
55 struct t_inputrec;
56 struct t_pbc;
57 struct t_commrec;
58 struct t_oriresdata;
59 class t_state;
60
61 namespace gmx
62 {
63 template<typename>
64 class ArrayRef;
65 } // namespace gmx
66
67 /*! \brief
68  * Decides whether orientation restraints can work, and initializes
69  * all the orientation restraint stuff in *od (and assumes *od is
70  * already allocated.
71  * If orientation restraint are used, globalState is read and modified
72  * on the master rank (which is the only rank, since orientation
73  * restraints can not run in parallel).
74  */
75 void init_orires(FILE*                 fplog,
76                  const gmx_mtop_t*     mtop,
77                  const t_inputrec*     ir,
78                  const t_commrec*      cr,
79                  const gmx_multisim_t* ms,
80                  t_state*              globalState,
81                  t_oriresdata*         od);
82
83 /*! \brief
84  * Calculates the time averaged D matrices, the S matrix for each experiment.
85  *
86  * Returns the weighted RMS deviation of the orientation restraints.
87  */
88 real calc_orires_dev(const gmx_multisim_t*          ms,
89                      int                            nfa,
90                      const t_iatom                  fa[],
91                      const t_iparams                ip[],
92                      const t_mdatoms*               md,
93                      gmx::ArrayRef<const gmx::RVec> xWholeMolecules,
94                      const rvec                     x[],
95                      const t_pbc*                   pbc,
96                      t_oriresdata*                  oriresdata,
97                      history_t*                     hist);
98
99 /*! \brief
100  * Diagonalizes the order tensor(s) of the orienation restraints.
101  *
102  * For each experiment eig containts first 3 eigenvalues and then
103  * the 3 eigenvectors. The eigenvalues are ordered on magnitude.
104  */
105 void diagonalize_orires_tensors(t_oriresdata* od);
106
107 //! Prints order parameter, eigenvalues and eigenvectors to the log file.
108 void print_orires_log(FILE* log, t_oriresdata* od);
109
110 //! Calculates the orientation restraint forces.
111 real orires(int              nfa,
112             const t_iatom    forceatoms[],
113             const t_iparams  ip[],
114             const rvec       x[],
115             rvec4            f[],
116             rvec             fshift[],
117             const t_pbc*     pbc,
118             real             lambda,
119             real*            dvdlambda,
120             const t_mdatoms* md,
121             t_fcdata*        fcd,
122             int*             global_atom_index);
123
124 //! Copies the new time averages that have been calculated in calc_orires_dev.
125 void update_orires_history(const t_oriresdata& oriresdata, history_t* hist);
126
127 #endif