Update instructions in containers.rst
[gromacs.git] / src / gromacs / fileio / espio.cpp
bloba380ffafc598e31da5d89a0b47a22620eeb5fef0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005, The GROMACS development team.
5 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
6 * Copyright (c) 2018,2019,2020, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "espio.h"
41 #include <cstdio>
42 #include <cstdlib>
43 #include <cstring>
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/topology/atoms.h"
49 #include "gromacs/topology/symtab.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/smalloc.h"
55 static int get_espresso_word(FILE* fp, char word[])
57 int ret, nc, i;
59 ret = 0;
60 nc = 0;
64 i = fgetc(fp);
65 if (i != EOF)
67 if (i == ' ' || i == '\n' || i == '\t')
69 if (nc > 0)
71 ret = 1;
74 else if (i == '{')
76 if (nc == 0)
78 word[nc++] = '{';
80 ret = 2;
82 else if (i == '}')
84 if (nc == 0)
86 word[nc++] = '}';
88 ret = 3;
90 else
92 word[nc++] = static_cast<char>(i);
95 } while (i != EOF && ret == 0);
97 word[nc] = '\0';
99 return ret;
102 static int check_open_parenthesis(FILE* fp, int r, const char* infile, const char* keyword)
104 int level_inc;
105 char word[STRLEN];
107 level_inc = 0;
108 if (r == 2)
110 level_inc++;
112 else
114 r = get_espresso_word(fp, word);
115 if (r == 2)
117 level_inc++;
119 else
121 gmx_fatal(FARGS, "Expected '{' after '%s' in file '%s'", keyword, infile);
125 return level_inc;
128 static int check_close_parenthesis(FILE* fp, int r, const char* infile, const char* keyword)
130 int level_inc;
131 char word[STRLEN];
133 level_inc = 0;
134 if (r == 3)
136 level_inc--;
138 else
140 r = get_espresso_word(fp, word);
141 if (r == 3)
143 level_inc--;
145 else
147 gmx_fatal(FARGS, "Expected '}' after section '%s' in file '%s'", keyword, infile);
151 return level_inc;
154 enum
156 espID,
157 espPOS,
158 espTYPE,
159 espQ,
160 espV,
161 espF,
162 espMOLECULE,
163 espNR
165 static const char* const esp_prop[espNR] = { "id", "pos", "type", "q", "v", "f", "molecule" };
167 void gmx_espresso_read_conf(const char* infile, t_symtab* symtab, char** name, t_atoms* atoms, rvec x[], rvec* v, matrix box)
169 FILE* fp;
170 char word[STRLEN], buf[STRLEN];
171 int level, r, nprop, p, i, m, molnr;
172 int prop[32];
173 double d;
174 gmx_bool bFoundParticles, bFoundProp, bFoundVariable, bMol;
176 if (name != nullptr)
178 // No title reading implemented for espresso files
179 *name = gmx_strdup("");
182 clear_mat(box);
184 atoms->haveMass = FALSE;
185 atoms->haveCharge = FALSE;
186 atoms->haveType = FALSE;
187 atoms->haveBState = FALSE;
188 atoms->havePdbInfo = FALSE;
190 fp = gmx_fio_fopen(infile, "r");
192 bFoundParticles = FALSE;
193 bFoundVariable = FALSE;
194 bMol = FALSE;
195 level = 0;
196 while ((r = get_espresso_word(fp, word)))
198 if (level == 1 && std::strcmp(word, "particles") == 0 && !bFoundParticles)
200 bFoundParticles = TRUE;
201 level += check_open_parenthesis(fp, r, infile, "particles");
202 nprop = 0;
203 while (level == 2 && (r = get_espresso_word(fp, word)))
205 bFoundProp = FALSE;
206 for (p = 0; p < espNR; p++)
208 if (strcmp(word, esp_prop[p]) == 0)
210 bFoundProp = TRUE;
211 prop[nprop++] = p;
212 if (p == espQ)
214 atoms->haveCharge = TRUE;
217 if (debug)
219 fprintf(debug, " prop[%d] = %s\n", nprop - 1, esp_prop[prop[nprop - 1]]);
223 if (!bFoundProp && word[0] != '}')
225 gmx_fatal(FARGS, "Can not read Espresso files with particle property '%s'", word);
227 if (bFoundProp && p == espMOLECULE)
229 bMol = TRUE;
231 if (r == 3)
233 level--;
237 i = 0;
238 while (level > 0 && (r = get_espresso_word(fp, word)))
240 if (r == 2)
242 level++;
244 else if (r == 3)
246 level--;
248 if (level == 2)
250 for (p = 0; p < nprop; p++)
252 switch (prop[p])
254 case espID:
255 r = get_espresso_word(fp, word);
256 /* Not used */
257 break;
258 case espPOS:
259 for (m = 0; m < 3; m++)
261 r = get_espresso_word(fp, word);
262 sscanf(word, "%lf", &d);
263 x[i][m] = d;
265 break;
266 case espTYPE:
267 r = get_espresso_word(fp, word);
268 atoms->atom[i].type = std::strtol(word, nullptr, 10);
269 break;
270 case espQ:
271 r = get_espresso_word(fp, word);
272 sscanf(word, "%lf", &d);
273 atoms->atom[i].q = d;
274 break;
275 case espV:
276 for (m = 0; m < 3; m++)
278 r = get_espresso_word(fp, word);
279 sscanf(word, "%lf", &d);
280 v[i][m] = d;
282 break;
283 case espF:
284 for (m = 0; m < 3; m++)
286 r = get_espresso_word(fp, word);
287 /* not used */
289 break;
290 case espMOLECULE:
291 r = get_espresso_word(fp, word);
292 molnr = std::strtol(word, nullptr, 10);
293 if (i == 0 || atoms->resinfo[atoms->atom[i - 1].resind].nr != molnr)
295 atoms->atom[i].resind = (i == 0 ? 0 : atoms->atom[i - 1].resind + 1);
296 atoms->resinfo[atoms->atom[i].resind].nr = molnr;
297 atoms->resinfo[atoms->atom[i].resind].ic = ' ';
298 atoms->resinfo[atoms->atom[i].resind].chainid = ' ';
299 atoms->resinfo[atoms->atom[i].resind].chainnum =
300 molnr; /* Not sure if this is right? */
302 else
304 atoms->atom[i].resind = atoms->atom[i - 1].resind;
306 break;
309 /* Generate an atom name from the particle type */
310 sprintf(buf, "T%hu", atoms->atom[i].type);
311 atoms->atomname[i] = put_symtab(symtab, buf);
312 if (bMol)
314 if (i == 0 || atoms->atom[i].resind != atoms->atom[i - 1].resind)
316 atoms->resinfo[atoms->atom[i].resind].name = put_symtab(symtab, "MOL");
319 else
321 /* Residue number is the atom number */
322 atoms->atom[i].resind = i;
323 /* Generate an residue name from the particle type */
324 if (atoms->atom[i].type < 26)
326 sprintf(buf, "T%c", 'A' + atoms->atom[i].type);
328 else
330 sprintf(buf, "T%c%c", 'A' + atoms->atom[i].type / 26,
331 'A' + atoms->atom[i].type % 26);
333 t_atoms_set_resinfo(atoms, i, symtab, buf, i, ' ', 0, ' ');
336 if (r == 3)
338 level--;
340 i++;
343 atoms->nres = atoms->nr;
345 if (i != atoms->nr)
347 gmx_fatal(FARGS,
348 "Internal inconsistency in Espresso routines, read %d atoms, expected %d "
349 "atoms",
350 i, atoms->nr);
353 else if (level == 1 && std::strcmp(word, "variable") == 0 && !bFoundVariable)
355 bFoundVariable = TRUE;
356 level += check_open_parenthesis(fp, r, infile, "variable");
357 while (level == 2 && (r = get_espresso_word(fp, word)))
359 if (level == 2 && std::strcmp(word, "box_l") == 0)
361 for (m = 0; m < 3; m++)
363 r = get_espresso_word(fp, word);
364 sscanf(word, "%lf", &d);
365 box[m][m] = d;
367 level += check_close_parenthesis(fp, r, infile, "box_l");
371 else if (r == 2)
373 level++;
375 else if (r == 3)
377 level--;
381 if (!bFoundParticles)
383 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n", infile);
386 gmx_fio_fclose(fp);
389 int get_espresso_coordnum(const char* infile)
391 FILE* fp;
392 char word[STRLEN];
393 int natoms, level, r;
394 gmx_bool bFoundParticles;
396 natoms = 0;
398 fp = gmx_fio_fopen(infile, "r");
400 bFoundParticles = FALSE;
401 level = 0;
402 while ((r = get_espresso_word(fp, word)) && !bFoundParticles)
404 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
406 bFoundParticles = TRUE;
407 level += check_open_parenthesis(fp, r, infile, "particles");
408 while (level > 0 && (r = get_espresso_word(fp, word)))
410 if (r == 2)
412 level++;
413 if (level == 2)
415 natoms++;
418 else if (r == 3)
420 level--;
424 else if (r == 2)
426 level++;
428 else if (r == 3)
430 level--;
433 if (!bFoundParticles)
435 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n", infile);
438 gmx_fio_fclose(fp);
440 return natoms;
443 void write_espresso_conf_indexed(FILE* out,
444 const char* title,
445 const t_atoms* atoms,
446 int nx,
447 const int* index,
448 const rvec* x,
449 const rvec* v,
450 const matrix box)
452 int i, j;
454 fprintf(out, "# %s\n", title);
455 if (TRICLINIC(box))
457 gmx_warning("The Espresso format does not support triclinic unit-cells");
459 fprintf(out, "{variable {box_l %f %f %f}}\n", box[0][0], box[1][1], box[2][2]);
461 fprintf(out, "{particles {id pos type q%s}\n", v ? " v" : "");
462 for (i = 0; i < nx; i++)
464 if (index)
466 j = index[i];
468 else
470 j = i;
472 fprintf(out, "\t{%d %f %f %f %hu %g", j, x[j][XX], x[j][YY], x[j][ZZ], atoms->atom[j].type,
473 atoms->atom[j].q);
474 if (v)
476 fprintf(out, " %f %f %f", v[j][XX], v[j][YY], v[j][ZZ]);
478 fprintf(out, "}\n");
480 fprintf(out, "}\n");