a380ffafc598e31da5d89a0b47a22620eeb5fef0
[alexxy/gromacs.git] / src / gromacs / fileio / espio.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 "espio.h"
40
41 #include <cstdio>
42 #include <cstdlib>
43 #include <cstring>
44
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"
54
55 static int get_espresso_word(FILE* fp, char word[])
56 {
57     int ret, nc, i;
58
59     ret = 0;
60     nc  = 0;
61
62     do
63     {
64         i = fgetc(fp);
65         if (i != EOF)
66         {
67             if (i == ' ' || i == '\n' || i == '\t')
68             {
69                 if (nc > 0)
70                 {
71                     ret = 1;
72                 }
73             }
74             else if (i == '{')
75             {
76                 if (nc == 0)
77                 {
78                     word[nc++] = '{';
79                 }
80                 ret = 2;
81             }
82             else if (i == '}')
83             {
84                 if (nc == 0)
85                 {
86                     word[nc++] = '}';
87                 }
88                 ret = 3;
89             }
90             else
91             {
92                 word[nc++] = static_cast<char>(i);
93             }
94         }
95     } while (i != EOF && ret == 0);
96
97     word[nc] = '\0';
98
99     return ret;
100 }
101
102 static int check_open_parenthesis(FILE* fp, int r, const char* infile, const char* keyword)
103 {
104     int  level_inc;
105     char word[STRLEN];
106
107     level_inc = 0;
108     if (r == 2)
109     {
110         level_inc++;
111     }
112     else
113     {
114         r = get_espresso_word(fp, word);
115         if (r == 2)
116         {
117             level_inc++;
118         }
119         else
120         {
121             gmx_fatal(FARGS, "Expected '{' after '%s' in file '%s'", keyword, infile);
122         }
123     }
124
125     return level_inc;
126 }
127
128 static int check_close_parenthesis(FILE* fp, int r, const char* infile, const char* keyword)
129 {
130     int  level_inc;
131     char word[STRLEN];
132
133     level_inc = 0;
134     if (r == 3)
135     {
136         level_inc--;
137     }
138     else
139     {
140         r = get_espresso_word(fp, word);
141         if (r == 3)
142         {
143             level_inc--;
144         }
145         else
146         {
147             gmx_fatal(FARGS, "Expected '}' after section '%s' in file '%s'", keyword, infile);
148         }
149     }
150
151     return level_inc;
152 }
153
154 enum
155 {
156     espID,
157     espPOS,
158     espTYPE,
159     espQ,
160     espV,
161     espF,
162     espMOLECULE,
163     espNR
164 };
165 static const char* const esp_prop[espNR] = { "id", "pos", "type", "q", "v", "f", "molecule" };
166
167 void gmx_espresso_read_conf(const char* infile, t_symtab* symtab, char** name, t_atoms* atoms, rvec x[], rvec* v, matrix box)
168 {
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;
175
176     if (name != nullptr)
177     {
178         // No title reading implemented for espresso files
179         *name = gmx_strdup("");
180     }
181
182     clear_mat(box);
183
184     atoms->haveMass    = FALSE;
185     atoms->haveCharge  = FALSE;
186     atoms->haveType    = FALSE;
187     atoms->haveBState  = FALSE;
188     atoms->havePdbInfo = FALSE;
189
190     fp = gmx_fio_fopen(infile, "r");
191
192     bFoundParticles = FALSE;
193     bFoundVariable  = FALSE;
194     bMol            = FALSE;
195     level           = 0;
196     while ((r = get_espresso_word(fp, word)))
197     {
198         if (level == 1 && std::strcmp(word, "particles") == 0 && !bFoundParticles)
199         {
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)))
204             {
205                 bFoundProp = FALSE;
206                 for (p = 0; p < espNR; p++)
207                 {
208                     if (strcmp(word, esp_prop[p]) == 0)
209                     {
210                         bFoundProp    = TRUE;
211                         prop[nprop++] = p;
212                         if (p == espQ)
213                         {
214                             atoms->haveCharge = TRUE;
215                         }
216
217                         if (debug)
218                         {
219                             fprintf(debug, "  prop[%d] = %s\n", nprop - 1, esp_prop[prop[nprop - 1]]);
220                         }
221                     }
222                 }
223                 if (!bFoundProp && word[0] != '}')
224                 {
225                     gmx_fatal(FARGS, "Can not read Espresso files with particle property '%s'", word);
226                 }
227                 if (bFoundProp && p == espMOLECULE)
228                 {
229                     bMol = TRUE;
230                 }
231                 if (r == 3)
232                 {
233                     level--;
234                 }
235             }
236
237             i = 0;
238             while (level > 0 && (r = get_espresso_word(fp, word)))
239             {
240                 if (r == 2)
241                 {
242                     level++;
243                 }
244                 else if (r == 3)
245                 {
246                     level--;
247                 }
248                 if (level == 2)
249                 {
250                     for (p = 0; p < nprop; p++)
251                     {
252                         switch (prop[p])
253                         {
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++)
260                                 {
261                                     r = get_espresso_word(fp, word);
262                                     sscanf(word, "%lf", &d);
263                                     x[i][m] = d;
264                                 }
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++)
277                                 {
278                                     r = get_espresso_word(fp, word);
279                                     sscanf(word, "%lf", &d);
280                                     v[i][m] = d;
281                                 }
282                                 break;
283                             case espF:
284                                 for (m = 0; m < 3; m++)
285                                 {
286                                     r = get_espresso_word(fp, word);
287                                     /* not used */
288                                 }
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)
294                                 {
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? */
301                                 }
302                                 else
303                                 {
304                                     atoms->atom[i].resind = atoms->atom[i - 1].resind;
305                                 }
306                                 break;
307                         }
308                     }
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)
313                     {
314                         if (i == 0 || atoms->atom[i].resind != atoms->atom[i - 1].resind)
315                         {
316                             atoms->resinfo[atoms->atom[i].resind].name = put_symtab(symtab, "MOL");
317                         }
318                     }
319                     else
320                     {
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)
325                         {
326                             sprintf(buf, "T%c", 'A' + atoms->atom[i].type);
327                         }
328                         else
329                         {
330                             sprintf(buf, "T%c%c", 'A' + atoms->atom[i].type / 26,
331                                     'A' + atoms->atom[i].type % 26);
332                         }
333                         t_atoms_set_resinfo(atoms, i, symtab, buf, i, ' ', 0, ' ');
334                     }
335
336                     if (r == 3)
337                     {
338                         level--;
339                     }
340                     i++;
341                 }
342             }
343             atoms->nres = atoms->nr;
344
345             if (i != atoms->nr)
346             {
347                 gmx_fatal(FARGS,
348                           "Internal inconsistency in Espresso routines, read %d atoms, expected %d "
349                           "atoms",
350                           i, atoms->nr);
351             }
352         }
353         else if (level == 1 && std::strcmp(word, "variable") == 0 && !bFoundVariable)
354         {
355             bFoundVariable = TRUE;
356             level += check_open_parenthesis(fp, r, infile, "variable");
357             while (level == 2 && (r = get_espresso_word(fp, word)))
358             {
359                 if (level == 2 && std::strcmp(word, "box_l") == 0)
360                 {
361                     for (m = 0; m < 3; m++)
362                     {
363                         r = get_espresso_word(fp, word);
364                         sscanf(word, "%lf", &d);
365                         box[m][m] = d;
366                     }
367                     level += check_close_parenthesis(fp, r, infile, "box_l");
368                 }
369             }
370         }
371         else if (r == 2)
372         {
373             level++;
374         }
375         else if (r == 3)
376         {
377             level--;
378         }
379     }
380
381     if (!bFoundParticles)
382     {
383         fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n", infile);
384     }
385
386     gmx_fio_fclose(fp);
387 }
388
389 int get_espresso_coordnum(const char* infile)
390 {
391     FILE*    fp;
392     char     word[STRLEN];
393     int      natoms, level, r;
394     gmx_bool bFoundParticles;
395
396     natoms = 0;
397
398     fp = gmx_fio_fopen(infile, "r");
399
400     bFoundParticles = FALSE;
401     level           = 0;
402     while ((r = get_espresso_word(fp, word)) && !bFoundParticles)
403     {
404         if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
405         {
406             bFoundParticles = TRUE;
407             level += check_open_parenthesis(fp, r, infile, "particles");
408             while (level > 0 && (r = get_espresso_word(fp, word)))
409             {
410                 if (r == 2)
411                 {
412                     level++;
413                     if (level == 2)
414                     {
415                         natoms++;
416                     }
417                 }
418                 else if (r == 3)
419                 {
420                     level--;
421                 }
422             }
423         }
424         else if (r == 2)
425         {
426             level++;
427         }
428         else if (r == 3)
429         {
430             level--;
431         }
432     }
433     if (!bFoundParticles)
434     {
435         fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n", infile);
436     }
437
438     gmx_fio_fclose(fp);
439
440     return natoms;
441 }
442
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)
451 {
452     int i, j;
453
454     fprintf(out, "# %s\n", title);
455     if (TRICLINIC(box))
456     {
457         gmx_warning("The Espresso format does not support triclinic unit-cells");
458     }
459     fprintf(out, "{variable {box_l %f %f %f}}\n", box[0][0], box[1][1], box[2][2]);
460
461     fprintf(out, "{particles {id pos type q%s}\n", v ? " v" : "");
462     for (i = 0; i < nx; i++)
463     {
464         if (index)
465         {
466             j = index[i];
467         }
468         else
469         {
470             j = i;
471         }
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)
475         {
476             fprintf(out, " %f %f %f", v[j][XX], v[j][YY], v[j][ZZ]);
477         }
478         fprintf(out, "}\n");
479     }
480     fprintf(out, "}\n");
481 }