Merge branch release-2018
[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, 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
36 #ifndef GMX_MDLIB_SHAKE_H
37 #define GMX_MDLIB_SHAKE_H
38
39 #include "gromacs/mdlib/constr.h"
40 #include "gromacs/topology/idef.h"
41
42 struct t_inputrec;
43
44 /* Abstract type for SHAKE that is defined only in the file that uses it */
45 typedef struct gmx_shakedata *gmx_shakedata_t;
46
47 gmx_shakedata_t shake_init();
48 /* Initializes and return the SHAKE data structure */
49
50 gmx_bool bshakef(FILE           *log,          /* Log file                      */
51                  gmx_shakedata_t shaked,       /* Total number of atoms */
52                  real            invmass[],    /* Atomic masses         */
53                  int             nblocks,      /* The number of shake blocks    */
54                  int             sblock[],     /* The shake blocks             */
55                  t_idef         *idef,         /* The interaction def           */
56                  t_inputrec     *ir,           /* Input record                  */
57                  rvec            x_s[],        /* Coords before update          */
58                  rvec            prime[],      /* Output coords         */
59                  t_nrnb         *nrnb,         /* Performance measure          */
60                  real * const    lagr,         /* The Lagrange multipliers     */
61                  real            lambda,       /* FEP lambda                   */
62                  real           *dvdlambda,    /* FEP force                    */
63                  real            invdt,        /* 1/delta_t                    */
64                  rvec           *v,            /* Also constrain v if v!=NULL  */
65                  gmx_bool        bCalcVir,     /* Calculate r x m delta_r      */
66                  tensor          vir_r_m_dr,   /* sum r x m delta_r            */
67                  gmx_bool        bDumpOnError, /* Dump debugging stuff on error*/
68                  int             econq);       /* which type of constraint is occurring */
69 /* Shake all the atoms blockwise. It is assumed that all the constraints
70  * in the idef->shakes field are sorted, to ascending block nr. The
71  * sblock array points into the idef->shakes.iatoms field, with block 0
72  * starting
73  * at sblock[0] and running to ( < ) sblock[1], block n running from
74  * sblock[n] to sblock[n+1]. Array sblock should be large enough.
75  * Return TRUE when OK, FALSE when shake-error
76  */
77
78 void cshake(const int iatom[], int ncon, int *nnit, int maxnit,
79             const real dist2[], real xp[], const real rij[], const real m2[], real omega,
80             const real invmass[], const real tt[], real lagr[], int *nerror);
81 /* Regular iterative shake */
82
83 #endif