Add cross-correlation as density simlarity measure
[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)",
63                   __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)
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                 printf("turning H bonds into constraints...\n");
88                 break;
89             case eshALLBONDS:
90                 printf("turning all bonds into constraints...\n");
91                 break;
92             case eshHANGLES:
93                 printf("turning all bonds and H angles into constraints...\n");
94                 break;
95             case eshALLANGLES:
96                 printf("turning all bonds and angles into constraints...\n");
97                 break;
98             default:
99                 gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
100         }
101
102         if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
103         {
104             /* Add all the angles with hydrogens to the shake list
105              * and remove them from the bond list
106              */
107             for (int ftype = 0; (ftype < F_NRE); ftype++)
108             {
109                 if (interaction_function[ftype].flags & IF_BTYPE)
110                 {
111                     InteractionsOfType *bonds = &(plist[ftype]);
112
113                     for (int ftype_a = 0; (gmx::ssize(*bonds) > 0 && ftype_a < F_NRE); ftype_a++)
114                     {
115                         if (interaction_function[ftype_a].flags & IF_ATYPE)
116                         {
117                             InteractionsOfType *pr = &(plist[ftype_a]);
118
119                             for (auto parm = pr->interactionTypes.begin(); parm != pr->interactionTypes.end(); )
120                             {
121                                 const InteractionOfType *ang = &(*parm);
122 #ifdef DEBUG
123                                 printf("Angle: %d-%d-%d\n", ang->ai(), ang->aj(), ang->ak());
124 #endif
125                                 int numhydrogens = count_hydrogens(info, 3, ang->atoms());
126                                 if ((nshake == eshALLANGLES) ||
127                                     (numhydrogens > 1) ||
128                                     (numhydrogens == 1 && toupper(**(info[ang->aj()])) == 'O'))
129                                 {
130                                     /* Can only add hydrogen angle shake, if the two bonds
131                                      * are constrained.
132                                      * append this angle to the shake list
133                                      */
134                                     std::vector<int> atomNumbers = {ang->ai(), ang->ak()};
135
136                                     /* Calculate length of constraint */
137                                     bool bFound = false;
138                                     b_ij   = b_jk = 0.0;
139                                     for (const auto &bond : bonds->interactionTypes)
140                                     {
141                                         if (((bond.ai() == ang->ai()) &&
142                                              (bond.aj() == ang->aj())) ||
143                                             ((bond.ai() == ang->aj()) &&
144                                              (bond.aj() == ang->ai())))
145                                         {
146                                             b_ij = bond.c0();
147                                         }
148                                         if (((bond.ai() == ang->ak()) &&
149                                              (bond.aj() == ang->aj())) ||
150                                             ((bond.ai() == ang->aj()) &&
151                                              (bond.aj() == ang->ak())))
152                                         {
153                                             b_jk = bond.c0();
154                                         }
155                                         bFound = (b_ij != 0.0) && (b_jk != 0.0);
156                                     }
157                                     if (bFound)
158                                     {
159                                         real              param = std::sqrt( b_ij*b_ij + b_jk*b_jk -
160                                                                              2.0*b_ij*b_jk*cos(DEG2RAD*ang->c0()));
161                                         std::vector<real> forceParm = {param, param};
162                                         /* apply law of cosines */
163 #ifdef DEBUG
164                                         printf("p: %d, q: %d, dist: %12.5e\n", atomNumbers[0],
165                                                atomNumbers[1], forceParm[0]);
166 #endif
167                                         add_param_to_list (&(plist[F_CONSTR]), InteractionOfType(atomNumbers, forceParm));
168                                         /* move the last bond to this position */
169                                         *parm = *(pr->interactionTypes.end() - 1);
170                                         pr->interactionTypes.erase(pr->interactionTypes.end() - 1);
171                                     }
172                                 }
173                                 else
174                                 {
175                                     ++parm;
176                                 }
177                             }
178                         } /* if IF_ATYPE */
179                     }     /* for ftype_A */
180                 }         /* if IF_BTYPE */
181             }             /* for ftype */
182         }                 /* if shake angles */
183
184         /* Add all the bonds with hydrogens to the shake list
185          * and remove them from the bond list
186          */
187         for (int ftype = 0; (ftype < F_NRE); ftype++)
188         {
189             if (interaction_function[ftype].flags & IF_BTYPE)
190             {
191                 InteractionsOfType *pr = &(plist[ftype]);
192                 for (auto parm = pr->interactionTypes.begin(); parm != pr->interactionTypes.end(); )
193                 {
194                     if ( (nshake != eshHBONDS) ||
195                          (count_hydrogens (info, 2, parm->atoms()) > 0) )
196                     {
197                         /* append this bond to the shake list */
198                         std::vector<int>  atomNumbers = {parm->ai(), parm->aj()};
199                         std::vector<real> forceParm   = { parm->c0(), parm->c2()};
200                         add_param_to_list (&(plist[F_CONSTR]), InteractionOfType(atomNumbers, forceParm));
201                         parm = pr->interactionTypes.erase(parm);
202                     }
203                     else
204                     {
205                         ++parm;
206                     }
207                 }
208             }
209         }
210     }
211 }