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