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