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,2021, 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/math/utilities.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/logger.h"
56 #include "gromacs/utility/smalloc.h"
58 static int count_hydrogens(char*** atomname, int nra, gmx::ArrayRef<const int> a)
64 gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)", __FILE__, __LINE__);
68 for (int i = 0; (i < nra); i++)
70 if (toupper(**(atomname[a[i]])) == 'H')
79 void make_shake(gmx::ArrayRef<InteractionsOfType> plist, t_atoms* atoms, int nshake, const gmx::MDLogger& logger)
81 char*** info = atoms->atomname;
83 if (nshake != eshNONE)
90 .appendTextFormatted("turning H bonds into constraints...");
95 .appendTextFormatted("turning all bonds into constraints...");
100 .appendTextFormatted("turning all bonds and H angles into constraints...");
105 .appendTextFormatted("turning all bonds and angles into constraints...");
107 default: gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
110 if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
112 /* Add all the angles with hydrogens to the shake list
113 * and remove them from the bond list
115 for (int ftype = 0; (ftype < F_NRE); ftype++)
117 if (interaction_function[ftype].flags & IF_CHEMBOND)
119 InteractionsOfType* bonds = &(plist[ftype]);
121 for (int ftype_a = 0; (gmx::ssize(*bonds) > 0 && ftype_a < F_NRE); ftype_a++)
123 if (interaction_function[ftype_a].flags & IF_ATYPE)
125 InteractionsOfType* pr = &(plist[ftype_a]);
127 for (auto parm = pr->interactionTypes.begin();
128 parm != pr->interactionTypes.end();)
130 const InteractionOfType* ang = &(*parm);
134 .appendTextFormatted(
135 "Angle: %d-%d-%d", ang->ai(), ang->aj(), ang->ak());
137 int numhydrogens = count_hydrogens(info, 3, ang->atoms());
138 if ((nshake == eshALLANGLES) || (numhydrogens > 1)
139 || (numhydrogens == 1 && toupper(**(info[ang->aj()])) == 'O'))
141 /* Can only add hydrogen angle shake, if the two bonds
143 * append this angle to the shake list
145 std::vector<int> atomNumbers = { ang->ai(), ang->ak() };
147 /* Calculate length of constraint */
150 for (const auto& bond : bonds->interactionTypes)
152 if (((bond.ai() == ang->ai()) && (bond.aj() == ang->aj()))
153 || ((bond.ai() == ang->aj()) && (bond.aj() == ang->ai())))
157 if (((bond.ai() == ang->ak()) && (bond.aj() == ang->aj()))
158 || ((bond.ai() == ang->aj()) && (bond.aj() == ang->ak())))
162 bFound = (b_ij != 0.0) && (b_jk != 0.0);
166 real param = std::sqrt(
167 b_ij * b_ij + b_jk * b_jk
168 - 2.0 * b_ij * b_jk * cos(gmx::c_deg2Rad * ang->c0()));
169 std::vector<real> forceParm = { param, param };
170 if (ftype == F_CONNBONDS || ftype_a == F_CONNBONDS)
173 "Can not constrain all angles when they "
174 "involved bonds of type %s",
175 interaction_function[F_CONNBONDS].longname);
177 /* apply law of cosines */
181 .appendTextFormatted("p: %d, q: %d, dist: %12.5e",
186 add_param_to_list(&(plist[F_CONSTR]),
187 InteractionOfType(atomNumbers, forceParm));
188 /* move the last bond to this position */
189 *parm = *(pr->interactionTypes.end() - 1);
190 pr->interactionTypes.erase(pr->interactionTypes.end() - 1);
202 } /* if shake angles */
204 /* Add all the bonds with hydrogens to the shake list
205 * and remove them from the bond list
207 for (int ftype = 0; (ftype < F_NRE); ftype++)
209 if (interaction_function[ftype].flags & IF_BTYPE)
211 InteractionsOfType* pr = &(plist[ftype]);
212 for (auto parm = pr->interactionTypes.begin(); parm != pr->interactionTypes.end();)
214 if ((nshake != eshHBONDS) || (count_hydrogens(info, 2, parm->atoms()) > 0))
216 /* append this bond to the shake list */
217 std::vector<int> atomNumbers = { parm->ai(), parm->aj() };
218 std::vector<real> forceParm = { parm->c0(), parm->c2() };
219 add_param_to_list(&(plist[F_CONSTR]), InteractionOfType(atomNumbers, forceParm));
220 parm = pr->interactionTypes.erase(parm);