Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / pgutil.c
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, 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
39 #include "gmxpre.h"
40
41 #include "pgutil.h"
42
43 #include <string.h>
44
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/fatalerror.h"
47
48 #define BUFSIZE 1024
49 static void atom_not_found(int fatal_errno, const char *file, int line,
50                            const char *atomname, int resind,
51                            const char *resname,
52                            const char *bondtype, gmx_bool bAllowMissing)
53 {
54     char message_buffer[BUFSIZE];
55     if (strcmp(bondtype, "check") != 0)
56     {
57         if (0 != strcmp(bondtype, "atom"))
58         {
59             snprintf(message_buffer, 1024,
60                      "Residue %d named %s of a molecule in the input file was mapped\n"
61                      "to an entry in the topology database, but the atom %s used in\n"
62                      "an interaction of type %s in that entry is not found in the\n"
63                      "input file. Perhaps your atom and/or residue naming needs to be\n"
64                      "fixed.\n",
65                      resind+1, resname, atomname, bondtype);
66         }
67         else
68         {
69             snprintf(message_buffer, 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                      "that entry is not found in the input file. Perhaps your atom\n"
73                      "and/or residue naming needs to be fixed.\n",
74                      resind+1, resname, atomname);
75         }
76         if (bAllowMissing)
77         {
78             gmx_warning("WARNING: %s", message_buffer);
79         }
80         else
81         {
82             gmx_fatal(fatal_errno, file, line, message_buffer);
83         }
84     }
85 }
86
87 atom_id search_atom(const char *type, int start,
88                     t_atoms *atoms,
89                     const char *bondtype, gmx_bool bAllowMissing)
90 {
91     int             i, resind = -1;
92     gmx_bool        bPrevious, bNext;
93     int             natoms = atoms->nr;
94     t_atom         *at     = atoms->atom;
95     char ** const * anm    = atoms->atomname;
96
97     bPrevious = (strchr(type, '-') != NULL);
98     bNext     = (strchr(type, '+') != NULL);
99
100     if (!bPrevious)
101     {
102         resind = at[start].resind;
103         if (bNext)
104         {
105             /* The next residue */
106             type++;
107             while ((start < natoms) && (at[start].resind == resind))
108             {
109                 start++;
110             }
111             if (start < natoms)
112             {
113                 resind = at[start].resind;
114             }
115         }
116
117         for (i = start; (i < natoms) && (bNext || (at[i].resind == resind)); i++)
118         {
119             if (anm[i] && gmx_strcasecmp(type, *(anm[i])) == 0)
120             {
121                 return (atom_id) i;
122             }
123         }
124         if (!(bNext && at[start].resind == at[natoms-1].resind))
125         {
126             atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype, bAllowMissing);
127         }
128     }
129     else
130     {
131         /* The previous residue */
132         type++;
133         if (start > 0)
134         {
135             resind = at[start-1].resind;
136         }
137         for (i = start-1; (i >= 0) /*&& (at[i].resind == resind)*/; i--)
138         {
139             if (gmx_strcasecmp(type, *(anm[i])) == 0)
140             {
141                 return (atom_id) i;
142             }
143         }
144         if (start > 0)
145         {
146             atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype, bAllowMissing);
147         }
148     }
149     return NO_ATID;
150 }
151
152 atom_id
153 search_res_atom(const char *type, int resind,
154                 t_atoms *atoms,
155                 const char *bondtype, gmx_bool bAllowMissing)
156 {
157     int i;
158
159     for (i = 0; (i < atoms->nr); i++)
160     {
161         if (atoms->atom[i].resind == resind)
162         {
163             return search_atom(type, i, atoms, bondtype, bAllowMissing);
164         }
165     }
166
167     return NO_ATID;
168 }
169
170
171 void set_at(t_atom *at, real m, real q, int type, int resind)
172 {
173     at->m      = m;
174     at->q      = q;
175     at->type   = type;
176     at->resind = resind;
177 }