2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \libinternal \file
36 * \brief Declares interface to SHAKE code.
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
45 #ifndef GMX_MDLIB_SHAKE_H
46 #define GMX_MDLIB_SHAKE_H
48 #include "gromacs/topology/block.h"
49 #include "gromacs/topology/idef.h"
59 enum class ConstraintVariable
: int;
61 /* Abstract type for SHAKE that is defined only in the file that uses it */
64 /*! \brief Initializes and return the SHAKE data structure */
65 shakedata
*shake_init();
67 //! Make SHAKE blocks when not using DD.
69 make_shake_sblock_serial(shakedata
*shaked
,
70 const t_idef
*idef
, const t_mdatoms
&md
);
72 //! Make SHAKE blocks when using DD.
74 make_shake_sblock_dd(shakedata
*shaked
,
75 const t_ilist
*ilcon
, const t_block
*cgs
,
76 const gmx_domdec_t
*dd
);
78 /*! \brief Shake all the atoms blockwise. It is assumed that all the constraints
79 * in the idef->shakes field are sorted, to ascending block nr. The
80 * sblock array points into the idef->shakes.iatoms field, with block 0
82 * at sblock[0] and running to ( < ) sblock[1], block n running from
83 * sblock[n] to sblock[n+1]. Array sblock should be large enough.
84 * Return TRUE when OK, FALSE when shake-error
87 constrain_shake(FILE *log
, /* Log file */
88 shakedata
*shaked
, /* Total number of atoms */
89 const real invmass
[], /* Atomic masses */
90 const t_idef
&idef
, /* The interaction def */
91 const t_inputrec
&ir
, /* Input record */
92 const rvec x_s
[], /* Coords before update */
93 rvec xprime
[], /* Output coords when constraining x */
94 rvec vprime
[], /* Output coords when constraining v */
95 t_nrnb
*nrnb
, /* Performance measure */
96 real lambda
, /* FEP lambda */
97 real
*dvdlambda
, /* FEP force */
98 real invdt
, /* 1/delta_t */
99 rvec
*v
, /* Also constrain v if v!=NULL */
100 bool bCalcVir
, /* Calculate r x m delta_r */
101 tensor vir_r_m_dr
, /* sum r x m delta_r */
102 bool bDumpOnError
, /* Dump debugging stuff on error*/
103 ConstraintVariable econq
); /* which type of constraint is occurring */
105 /*! \brief Regular iterative shake */
106 void cshake(const int iatom
[], int ncon
, int *nnit
, int maxnit
,
107 const real dist2
[], real xp
[], const real rij
[], const real m2
[], real omega
,
108 const real invmass
[], const real tt
[], real lagr
[], int *nerror
);