Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / topshake.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "topshake.h"
41
42 #include <cctype>
43 #include <cmath>
44
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"
56
57 static int count_hydrogens(char*** atomname, int nra, gmx::ArrayRef<const int> a)
58 {
59     int nh;
60
61     if (!atomname)
62     {
63         gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)", __FILE__, __LINE__);
64     }
65
66     nh = 0;
67     for (int i = 0; (i < nra); i++)
68     {
69         if (toupper(**(atomname[a[i]])) == 'H')
70         {
71             nh++;
72         }
73     }
74
75     return nh;
76 }
77
78 void make_shake(gmx::ArrayRef<InteractionsOfType> plist, t_atoms* atoms, int nshake, const gmx::MDLogger& logger)
79 {
80     char*** info = atoms->atomname;
81     real    b_ij, b_jk;
82     if (nshake != eshNONE)
83     {
84         switch (nshake)
85         {
86             case eshHBONDS:
87                 GMX_LOG(logger.info)
88                         .asParagraph()
89                         .appendTextFormatted("turning H bonds into constraints...");
90                 break;
91             case eshALLBONDS:
92                 GMX_LOG(logger.info)
93                         .asParagraph()
94                         .appendTextFormatted("turning all bonds into constraints...");
95                 break;
96             case eshHANGLES:
97                 GMX_LOG(logger.info)
98                         .asParagraph()
99                         .appendTextFormatted("turning all bonds and H angles into constraints...");
100                 break;
101             case eshALLANGLES:
102                 GMX_LOG(logger.info)
103                         .asParagraph()
104                         .appendTextFormatted("turning all bonds and angles into constraints...");
105                 break;
106             default: gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
107         }
108
109         if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
110         {
111             /* Add all the angles with hydrogens to the shake list
112              * and remove them from the bond list
113              */
114             for (int ftype = 0; (ftype < F_NRE); ftype++)
115             {
116                 if (interaction_function[ftype].flags & IF_CHEMBOND)
117                 {
118                     InteractionsOfType* bonds = &(plist[ftype]);
119
120                     for (int ftype_a = 0; (gmx::ssize(*bonds) > 0 && ftype_a < F_NRE); ftype_a++)
121                     {
122                         if (interaction_function[ftype_a].flags & IF_ATYPE)
123                         {
124                             InteractionsOfType* pr = &(plist[ftype_a]);
125
126                             for (auto parm = pr->interactionTypes.begin();
127                                  parm != pr->interactionTypes.end();)
128                             {
129                                 const InteractionOfType* ang = &(*parm);
130 #ifdef DEBUG
131                                 GMX_LOG(logger.info)
132                                         .asParagraph()
133                                         .appendTextFormatted(
134                                                 "Angle: %d-%d-%d", ang->ai(), ang->aj(), ang->ak());
135 #endif
136                                 int numhydrogens = count_hydrogens(info, 3, ang->atoms());
137                                 if ((nshake == eshALLANGLES) || (numhydrogens > 1)
138                                     || (numhydrogens == 1 && toupper(**(info[ang->aj()])) == 'O'))
139                                 {
140                                     /* Can only add hydrogen angle shake, if the two bonds
141                                      * are constrained.
142                                      * append this angle to the shake list
143                                      */
144                                     std::vector<int> atomNumbers = { ang->ai(), ang->ak() };
145
146                                     /* Calculate length of constraint */
147                                     bool bFound = false;
148                                     b_ij = b_jk = 0.0;
149                                     for (const auto& bond : bonds->interactionTypes)
150                                     {
151                                         if (((bond.ai() == ang->ai()) && (bond.aj() == ang->aj()))
152                                             || ((bond.ai() == ang->aj()) && (bond.aj() == ang->ai())))
153                                         {
154                                             b_ij = bond.c0();
155                                         }
156                                         if (((bond.ai() == ang->ak()) && (bond.aj() == ang->aj()))
157                                             || ((bond.ai() == ang->aj()) && (bond.aj() == ang->ak())))
158                                         {
159                                             b_jk = bond.c0();
160                                         }
161                                         bFound = (b_ij != 0.0) && (b_jk != 0.0);
162                                     }
163                                     if (bFound)
164                                     {
165                                         real param = std::sqrt(b_ij * b_ij + b_jk * b_jk
166                                                                - 2.0 * b_ij * b_jk
167                                                                          * cos(DEG2RAD * ang->c0()));
168                                         std::vector<real> forceParm = { param, param };
169                                         if (ftype == F_CONNBONDS || ftype_a == F_CONNBONDS)
170                                         {
171                                             gmx_fatal(FARGS,
172                                                       "Can not constrain all angles when they "
173                                                       "involved bonds of type %s",
174                                                       interaction_function[F_CONNBONDS].longname);
175                                         }
176                                         /* apply law of cosines */
177 #ifdef DEBUG
178                                         GMX_LOG(logger.info)
179                                                 .asParagraph()
180                                                 .appendTextFormatted("p: %d, q: %d, dist: %12.5e",
181                                                                      atomNumbers[0],
182                                                                      atomNumbers[1],
183                                                                      forceParm[0]);
184 #endif
185                                         add_param_to_list(&(plist[F_CONSTR]),
186                                                           InteractionOfType(atomNumbers, forceParm));
187                                         /* move the last bond to this position */
188                                         *parm = *(pr->interactionTypes.end() - 1);
189                                         pr->interactionTypes.erase(pr->interactionTypes.end() - 1);
190                                     }
191                                 }
192                                 else
193                                 {
194                                     ++parm;
195                                 }
196                             }
197                         } /* if IF_ATYPE */
198                     }     /* for ftype_A */
199                 }         /* if IF_BTYPE */
200             }             /* for ftype */
201         }                 /* if shake angles */
202
203         /* Add all the bonds with hydrogens to the shake list
204          * and remove them from the bond list
205          */
206         for (int ftype = 0; (ftype < F_NRE); ftype++)
207         {
208             if (interaction_function[ftype].flags & IF_BTYPE)
209             {
210                 InteractionsOfType* pr = &(plist[ftype]);
211                 for (auto parm = pr->interactionTypes.begin(); parm != pr->interactionTypes.end();)
212                 {
213                     if ((nshake != eshHBONDS) || (count_hydrogens(info, 2, parm->atoms()) > 0))
214                     {
215                         /* append this bond to the shake list */
216                         std::vector<int>  atomNumbers = { parm->ai(), parm->aj() };
217                         std::vector<real> forceParm   = { parm->c0(), parm->c2() };
218                         add_param_to_list(&(plist[F_CONSTR]), InteractionOfType(atomNumbers, forceParm));
219                         parm = pr->interactionTypes.erase(parm);
220                     }
221                     else
222                     {
223                         ++parm;
224                     }
225                 }
226             }
227         }
228     }
229 }