bffcbdf8755ab2605a18c9a48086f3e944d06bc8
[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(InteractionsOfType* ps, int ai, int aj, gmx::ArrayRef<const real> c, const char* s)
56 {
57     if ((ai < 0) || (aj < 0))
58     {
59         gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
60     }
61     std::vector<int>  atoms = { ai, aj };
62     std::vector<real> forceParm(c.begin(), c.end());
63
64     ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm, s ? s : ""));
65 }
66
67 void add_cmap_param(InteractionsOfType* ps, int ai, int aj, int ak, int al, int am, const char* s)
68 {
69     std::vector<int> atoms = { ai, aj, ak, al, am };
70     ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}, s ? s : ""));
71 }
72
73 void add_vsite2_param(InteractionsOfType* ps, int ai, int aj, int ak, real c0)
74 {
75     std::vector<int>  atoms     = { ai, aj, ak };
76     std::vector<real> forceParm = { c0 };
77     ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm));
78 }
79
80 void add_vsite3_param(InteractionsOfType* ps, int ai, int aj, int ak, int al, real c0, real c1)
81 {
82     std::vector<int>  atoms     = { ai, aj, ak, al };
83     std::vector<real> forceParm = { c0, c1 };
84     ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm));
85 }
86
87 void add_vsite3_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, bool bSwapParity)
88 {
89     std::vector<int> atoms = { ai, aj, ak, al };
90     ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
91
92     if (bSwapParity)
93     {
94         ps->interactionTypes.back().setForceParameter(1, -1);
95     }
96 }
97
98 void add_vsite4_atoms(InteractionsOfType* ps, int ai, int aj, int ak, int al, int am)
99 {
100     std::vector<int> atoms = { ai, aj, ak, al, am };
101     ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
102 }
103
104 int search_jtype(const PreprocessResidue& localPpResidue, const char* name, bool bNterm)
105 {
106     int    niter, jmax;
107     size_t k, kmax, minstrlen;
108     char * rtpname, searchname[12];
109
110     strcpy(searchname, name);
111
112     /* Do a best match comparison */
113     /* for protein N-terminus, allow renaming of H1, H2 and H3 to H */
114     if (bNterm && (strlen(searchname) == 2) && (searchname[0] == 'H')
115         && ((searchname[1] == '1') || (searchname[1] == '2') || (searchname[1] == '3')))
116     {
117         niter = 2;
118     }
119     else
120     {
121         niter = 1;
122     }
123     kmax = 0;
124     jmax = -1;
125     for (int iter = 0; (iter < niter && jmax == -1); iter++)
126     {
127         if (iter == 1)
128         {
129             /* Try without the hydrogen number in the N-terminus */
130             searchname[1] = '\0';
131         }
132         for (int j = 0; (j < localPpResidue.natom()); j++)
133         {
134             rtpname = *(localPpResidue.atomname[j]);
135             if (gmx_strcasecmp(searchname, rtpname) == 0)
136             {
137                 jmax = j;
138                 kmax = strlen(searchname);
139                 break;
140             }
141             if (iter == niter - 1)
142             {
143                 minstrlen = std::min(strlen(searchname), strlen(rtpname));
144                 for (k = 0; k < minstrlen; k++)
145                 {
146                     if (searchname[k] != rtpname[k])
147                     {
148                         break;
149                     }
150                 }
151                 if (k > kmax)
152                 {
153                     kmax = k;
154                     jmax = j;
155                 }
156             }
157         }
158     }
159     if (jmax == -1)
160     {
161         gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s", searchname,
162                   localPpResidue.resname.c_str());
163     }
164     if (kmax != strlen(searchname))
165     {
166         gmx_fatal(FARGS,
167                   "Atom %s not found in rtp database in residue %s, "
168                   "it looks a bit like %s",
169                   searchname, localPpResidue.resname.c_str(), *(localPpResidue.atomname[jmax]));
170     }
171     return jmax;
172 }