Remove obsolete mdp option ns-type
[gromacs.git] / src / gromacs / gmxpreprocess / pgutil.cpp
blob2f83195e97fdcdb7d1b2912d9cdab96766c1711d
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) 2012,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! */
39 #include "gmxpre.h"
41 #include "pgutil.h"
43 #include <cstring>
45 #include "gromacs/topology/atoms.h"
46 #include "gromacs/utility/cstringutil.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/snprintf.h"
50 #define BUFSIZE 1024
51 static void atom_not_found(int fatal_errno, const char *file, int line,
52 const char *atomname, int resind,
53 const char *resname,
54 const char *bondtype, bool bAllowMissing)
56 char message_buffer[BUFSIZE];
57 if (strcmp(bondtype, "check") != 0)
59 if (0 != strcmp(bondtype, "atom"))
61 snprintf(message_buffer, 1024,
62 "Residue %d named %s of a molecule in the input file was mapped\n"
63 "to an entry in the topology database, but the atom %s used in\n"
64 "an interaction of type %s in that entry is not found in the\n"
65 "input file. Perhaps your atom and/or residue naming needs to be\n"
66 "fixed.\n",
67 resind+1, resname, atomname, bondtype);
69 else
71 snprintf(message_buffer, 1024,
72 "Residue %d named %s of a molecule in the input file was mapped\n"
73 "to an entry in the topology database, but the atom %s used in\n"
74 "that entry is not found in the input file. Perhaps your atom\n"
75 "and/or residue naming needs to be fixed.\n",
76 resind+1, resname, atomname);
78 if (bAllowMissing)
80 gmx_warning("WARNING: %s", message_buffer);
82 else
84 gmx_fatal(fatal_errno, file, line, "%s", message_buffer);
89 int search_atom(const char *type, int start,
90 const t_atoms *atoms,
91 const char *bondtype, bool bAllowMissing)
93 int i, resind = -1;
94 bool bPrevious, bNext;
95 int natoms = atoms->nr;
96 t_atom *at = atoms->atom;
97 char ** const * anm = atoms->atomname;
99 bPrevious = (strchr(type, '-') != nullptr);
100 bNext = (strchr(type, '+') != nullptr);
102 if (!bPrevious)
104 resind = at[start].resind;
105 if (bNext)
107 /* The next residue */
108 type++;
109 while ((start < natoms) && (at[start].resind == resind))
111 start++;
113 if (start < natoms)
115 resind = at[start].resind;
119 for (i = start; (i < natoms) && (bNext || (at[i].resind == resind)); i++)
121 if (anm[i] && gmx_strcasecmp(type, *(anm[i])) == 0)
123 return i;
126 if (!(bNext && at[start].resind == at[natoms-1].resind))
128 atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype, bAllowMissing);
131 else
133 /* The previous residue */
134 type++;
135 if (start > 0)
137 resind = at[start-1].resind;
139 for (i = start-1; (i >= 0) /*&& (at[i].resind == resind)*/; i--)
141 if (gmx_strcasecmp(type, *(anm[i])) == 0)
143 return i;
146 if (start > 0)
148 atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype, bAllowMissing);
151 return -1;
155 search_res_atom(const char *type, int resind,
156 const t_atoms *atoms,
157 const char *bondtype, bool bAllowMissing)
159 int i;
161 for (i = 0; (i < atoms->nr); i++)
163 if (atoms->atom[i].resind == resind)
165 return search_atom(type, i, atoms, bondtype, bAllowMissing);
169 return -1;
173 void set_at(t_atom *at, real m, real q, int type, int resind)
175 at->m = m;
176 at->q = q;
177 at->type = type;
178 at->resind = resind;