Fix random typos
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / add_par.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) 2011,2014,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "add_par.h"
42
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/gmxpreprocess/grompp_impl.h"
48 #include "gromacs/gmxpreprocess/notset.h"
49 #include "gromacs/gmxpreprocess/toputil.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
53
54 #include "hackblock.h"
55
56 void add_param(InteractionsOfType* ps, int ai, int aj, gmx::ArrayRef<const real> c, const char* s)
57 {
58     if ((ai < 0) || (aj < 0))
59     {
60         gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
61     }
62     std::vector<int>  atoms = { ai, aj };
63     std::vector<real> forceParm(c.begin(), c.end());
64
65     ps->interactionTypes.emplace_back(atoms, forceParm, s ? s : "");
66 }
67
68 void add_cmap_param(InteractionsOfType* ps, int ai, int aj, int ak, int al, int am, const char* s)
69 {
70     std::vector<int> atoms = { ai, aj, ak, al, am };
71     ps->interactionTypes.emplace_back(atoms, gmx::ArrayRef<const real>{}, s ? s : "");
72 }
73
74 void add_vsite2_param(InteractionsOfType* ps, int ai, int aj, int ak, real c0)
75 {
76     std::vector<int>  atoms     = { ai, aj, ak };
77     std::vector<real> forceParm = { c0 };
78     ps->interactionTypes.emplace_back(atoms, forceParm);
79 }
80
81 void add_vsite3_param(InteractionsOfType* ps, int ai, int aj, int ak, int al, real c0, real c1)
82 {
83     std::vector<int>  atoms     = { ai, aj, ak, al };
84     std::vector<real> forceParm = { c0, c1 };
85     ps->interactionTypes.emplace_back(atoms, forceParm);
86 }
87
88 void add_vsite3_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, bool bSwapParity)
89 {
90     std::vector<int> atoms = { ai, aj, ak, al };
91     ps->interactionTypes.emplace_back(atoms, gmx::ArrayRef<const real>{});
92
93     if (bSwapParity)
94     {
95         ps->interactionTypes.back().setForceParameter(1, -1);
96     }
97 }
98
99 void add_vsite4_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, int am)
100 {
101     std::vector<int> atoms = { ai, aj, ak, al, am };
102     ps->interactionTypes.emplace_back(atoms, gmx::ArrayRef<const real>{});
103 }
104
105 int search_jtype(const PreprocessResidue& localPpResidue, const char* name, bool bNterm)
106 {
107     int    niter, jmax;
108     size_t k, kmax, minstrlen;
109     char * rtpname, searchname[12];
110
111     strcpy(searchname, name);
112
113     /* Do a best match comparison */
114     /* for protein N-terminus, allow renaming of H1, H2 and H3 to H */
115     if (bNterm && (strlen(searchname) == 2) && (searchname[0] == 'H')
116         && ((searchname[1] == '1') || (searchname[1] == '2') || (searchname[1] == '3')))
117     {
118         niter = 2;
119     }
120     else
121     {
122         niter = 1;
123     }
124     kmax = 0;
125     jmax = -1;
126     for (int iter = 0; (iter < niter && jmax == -1); iter++)
127     {
128         if (iter == 1)
129         {
130             /* Try without the hydrogen number in the N-terminus */
131             searchname[1] = '\0';
132         }
133         for (int j = 0; (j < localPpResidue.natom()); j++)
134         {
135             rtpname = *(localPpResidue.atomname[j]);
136             if (gmx_strcasecmp(searchname, rtpname) == 0)
137             {
138                 jmax = j;
139                 kmax = strlen(searchname);
140                 break;
141             }
142             if (iter == niter - 1)
143             {
144                 minstrlen = std::min(strlen(searchname), strlen(rtpname));
145                 for (k = 0; k < minstrlen; k++)
146                 {
147                     if (searchname[k] != rtpname[k])
148                     {
149                         break;
150                     }
151                 }
152                 if (k > kmax)
153                 {
154                     kmax = k;
155                     jmax = j;
156                 }
157             }
158         }
159     }
160     if (jmax == -1)
161     {
162         gmx_fatal(FARGS,
163                   "Atom %s not found in rtp database in residue %s",
164                   searchname,
165                   localPpResidue.resname.c_str());
166     }
167     if (kmax != strlen(searchname))
168     {
169         gmx_fatal(FARGS,
170                   "Atom %s not found in rtp database in residue %s, "
171                   "it looks a bit like %s",
172                   searchname,
173                   localPpResidue.resname.c_str(),
174                   *(localPpResidue.atomname[jmax]));
175     }
176     return jmax;
177 }