Tidy: modernize-use-nullptr
[alexxy/gromacs.git] / src / gromacs / fileio / confio.cpp
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,2015,2016,2017, 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 "confio.h"
40
41 #include <cstdio>
42 #include <cstring>
43
44 #include "gromacs/fileio/espio.h"
45 #include "gromacs/fileio/filetypes.h"
46 #include "gromacs/fileio/g96io.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/groio.h"
49 #include "gromacs/fileio/pdbio.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/topology/atoms.h"
54 #include "gromacs/topology/block.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/symtab.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/smalloc.h"
63
64 void write_sto_conf_indexed(const char *outfile, const char *title,
65                             const t_atoms *atoms,
66                             const rvec x[], const rvec *v, int ePBC, const matrix box,
67                             int nindex, int index[])
68 {
69     FILE       *out;
70     int         ftp;
71     t_trxframe  fr;
72
73     ftp = fn2ftp(outfile);
74     switch (ftp)
75     {
76         case efGRO:
77             out = gmx_fio_fopen(outfile, "w");
78             write_hconf_indexed_p(out, title, atoms, nindex, index, x, v, box);
79             gmx_fio_fclose(out);
80             break;
81         case efG96:
82             clear_trxframe(&fr, TRUE);
83             fr.bTitle = TRUE;
84             fr.title  = title;
85             fr.natoms = atoms->nr;
86             fr.bAtoms = TRUE;
87             fr.atoms  = const_cast<t_atoms *>(atoms);
88             fr.bX     = TRUE;
89             fr.x      = const_cast<rvec *>(x);
90             if (v)
91             {
92                 fr.bV = TRUE;
93                 fr.v  = const_cast<rvec *>(v);
94             }
95             fr.bBox = TRUE;
96             copy_mat(box, fr.box);
97             out = gmx_fio_fopen(outfile, "w");
98             write_g96_conf(out, &fr, nindex, index);
99             gmx_fio_fclose(out);
100             break;
101         case efPDB:
102         case efBRK:
103         case efENT:
104         case efPQR:
105             out = gmx_fio_fopen(outfile, "w");
106             write_pdbfile_indexed(out, title, atoms, x, ePBC, box, ' ', -1, nindex, index, nullptr, TRUE);
107             gmx_fio_fclose(out);
108             break;
109         case efESP:
110             out = gmx_fio_fopen(outfile, "w");
111             write_espresso_conf_indexed(out, title, atoms, nindex, index, x, v, box);
112             gmx_fio_fclose(out);
113             break;
114         case efTPR:
115             gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
116             break;
117         default:
118             gmx_incons("Not supported in write_sto_conf_indexed");
119     }
120 }
121
122 void write_sto_conf(const char *outfile, const char *title, const t_atoms *atoms,
123                     const rvec x[], const rvec *v, int ePBC, const matrix box)
124 {
125     FILE       *out;
126     int         ftp;
127     t_trxframe  fr;
128
129     ftp = fn2ftp(outfile);
130     switch (ftp)
131     {
132         case efGRO:
133             write_conf_p(outfile, title, atoms, x, v, box);
134             break;
135         case efG96:
136             clear_trxframe(&fr, TRUE);
137             fr.bTitle = TRUE;
138             fr.title  = title;
139             fr.natoms = atoms->nr;
140             fr.bAtoms = TRUE;
141             fr.atoms  = const_cast<t_atoms *>(atoms); // TODO check
142             fr.bX     = TRUE;
143             fr.x      = const_cast<rvec *>(x);
144             if (v)
145             {
146                 fr.bV = TRUE;
147                 fr.v  = const_cast<rvec *>(v);
148             }
149             fr.bBox = TRUE;
150             copy_mat(box, fr.box);
151             out = gmx_fio_fopen(outfile, "w");
152             write_g96_conf(out, &fr, -1, nullptr);
153             gmx_fio_fclose(out);
154             break;
155         case efPDB:
156         case efBRK:
157         case efENT:
158             out = gmx_fio_fopen(outfile, "w");
159             write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, nullptr, TRUE);
160             gmx_fio_fclose(out);
161             break;
162         case efESP:
163             out = gmx_fio_fopen(outfile, "w");
164             write_espresso_conf_indexed(out, title, atoms, atoms->nr, nullptr, x, v, box);
165             gmx_fio_fclose(out);
166             break;
167         case efTPR:
168             gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
169             break;
170         default:
171             gmx_incons("Not supported in write_sto_conf");
172     }
173 }
174
175 void write_sto_conf_mtop(const char *outfile, const char *title,
176                          gmx_mtop_t *mtop,
177                          const rvec x[], const rvec *v, int ePBC, const matrix box)
178 {
179     int     ftp;
180     FILE   *out;
181     t_atoms atoms;
182
183     ftp = fn2ftp(outfile);
184     switch (ftp)
185     {
186         case efGRO:
187             out = gmx_fio_fopen(outfile, "w");
188             write_hconf_mtop(out, title, mtop, x, v, box);
189             gmx_fio_fclose(out);
190             break;
191         default:
192             /* This is a brute force approach which requires a lot of memory.
193              * We should implement mtop versions of all writing routines.
194              */
195             atoms = gmx_mtop_global_atoms(mtop);
196
197             write_sto_conf(outfile, title, &atoms, x, v, ePBC, box);
198
199             done_atom(&atoms);
200             break;
201     }
202 }
203
204 static void get_stx_coordnum(const char *infile, int *natoms)
205 {
206     FILE      *in;
207     int        ftp;
208     t_trxframe fr;
209     char       g96_line[STRLEN+1];
210
211     ftp = fn2ftp(infile);
212     range_check(ftp, 0, efNR);
213     switch (ftp)
214     {
215         case efGRO:
216             get_coordnum(infile, natoms);
217             break;
218         case efG96:
219         {
220             in        = gmx_fio_fopen(infile, "r");
221             fr.title  = nullptr;
222             fr.natoms = -1;
223             fr.atoms  = nullptr;
224             fr.x      = nullptr;
225             fr.v      = nullptr;
226             fr.f      = nullptr;
227             *natoms   = read_g96_conf(in, infile, &fr, nullptr, g96_line);
228             sfree(const_cast<char *>(fr.title));
229             gmx_fio_fclose(in);
230             break;
231         }
232         case efPDB:
233         case efBRK:
234         case efENT:
235             in = gmx_fio_fopen(infile, "r");
236             get_pdb_coordnum(in, natoms);
237             gmx_fio_fclose(in);
238             break;
239         case efESP:
240             *natoms = get_espresso_coordnum(infile);
241             break;
242         default:
243             gmx_fatal(FARGS, "File type %s not supported in get_stx_coordnum",
244                       ftp2ext(ftp));
245     }
246 }
247
248 static void tpx_make_chain_identifiers(t_atoms *atoms, t_block *mols)
249 {
250     /* We always assign a new chain number, but save the chain id characters
251      * for larger molecules.
252      */
253     const int chainMinAtoms = 15;
254
255     int       chainnum = 0;
256     char      chainid  = 'A';
257     bool      outOfIds = false;
258     for (int m = 0; m < mols->nr; m++)
259     {
260         int a0 = mols->index[m];
261         int a1 = mols->index[m+1];
262         int c;
263         if (a1 - a0 >= chainMinAtoms && !outOfIds)
264         {
265             /* Set the chain id for the output */
266             c = chainid;
267             /* Here we allow for the max possible 2*26+10=62 chain ids */
268             if (chainid == 'Z')
269             {
270                 chainid = 'a';
271             }
272             else if (chainid == 'z')
273             {
274                 chainid = '0';
275             }
276             else if (chainid == '9')
277             {
278                 outOfIds = true;
279             }
280             else
281             {
282                 chainid++;
283             }
284         }
285         else
286         {
287             c = ' ';
288         }
289         for (int a = a0; a < a1; a++)
290         {
291             atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
292             atoms->resinfo[atoms->atom[a].resind].chainid  = c;
293         }
294         chainnum++;
295     }
296
297     /* Blank out the chain id if there was only one chain */
298     if (chainid == 'B')
299     {
300         for (int r = 0; r < atoms->nres; r++)
301         {
302             atoms->resinfo[r].chainid = ' ';
303         }
304     }
305 }
306
307 static void read_stx_conf(const char *infile,
308                           t_symtab *symtab, char ***name, t_atoms *atoms,
309                           rvec x[], rvec *v, int *ePBC, matrix box)
310 {
311     FILE       *in;
312     t_trxframe  fr;
313     int         ftp;
314     char        g96_line[STRLEN+1];
315
316     if (atoms->nr == 0)
317     {
318         fprintf(stderr, "Warning: Number of atoms in %s is 0\n", infile);
319     }
320     else if (atoms->atom == nullptr)
321     {
322         gmx_mem("Uninitialized array atom");
323     }
324
325     if (ePBC)
326     {
327         *ePBC = -1;
328     }
329
330     ftp = fn2ftp(infile);
331     switch (ftp)
332     {
333         case efGRO:
334             gmx_gro_read_conf(infile, symtab, name, atoms, x, v, box);
335             break;
336         case efG96:
337             fr.title  = nullptr;
338             fr.natoms = atoms->nr;
339             fr.atoms  = atoms;
340             fr.x      = x;
341             fr.v      = v;
342             fr.f      = nullptr;
343             in        = gmx_fio_fopen(infile, "r");
344             read_g96_conf(in, infile, &fr, symtab, g96_line);
345             gmx_fio_fclose(in);
346             copy_mat(fr.box, box);
347             *name     = put_symtab(symtab, fr.title);
348             sfree(const_cast<char *>(fr.title));
349             break;
350         case efPDB:
351         case efBRK:
352         case efENT:
353             gmx_pdb_read_conf(infile, symtab, name, atoms, x, ePBC, box);
354             break;
355         case efESP:
356             gmx_espresso_read_conf(infile, symtab, name, atoms, x, v, box);
357             break;
358         default:
359             gmx_incons("Not supported in read_stx_conf");
360     }
361 }
362
363 static void readConfAndAtoms(const char *infile,
364                              t_symtab *symtab, char ***name, t_atoms *atoms,
365                              int *ePBC,
366                              rvec **x, rvec **v, matrix box)
367 {
368     int natoms;
369     get_stx_coordnum(infile, &natoms);
370
371     init_t_atoms(atoms, natoms, (fn2ftp(infile) == efPDB));
372
373     bool xIsNull = false;
374     if (x == nullptr)
375     {
376         snew(x, 1);
377         xIsNull = true;
378     }
379     snew(*x, natoms);
380     if (v)
381     {
382         snew(*v, natoms);
383     }
384     read_stx_conf(infile,
385                   symtab, name, atoms,
386                   *x, (v == nullptr) ? nullptr : *v, ePBC, box);
387     if (xIsNull)
388     {
389         sfree(*x);
390         sfree(x);
391     }
392 }
393
394 void readConfAndTopology(const char *infile,
395                          bool *haveTopology, gmx_mtop_t *mtop,
396                          int *ePBC,
397                          rvec **x, rvec **v, matrix box)
398 {
399     GMX_RELEASE_ASSERT(mtop != nullptr, "readConfAndTopology requires mtop!=NULL");
400
401     if (ePBC != nullptr)
402     {
403         *ePBC = -1;
404     }
405
406     *haveTopology = fn2bTPX(infile);
407     if (*haveTopology)
408     {
409         t_tpxheader header;
410         read_tpxheader(infile, &header, TRUE);
411         if (x)
412         {
413             snew(*x, header.natoms);
414         }
415         if (v)
416         {
417             snew(*v, header.natoms);
418         }
419         int natoms;
420         int ePBC_tmp
421             = read_tpx(infile, nullptr, box, &natoms,
422                        (x == nullptr) ? nullptr : *x, (v == nullptr) ? nullptr : *v, mtop);
423         if (ePBC != nullptr)
424         {
425             *ePBC = ePBC_tmp;
426         }
427     }
428     else
429     {
430         t_symtab   symtab;
431         char     **name;
432         t_atoms    atoms;
433
434         open_symtab(&symtab);
435
436         readConfAndAtoms(infile, &symtab, &name, &atoms, ePBC, x, v, box);
437
438         init_mtop(mtop);
439         convertAtomsToMtop(&symtab, name, &atoms, mtop);
440     }
441 }
442
443 gmx_bool read_tps_conf(const char *infile, t_topology *top, int *ePBC,
444                        rvec **x, rvec **v, matrix box, gmx_bool requireMasses)
445 {
446     bool        haveTopology;
447     gmx_mtop_t *mtop;
448
449     // Note: We should have an initializer instead of relying on snew
450     snew(mtop, 1);
451     readConfAndTopology(infile, &haveTopology, mtop, ePBC, x, v, box);
452
453     *top = gmx_mtop_t_to_t_topology(mtop, true);
454     sfree(mtop);
455
456     tpx_make_chain_identifiers(&top->atoms, &top->mols);
457
458     if (requireMasses && !top->atoms.haveMass)
459     {
460         atomsSetMassesBasedOnNames(&top->atoms, TRUE);
461
462         if (!top->atoms.haveMass)
463         {
464             gmx_fatal(FARGS, "Masses were requested, but for some atom(s) masses could not be found in the database. Use a tpr file as input, if possible, or add these atoms to the mass database.");
465         }
466     }
467
468     return haveTopology;
469 }