Add second LINCS atom update task
[alexxy/gromacs.git] / src / gromacs / mdlib / settle.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020,2021, 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 SETTLE code.
37  *
38  * \author Berk Hess <hess@kth.se>
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  * \ingroup module_mdlib
41  * \inlibraryapi
42  */
43
44 #ifndef GMX_MDLIB_SETTLE_H
45 #define GMX_MDLIB_SETTLE_H
46
47 #include "gromacs/topology/idef.h"
48 #include "gromacs/utility/alignedallocator.h"
49
50 struct gmx_mtop_t;
51 struct InteractionList;
52 struct t_inputrec;
53 struct t_pbc;
54
55 namespace gmx
56 {
57
58 template<typename>
59 class ArrayRef;
60 template<typename>
61 class ArrayRefWithPadding;
62 enum class ConstraintVariable : int;
63
64 /* \brief Parameters for SETTLE algorithm.
65  *
66  * \todo Remove duplicates, check if recomputing makes more sense in some cases.
67  * \todo Move the projection parameters into separate structure.
68  */
69 struct SettleParameters
70 {
71     //! Mass of oxygen atom
72     real mO;
73     //! Mass of hydrogen atom
74     real mH;
75     //! Relative hydrogen mass (i.e. mH/(mO+2*mH))
76     real wh;
77     //! Target distance between oxygen and hydrogen
78     real dOH;
79     //! Target distance between hydrogen atoms
80     real dHH;
81     //! Double relative H mass times height of the H-O-H triangle.
82     real ra;
83     //! Height of the H-O-H triangle minus ra.
84     real rb;
85     //! Half the H-H distance.
86     real rc;
87     //! Reciprocal H-H distance
88     real irc2;
89
90     /* Parameters below are used for projection */
91     //! Reciprocal oxygen mass
92     real imO;
93     //! Reciprocal hydrogen mass
94     real imH;
95     //! Reciprocal O-H distance
96     real invdOH;
97     //! Reciprocal H-H distance (again)
98     real invdHH;
99     //! Reciprocal projection matrix
100     matrix invmat;
101 };
102
103 /*! \brief Computes and returns settle parameters.
104  *
105  * \param[in]  mO     Mass of oxygen atom
106  * \param[in]  mH     Mass of hydrogen atom
107  * \param[in]  invmO  Reciprocal mass of oxygen atom
108  * \param[in]  invmH  Reciprocal mass of hydrogen atom
109  * \param[in]  dOH    Target O-H bond length
110  * \param[in]  dHH    Target H-H bond length
111  */
112 SettleParameters settleParameters(real mO, real mH, real dOH, real invmO, real invmH, real dHH);
113
114 /*! \libinternal
115  * \brief Data for executing SETTLE constraining
116  */
117 class SettleData
118 {
119 public:
120     //! Constructor
121     SettleData(const gmx_mtop_t& mtop);
122
123     //! Sets the constraints from the interaction list and the masses
124     void setConstraints(const InteractionList&    il_settle,
125                         int                       numHomeAtoms,
126                         gmx::ArrayRef<const real> masses,
127                         gmx::ArrayRef<const real> inverseMasses);
128
129     //! Returns settle parameters for constraining coordinates and forces
130     const SettleParameters& parametersMassWeighted() const { return parametersMassWeighted_; }
131
132     //! Returns settle parameters for constraining forces: all masses are set to 1
133     const SettleParameters& parametersAllMasses1() const { return parametersAllMasses1_; }
134
135     //! Returns the number of SETTLEs
136     int numSettles() const { return numSettles_; }
137
138     //! Returns a pointer to the indices of oxygens for each SETTLE
139     const int* ow1() const { return ow1_.data(); }
140     //! Returns a pointer to the indices of the first hydrogen for each SETTLE
141     const int* hw2() const { return hw2_.data(); }
142     //! Returns a pointer to the indices of the second hydrogen for each SETTLE
143     const int* hw3() const { return hw3_.data(); }
144     //! Returns a pointer to the virial contribution factor (either 1 or 0) for each SETTLE
145     const real* virfac() const { return virfac_.data(); }
146
147     //! Returns whether we should use SIMD intrinsics code
148     bool useSimd() const { return useSimd_; }
149
150 private:
151     //! Parameters for SETTLE for coordinates
152     SettleParameters parametersMassWeighted_;
153     //! Parameters with all masses 1, for forces
154     SettleParameters parametersAllMasses1_;
155
156     //! The number of settles on our rank
157     int numSettles_;
158
159     //! Index to OW1 atoms, size numSettles_ + SIMD padding
160     std::vector<int, gmx::AlignedAllocator<int>> ow1_;
161     //! Index to HW2 atoms, size numSettles_ + SIMD padding
162     std::vector<int, gmx::AlignedAllocator<int>> hw2_;
163     //! Index to HW3 atoms, size numSettles_ + SIMD padding
164     std::vector<int, gmx::AlignedAllocator<int>> hw3_;
165     //! Virial factor 0 or 1, size numSettles_ + SIMD padding
166     std::vector<real, gmx::AlignedAllocator<real>> virfac_;
167
168     //! Tells whether we will use SIMD intrinsics code
169     bool useSimd_;
170 };
171
172 /*! \brief Constrain coordinates using SETTLE.
173  * Can be called on any number of threads.
174  */
175 void csettle(const SettleData&               settled,     /* The SETTLE structure */
176              int                             nthread,     /* The number of threads used */
177              int                             thread,      /* Our thread index */
178              const t_pbc*                    pbc,         /* PBC data pointer, can be NULL */
179              ArrayRefWithPadding<const RVec> x,           /* Reference coordinates */
180              ArrayRefWithPadding<RVec>       xprime,      /* New coords, to be settled */
181              real                            invdt,       /* 1/delta_t */
182              ArrayRefWithPadding<RVec>       v,           /* Also constrain v if v!=NULL */
183              bool                            bCalcVirial, /* Calculate the virial contribution */
184              tensor                          vir_r_m_dr,  /* sum r x m delta_r */
185              bool*                           bErrorHasOccurred /* True if a settle error occurred */
186 );
187
188 /*! \brief Analytical algorithm to subtract the components of derivatives
189  * of coordinates working on settle type constraint.
190  */
191 void settle_proj(const SettleData&    settled,
192                  ConstraintVariable   econq,
193                  int                  nsettle,
194                  const int            iatoms[],
195                  const t_pbc*         pbc, /* PBC data pointer, can be NULL  */
196                  ArrayRef<const RVec> x,
197                  ArrayRef<RVec>       der,
198                  ArrayRef<RVec>       derp,
199                  int                  CalcVirAtomEnd,
200                  tensor               vir_r_m_dder);
201
202 } // namespace gmx
203
204 #endif