Remove obsolete mdp option ns-type
[gromacs.git] / src / gromacs / gmxpreprocess / add_par.cpp
blob7d220c170764afe41f2a13afd07749b7f49bff08
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "add_par.h"
42 #include <cstring>
44 #include <algorithm>
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"
53 #include "hackblock.h"
55 void add_param(InteractionsOfType *ps,
56 int ai,
57 int aj,
58 gmx::ArrayRef<const real> c,
59 const char *s)
61 if ((ai < 0) || (aj < 0))
63 gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
65 std::vector<int> atoms = {ai, aj};
66 std::vector<real> forceParm(c.begin(), c.end());
68 ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm, s ? s : ""));
71 void add_imp_param(InteractionsOfType *ps, int ai, int aj, int ak, int al, real c0, real c1,
72 const char *s)
74 std::vector<int> atoms = {ai, aj, ak, al};
75 std::vector<real> forceParm = {c0, c1};
76 ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm, s ? s : ""));
79 void add_dih_param(InteractionsOfType *ps, int ai, int aj, int ak, int al, real c0, real c1,
80 real c2, const char *s)
82 std::vector<int> atoms = {ai, aj, ak, al};
83 std::vector<real> forceParm = {c0, c1, c2};
84 ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm, s ? s : ""));
87 void add_cmap_param(InteractionsOfType *ps, int ai, int aj, int ak, int al, int am, const char *s)
89 std::vector<int> atoms = {ai, aj, ak, al, am};
90 ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}, s ? s : ""));
93 void add_vsite2_atoms(InteractionsOfType *ps, int ai, int aj, int ak)
95 std::vector<int> atoms = {ai, aj, ak};
96 ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
99 void add_vsite2_param(InteractionsOfType *ps, int ai, int aj, int ak, real c0)
101 std::vector<int> atoms = {ai, aj, ak};
102 std::vector<real> forceParm = {c0};
103 ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm));
106 void add_vsite3_param(InteractionsOfType *ps, int ai, int aj, int ak, int al,
107 real c0, real c1)
109 std::vector<int> atoms = {ai, aj, ak, al};
110 std::vector<real> forceParm = {c0, c1};
111 ps->interactionTypes.emplace_back(InteractionOfType(atoms, forceParm));
114 void add_vsite3_atoms(InteractionsOfType *ps, int ai, int aj, int ak, int al, bool bSwapParity)
116 std::vector<int> atoms = {ai, aj, ak, al};
117 ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
119 if (bSwapParity)
121 ps->interactionTypes.back().setForceParameter(1, -1);
125 void add_vsite4_atoms(InteractionsOfType *ps, int ai, int aj, int ak, int al, int am)
127 std::vector<int> atoms = {ai, aj, ak, al, am};
128 ps->interactionTypes.emplace_back(InteractionOfType(atoms, {}));
131 int search_jtype(const PreprocessResidue &localPpResidue, const char *name, bool bNterm)
133 int niter, jmax;
134 size_t k, kmax, minstrlen;
135 char *rtpname, searchname[12];
137 strcpy(searchname, name);
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') ) )
145 niter = 2;
147 else
149 niter = 1;
151 kmax = 0;
152 jmax = -1;
153 for (int iter = 0; (iter < niter && jmax == -1); iter++)
155 if (iter == 1)
157 /* Try without the hydrogen number in the N-terminus */
158 searchname[1] = '\0';
160 for (int j = 0; (j < localPpResidue.natom()); j++)
162 rtpname = *(localPpResidue.atomname[j]);
163 if (gmx_strcasecmp(searchname, rtpname) == 0)
165 jmax = j;
166 kmax = strlen(searchname);
167 break;
169 if (iter == niter - 1)
171 minstrlen = std::min(strlen(searchname), strlen(rtpname));
172 for (k = 0; k < minstrlen; k++)
174 if (searchname[k] != rtpname[k])
176 break;
179 if (k > kmax)
181 kmax = k;
182 jmax = j;
187 if (jmax == -1)
189 gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s",
190 searchname, localPpResidue.resname.c_str());
192 if (kmax != strlen(searchname))
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]));
198 return jmax;