21f70b014f3f0f8a65529eeeb8789cc151e8410d
[alexxy/gromacs.git] / src / gromacs / fileio / strdb.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) 2013,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 #include "gmxpre.h"
38
39 #include "strdb.h"
40
41 #include <stdio.h>
42 #include <stdlib.h>
43 #include <string.h>
44
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/utility/cstringutil.h"
48
49 #include "gromacs/utility/futil.h"
50
51 gmx_bool get_a_line(FILE *fp, char line[], int n)
52 {
53     char *line0;
54     char *dum;
55
56     snew(line0, n+1);
57
58     do
59     {
60         if (!fgets(line0, n+1, fp))
61         {
62             sfree(line0);
63             return FALSE;
64         }
65         dum = strchr(line0, '\n');
66         if (dum)
67         {
68             dum[0] = '\0';
69         }
70         else if (strlen(line0) == n)
71         {
72             fprintf(stderr, "Warning: line length exceeds buffer length (%d), data might be corrupted\n", n);
73             line0[n-1] = '\0';
74         }
75         else
76         {
77             fprintf(stderr, "Warning: file does not end with a newline, last line:\n%s\n",
78                     line0);
79         }
80         dum = strchr(line0, ';');
81         if (dum)
82         {
83             dum[0] = '\0';
84         }
85         strncpy(line, line0, n);
86         dum = line0;
87         ltrim(dum);
88     }
89     while (dum[0] == '\0');
90
91     sfree(line0);
92     return TRUE;
93 }
94
95 gmx_bool get_header(char line[], char *header)
96 {
97     char temp[STRLEN], *dum;
98
99     strcpy(temp, line);
100     dum = strchr(temp, '[');
101     if (dum == NULL)
102     {
103         return FALSE;
104     }
105     dum[0] = ' ';
106     dum    = strchr(temp, ']');
107     if (dum == NULL)
108     {
109         gmx_fatal(FARGS, "header is not terminated on line:\n'%s'\n", line);
110         return FALSE;
111     }
112     dum[0] = '\0';
113     if (sscanf(temp, "%s%*s", header) != 1)
114     {
115         return FALSE;
116     }
117
118     return TRUE;
119 }
120
121 int get_strings(const char *db, char ***strings)
122 {
123     FILE  *in;
124     char **ptr;
125     char   buf[256];
126     int    i, nstr;
127
128     in = libopen(db);
129
130     if (fscanf(in, "%d", &nstr) != 1)
131     {
132         gmx_warning("File %s is empty", db);
133         gmx_ffclose(in);
134         return 0;
135     }
136     snew(ptr, nstr);
137     for (i = 0; (i < nstr); i++)
138     {
139         if (1 != fscanf(in, "%s", buf))
140         {
141             gmx_fatal(FARGS, "Cannot read string from buffer");
142         }
143 #ifdef DEBUG
144         fprintf(stderr, "Have read: %s\n", buf);
145 #endif
146         ptr[i] = gmx_strdup(buf);
147     }
148     gmx_ffclose(in);
149
150     *strings = ptr;
151
152     return nstr;
153 }
154
155 int search_str(int nstr, char **str, char *key)
156 {
157     int i;
158
159     /* Linear search */
160     for (i = 0; (i < nstr); i++)
161     {
162         if (gmx_strcasecmp(str[i], key) == 0)
163         {
164             return i;
165         }
166     }
167
168     return -1;
169 }
170
171 int fget_lines(FILE *in, char ***strings)
172 {
173     char **ptr;
174     char   buf[STRLEN];
175     int    i, nstr;
176     char  *pret;
177
178     pret = fgets(buf, STRLEN-1, in);
179     if (pret == NULL  || sscanf(buf, "%d", &nstr) != 1)
180     {
181         gmx_warning("File is empty");
182         gmx_ffclose(in);
183
184         return 0;
185     }
186     snew(ptr, nstr);
187     for (i = 0; (i < nstr); i++)
188     {
189         fgets2(buf, STRLEN-1, in);
190         ptr[i] = gmx_strdup(buf);
191     }
192
193     (*strings) = ptr;
194
195     return nstr;
196 }
197
198 int get_lines(const char *db, char ***strings)
199 {
200     FILE *in;
201     int   nstr;
202
203     in   = libopen(db);
204     nstr = fget_lines(in, strings);
205     gmx_ffclose(in);
206
207     return nstr;
208 }
209
210 int get_file(const char *db, char ***strings)
211 {
212     FILE  *in;
213     char **ptr = NULL;
214     char   buf[STRLEN];
215     int    i, nstr, maxi;
216
217     in = libopen(db);
218
219     i = maxi = 0;
220     while (fgets2(buf, STRLEN-1, in))
221     {
222         if (i >= maxi)
223         {
224             maxi += 50;
225             srenew(ptr, maxi);
226         }
227         ptr[i] = gmx_strdup(buf);
228         i++;
229     }
230     nstr = i;
231     gmx_ffclose(in);
232     srenew(ptr, nstr);
233     *strings = ptr;
234
235     return nstr;
236 }