2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2018,2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /* This file is completely threadsafe - keep it that way! */
45 #include "gromacs/gmxpreprocess/grompp_impl.h"
46 #include "gromacs/gmxpreprocess/notset.h"
47 #include "gromacs/gmxpreprocess/readir.h"
48 #include "gromacs/gmxpreprocess/topdirs.h"
49 #include "gromacs/gmxpreprocess/toppush.h"
50 #include "gromacs/gmxpreprocess/toputil.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/logger.h"
55 #include "gromacs/utility/smalloc.h"
57 static int count_hydrogens(char*** atomname
, int nra
, gmx::ArrayRef
<const int> a
)
63 gmx_fatal(FARGS
, "Cannot call count_hydrogens with no atomname (%s %d)", __FILE__
, __LINE__
);
67 for (int i
= 0; (i
< nra
); i
++)
69 if (toupper(**(atomname
[a
[i
]])) == 'H')
78 void make_shake(gmx::ArrayRef
<InteractionsOfType
> plist
, t_atoms
* atoms
, int nshake
, const gmx::MDLogger
& logger
)
80 char*** info
= atoms
->atomname
;
82 if (nshake
!= eshNONE
)
89 .appendTextFormatted("turning H bonds into constraints...");
94 .appendTextFormatted("turning all bonds into constraints...");
99 .appendTextFormatted("turning all bonds and H angles into constraints...");
104 .appendTextFormatted("turning all bonds and angles into constraints...");
106 default: gmx_fatal(FARGS
, "Invalid option for make_shake (%d)", nshake
);
109 if ((nshake
== eshHANGLES
) || (nshake
== eshALLANGLES
))
111 /* Add all the angles with hydrogens to the shake list
112 * and remove them from the bond list
114 for (int ftype
= 0; (ftype
< F_NRE
); ftype
++)
116 if (interaction_function
[ftype
].flags
& IF_CHEMBOND
)
118 InteractionsOfType
* bonds
= &(plist
[ftype
]);
120 for (int ftype_a
= 0; (gmx::ssize(*bonds
) > 0 && ftype_a
< F_NRE
); ftype_a
++)
122 if (interaction_function
[ftype_a
].flags
& IF_ATYPE
)
124 InteractionsOfType
* pr
= &(plist
[ftype_a
]);
126 for (auto parm
= pr
->interactionTypes
.begin();
127 parm
!= pr
->interactionTypes
.end();)
129 const InteractionOfType
* ang
= &(*parm
);
133 .appendTextFormatted("Angle: %d-%d-%d", ang
->ai(),
134 ang
->aj(), ang
->ak());
136 int numhydrogens
= count_hydrogens(info
, 3, ang
->atoms());
137 if ((nshake
== eshALLANGLES
) || (numhydrogens
> 1)
138 || (numhydrogens
== 1 && toupper(**(info
[ang
->aj()])) == 'O'))
140 /* Can only add hydrogen angle shake, if the two bonds
142 * append this angle to the shake list
144 std::vector
<int> atomNumbers
= { ang
->ai(), ang
->ak() };
146 /* Calculate length of constraint */
149 for (const auto& bond
: bonds
->interactionTypes
)
151 if (((bond
.ai() == ang
->ai()) && (bond
.aj() == ang
->aj()))
152 || ((bond
.ai() == ang
->aj()) && (bond
.aj() == ang
->ai())))
156 if (((bond
.ai() == ang
->ak()) && (bond
.aj() == ang
->aj()))
157 || ((bond
.ai() == ang
->aj()) && (bond
.aj() == ang
->ak())))
161 bFound
= (b_ij
!= 0.0) && (b_jk
!= 0.0);
165 real param
= std::sqrt(b_ij
* b_ij
+ b_jk
* b_jk
167 * cos(DEG2RAD
* ang
->c0()));
168 std::vector
<real
> forceParm
= { param
, param
};
169 if (ftype
== F_CONNBONDS
|| ftype_a
== F_CONNBONDS
)
172 "Can not constrain all angles when they "
173 "involved bonds of type %s",
174 interaction_function
[F_CONNBONDS
].longname
);
176 /* apply law of cosines */
180 .appendTextFormatted("p: %d, q: %d, dist: %12.5e",
181 atomNumbers
[0], atomNumbers
[1],
184 add_param_to_list(&(plist
[F_CONSTR
]),
185 InteractionOfType(atomNumbers
, forceParm
));
186 /* move the last bond to this position */
187 *parm
= *(pr
->interactionTypes
.end() - 1);
188 pr
->interactionTypes
.erase(pr
->interactionTypes
.end() - 1);
200 } /* if shake angles */
202 /* Add all the bonds with hydrogens to the shake list
203 * and remove them from the bond list
205 for (int ftype
= 0; (ftype
< F_NRE
); ftype
++)
207 if (interaction_function
[ftype
].flags
& IF_BTYPE
)
209 InteractionsOfType
* pr
= &(plist
[ftype
]);
210 for (auto parm
= pr
->interactionTypes
.begin(); parm
!= pr
->interactionTypes
.end();)
212 if ((nshake
!= eshHBONDS
) || (count_hydrogens(info
, 2, parm
->atoms()) > 0))
214 /* append this bond to the shake list */
215 std::vector
<int> atomNumbers
= { parm
->ai(), parm
->aj() };
216 std::vector
<real
> forceParm
= { parm
->c0(), parm
->c2() };
217 add_param_to_list(&(plist
[F_CONSTR
]), InteractionOfType(atomNumbers
, forceParm
));
218 parm
= pr
->interactionTypes
.erase(parm
);