Remove dependence of constraints on t_mdatoms
[gromacs.git] / src / gromacs / mdlib / settle.h
blob85a499d112ee969e395d18d246834ed4cfb0447d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020, 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 SETTLE code.
38 * \author Berk Hess <hess@kth.se>
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_mdlib
41 * \inlibraryapi
44 #ifndef GMX_MDLIB_SETTLE_H
45 #define GMX_MDLIB_SETTLE_H
47 #include "gromacs/topology/idef.h"
49 struct gmx_mtop_t;
50 struct InteractionList;
51 struct t_inputrec;
52 struct t_pbc;
54 namespace gmx
57 template<typename>
58 class ArrayRef;
59 template<typename>
60 class ArrayRefWithPadding;
61 enum class ConstraintVariable : int;
63 /* Abstract type for SETTLE that is defined only in the file that uses it */
64 struct settledata;
66 /*! \brief Initializes and returns a structure with SETTLE parameters */
67 settledata* settle_init(const gmx_mtop_t& mtop);
69 //! Cleans up.
70 void settle_free(settledata* settled);
72 /*! \brief Set up the indices for the settle constraints */
73 void settle_set_constraints(settledata* settled,
74 const InteractionList& il_settle,
75 int numHomeAtoms,
76 const real* masses,
77 const real* inverseMasses);
79 /*! \brief Constrain coordinates using SETTLE.
80 * Can be called on any number of threads.
82 void csettle(settledata* settled, /* The SETTLE structure */
83 int nthread, /* The number of threads used */
84 int thread, /* Our thread index */
85 const t_pbc* pbc, /* PBC data pointer, can be NULL */
86 ArrayRefWithPadding<const RVec> x, /* Reference coordinates */
87 ArrayRefWithPadding<RVec> xprime, /* New coords, to be settled */
88 real invdt, /* 1/delta_t */
89 ArrayRefWithPadding<RVec> v, /* Also constrain v if v!=NULL */
90 bool bCalcVirial, /* Calculate the virial contribution */
91 tensor vir_r_m_dr, /* sum r x m delta_r */
92 bool* bErrorHasOccurred /* True if a settle error occurred */
95 /*! \brief Analytical algorithm to subtract the components of derivatives
96 * of coordinates working on settle type constraint.
98 void settle_proj(settledata* settled,
99 ConstraintVariable econq,
100 int nsettle,
101 const int iatoms[],
102 const t_pbc* pbc, /* PBC data pointer, can be NULL */
103 ArrayRef<const RVec> x,
104 ArrayRef<RVec> der,
105 ArrayRef<RVec> derp,
106 int CalcVirAtomEnd,
107 tensor vir_r_m_dder);
109 } // namespace gmx
111 #endif