Refactor preprocessing atom types.
[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 void copy_bond (InteractionTypeParameters *pr, int to, int from)
57 /* copies an entry in a bond list to another position.
58  * does no allocing or freeing of memory
59  */
60 {
61     /*memcpy((char*) &(pr->param[to]),(char*) &(pr->param[from]),
62        (size_t)sizeof(pr->param[from]));*/
63     int i;
64
65     if (to != from)
66     {
67         range_check(to, 0, pr->nr);
68         range_check(from, 0, pr->nr);
69         for (i = 0; (i < MAXATOMLIST); i++)
70         {
71             pr->param[to].a[i] = pr->param[from].a[i];
72         }
73         for (i = 0; (i < MAXFORCEPARAM); i++)
74         {
75             pr->param[to].c[i] = pr->param[from].c[i];
76         }
77         for (i = 0; (i < MAXSLEN); i++)
78         {
79             pr->param[to].s[i] = pr->param[from].s[i];
80         }
81     }
82 }
83
84 static int count_hydrogens (char ***atomname, int nra, const int a[])
85 {
86     int  nh;
87
88     if (!atomname)
89     {
90         gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)",
91                   __FILE__, __LINE__);
92     }
93
94     nh = 0;
95     for (int i = 0; (i < nra); i++)
96     {
97         if (toupper(**(atomname[a[i]])) == 'H')
98         {
99             nh++;
100         }
101     }
102
103     return nh;
104 }
105
106 void make_shake(gmx::ArrayRef<InteractionTypeParameters> plist, t_atoms *atoms, int nshake)
107 {
108     char                  ***info = atoms->atomname;
109     real                     b_ij, b_jk;
110     int                      i, j;
111
112     if (nshake != eshNONE)
113     {
114         switch (nshake)
115         {
116             case eshHBONDS:
117                 printf("turning H bonds into constraints...\n");
118                 break;
119             case eshALLBONDS:
120                 printf("turning all bonds into constraints...\n");
121                 break;
122             case eshHANGLES:
123                 printf("turning all bonds and H angles into constraints...\n");
124                 break;
125             case eshALLANGLES:
126                 printf("turning all bonds and angles into constraints...\n");
127                 break;
128             default:
129                 gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
130         }
131
132         if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
133         {
134             /* Add all the angles with hydrogens to the shake list
135              * and remove them from the bond list
136              */
137             for (int ftype = 0; (ftype < F_NRE); ftype++)
138             {
139                 if (interaction_function[ftype].flags & IF_BTYPE)
140                 {
141                     InteractionTypeParameters *bonds = &(plist[ftype]);
142
143                     for (int ftype_a = 0; (bonds->nr > 0 && ftype_a < F_NRE); ftype_a++)
144                     {
145                         if (interaction_function[ftype_a].flags & IF_ATYPE)
146                         {
147                             InteractionTypeParameters *pr = &(plist[ftype_a]);
148
149                             for (int i = 0; (i < pr->nr); )
150                             {
151                                 t_param *ang = &(pr->param[i]);
152 #ifdef DEBUG
153                                 printf("Angle: %d-%d-%d\n", ang->ai(), ang->aj(), ang->ak());
154 #endif
155                                 int numhydrogens = count_hydrogens(info, 3, ang->a);
156                                 if ((nshake == eshALLANGLES) ||
157                                     (numhydrogens > 1) ||
158                                     (numhydrogens == 1 && toupper(**(info[ang->a[1]])) == 'O'))
159                                 {
160                                     /* Can only add hydrogen angle shake, if the two bonds
161                                      * are constrained.
162                                      * append this angle to the shake list
163                                      */
164                                     t_param p;
165                                     p.ai() = ang->ai();
166                                     p.aj() = ang->ak();
167
168                                     /* Calculate length of constraint */
169                                     bool bFound = false;
170                                     b_ij   = b_jk = 0.0;
171                                     for (int j = 0; !bFound && (j < bonds->nr); j++)
172                                     {
173                                         t_param *bond = &(bonds->param[j]);
174                                         if (((bond->ai() == ang->ai()) &&
175                                              (bond->aj() == ang->aj())) ||
176                                             ((bond->ai() == ang->aj()) &&
177                                              (bond->aj() == ang->ai())))
178                                         {
179                                             b_ij = bond->c0();
180                                         }
181                                         if (((bond->ai() == ang->ak()) &&
182                                              (bond->aj() == ang->aj())) ||
183                                             ((bond->ai() == ang->aj()) &&
184                                              (bond->aj() == ang->ak())))
185                                         {
186                                             b_jk = bond->c0();
187                                         }
188                                         bFound = (b_ij != 0.0) && (b_jk != 0.0);
189                                     }
190                                     if (bFound)
191                                     {
192                                         /* apply law of cosines */
193                                         p.c0() = std::sqrt( b_ij*b_ij + b_jk*b_jk -
194                                                             2.0*b_ij*b_jk*cos(DEG2RAD*ang->c0()) );
195                                         p.c1() = p.c0();
196 #ifdef DEBUG
197                                         printf("p: %d, q: %d, dist: %12.5e\n", p.ai(), p.aj(), p.c0());
198 #endif
199                                         add_param_to_list (&(plist[F_CONSTR]), &p);
200                                         /* move the last bond to this position */
201                                         copy_bond (pr, i, pr->nr-1);
202                                         /* should free memory here!! */
203                                         pr->nr--;
204                                     }
205                                 }
206                                 else
207                                 {
208                                     i++;
209                                 }
210                             }
211                         } /* if IF_ATYPE */
212                     }     /* for ftype_A */
213                 }         /* if IF_BTYPE */
214             }             /* for ftype */
215         }                 /* if shake angles */
216
217         /* Add all the bonds with hydrogens to the shake list
218          * and remove them from the bond list
219          */
220         for (int ftype = 0; (ftype < F_NRE); ftype++)
221         {
222             if (interaction_function[ftype].flags & IF_BTYPE)
223             {
224                 InteractionTypeParameters *pr = &(plist[ftype]);
225                 j  = 0;
226                 for (i = 0; i < pr->nr; i++)
227                 {
228                     if ( (nshake != eshHBONDS) ||
229                          (count_hydrogens (info, 2, pr->param[i].a) > 0) )
230                     {
231                         /* append this bond to the shake list */
232                         t_param p;
233                         p.ai() = pr->param[i].ai();
234                         p.aj() = pr->param[i].aj();
235                         p.c0() = pr->param[i].c0();
236                         p.c1() = pr->param[i].c2();
237                         add_param_to_list (&(plist[F_CONSTR]), &p);
238                     }
239                     else
240                     {
241                         copy_bond(pr, j++, i);
242                     }
243                 }
244                 pr->nr = j;
245             }
246         }
247     }
248 }