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