1249d3c1f577094c7bba5e64d229b029d65bf627
[alexxy/gromacs.git] / src / gromacs / mdlib / shake.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \libinternal \file
36  * \brief Declares interface to SHAKE code.
37  *
38  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
39  * \author Berk Hess <hess@kth.se>
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_mdlib
42  * \inlibraryapi
43  */
44
45 #ifndef GMX_MDLIB_SHAKE_H
46 #define GMX_MDLIB_SHAKE_H
47
48 #include "gromacs/topology/block.h"
49 #include "gromacs/topology/idef.h"
50
51 struct gmx_domdec_t;
52 struct t_inputrec;
53 struct t_mdatoms;
54 struct t_nrnb;
55
56 namespace gmx
57 {
58
59 enum class ConstraintVariable : int;
60
61 /* Abstract type for SHAKE that is defined only in the file that uses it */
62 struct shakedata;
63
64 /*! \brief Initializes and return the SHAKE data structure */
65 shakedata* shake_init();
66
67 //! Destroy SHAKE. Needed to solve memory leaks.
68 void done_shake(shakedata* d);
69
70 //! Make SHAKE blocks when not using DD.
71 void make_shake_sblock_serial(shakedata* shaked, const t_idef* idef, const t_mdatoms& md);
72
73 //! Make SHAKE blocks when using DD.
74 void make_shake_sblock_dd(shakedata* shaked, const t_ilist* ilcon, const gmx_domdec_t* dd);
75
76 /*! \brief Shake all the atoms blockwise. It is assumed that all the constraints
77  * in the idef->shakes field are sorted, to ascending block nr. The
78  * sblock array points into the idef->shakes.iatoms field, with block 0
79  * starting
80  * at sblock[0] and running to ( < ) sblock[1], block n running from
81  * sblock[n] to sblock[n+1]. Array sblock should be large enough.
82  * Return TRUE when OK, FALSE when shake-error
83  */
84 bool constrain_shake(FILE*              log,          /* Log file                       */
85                      shakedata*         shaked,       /* Total number of atoms  */
86                      const real         invmass[],    /* Atomic masses          */
87                      const t_idef&      idef,         /* The interaction def            */
88                      const t_inputrec&  ir,           /* Input record                   */
89                      const rvec         x_s[],        /* Coords before update           */
90                      rvec               xprime[],     /* Output coords when constraining x */
91                      rvec               vprime[],     /* Output coords when constraining v */
92                      t_nrnb*            nrnb,         /* Performance measure          */
93                      real               lambda,       /* FEP lambda                   */
94                      real*              dvdlambda,    /* FEP force                    */
95                      real               invdt,        /* 1/delta_t                    */
96                      rvec*              v,            /* Also constrain v if v!=NULL  */
97                      bool               bCalcVir,     /* Calculate r x m delta_r      */
98                      tensor             vir_r_m_dr,   /* sum r x m delta_r            */
99                      bool               bDumpOnError, /* Dump debugging stuff on error*/
100                      ConstraintVariable econq);       /* which type of constraint is occurring */
101
102 /*! \brief Regular iterative shake */
103 void cshake(const int  iatom[],
104             int        ncon,
105             int*       nnit,
106             int        maxnit,
107             const real dist2[],
108             real       xp[],
109             const real rij[],
110             const real m2[],
111             real       omega,
112             const real invmass[],
113             const real tt[],
114             real       lagr[],
115             int*       nerror);
116
117 } // namespace gmx
118
119 #endif