Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / pgutil.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) 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39
40 #include "gmxpre.h"
41
42 #include "pgutil.h"
43
44 #include <cstring>
45
46 #include <algorithm>
47
48 #include "gromacs/topology/atoms.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/snprintf.h"
52
53 #define BUFSIZE 1024
54 static void atom_not_found(int         fatal_errno,
55                            const char* file,
56                            int         line,
57                            const char* atomname,
58                            int         resind,
59                            const char* resname,
60                            const char* bondtype,
61                            bool        bAllowMissing)
62 {
63     char message_buffer[BUFSIZE];
64     if (strcmp(bondtype, "check") != 0)
65     {
66         if (0 != strcmp(bondtype, "atom"))
67         {
68             snprintf(message_buffer,
69                      1024,
70                      "Residue %d named %s of a molecule in the input file was mapped\n"
71                      "to an entry in the topology database, but the atom %s used in\n"
72                      "an interaction of type %s in that entry is not found in the\n"
73                      "input file. Perhaps your atom and/or residue naming needs to be\n"
74                      "fixed.\n",
75                      resind + 1,
76                      resname,
77                      atomname,
78                      bondtype);
79         }
80         else
81         {
82             snprintf(message_buffer,
83                      1024,
84                      "Residue %d named %s of a molecule in the input file was mapped\n"
85                      "to an entry in the topology database, but the atom %s used in\n"
86                      "that entry is not found in the input file. Perhaps your atom\n"
87                      "and/or residue naming needs to be fixed.\n",
88                      resind + 1,
89                      resname,
90                      atomname);
91         }
92         if (bAllowMissing)
93         {
94             gmx_warning("WARNING: %s", message_buffer);
95         }
96         else
97         {
98             gmx_fatal(fatal_errno, file, line, "%s", message_buffer);
99         }
100     }
101 }
102
103 int search_atom(const char*              type,
104                 int                      start,
105                 const t_atoms*           atoms,
106                 const char*              bondtype,
107                 bool                     bAllowMissing,
108                 gmx::ArrayRef<const int> cyclicBondsIndex)
109 {
110     int                                i, resind = -1;
111     bool                               bPrevious, bNext, bOverring;
112     int                                natoms = atoms->nr;
113     t_atom*                            at     = atoms->atom;
114     char** const*                      anm    = atoms->atomname;
115     gmx::ArrayRef<const int>::iterator cyclicBondsIterator;
116
117     bPrevious = (strchr(type, '-') != nullptr);
118     bNext     = (strchr(type, '+') != nullptr);
119
120     if (!bPrevious)
121     {
122         resind = at[start].resind;
123         if (bNext)
124         {
125             /* The next residue */
126             type++;
127             bOverring = !cyclicBondsIndex.empty()
128                         && (cyclicBondsIterator =
129                                     std::find(cyclicBondsIndex.begin(), cyclicBondsIndex.end(), resind))
130                                    != cyclicBondsIndex.end();
131             if (bOverring && ((cyclicBondsIterator - cyclicBondsIndex.begin()) & 1))
132             {
133                 resind = *(--cyclicBondsIterator);
134                 return search_res_atom(type, resind, atoms, bondtype, false);
135             }
136             else
137             {
138                 while ((start < natoms) && (at[start].resind == resind))
139                 {
140                     start++;
141                 }
142                 if (start < natoms)
143                 {
144                     resind = at[start].resind;
145                 }
146             }
147         }
148
149         for (i = start; (i < natoms) && (bNext || (at[i].resind == resind)); i++)
150         {
151             if (anm[i] && gmx_strcasecmp(type, *(anm[i])) == 0)
152             {
153                 return i;
154             }
155         }
156         if (!(bNext && at[start].resind == at[natoms - 1].resind))
157         {
158             atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype, bAllowMissing);
159         }
160     }
161     else
162     {
163         /* The previous residue */
164         type++;
165         resind    = at[start].resind;
166         bOverring = !cyclicBondsIndex.empty()
167                     && (cyclicBondsIterator =
168                                 std::find(cyclicBondsIndex.begin(), cyclicBondsIndex.end(), resind))
169                                != cyclicBondsIndex.end();
170
171         if (bOverring && !((cyclicBondsIterator - cyclicBondsIndex.begin()) & 1))
172         {
173             resind = *(++cyclicBondsIterator);
174             return search_res_atom(type, resind, atoms, bondtype, false);
175         }
176         else
177         {
178             while ((start >= 0) && (at[start].resind == resind))
179             {
180                 start--;
181             }
182             if (start >= 0)
183             {
184                 resind = at[start].resind;
185                 start++;
186             }
187         }
188         for (i = start - 1; (i >= 0) && (at[i].resind == resind); i--)
189         {
190             if (gmx_strcasecmp(type, *(anm[i])) == 0)
191             {
192                 return i;
193             }
194         }
195         if (start > 0)
196         {
197             atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype, bAllowMissing);
198         }
199     }
200     return -1;
201 }
202
203 int search_res_atom(const char* type, int resind, const t_atoms* atoms, const char* bondtype, bool bAllowMissing)
204 {
205     int i;
206
207     for (i = 0; (i < atoms->nr); i++)
208     {
209         if (atoms->atom[i].resind == resind)
210         {
211             return search_atom(type, i, atoms, bondtype, bAllowMissing, gmx::ArrayRef<const int>());
212         }
213     }
214
215     return -1;
216 }