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