Clean up grompp memory leaks in rotation and pull structures.
[gromacs.git] / src / gromacs / gmxpreprocess / pgutil.cpp
blobcced555d25f0f2e14b653d53f057a5a5203b9970
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 by the GROMACS development team.
7 * Copyright (c) 2019,2020, 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.
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.
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.
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.
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.
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.
38 /* This file is completely threadsafe - keep it that way! */
40 #include "gmxpre.h"
42 #include "pgutil.h"
44 #include <cstring>
46 #include "gromacs/topology/atoms.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/snprintf.h"
51 #define BUFSIZE 1024
52 static void atom_not_found(int fatal_errno,
53 const char* file,
54 int line,
55 const char* atomname,
56 int resind,
57 const char* resname,
58 const char* bondtype,
59 bool bAllowMissing)
61 char message_buffer[BUFSIZE];
62 if (strcmp(bondtype, "check") != 0)
64 if (0 != strcmp(bondtype, "atom"))
66 snprintf(message_buffer, 1024,
67 "Residue %d named %s of a molecule in the input file was mapped\n"
68 "to an entry in the topology database, but the atom %s used in\n"
69 "an interaction of type %s in that entry is not found in the\n"
70 "input file. Perhaps your atom and/or residue naming needs to be\n"
71 "fixed.\n",
72 resind + 1, resname, atomname, bondtype);
74 else
76 snprintf(message_buffer, 1024,
77 "Residue %d named %s of a molecule in the input file was mapped\n"
78 "to an entry in the topology database, but the atom %s used in\n"
79 "that entry is not found in the input file. Perhaps your atom\n"
80 "and/or residue naming needs to be fixed.\n",
81 resind + 1, resname, atomname);
83 if (bAllowMissing)
85 gmx_warning("WARNING: %s", message_buffer);
87 else
89 gmx_fatal(fatal_errno, file, line, "%s", message_buffer);
94 int search_atom(const char* type, int start, const t_atoms* atoms, const char* bondtype, bool bAllowMissing)
96 int i, resind = -1;
97 bool bPrevious, bNext;
98 int natoms = atoms->nr;
99 t_atom* at = atoms->atom;
100 char** const* anm = atoms->atomname;
102 bPrevious = (strchr(type, '-') != nullptr);
103 bNext = (strchr(type, '+') != nullptr);
105 if (!bPrevious)
107 resind = at[start].resind;
108 if (bNext)
110 /* The next residue */
111 type++;
112 while ((start < natoms) && (at[start].resind == resind))
114 start++;
116 if (start < natoms)
118 resind = at[start].resind;
122 for (i = start; (i < natoms) && (bNext || (at[i].resind == resind)); i++)
124 if (anm[i] && gmx_strcasecmp(type, *(anm[i])) == 0)
126 return i;
129 if (!(bNext && at[start].resind == at[natoms - 1].resind))
131 atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype,
132 bAllowMissing);
135 else
137 /* The previous residue */
138 type++;
139 if (start > 0)
141 resind = at[start - 1].resind;
143 for (i = start - 1; (i >= 0) /*&& (at[i].resind == resind)*/; i--)
145 if (gmx_strcasecmp(type, *(anm[i])) == 0)
147 return i;
150 if (start > 0)
152 atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype,
153 bAllowMissing);
156 return -1;
159 int search_res_atom(const char* type, int resind, const t_atoms* atoms, const char* bondtype, bool bAllowMissing)
161 int i;
163 for (i = 0; (i < atoms->nr); i++)
165 if (atoms->atom[i].resind == resind)
167 return search_atom(type, i, atoms, bondtype, bAllowMissing);
171 return -1;