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