Refactor t_param to InteractionType
[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,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 "add_par.h"
41
42 #include <cstring>
43
44 #include <algorithm>
45
46 #include "gromacs/gmxpreprocess/grompp_impl.h"
47 #include "gromacs/gmxpreprocess/notset.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
52
53 #include "hackblock.h"
54
55 void add_param(InteractionTypeParameters *ps,
56                int                        ai,
57                int                        aj,
58                gmx::ArrayRef<const real>  c,
59                const char                *s)
60 {
61     if ((ai < 0) || (aj < 0))
62     {
63         gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
64     }
65     std::vector<int>  atoms = {ai, aj};
66     std::vector<real> forceParm(c.begin(), c.end());
67
68     ps->interactionTypes.emplace_back(InteractionType(atoms, forceParm, s ? s : ""));
69 }
70
71 void add_imp_param(InteractionTypeParameters *ps, int ai, int aj, int ak, int al, real c0, real c1,
72                    const char *s)
73 {
74     std::vector<int>  atoms     = {ai, aj, ak, al};
75     std::vector<real> forceParm = {c0, c1};
76     ps->interactionTypes.emplace_back(InteractionType(atoms, forceParm, s ? s : ""));
77 }
78
79 void add_dih_param(InteractionTypeParameters *ps, int ai, int aj, int ak, int al, real c0, real c1,
80                    real c2, const char *s)
81 {
82     std::vector<int>  atoms     = {ai, aj, ak, al};
83     std::vector<real> forceParm = {c0, c1, c2};
84     ps->interactionTypes.emplace_back(InteractionType(atoms, forceParm, s ? s : ""));
85 }
86
87 void add_cmap_param(InteractionTypeParameters *ps, int ai, int aj, int ak, int al, int am, const char *s)
88 {
89     std::vector<int> atoms = {ai, aj, ak, al, am};
90     ps->interactionTypes.emplace_back(InteractionType(atoms, {}, s ? s : ""));
91 }
92
93 void add_vsite2_atoms(InteractionTypeParameters *ps, int ai, int aj, int ak)
94 {
95     std::vector<int> atoms = {ai, aj, ak};
96     ps->interactionTypes.emplace_back(InteractionType(atoms, {}));
97 }
98
99 void add_vsite2_param(InteractionTypeParameters *ps, int ai, int aj, int ak, real c0)
100 {
101     std::vector<int>  atoms     = {ai, aj, ak};
102     std::vector<real> forceParm = {c0};
103     ps->interactionTypes.emplace_back(InteractionType(atoms, forceParm));
104 }
105
106 void add_vsite3_param(InteractionTypeParameters *ps, int ai, int aj, int ak, int al,
107                       real c0, real c1)
108 {
109     std::vector<int>  atoms     = {ai, aj, ak, al};
110     std::vector<real> forceParm = {c0, c1};
111     ps->interactionTypes.emplace_back(InteractionType(atoms, forceParm));
112 }
113
114 void add_vsite3_atoms(InteractionTypeParameters *ps, int ai, int aj, int ak, int al, bool bSwapParity)
115 {
116     std::vector<int>  atoms = {ai, aj, ak, al};
117     ps->interactionTypes.emplace_back(InteractionType(atoms, {}));
118
119     if (bSwapParity)
120     {
121         ps->interactionTypes.back().setForceParameter(1, -1);
122     }
123 }
124
125 void add_vsite4_atoms(InteractionTypeParameters *ps, int ai, int aj, int ak, int al, int am)
126 {
127     std::vector<int> atoms = {ai, aj, ak, al, am};
128     ps->interactionTypes.emplace_back(InteractionType(atoms, {}));
129 }
130
131 int search_jtype(const PreprocessResidue &localPpResidue, const char *name, bool bNterm)
132 {
133     int    niter, jmax;
134     size_t k, kmax, minstrlen;
135     char  *rtpname, searchname[12];
136
137     strcpy(searchname, name);
138
139     /* Do a best match comparison */
140     /* for protein N-terminus, allow renaming of H1, H2 and H3 to H */
141     if (bNterm && (strlen(searchname) == 2) && (searchname[0] == 'H') &&
142         ( (searchname[1] == '1') || (searchname[1] == '2') ||
143           (searchname[1] == '3') ) )
144     {
145         niter = 2;
146     }
147     else
148     {
149         niter = 1;
150     }
151     kmax = 0;
152     jmax = -1;
153     for (int iter = 0; (iter < niter && jmax == -1); iter++)
154     {
155         if (iter == 1)
156         {
157             /* Try without the hydrogen number in the N-terminus */
158             searchname[1] = '\0';
159         }
160         for (int j = 0; (j < localPpResidue.natom()); j++)
161         {
162             rtpname = *(localPpResidue.atomname[j]);
163             if (gmx_strcasecmp(searchname, rtpname) == 0)
164             {
165                 jmax = j;
166                 kmax = strlen(searchname);
167                 break;
168             }
169             if (iter == niter - 1)
170             {
171                 minstrlen = std::min(strlen(searchname), strlen(rtpname));
172                 for (k = 0; k < minstrlen; k++)
173                 {
174                     if (searchname[k] != rtpname[k])
175                     {
176                         break;
177                     }
178                 }
179                 if (k > kmax)
180                 {
181                     kmax = k;
182                     jmax = j;
183                 }
184             }
185         }
186     }
187     if (jmax == -1)
188     {
189         gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s",
190                   searchname, localPpResidue.resname.c_str());
191     }
192     if (kmax != strlen(searchname))
193     {
194         gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s, "
195                   "it looks a bit like %s",
196                   searchname, localPpResidue.resname.c_str(), *(localPpResidue.atomname[jmax]));
197     }
198     return jmax;
199 }