Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / editconf.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.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "editconf.h"
41
42 #include <cmath>
43 #include <cstring>
44
45 #include <algorithm>
46 #include <string>
47
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/commandline/viewit.h"
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/fileio/pdbio.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/gmxana/princ.h"
55 #include "gromacs/gmxlib/conformation_utilities.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/pbcutil/rmpbc.h"
61 #include "gromacs/topology/atomprop.h"
62 #include "gromacs/topology/index.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/utility/arraysize.h"
65 #include "gromacs/utility/cstringutil.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/futil.h"
68 #include "gromacs/utility/gmxassert.h"
69 #include "gromacs/utility/smalloc.h"
70 #include "gromacs/utility/strdb.h"
71
72 static real calc_mass(t_atoms* atoms, gmx_bool bGetMass, AtomProperties* aps)
73 {
74     real tmass;
75     int  i;
76
77     tmass = 0;
78     for (i = 0; (i < atoms->nr); i++)
79     {
80         if (bGetMass)
81         {
82             aps->setAtomProperty(epropMass,
83                                  std::string(*atoms->resinfo[atoms->atom[i].resind].name),
84                                  std::string(*atoms->atomname[i]),
85                                  &(atoms->atom[i].m));
86         }
87         tmass += atoms->atom[i].m;
88     }
89
90     return tmass;
91 }
92
93 static real calc_geom(int isize, const int* index, rvec* x, rvec geom_center, rvec minval, rvec maxval, gmx_bool bDiam)
94 {
95     real diam2, d;
96     int  ii, i, j;
97
98     clear_rvec(geom_center);
99     diam2 = 0;
100     if (isize == 0)
101     {
102         clear_rvec(minval);
103         clear_rvec(maxval);
104     }
105     else
106     {
107         if (index)
108         {
109             ii = index[0];
110         }
111         else
112         {
113             ii = 0;
114         }
115         for (j = 0; j < DIM; j++)
116         {
117             minval[j] = maxval[j] = x[ii][j];
118         }
119         for (i = 0; i < isize; i++)
120         {
121             if (index)
122             {
123                 ii = index[i];
124             }
125             else
126             {
127                 ii = i;
128             }
129             rvec_inc(geom_center, x[ii]);
130             for (j = 0; j < DIM; j++)
131             {
132                 if (x[ii][j] < minval[j])
133                 {
134                     minval[j] = x[ii][j];
135                 }
136                 if (x[ii][j] > maxval[j])
137                 {
138                     maxval[j] = x[ii][j];
139                 }
140             }
141             if (bDiam)
142             {
143                 if (index)
144                 {
145                     for (j = i + 1; j < isize; j++)
146                     {
147                         d     = distance2(x[ii], x[index[j]]);
148                         diam2 = std::max(d, diam2);
149                     }
150                 }
151                 else
152                 {
153                     for (j = i + 1; j < isize; j++)
154                     {
155                         d     = distance2(x[i], x[j]);
156                         diam2 = std::max(d, diam2);
157                     }
158                 }
159             }
160         }
161         svmul(1.0 / isize, geom_center, geom_center);
162     }
163
164     return std::sqrt(diam2);
165 }
166
167 static void center_conf(int natom, rvec* x, rvec center, rvec geom_cent)
168 {
169     int  i;
170     rvec shift;
171
172     rvec_sub(center, geom_cent, shift);
173
174     printf("    shift       :%7.3f%7.3f%7.3f (nm)\n", shift[XX], shift[YY], shift[ZZ]);
175
176     for (i = 0; (i < natom); i++)
177     {
178         rvec_inc(x[i], shift);
179     }
180 }
181
182 static void scale_conf(int natom, rvec x[], matrix box, const rvec scale)
183 {
184     int i, j;
185
186     for (i = 0; i < natom; i++)
187     {
188         for (j = 0; j < DIM; j++)
189         {
190             x[i][j] *= scale[j];
191         }
192     }
193     for (i = 0; i < DIM; i++)
194     {
195         for (j = 0; j < DIM; j++)
196         {
197             box[i][j] *= scale[j];
198         }
199     }
200 }
201
202 static void read_bfac(const char* fn, int* n_bfac, double** bfac_val, int** bfac_nr)
203 {
204     int    i;
205     char** bfac_lines;
206
207     *n_bfac = get_lines(fn, &bfac_lines);
208     snew(*bfac_val, *n_bfac);
209     snew(*bfac_nr, *n_bfac);
210     fprintf(stderr, "Reading %d B-factors from %s\n", *n_bfac, fn);
211     for (i = 0; (i < *n_bfac); i++)
212     {
213         sscanf(bfac_lines[i], "%d %lf", &(*bfac_nr)[i], &(*bfac_val)[i]);
214     }
215 }
216
217 static void set_pdb_conf_bfac(int natoms, int nres, t_atoms* atoms, int n_bfac, double* bfac, int* bfac_nr, gmx_bool peratom)
218 {
219     real     bfac_min, bfac_max;
220     int      i, n;
221     gmx_bool found;
222
223     if (n_bfac > atoms->nres)
224     {
225         peratom = TRUE;
226     }
227
228     bfac_max = -1e10;
229     bfac_min = 1e10;
230     for (i = 0; (i < n_bfac); i++)
231     {
232         /*    if ((bfac_nr[i]-1<0) || (bfac_nr[i]-1>=atoms->nr))
233            gmx_fatal(FARGS,"Index of B-Factor %d is out of range: %d (%g)",
234            i+1,bfac_nr[i],bfac[i]); */
235         if (bfac[i] > bfac_max)
236         {
237             bfac_max = bfac[i];
238         }
239         if (bfac[i] < bfac_min)
240         {
241             bfac_min = bfac[i];
242         }
243     }
244     while ((bfac_max > 99.99) || (bfac_min < -99.99))
245     {
246         fprintf(stderr,
247                 "Range of values for B-factors too large (min %g, max %g) "
248                 "will scale down a factor 10\n",
249                 bfac_min,
250                 bfac_max);
251         for (i = 0; (i < n_bfac); i++)
252         {
253             bfac[i] /= 10;
254         }
255         bfac_max /= 10;
256         bfac_min /= 10;
257     }
258     while ((std::abs(bfac_max) < 0.5) && (std::abs(bfac_min) < 0.5))
259     {
260         fprintf(stderr,
261                 "Range of values for B-factors too small (min %g, max %g) "
262                 "will scale up a factor 10\n",
263                 bfac_min,
264                 bfac_max);
265         for (i = 0; (i < n_bfac); i++)
266         {
267             bfac[i] *= 10;
268         }
269         bfac_max *= 10;
270         bfac_min *= 10;
271     }
272
273     for (i = 0; (i < natoms); i++)
274     {
275         atoms->pdbinfo[i].bfac = 0;
276     }
277
278     if (!peratom)
279     {
280         fprintf(stderr, "Will attach %d B-factors to %d residues\n", n_bfac, nres);
281         for (i = 0; (i < n_bfac); i++)
282         {
283             found = FALSE;
284             for (n = 0; (n < natoms); n++)
285             {
286                 if (bfac_nr[i] == atoms->resinfo[atoms->atom[n].resind].nr)
287                 {
288                     atoms->pdbinfo[n].bfac = bfac[i];
289                     found                  = TRUE;
290                 }
291             }
292             if (!found)
293             {
294                 gmx_warning("Residue nr %d not found\n", bfac_nr[i]);
295             }
296         }
297     }
298     else
299     {
300         fprintf(stderr, "Will attach %d B-factors to %d atoms\n", n_bfac, natoms);
301         for (i = 0; (i < n_bfac); i++)
302         {
303             atoms->pdbinfo[bfac_nr[i] - 1].bfac = bfac[i];
304         }
305     }
306 }
307
308 static void pdb_legend(FILE* out, int natoms, int nres, t_atoms* atoms, rvec x[])
309 {
310     real bfac_min, bfac_max, xmin, ymin, zmin;
311     int  i;
312     int  space = ' ';
313
314     bfac_max = -1e10;
315     bfac_min = 1e10;
316     xmin     = 1e10;
317     ymin     = 1e10;
318     zmin     = 1e10;
319     for (i = 0; (i < natoms); i++)
320     {
321         xmin     = std::min(xmin, x[i][XX]);
322         ymin     = std::min(ymin, x[i][YY]);
323         zmin     = std::min(zmin, x[i][ZZ]);
324         bfac_min = std::min(bfac_min, atoms->pdbinfo[i].bfac);
325         bfac_max = std::max(bfac_max, atoms->pdbinfo[i].bfac);
326     }
327     fprintf(stderr, "B-factors range from %g to %g\n", bfac_min, bfac_max);
328     for (i = 1; (i < 12); i++)
329     {
330         fprintf(out,
331                 "%-6s%5d  %-4.4s%3.3s %c%4d%c   %8.3f%8.3f%8.3f%6.2f%6.2f\n",
332                 "ATOM  ",
333                 natoms + 1 + i,
334                 "CA",
335                 "LEG",
336                 space,
337                 nres + 1,
338                 space,
339                 (xmin + (i * 0.12)) * 10,
340                 ymin * 10,
341                 zmin * 10,
342                 1.0,
343                 bfac_min + ((i - 1.0) * (bfac_max - bfac_min) / 10));
344     }
345 }
346
347 static void visualize_images(const char* fn, PbcType pbcType, matrix box)
348 {
349     t_atoms atoms;
350     rvec*   img;
351     char *  c, *ala;
352     int     nat, i;
353
354     nat = NTRICIMG + 1;
355     init_t_atoms(&atoms, nat, FALSE);
356     atoms.nr = nat;
357     snew(img, nat);
358     /* FIXME: Constness should not be cast away */
359     c   = const_cast<char*>("C");
360     ala = const_cast<char*>("ALA");
361     for (i = 0; i < nat; i++)
362     {
363         atoms.atomname[i]        = &c;
364         atoms.atom[i].resind     = i;
365         atoms.resinfo[i].name    = &ala;
366         atoms.resinfo[i].nr      = i + 1;
367         atoms.resinfo[i].chainid = 'A' + i / NCUCVERT;
368     }
369     calc_triclinic_images(box, img + 1);
370
371     write_sto_conf(fn, "Images", &atoms, img, nullptr, pbcType, box);
372
373     done_atom(&atoms);
374     sfree(img);
375 }
376
377 static void visualize_box(FILE* out, int a0, int r0, matrix box, const rvec gridsize)
378 {
379     int*  edge;
380     rvec *vert, shift;
381     int   nx, ny, nz, nbox, nat;
382     int   i, j, x, y, z;
383     int   rectedge[24] = { 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6, 4 };
384
385     a0++;
386     r0++;
387
388     nx   = gmx::roundToInt(gridsize[XX]);
389     ny   = gmx::roundToInt(gridsize[YY]);
390     nz   = gmx::roundToInt(gridsize[ZZ]);
391     nbox = nx * ny * nz;
392     if (TRICLINIC(box))
393     {
394         nat = nbox * NCUCVERT;
395         snew(vert, nat);
396         calc_compact_unitcell_vertices(ecenterDEF, box, vert);
397         j = 0;
398         for (z = 0; z < nz; z++)
399         {
400             for (y = 0; y < ny; y++)
401             {
402                 for (x = 0; x < nx; x++)
403                 {
404                     for (i = 0; i < DIM; i++)
405                     {
406                         shift[i] = x * box[0][i] + y * box[1][i] + z * box[2][i];
407                     }
408                     for (i = 0; i < NCUCVERT; i++)
409                     {
410                         rvec_add(vert[i], shift, vert[j]);
411                         j++;
412                     }
413                 }
414             }
415         }
416
417         for (i = 0; i < nat; i++)
418         {
419             gmx_fprintf_pdb_atomline(out,
420                                      epdbATOM,
421                                      a0 + i,
422                                      "C",
423                                      ' ',
424                                      "BOX",
425                                      'K' + i / NCUCVERT,
426                                      r0 + i,
427                                      ' ',
428                                      10 * vert[i][XX],
429                                      10 * vert[i][YY],
430                                      10 * vert[i][ZZ],
431                                      1.0,
432                                      0.0,
433                                      "");
434         }
435
436         edge = compact_unitcell_edges();
437         for (j = 0; j < nbox; j++)
438         {
439             for (i = 0; i < NCUCEDGE; i++)
440             {
441                 fprintf(out,
442                         "CONECT%5d%5d\n",
443                         a0 + j * NCUCVERT + edge[2 * i],
444                         a0 + j * NCUCVERT + edge[2 * i + 1]);
445             }
446         }
447
448         sfree(vert);
449     }
450     else
451     {
452         i = 0;
453         for (z = 0; z <= 1; z++)
454         {
455             for (y = 0; y <= 1; y++)
456             {
457                 for (x = 0; x <= 1; x++)
458                 {
459                     gmx_fprintf_pdb_atomline(out,
460                                              epdbATOM,
461                                              a0 + i,
462                                              "C",
463                                              ' ',
464                                              "BOX",
465                                              'K' + i / 8,
466                                              r0 + i,
467                                              ' ',
468                                              x * 10 * box[XX][XX],
469                                              y * 10 * box[YY][YY],
470                                              z * 10 * box[ZZ][ZZ],
471                                              1.0,
472                                              0.0,
473                                              "");
474                     i++;
475                 }
476             }
477         }
478         for (i = 0; i < 24; i += 2)
479         {
480             fprintf(out, "CONECT%5d%5d\n", a0 + rectedge[i], a0 + rectedge[i + 1]);
481         }
482     }
483 }
484
485 static void calc_rotmatrix(rvec principal_axis, rvec targetvec, matrix rotmatrix)
486 {
487     rvec rotvec;
488     real ux, uy, uz, costheta, sintheta;
489
490     costheta = cos_angle(principal_axis, targetvec);
491     sintheta = std::sqrt(1.0 - costheta * costheta); /* sign is always positive since 0<theta<pi */
492
493     /* Determine rotation from cross product with target vector */
494     cprod(principal_axis, targetvec, rotvec);
495     unitv(rotvec, rotvec);
496     printf("Aligning %g %g %g to %g %g %g : xprod  %g %g %g\n",
497            principal_axis[XX],
498            principal_axis[YY],
499            principal_axis[ZZ],
500            targetvec[XX],
501            targetvec[YY],
502            targetvec[ZZ],
503            rotvec[XX],
504            rotvec[YY],
505            rotvec[ZZ]);
506
507     ux              = rotvec[XX];
508     uy              = rotvec[YY];
509     uz              = rotvec[ZZ];
510     rotmatrix[0][0] = ux * ux + (1.0 - ux * ux) * costheta;
511     rotmatrix[0][1] = ux * uy * (1 - costheta) - uz * sintheta;
512     rotmatrix[0][2] = ux * uz * (1 - costheta) + uy * sintheta;
513     rotmatrix[1][0] = ux * uy * (1 - costheta) + uz * sintheta;
514     rotmatrix[1][1] = uy * uy + (1.0 - uy * uy) * costheta;
515     rotmatrix[1][2] = uy * uz * (1 - costheta) - ux * sintheta;
516     rotmatrix[2][0] = ux * uz * (1 - costheta) - uy * sintheta;
517     rotmatrix[2][1] = uy * uz * (1 - costheta) + ux * sintheta;
518     rotmatrix[2][2] = uz * uz + (1.0 - uz * uz) * costheta;
519
520     printf("Rotation matrix: \n%g %g %g\n%g %g %g\n%g %g %g\n",
521            rotmatrix[0][0],
522            rotmatrix[0][1],
523            rotmatrix[0][2],
524            rotmatrix[1][0],
525            rotmatrix[1][1],
526            rotmatrix[1][2],
527            rotmatrix[2][0],
528            rotmatrix[2][1],
529            rotmatrix[2][2]);
530 }
531
532 static void renum_resnr(t_atoms* atoms, int isize, const int* index, int resnr_start)
533 {
534     int i, resind_prev, resind;
535
536     resind_prev = -1;
537     for (i = 0; i < isize; i++)
538     {
539         resind = atoms->atom[index == nullptr ? i : index[i]].resind;
540         if (resind != resind_prev)
541         {
542             atoms->resinfo[resind].nr = resnr_start;
543             resnr_start++;
544         }
545         resind_prev = resind;
546     }
547 }
548
549 int gmx_editconf(int argc, char* argv[])
550 {
551     const char* desc[] = {
552         "[THISMODULE] converts generic structure format to [REF].gro[ref], [TT].g96[tt]",
553         "or [REF].pdb[ref].",
554         "[PAR]",
555         "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and",
556         "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]",
557         "will center the system in the box, unless [TT]-noc[tt] is used.",
558         "The [TT]-center[tt] option can be used to shift the geometric center",
559         "of the system from the default of (x/2, y/2, z/2) implied by [TT]-c[tt]",
560         "to some other value.",
561         "[PAR]",
562         "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
563         "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
564         "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
565         "[TT]octahedron[tt] is a truncated octahedron.",
566         "The last two are special cases of a triclinic box.",
567         "The length of the three box vectors of the truncated octahedron is the",
568         "shortest distance between two opposite hexagons.",
569         "Relative to a cubic box with some periodic image distance, the volume of a ",
570         "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
571         "and that of a truncated octahedron is 0.77 times.",
572         "[PAR]",
573         "Option [TT]-box[tt] requires only",
574         "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
575         "[PAR]",
576         "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, ",
577         "[IT]y[it]-,",
578         "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
579         "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
580         "to the diameter of the system (largest distance between atoms) plus twice",
581         "the specified distance.",
582         "[PAR]",
583         "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
584         "a triclinic box and cannot be used with option [TT]-d[tt].",
585         "[PAR]",
586         "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
587         "can be selected for calculating the size and the geometric center,",
588         "otherwise the whole system is used.",
589         "[PAR]",
590         "[TT]-rotate[tt] rotates the coordinates and velocities.",
591         "[PAR]",
592         "[TT]-princ[tt] aligns the principal axes of the system along the",
593         "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
594         "This may allow you to decrease the box volume,",
595         "but beware that molecules can rotate significantly in a nanosecond.",
596         "[PAR]",
597         "Scaling is applied before any of the other operations are",
598         "performed. Boxes and coordinates can be scaled to give a certain density (option",
599         "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
600         "file is given as input. A special feature of the scaling option is that when the",
601         "factor -1 is given in one dimension, one obtains a mirror image,",
602         "mirrored in one of the planes. When one uses -1 in three dimensions, ",
603         "a point-mirror image is obtained.[PAR]",
604         "Groups are selected after all operations have been applied.[PAR]",
605         "Periodicity can be removed in a crude manner.",
606         "It is important that the box vectors at the bottom of your input file",
607         "are correct when the periodicity is to be removed.",
608         "[PAR]",
609         "When writing [REF].pdb[ref] files, B-factors can be",
610         "added with the [TT]-bf[tt] option. B-factors are read",
611         "from a file with with following format: first line states number of",
612         "entries in the file, next lines state an index",
613         "followed by a B-factor. The B-factors will be attached per residue",
614         "unless the number of B-factors is larger than the number of the residues or unless the",
615         "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
616         "be added instead of B-factors. [TT]-legend[tt] will produce",
617         "a row of CA atoms with B-factors ranging from the minimum to the",
618         "maximum value found, effectively making a legend for viewing.",
619         "[PAR]",
620         "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
621         "file for the MEAD electrostatics",
622         "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
623         "is that the input file is a run input file.",
624         "The B-factor field is then filled with the Van der Waals radius",
625         "of the atoms while the occupancy field will hold the charge.",
626         "[PAR]",
627         "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
628         "and the radius in the occupancy.",
629         "[PAR]",
630         "Option [TT]-align[tt] allows alignment",
631         "of the principal axis of a specified group against the given vector, ",
632         "with an optional center of rotation specified by [TT]-aligncenter[tt].",
633         "[PAR]",
634         "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
635         "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
636         "[PAR]",
637         "To convert a truncated octrahedron file produced by a package which uses",
638         "a cubic box with the corners cut off (such as GROMOS), use::",
639         "",
640         "  gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
641         "",
642         "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
643     };
644     const char* bugs[] = {
645         "For complex molecules, the periodicity removal routine may break down, ",
646         "in that case you can use [gmx-trjconv]."
647     };
648     static real     dist = 0.0;
649     static gmx_bool bNDEF = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW = FALSE, bCONECT = FALSE;
650     static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead = FALSE,
651                     bGrasp = FALSE, bSig56 = FALSE;
652     static rvec scale = { 1, 1, 1 }, newbox = { 0, 0, 0 }, newang = { 90, 90, 90 };
653     static real rho = 1000.0, rvdw = 0.12;
654     static rvec center = { 0, 0, 0 }, translation = { 0, 0, 0 }, rotangles = { 0, 0, 0 },
655                 aligncenter = { 0, 0, 0 }, targetvec = { 0, 0, 0 };
656     static const char *btype[] = { nullptr,        "triclinic",  "cubic",
657                                    "dodecahedron", "octahedron", nullptr },
658                       *label   = "A";
659     static rvec visbox         = { 0, 0, 0 };
660     static int  resnr_start    = -1;
661     t_pargs     pa[]           = {
662         { "-ndef", FALSE, etBOOL, { &bNDEF }, "Choose output from default index groups" },
663         { "-visbox",
664           FALSE,
665           etRVEC,
666           { visbox },
667           "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
668         { "-bt", FALSE, etENUM, { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
669         { "-box", FALSE, etRVEC, { newbox }, "Box vector lengths (a,b,c)" },
670         { "-angles", FALSE, etRVEC, { newang }, "Angles between the box vectors (bc,ac,ab)" },
671         { "-d", FALSE, etREAL, { &dist }, "Distance between the solute and the box" },
672         { "-c",
673           FALSE,
674           etBOOL,
675           { &bCenter },
676           "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
677         { "-center", FALSE, etRVEC, { center }, "Shift the geometrical center to (x,y,z)" },
678         { "-aligncenter", FALSE, etRVEC, { aligncenter }, "Center of rotation for alignment" },
679         { "-align", FALSE, etRVEC, { targetvec }, "Align to target vector" },
680         { "-translate", FALSE, etRVEC, { translation }, "Translation" },
681         { "-rotate",
682           FALSE,
683           etRVEC,
684           { rotangles },
685           "Rotation around the X, Y and Z axes in degrees" },
686         { "-princ", FALSE, etBOOL, { &bOrient }, "Orient molecule(s) along their principal axes" },
687         { "-scale", FALSE, etRVEC, { scale }, "Scaling factor" },
688         { "-density",
689           FALSE,
690           etREAL,
691           { &rho },
692           "Density (g/L) of the output box achieved by scaling" },
693         { "-pbc", FALSE, etBOOL, { &bRMPBC }, "Remove the periodicity (make molecule whole again)" },
694         { "-resnr", FALSE, etINT, { &resnr_start }, " Renumber residues starting from resnr" },
695         { "-grasp",
696           FALSE,
697           etBOOL,
698           { &bGrasp },
699           "Store the charge of the atom in the B-factor field and the radius of the atom in the "
700           "occupancy field" },
701         { "-rvdw",
702           FALSE,
703           etREAL,
704           { &rvdw },
705           "Default Van der Waals radius (in nm) if one can not be found in the database or if no "
706           "parameters are present in the topology file" },
707         { "-sig56",
708           FALSE,
709           etBOOL,
710           { &bSig56 },
711           "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
712         { "-vdwread",
713           FALSE,
714           etBOOL,
715           { &bReadVDW },
716           "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing "
717           "the radii based on the force field" },
718         { "-atom", FALSE, etBOOL, { &peratom }, "Force B-factor attachment per atom" },
719         { "-legend", FALSE, etBOOL, { &bLegend }, "Make B-factor legend" },
720         { "-label", FALSE, etSTR, { &label }, "Add chain label for all residues" },
721         { "-conect",
722           FALSE,
723           etBOOL,
724           { &bCONECT },
725           "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a "
726           "topology is present" }
727     };
728 #define NPA asize(pa)
729
730     FILE*             out;
731     const char *      infile, *outfile;
732     int               outftp, inftp, natom, i, j, n_bfac, itype, ntype;
733     double *          bfac    = nullptr, c6, c12;
734     int*              bfac_nr = nullptr;
735     t_topology*       top     = nullptr;
736     char *            grpname, *sgrpname, *agrpname;
737     int               isize, ssize, numAlignmentAtoms;
738     int *             index, *sindex, *aindex;
739     rvec *            x, *v, gc, rmin, rmax, size;
740     PbcType           pbcType;
741     matrix            box, rotmatrix, trans;
742     rvec              princd, tmpvec;
743     gmx_bool          bIndex, bSetSize, bSetAng, bDist, bSetCenter, bAlign;
744     gmx_bool          bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
745     real              diam = 0, mass = 0, d, vdw;
746     gmx_conect        conect;
747     gmx_output_env_t* oenv;
748     t_filenm          fnm[] = { { efSTX, "-f", nullptr, ffREAD },
749                        { efNDX, "-n", nullptr, ffOPTRD },
750                        { efSTO, nullptr, nullptr, ffOPTWR },
751                        { efPQR, "-mead", "mead", ffOPTWR },
752                        { efDAT, "-bf", "bfact", ffOPTRD } };
753 #define NFILE asize(fnm)
754
755     if (!parse_common_args(
756                 &argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa, asize(desc), desc, asize(bugs), bugs, &oenv))
757     {
758         return 0;
759     }
760     fprintf(stdout,
761             "Note that major changes are planned in future for "
762             "editconf, to improve usability and utility.\n");
763
764     bIndex     = opt2bSet("-n", NFILE, fnm) || bNDEF;
765     bMead      = opt2bSet("-mead", NFILE, fnm);
766     bSetSize   = opt2parg_bSet("-box", NPA, pa);
767     bSetAng    = opt2parg_bSet("-angles", NPA, pa);
768     bSetCenter = opt2parg_bSet("-center", NPA, pa);
769     bDist      = opt2parg_bSet("-d", NPA, pa);
770     bAlign     = opt2parg_bSet("-align", NPA, pa);
771     /* Only automatically turn on centering without -noc */
772     if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
773     {
774         bCenter = TRUE;
775     }
776     bScale     = opt2parg_bSet("-scale", NPA, pa);
777     bRho       = opt2parg_bSet("-density", NPA, pa);
778     bTranslate = opt2parg_bSet("-translate", NPA, pa);
779     bRotate    = opt2parg_bSet("-rotate", NPA, pa);
780     if (bScale && bRho)
781     {
782         fprintf(stderr, "WARNING: setting -density overrides -scale\n");
783     }
784     bScale    = bScale || bRho;
785     bCalcGeom = bCenter || bRotate || bOrient || bScale;
786
787     GMX_RELEASE_ASSERT(btype[0] != nullptr, "Option setting inconsistency; btype[0] is NULL");
788
789     bCalcDiam = (btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o');
790
791     infile = ftp2fn(efSTX, NFILE, fnm);
792     if (bMead)
793     {
794         outfile = ftp2fn(efPQR, NFILE, fnm);
795     }
796     else
797     {
798         outfile = ftp2fn(efSTO, NFILE, fnm);
799     }
800     outftp = fn2ftp(outfile);
801     inftp  = fn2ftp(infile);
802
803     AtomProperties aps;
804
805     if (bMead && bGrasp)
806     {
807         printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
808         bGrasp = FALSE;
809     }
810     if (bGrasp && (outftp != efPDB))
811     {
812         gmx_fatal(FARGS,
813                   "Output file should be a .pdb file"
814                   " when using the -grasp option\n");
815     }
816     if ((bMead || bGrasp) && (fn2ftp(infile) != efTPR))
817     {
818         gmx_fatal(FARGS,
819                   "Input file should be a .tpr file"
820                   " when using the -mead option\n");
821     }
822
823     t_symtab symtab;
824     char*    name;
825     t_atoms  atoms;
826     open_symtab(&symtab);
827     readConfAndAtoms(infile, &symtab, &name, &atoms, &pbcType, &x, &v, box);
828     natom = atoms.nr;
829     if (atoms.pdbinfo == nullptr)
830     {
831         snew(atoms.pdbinfo, atoms.nr);
832     }
833     atoms.havePdbInfo = TRUE;
834
835     if (fn2ftp(infile) == efPDB)
836     {
837         get_pdb_atomnumber(&atoms, &aps);
838     }
839     printf("Read %d atoms\n", atoms.nr);
840
841     /* Get the element numbers if available in a pdb file */
842     if (fn2ftp(infile) == efPDB)
843     {
844         get_pdb_atomnumber(&atoms, &aps);
845     }
846
847     if (pbcType != PbcType::No)
848     {
849         real vol = det(box);
850         printf("Volume: %g nm^3, corresponds to roughly %d electrons\n",
851                vol,
852                100 * (static_cast<int>(vol * 4.5)));
853     }
854
855     if (bMead || bGrasp || bCONECT)
856     {
857         top = read_top(infile, nullptr);
858     }
859
860     if (bMead || bGrasp)
861     {
862         if (atoms.nr != top->atoms.nr)
863         {
864             gmx_fatal(FARGS, "Atom numbers don't match (%d vs. %d)", atoms.nr, top->atoms.nr);
865         }
866         snew(atoms.pdbinfo, top->atoms.nr);
867         ntype = top->idef.atnr;
868         for (i = 0; (i < atoms.nr); i++)
869         {
870             /* Determine the Van der Waals radius from the force field */
871             if (bReadVDW)
872             {
873                 if (!aps.setAtomProperty(epropVDW,
874                                          *top->atoms.resinfo[top->atoms.atom[i].resind].name,
875                                          *top->atoms.atomname[i],
876                                          &vdw))
877                 {
878                     vdw = rvdw;
879                 }
880             }
881             else
882             {
883                 itype = top->atoms.atom[i].type;
884                 c12   = top->idef.iparams[itype * ntype + itype].lj.c12;
885                 c6    = top->idef.iparams[itype * ntype + itype].lj.c6;
886                 if ((c6 != 0) && (c12 != 0))
887                 {
888                     real sig6;
889                     if (bSig56)
890                     {
891                         sig6 = 2 * c12 / c6;
892                     }
893                     else
894                     {
895                         sig6 = c12 / c6;
896                     }
897                     vdw = 0.5 * gmx::sixthroot(sig6);
898                 }
899                 else
900                 {
901                     vdw = rvdw;
902                 }
903             }
904             /* Factor of 10 for nm -> Angstroms */
905             vdw *= 10;
906
907             if (bMead)
908             {
909                 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
910                 atoms.pdbinfo[i].bfac  = vdw;
911             }
912             else
913             {
914                 atoms.pdbinfo[i].occup = vdw;
915                 atoms.pdbinfo[i].bfac  = top->atoms.atom[i].q;
916             }
917         }
918     }
919     bHaveV = FALSE;
920     for (i = 0; (i < natom) && !bHaveV; i++)
921     {
922         for (j = 0; (j < DIM) && !bHaveV; j++)
923         {
924             bHaveV = bHaveV || (v[i][j] != 0);
925         }
926     }
927     printf("%selocities found\n", bHaveV ? "V" : "No v");
928
929     if (visbox[0] > 0)
930     {
931         if (bIndex)
932         {
933             gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
934         }
935         if (outftp != efPDB)
936         {
937             gmx_fatal(FARGS, "Sorry, can only visualize box with a pdb file");
938         }
939     }
940     else if (visbox[0] == -1)
941     {
942         visualize_images("images.pdb", pbcType, box);
943     }
944
945     /* remove pbc */
946     if (bRMPBC)
947     {
948         rm_gropbc(&atoms, x, box);
949     }
950
951     if (bCalcGeom)
952     {
953         if (bIndex)
954         {
955             fprintf(stderr, "\nSelect a group for determining the system size:\n");
956             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ssize, &sindex, &sgrpname);
957         }
958         else
959         {
960             ssize  = atoms.nr;
961             sindex = nullptr;
962         }
963         diam = calc_geom(ssize, sindex, x, gc, rmin, rmax, bCalcDiam);
964         rvec_sub(rmax, rmin, size);
965         printf("    system size :%7.3f%7.3f%7.3f (nm)\n", size[XX], size[YY], size[ZZ]);
966         if (bCalcDiam)
967         {
968             printf("    diameter    :%7.3f               (nm)\n", diam);
969         }
970         printf("    center      :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
971         printf("    box vectors :%7.3f%7.3f%7.3f (nm)\n", norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
972         printf("    box angles  :%7.2f%7.2f%7.2f (degrees)\n",
973                norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[YY], box[ZZ]),
974                norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[ZZ]),
975                norm2(box[YY]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[YY]));
976         printf("    box volume  :%7.2f               (nm^3)\n", det(box));
977     }
978
979     if (bRho || bOrient || bAlign)
980     {
981         mass = calc_mass(&atoms, !fn2bTPX(infile), &aps);
982     }
983
984     if (bOrient)
985     {
986         int*  index;
987         char* grpnames;
988
989         /* Get a group for principal component analysis */
990         fprintf(stderr, "\nSelect group for the determining the orientation\n");
991         get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
992
993         /* Orient the principal axes along the coordinate axes */
994         orient_princ(&atoms, isize, index, natom, x, bHaveV ? v : nullptr, nullptr);
995         sfree(index);
996         sfree(grpnames);
997     }
998
999     if (bScale)
1000     {
1001         /* scale coordinates and box */
1002         if (bRho)
1003         {
1004             /* Compute scaling constant */
1005             real vol, dens;
1006
1007             vol  = det(box);
1008             dens = (mass * AMU) / (vol * NANO * NANO * NANO);
1009             fprintf(stderr, "Volume  of input %g (nm^3)\n", vol);
1010             fprintf(stderr, "Mass    of input %g (a.m.u.)\n", mass);
1011             fprintf(stderr, "Density of input %g (g/l)\n", dens);
1012             if (vol == 0 || mass == 0)
1013             {
1014                 gmx_fatal(FARGS,
1015                           "Cannot scale density with "
1016                           "zero mass (%g) or volume (%g)\n",
1017                           mass,
1018                           vol);
1019             }
1020
1021             scale[XX] = scale[YY] = scale[ZZ] = std::cbrt(dens / rho);
1022             fprintf(stderr, "Scaling all box vectors by %g\n", scale[XX]);
1023         }
1024         scale_conf(atoms.nr, x, box, scale);
1025     }
1026
1027     if (bAlign)
1028     {
1029         if (bIndex)
1030         {
1031             fprintf(stderr, "\nSelect a group that you want to align:\n");
1032             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &numAlignmentAtoms, &aindex, &agrpname);
1033         }
1034         else
1035         {
1036             numAlignmentAtoms = atoms.nr;
1037             snew(aindex, numAlignmentAtoms);
1038             for (i = 0; i < numAlignmentAtoms; i++)
1039             {
1040                 aindex[i] = i;
1041             }
1042         }
1043         printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n",
1044                numAlignmentAtoms,
1045                natom,
1046                targetvec[XX],
1047                targetvec[YY],
1048                targetvec[ZZ],
1049                aligncenter[XX],
1050                aligncenter[YY],
1051                aligncenter[ZZ]);
1052         /*subtract out pivot point*/
1053         for (i = 0; i < numAlignmentAtoms; i++)
1054         {
1055             rvec_dec(x[aindex[i]], aligncenter);
1056         }
1057         /*now determine transform and rotate*/
1058         /*will this work?*/
1059         principal_comp(numAlignmentAtoms, aindex, atoms.atom, x, trans, princd);
1060
1061         unitv(targetvec, targetvec);
1062         printf("Using %g %g %g as principal axis\n", trans[0][2], trans[1][2], trans[2][2]);
1063         tmpvec[XX] = trans[0][2];
1064         tmpvec[YY] = trans[1][2];
1065         tmpvec[ZZ] = trans[2][2];
1066         calc_rotmatrix(tmpvec, targetvec, rotmatrix);
1067         /* rotmatrix finished */
1068
1069         for (i = 0; i < numAlignmentAtoms; ++i)
1070         {
1071             mvmul(rotmatrix, x[aindex[i]], tmpvec);
1072             copy_rvec(tmpvec, x[aindex[i]]);
1073         }
1074
1075         /*add pivot point back*/
1076         for (i = 0; i < numAlignmentAtoms; i++)
1077         {
1078             rvec_inc(x[aindex[i]], aligncenter);
1079         }
1080         if (!bIndex)
1081         {
1082             sfree(aindex);
1083         }
1084     }
1085
1086     if (bTranslate)
1087     {
1088         if (bIndex)
1089         {
1090             fprintf(stderr, "\nSelect a group that you want to translate:\n");
1091             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ssize, &sindex, &sgrpname);
1092         }
1093         else
1094         {
1095             ssize  = atoms.nr;
1096             sindex = nullptr;
1097         }
1098         printf("Translating %d atoms (out of %d) by %g %g %g nm\n",
1099                ssize,
1100                natom,
1101                translation[XX],
1102                translation[YY],
1103                translation[ZZ]);
1104         if (sindex)
1105         {
1106             for (i = 0; i < ssize; i++)
1107             {
1108                 rvec_inc(x[sindex[i]], translation);
1109             }
1110         }
1111         else
1112         {
1113             for (i = 0; i < natom; i++)
1114             {
1115                 rvec_inc(x[i], translation);
1116             }
1117         }
1118     }
1119     if (bRotate)
1120     {
1121         /* Rotate */
1122         printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n",
1123                rotangles[XX],
1124                rotangles[YY],
1125                rotangles[ZZ]);
1126         for (i = 0; i < DIM; i++)
1127         {
1128             rotangles[i] *= DEG2RAD;
1129         }
1130         rotate_conf(natom, x, v, rotangles[XX], rotangles[YY], rotangles[ZZ]);
1131     }
1132
1133     if (bCalcGeom)
1134     {
1135         /* recalc geometrical center and max and min coordinates and size */
1136         calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1137         rvec_sub(rmax, rmin, size);
1138         if (bScale || bOrient || bRotate)
1139         {
1140             printf("new system size : %6.3f %6.3f %6.3f\n", size[XX], size[YY], size[ZZ]);
1141         }
1142     }
1143
1144     if ((btype[0] != nullptr) && (bSetSize || bDist || (btype[0][0] == 't' && bSetAng)))
1145     {
1146         pbcType = PbcType::Xyz;
1147         if (!(bSetSize || bDist))
1148         {
1149             for (i = 0; i < DIM; i++)
1150             {
1151                 newbox[i] = norm(box[i]);
1152             }
1153         }
1154         clear_mat(box);
1155         /* calculate new boxsize */
1156         switch (btype[0][0])
1157         {
1158             case 't':
1159                 if (bDist)
1160                 {
1161                     for (i = 0; i < DIM; i++)
1162                     {
1163                         newbox[i] = size[i] + 2 * dist;
1164                     }
1165                 }
1166                 if (!bSetAng)
1167                 {
1168                     box[XX][XX] = newbox[XX];
1169                     box[YY][YY] = newbox[YY];
1170                     box[ZZ][ZZ] = newbox[ZZ];
1171                 }
1172                 else
1173                 {
1174                     matrix_convert(box, newbox, newang);
1175                 }
1176                 break;
1177             case 'c':
1178             case 'd':
1179             case 'o':
1180                 if (bSetSize)
1181                 {
1182                     d = newbox[0];
1183                 }
1184                 else
1185                 {
1186                     d = diam + 2 * dist;
1187                 }
1188                 if (btype[0][0] == 'c')
1189                 {
1190                     for (i = 0; i < DIM; i++)
1191                     {
1192                         box[i][i] = d;
1193                     }
1194                 }
1195                 else if (btype[0][0] == 'd')
1196                 {
1197                     box[XX][XX] = d;
1198                     box[YY][YY] = d;
1199                     box[ZZ][XX] = d / 2;
1200                     box[ZZ][YY] = d / 2;
1201                     box[ZZ][ZZ] = d * std::sqrt(2.0) / 2.0;
1202                 }
1203                 else
1204                 {
1205                     box[XX][XX] = d;
1206                     box[YY][XX] = d / 3;
1207                     box[YY][YY] = d * std::sqrt(2.0) * 2.0 / 3.0;
1208                     box[ZZ][XX] = -d / 3;
1209                     box[ZZ][YY] = d * std::sqrt(2.0) / 3.0;
1210                     box[ZZ][ZZ] = d * std::sqrt(6.0) / 3.0;
1211                 }
1212                 break;
1213         }
1214     }
1215
1216     /* calculate new coords for geometrical center */
1217     if (!bSetCenter)
1218     {
1219         calc_box_center(ecenterDEF, box, center);
1220     }
1221
1222     /* center molecule on 'center' */
1223     if (bCenter)
1224     {
1225         center_conf(natom, x, center, gc);
1226     }
1227
1228     /* print some */
1229     if (bCalcGeom)
1230     {
1231         calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1232         printf("new center      :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
1233     }
1234     if (bOrient || bScale || bDist || bSetSize)
1235     {
1236         printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n", norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
1237         printf("new box angles  :%7.2f%7.2f%7.2f (degrees)\n",
1238                norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[YY], box[ZZ]),
1239                norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[ZZ]),
1240                norm2(box[YY]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[YY]));
1241         printf("new box volume  :%7.2f               (nm^3)\n", det(box));
1242     }
1243
1244     if (check_box(PbcType::Xyz, box))
1245     {
1246         printf("\nWARNING: %s\n"
1247                "See the GROMACS manual for a description of the requirements that\n"
1248                "must be satisfied by descriptions of simulation cells.\n",
1249                check_box(PbcType::Xyz, box));
1250     }
1251
1252     if (bDist && btype[0][0] == 't')
1253     {
1254         if (TRICLINIC(box))
1255         {
1256             printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1257                    "distance from the solute to a box surface along the corresponding normal\n"
1258                    "vector might be somewhat smaller than your specified value %f.\n"
1259                    "You can check the actual value with g_mindist -pi\n",
1260                    dist);
1261         }
1262         else if (!opt2parg_bSet("-bt", NPA, pa))
1263         {
1264             printf("\nWARNING: No boxtype specified - distance condition applied in each "
1265                    "dimension.\n"
1266                    "If the molecule rotates the actual distance will be smaller. You might want\n"
1267                    "to use a cubic box instead, or why not try a dodecahedron today?\n");
1268         }
1269     }
1270     if (bCONECT && (outftp == efPDB) && (inftp == efTPR))
1271     {
1272         conect = gmx_conect_generate(top);
1273     }
1274     else
1275     {
1276         conect = nullptr;
1277     }
1278
1279     if (bIndex)
1280     {
1281         fprintf(stderr, "\nSelect a group for output:\n");
1282         get_index(&atoms, opt2fn_null("-n", NFILE, fnm), 1, &isize, &index, &grpname);
1283
1284         if (resnr_start >= 0)
1285         {
1286             renum_resnr(&atoms, isize, index, resnr_start);
1287         }
1288
1289         if (opt2parg_bSet("-label", NPA, pa))
1290         {
1291             for (i = 0; (i < atoms.nr); i++)
1292             {
1293                 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1294             }
1295         }
1296
1297         if (opt2bSet("-bf", NFILE, fnm) || bLegend)
1298         {
1299             gmx_fatal(FARGS, "Sorry, cannot do bfactors with an index group.");
1300         }
1301
1302         if (outftp == efPDB)
1303         {
1304             out = gmx_ffopen(outfile, "w");
1305             write_pdbfile_indexed(out, name, &atoms, x, pbcType, box, ' ', 1, isize, index, conect, FALSE);
1306             gmx_ffclose(out);
1307         }
1308         else
1309         {
1310             write_sto_conf_indexed(
1311                     outfile, name, &atoms, x, bHaveV ? v : nullptr, pbcType, box, isize, index);
1312         }
1313         sfree(grpname);
1314         sfree(index);
1315     }
1316     else
1317     {
1318         if (resnr_start >= 0)
1319         {
1320             renum_resnr(&atoms, atoms.nr, nullptr, resnr_start);
1321         }
1322
1323         if ((outftp == efPDB) || (outftp == efPQR))
1324         {
1325             out = gmx_ffopen(outfile, "w");
1326             if (bMead)
1327             {
1328                 fprintf(out,
1329                         "REMARK    "
1330                         "The B-factors in this file hold atomic radii\n");
1331                 fprintf(out,
1332                         "REMARK    "
1333                         "The occupancy in this file hold atomic charges\n");
1334             }
1335             else if (bGrasp)
1336             {
1337                 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1338                 fprintf(out,
1339                         "REMARK    "
1340                         "The B-factors in this file hold atomic charges\n");
1341                 fprintf(out,
1342                         "REMARK    "
1343                         "The occupancy in this file hold atomic radii\n");
1344             }
1345             else if (opt2bSet("-bf", NFILE, fnm))
1346             {
1347                 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1348                 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms, n_bfac, bfac, bfac_nr, peratom);
1349             }
1350             if (opt2parg_bSet("-label", NPA, pa))
1351             {
1352                 for (i = 0; (i < atoms.nr); i++)
1353                 {
1354                     atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1355                 }
1356             }
1357             /* Need to bypass the regular write_pdbfile because I don't want to change
1358              * all instances to include the boolean flag for writing out PQR files.
1359              */
1360             int* index;
1361             snew(index, atoms.nr);
1362             for (int i = 0; i < atoms.nr; i++)
1363             {
1364                 index[i] = i;
1365             }
1366             write_pdbfile_indexed(
1367                     out, name, &atoms, x, pbcType, box, ' ', -1, atoms.nr, index, conect, outftp == efPQR);
1368             sfree(index);
1369             if (bLegend)
1370             {
1371                 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1372             }
1373             if (visbox[0] > 0)
1374             {
1375                 visualize_box(out,
1376                               bLegend ? atoms.nr + 12 : atoms.nr,
1377                               bLegend ? atoms.nres = 12 : atoms.nres,
1378                               box,
1379                               visbox);
1380             }
1381             gmx_ffclose(out);
1382         }
1383         else
1384         {
1385             write_sto_conf(outfile, name, &atoms, x, bHaveV ? v : nullptr, pbcType, box);
1386         }
1387     }
1388     done_atom(&atoms);
1389     done_symtab(&symtab);
1390     sfree(name);
1391     if (x)
1392     {
1393         sfree(x);
1394     }
1395     if (v)
1396     {
1397         sfree(v);
1398     }
1399     do_view(oenv, outfile, nullptr);
1400     output_env_done(oenv);
1401
1402     return 0;
1403 }