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