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