Merge branch release-5-1 into release-2016
[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, 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 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 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 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 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 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 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 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 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, NULL, ePBC, box);
364
365     done_atom(&atoms);
366     sfree(img);
367 }
368
369 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 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 == NULL ? 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         "[PAR]",
517         "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
518         "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
519         "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
520         "[TT]octahedron[tt] is a truncated octahedron.",
521         "The last two are special cases of a triclinic box.",
522         "The length of the three box vectors of the truncated octahedron is the",
523         "shortest distance between two opposite hexagons.",
524         "Relative to a cubic box with some periodic image distance, the volume of a ",
525         "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
526         "and that of a truncated octahedron is 0.77 times.",
527         "[PAR]",
528         "Option [TT]-box[tt] requires only",
529         "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
530         "[PAR]",
531         "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, [IT]y[it]-,",
532         "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
533         "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
534         "to the diameter of the system (largest distance between atoms) plus twice",
535         "the specified distance.",
536         "[PAR]",
537         "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
538         "a triclinic box and cannot be used with option [TT]-d[tt].",
539         "[PAR]",
540         "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
541         "can be selected for calculating the size and the geometric center,",
542         "otherwise the whole system is used.",
543         "[PAR]",
544         "[TT]-rotate[tt] rotates the coordinates and velocities.",
545         "[PAR]",
546         "[TT]-princ[tt] aligns the principal axes of the system along the",
547         "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
548         "This may allow you to decrease the box volume,",
549         "but beware that molecules can rotate significantly in a nanosecond.",
550         "[PAR]",
551         "Scaling is applied before any of the other operations are",
552         "performed. Boxes and coordinates can be scaled to give a certain density (option",
553         "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
554         "file is given as input. A special feature of the scaling option is that when the",
555         "factor -1 is given in one dimension, one obtains a mirror image,",
556         "mirrored in one of the planes. When one uses -1 in three dimensions, ",
557         "a point-mirror image is obtained.[PAR]",
558         "Groups are selected after all operations have been applied.[PAR]",
559         "Periodicity can be removed in a crude manner.",
560         "It is important that the box vectors at the bottom of your input file",
561         "are correct when the periodicity is to be removed.",
562         "[PAR]",
563         "When writing [REF].pdb[ref] files, B-factors can be",
564         "added with the [TT]-bf[tt] option. B-factors are read",
565         "from a file with with following format: first line states number of",
566         "entries in the file, next lines state an index",
567         "followed by a B-factor. The B-factors will be attached per residue",
568         "unless the number of B-factors is larger than the number of the residues or unless the",
569         "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
570         "be added instead of B-factors. [TT]-legend[tt] will produce",
571         "a row of CA atoms with B-factors ranging from the minimum to the",
572         "maximum value found, effectively making a legend for viewing.",
573         "[PAR]",
574         "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
575         "file for the MEAD electrostatics",
576         "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
577         "is that the input file is a run input file.",
578         "The B-factor field is then filled with the Van der Waals radius",
579         "of the atoms while the occupancy field will hold the charge.",
580         "[PAR]",
581         "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
582         "and the radius in the occupancy.",
583         "[PAR]",
584         "Option [TT]-align[tt] allows alignment",
585         "of the principal axis of a specified group against the given vector, ",
586         "with an optional center of rotation specified by [TT]-aligncenter[tt].",
587         "[PAR]",
588         "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
589         "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
590         "[PAR]",
591         "To convert a truncated octrahedron file produced by a package which uses",
592         "a cubic box with the corners cut off (such as GROMOS), use::",
593         "",
594         "  gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
595         "",
596         "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
597     };
598     const char     *bugs[] =
599     {
600         "For complex molecules, the periodicity removal routine may break down, "
601         "in that case you can use [gmx-trjconv]."
602     };
603     static real     dist    = 0.0;
604     static gmx_bool bNDEF   = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW =
605         FALSE, bCONECT      = FALSE;
606     static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead =
607         FALSE, bGrasp       = FALSE, bSig56 = FALSE;
608     static rvec     scale   =
609     { 1, 1, 1 }, newbox     =
610     { 0, 0, 0 }, newang     =
611     { 90, 90, 90 };
612     static real rho          = 1000.0, rvdw = 0.12;
613     static rvec center       =
614     { 0, 0, 0 }, translation =
615     { 0, 0, 0 }, rotangles   =
616     { 0, 0, 0 }, aligncenter =
617     { 0, 0, 0 }, targetvec   =
618     { 0, 0, 0 };
619     static const char *btype[] =
620     { NULL, "triclinic", "cubic", "dodecahedron", "octahedron", NULL },
621     *label             = "A";
622     static rvec visbox =
623     { 0, 0, 0 };
624     static int  resnr_start = -1;
625     t_pargs
626                 pa[] =
627     {
628         { "-ndef", FALSE, etBOOL,
629           { &bNDEF }, "Choose output from default index groups" },
630         { "-visbox", FALSE, etRVEC,
631           { visbox },
632           "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
633         { "-bt", FALSE, etENUM,
634           { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
635         { "-box", FALSE, etRVEC,
636           { newbox }, "Box vector lengths (a,b,c)" },
637         { "-angles", FALSE, etRVEC,
638           { newang }, "Angles between the box vectors (bc,ac,ab)" },
639         { "-d", FALSE, etREAL,
640           { &dist }, "Distance between the solute and the box" },
641         { "-c", FALSE, etBOOL,
642           { &bCenter },
643           "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
644         { "-center", FALSE, etRVEC,
645           { center }, "Coordinates of geometrical center" },
646         { "-aligncenter", FALSE, etRVEC,
647           { aligncenter }, "Center of rotation for alignment" },
648         { "-align", FALSE, etRVEC,
649           { targetvec },
650           "Align to target vector" },
651         { "-translate", FALSE, etRVEC,
652           { translation }, "Translation" },
653         { "-rotate", FALSE, etRVEC,
654           { rotangles },
655           "Rotation around the X, Y and Z axes in degrees" },
656         { "-princ", FALSE, etBOOL,
657           { &bOrient },
658           "Orient molecule(s) along their principal axes" },
659         { "-scale", FALSE, etRVEC,
660           { scale }, "Scaling factor" },
661         { "-density", FALSE, etREAL,
662           { &rho },
663           "Density (g/L) of the output box achieved by scaling" },
664         { "-pbc", FALSE, etBOOL,
665           { &bRMPBC },
666           "Remove the periodicity (make molecule whole again)" },
667         { "-resnr", FALSE, etINT,
668           { &resnr_start },
669           " Renumber residues starting from resnr" },
670         { "-grasp", FALSE, etBOOL,
671           { &bGrasp },
672           "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" },
673         {
674             "-rvdw", FALSE, etREAL,
675             { &rvdw },
676             "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"
677         },
678         { "-sig56", FALSE, etBOOL,
679           { &bSig56 },
680           "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
681         {
682             "-vdwread", FALSE, etBOOL,
683             { &bReadVDW },
684             "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing the radii based on the force field"
685         },
686         { "-atom", FALSE, etBOOL,
687           { &peratom }, "Force B-factor attachment per atom" },
688         { "-legend", FALSE, etBOOL,
689           { &bLegend }, "Make B-factor legend" },
690         { "-label", FALSE, etSTR,
691           { &label }, "Add chain label for all residues" },
692         {
693             "-conect", FALSE, etBOOL,
694             { &bCONECT },
695             "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a topology is present"
696         }
697     };
698 #define NPA asize(pa)
699
700     FILE             *out;
701     const char       *infile, *outfile;
702     int               outftp, inftp, natom, i, j, n_bfac, itype, ntype;
703     double           *bfac    = NULL, c6, c12;
704     int              *bfac_nr = NULL;
705     t_topology       *top     = NULL;
706     char             *grpname, *sgrpname, *agrpname;
707     int               isize, ssize, asize;
708     int              *index, *sindex, *aindex;
709     rvec             *x, *v, gc, rmin, rmax, size;
710     int               ePBC;
711     matrix            box, rotmatrix, trans;
712     rvec              princd, tmpvec;
713     gmx_bool          bIndex, bSetSize, bSetAng, bDist, bSetCenter, bAlign;
714     gmx_bool          bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
715     real              diam = 0, mass = 0, d, vdw;
716     gmx_atomprop_t    aps;
717     gmx_conect        conect;
718     gmx_output_env_t *oenv;
719     t_filenm          fnm[] =
720     {
721         { efSTX, "-f", NULL, ffREAD },
722         { efNDX, "-n", NULL, ffOPTRD },
723         { efSTO, NULL, NULL, ffOPTWR },
724         { efPQR, "-mead", "mead", ffOPTWR },
725         { efDAT, "-bf", "bfact", ffOPTRD }
726     };
727 #define NFILE asize(fnm)
728
729     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa,
730                            asize(desc), desc, asize(bugs), bugs, &oenv))
731     {
732         return 0;
733     }
734
735     bIndex     = opt2bSet("-n", NFILE, fnm) || bNDEF;
736     bMead      = opt2bSet("-mead", NFILE, fnm);
737     bSetSize   = opt2parg_bSet("-box", NPA, pa);
738     bSetAng    = opt2parg_bSet("-angles", NPA, pa);
739     bSetCenter = opt2parg_bSet("-center", NPA, pa);
740     bDist      = opt2parg_bSet("-d", NPA, pa);
741     bAlign     = opt2parg_bSet("-align", NPA, pa);
742     /* Only automatically turn on centering without -noc */
743     if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
744     {
745         bCenter = TRUE;
746     }
747     bScale     = opt2parg_bSet("-scale", NPA, pa);
748     bRho       = opt2parg_bSet("-density", NPA, pa);
749     bTranslate = opt2parg_bSet("-translate", NPA, pa);
750     bRotate    = opt2parg_bSet("-rotate", NPA, pa);
751     if (bScale && bRho)
752     {
753         fprintf(stderr, "WARNING: setting -density overrides -scale\n");
754     }
755     bScale    = bScale || bRho;
756     bCalcGeom = bCenter || bRotate || bOrient || bScale;
757
758     GMX_RELEASE_ASSERT(btype[0] != NULL, "Option setting inconsistency; btype[0] is NULL");
759
760     bCalcDiam = (btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o');
761
762     infile = ftp2fn(efSTX, NFILE, fnm);
763     if (bMead)
764     {
765         outfile = ftp2fn(efPQR, NFILE, fnm);
766     }
767     else
768     {
769         outfile = ftp2fn(efSTO, NFILE, fnm);
770     }
771     outftp = fn2ftp(outfile);
772     inftp  = fn2ftp(infile);
773
774     aps = gmx_atomprop_init();
775
776     if (bMead && bGrasp)
777     {
778         printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
779         bGrasp = FALSE;
780     }
781     if (bGrasp && (outftp != efPDB))
782     {
783         gmx_fatal(FARGS, "Output file should be a .pdb file"
784                   " when using the -grasp option\n");
785     }
786     if ((bMead || bGrasp) && (fn2ftp(infile) != efTPR))
787     {
788         gmx_fatal(FARGS, "Input file should be a .tpr file"
789                   " when using the -mead option\n");
790     }
791
792     t_topology *top_tmp;
793     snew(top_tmp, 1);
794     read_tps_conf(infile, top_tmp, &ePBC, &x, &v, box, FALSE);
795     t_atoms  &atoms = top_tmp->atoms;
796     natom = atoms.nr;
797     if (atoms.pdbinfo == NULL)
798     {
799         snew(atoms.pdbinfo, atoms.nr);
800     }
801     if (fn2ftp(infile) == efPDB)
802     {
803         get_pdb_atomnumber(&atoms, aps);
804     }
805     printf("Read %d atoms\n", atoms.nr);
806
807     /* Get the element numbers if available in a pdb file */
808     if (fn2ftp(infile) == efPDB)
809     {
810         get_pdb_atomnumber(&atoms, aps);
811     }
812
813     if (ePBC != epbcNONE)
814     {
815         real vol = det(box);
816         printf("Volume: %g nm^3, corresponds to roughly %d electrons\n",
817                vol, 100*(static_cast<int>(vol*4.5)));
818     }
819
820     if (bMead || bGrasp || bCONECT)
821     {
822         top = read_top(infile, NULL);
823     }
824
825     if (bMead || bGrasp)
826     {
827         if (atoms.nr != top->atoms.nr)
828         {
829             gmx_fatal(FARGS, "Atom numbers don't match (%d vs. %d)", atoms.nr, top->atoms.nr);
830         }
831         snew(atoms.pdbinfo, top->atoms.nr);
832         ntype = top->idef.atnr;
833         for (i = 0; (i < atoms.nr); i++)
834         {
835             /* Determine the Van der Waals radius from the force field */
836             if (bReadVDW)
837             {
838                 if (!gmx_atomprop_query(aps, epropVDW,
839                                         *top->atoms.resinfo[top->atoms.atom[i].resind].name,
840                                         *top->atoms.atomname[i], &vdw))
841                 {
842                     vdw = rvdw;
843                 }
844             }
845             else
846             {
847                 itype = top->atoms.atom[i].type;
848                 c12   = top->idef.iparams[itype*ntype+itype].lj.c12;
849                 c6    = top->idef.iparams[itype*ntype+itype].lj.c6;
850                 if ((c6 != 0) && (c12 != 0))
851                 {
852                     real sig6;
853                     if (bSig56)
854                     {
855                         sig6 = 2*c12/c6;
856                     }
857                     else
858                     {
859                         sig6 = c12/c6;
860                     }
861                     vdw   = 0.5*gmx::sixthroot(sig6);
862                 }
863                 else
864                 {
865                     vdw = rvdw;
866                 }
867             }
868             /* Factor of 10 for nm -> Angstroms */
869             vdw *= 10;
870
871             if (bMead)
872             {
873                 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
874                 atoms.pdbinfo[i].bfac  = vdw;
875             }
876             else
877             {
878                 atoms.pdbinfo[i].occup = vdw;
879                 atoms.pdbinfo[i].bfac  = top->atoms.atom[i].q;
880             }
881         }
882     }
883     bHaveV = FALSE;
884     for (i = 0; (i < natom) && !bHaveV; i++)
885     {
886         for (j = 0; (j < DIM) && !bHaveV; j++)
887         {
888             bHaveV = bHaveV || (v[i][j] != 0);
889         }
890     }
891     printf("%selocities found\n", bHaveV ? "V" : "No v");
892
893     if (visbox[0] > 0)
894     {
895         if (bIndex)
896         {
897             gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
898         }
899         if (outftp != efPDB)
900         {
901             gmx_fatal(FARGS, "Sorry, can only visualize box with a pdb file");
902         }
903     }
904     else if (visbox[0] == -1)
905     {
906         visualize_images("images.pdb", ePBC, box);
907     }
908
909     /* remove pbc */
910     if (bRMPBC)
911     {
912         rm_gropbc(&atoms, x, box);
913     }
914
915     if (bCalcGeom)
916     {
917         if (bIndex)
918         {
919             fprintf(stderr, "\nSelect a group for determining the system size:\n");
920             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
921                       1, &ssize, &sindex, &sgrpname);
922         }
923         else
924         {
925             ssize  = atoms.nr;
926             sindex = NULL;
927         }
928         diam = calc_geom(ssize, sindex, x, gc, rmin, rmax, bCalcDiam);
929         rvec_sub(rmax, rmin, size);
930         printf("    system size :%7.3f%7.3f%7.3f (nm)\n",
931                size[XX], size[YY], size[ZZ]);
932         if (bCalcDiam)
933         {
934             printf("    diameter    :%7.3f               (nm)\n", diam);
935         }
936         printf("    center      :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
937         printf("    box vectors :%7.3f%7.3f%7.3f (nm)\n",
938                norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
939         printf("    box angles  :%7.2f%7.2f%7.2f (degrees)\n",
940                norm2(box[ZZ]) == 0 ? 0 :
941                RAD2DEG*gmx_angle(box[YY], box[ZZ]),
942                norm2(box[ZZ]) == 0 ? 0 :
943                RAD2DEG*gmx_angle(box[XX], box[ZZ]),
944                norm2(box[YY]) == 0 ? 0 :
945                RAD2DEG*gmx_angle(box[XX], box[YY]));
946         printf("    box volume  :%7.2f               (nm^3)\n", det(box));
947     }
948
949     if (bRho || bOrient || bAlign)
950     {
951         mass = calc_mass(&atoms, !fn2bTPX(infile), aps);
952     }
953
954     if (bOrient)
955     {
956         int     *index;
957         char    *grpnames;
958
959         /* Get a group for principal component analysis */
960         fprintf(stderr, "\nSelect group for the determining the orientation\n");
961         get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
962
963         /* Orient the principal axes along the coordinate axes */
964         orient_princ(&atoms, isize, index, natom, x, bHaveV ? v : NULL, NULL);
965         sfree(index);
966         sfree(grpnames);
967     }
968
969     if (bScale)
970     {
971         /* scale coordinates and box */
972         if (bRho)
973         {
974             /* Compute scaling constant */
975             real vol, dens;
976
977             vol  = det(box);
978             dens = (mass*AMU)/(vol*NANO*NANO*NANO);
979             fprintf(stderr, "Volume  of input %g (nm^3)\n", vol);
980             fprintf(stderr, "Mass    of input %g (a.m.u.)\n", mass);
981             fprintf(stderr, "Density of input %g (g/l)\n", dens);
982             if (vol == 0 || mass == 0)
983             {
984                 gmx_fatal(FARGS, "Cannot scale density with "
985                           "zero mass (%g) or volume (%g)\n", mass, vol);
986             }
987
988             scale[XX] = scale[YY] = scale[ZZ] = std::cbrt(dens/rho);
989             fprintf(stderr, "Scaling all box vectors by %g\n", scale[XX]);
990         }
991         scale_conf(atoms.nr, x, box, scale);
992     }
993
994     if (bAlign)
995     {
996         if (bIndex)
997         {
998             fprintf(stderr, "\nSelect a group that you want to align:\n");
999             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
1000                       1, &asize, &aindex, &agrpname);
1001         }
1002         else
1003         {
1004             asize = atoms.nr;
1005             snew(aindex, asize);
1006             for (i = 0; i < asize; i++)
1007             {
1008                 aindex[i] = i;
1009             }
1010         }
1011         printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n", asize, natom,
1012                targetvec[XX], targetvec[YY], targetvec[ZZ],
1013                aligncenter[XX], aligncenter[YY], aligncenter[ZZ]);
1014         /*subtract out pivot point*/
1015         for (i = 0; i < asize; i++)
1016         {
1017             rvec_dec(x[aindex[i]], aligncenter);
1018         }
1019         /*now determine transform and rotate*/
1020         /*will this work?*/
1021         principal_comp(asize, aindex, atoms.atom, x, trans, princd);
1022
1023         unitv(targetvec, targetvec);
1024         printf("Using %g %g %g as principal axis\n", trans[0][2], trans[1][2], trans[2][2]);
1025         tmpvec[XX] = trans[0][2]; tmpvec[YY] = trans[1][2]; tmpvec[ZZ] = trans[2][2];
1026         calc_rotmatrix(tmpvec, targetvec, rotmatrix);
1027         /* rotmatrix finished */
1028
1029         for (i = 0; i < asize; ++i)
1030         {
1031             mvmul(rotmatrix, x[aindex[i]], tmpvec);
1032             copy_rvec(tmpvec, x[aindex[i]]);
1033         }
1034
1035         /*add pivot point back*/
1036         for (i = 0; i < asize; i++)
1037         {
1038             rvec_inc(x[aindex[i]], aligncenter);
1039         }
1040         if (!bIndex)
1041         {
1042             sfree(aindex);
1043         }
1044     }
1045
1046     if (bTranslate)
1047     {
1048         if (bIndex)
1049         {
1050             fprintf(stderr, "\nSelect a group that you want to translate:\n");
1051             get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
1052                       1, &ssize, &sindex, &sgrpname);
1053         }
1054         else
1055         {
1056             ssize  = atoms.nr;
1057             sindex = NULL;
1058         }
1059         printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize, natom,
1060                translation[XX], translation[YY], translation[ZZ]);
1061         if (sindex)
1062         {
1063             for (i = 0; i < ssize; i++)
1064             {
1065                 rvec_inc(x[sindex[i]], translation);
1066             }
1067         }
1068         else
1069         {
1070             for (i = 0; i < natom; i++)
1071             {
1072                 rvec_inc(x[i], translation);
1073             }
1074         }
1075     }
1076     if (bRotate)
1077     {
1078         /* Rotate */
1079         printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n", rotangles[XX], rotangles[YY], rotangles[ZZ]);
1080         for (i = 0; i < DIM; i++)
1081         {
1082             rotangles[i] *= DEG2RAD;
1083         }
1084         rotate_conf(natom, x, v, rotangles[XX], rotangles[YY], rotangles[ZZ]);
1085     }
1086
1087     if (bCalcGeom)
1088     {
1089         /* recalc geometrical center and max and min coordinates and size */
1090         calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1091         rvec_sub(rmax, rmin, size);
1092         if (bScale || bOrient || bRotate)
1093         {
1094             printf("new system size : %6.3f %6.3f %6.3f\n",
1095                    size[XX], size[YY], size[ZZ]);
1096         }
1097     }
1098
1099     if ((btype[0] != NULL) && (bSetSize || bDist || (btype[0][0] == 't' && bSetAng)))
1100     {
1101         ePBC = epbcXYZ;
1102         if (!(bSetSize || bDist))
1103         {
1104             for (i = 0; i < DIM; i++)
1105             {
1106                 newbox[i] = norm(box[i]);
1107             }
1108         }
1109         clear_mat(box);
1110         /* calculate new boxsize */
1111         switch (btype[0][0])
1112         {
1113             case 't':
1114                 if (bDist)
1115                 {
1116                     for (i = 0; i < DIM; i++)
1117                     {
1118                         newbox[i] = size[i]+2*dist;
1119                     }
1120                 }
1121                 if (!bSetAng)
1122                 {
1123                     box[XX][XX] = newbox[XX];
1124                     box[YY][YY] = newbox[YY];
1125                     box[ZZ][ZZ] = newbox[ZZ];
1126                 }
1127                 else
1128                 {
1129                     matrix_convert(box, newbox, newang);
1130                 }
1131                 break;
1132             case 'c':
1133             case 'd':
1134             case 'o':
1135                 if (bSetSize)
1136                 {
1137                     d = newbox[0];
1138                 }
1139                 else
1140                 {
1141                     d = diam+2*dist;
1142                 }
1143                 if (btype[0][0] == 'c')
1144                 {
1145                     for (i = 0; i < DIM; i++)
1146                     {
1147                         box[i][i] = d;
1148                     }
1149                 }
1150                 else if (btype[0][0] == 'd')
1151                 {
1152                     box[XX][XX] = d;
1153                     box[YY][YY] = d;
1154                     box[ZZ][XX] = d/2;
1155                     box[ZZ][YY] = d/2;
1156                     box[ZZ][ZZ] = d*std::sqrt(2.0)/2.0;
1157                 }
1158                 else
1159                 {
1160                     box[XX][XX] = d;
1161                     box[YY][XX] = d/3;
1162                     box[YY][YY] = d*std::sqrt(2.0)*2.0/3.0;
1163                     box[ZZ][XX] = -d/3;
1164                     box[ZZ][YY] = d*std::sqrt(2.0)/3.0;
1165                     box[ZZ][ZZ] = d*std::sqrt(6.0)/3.0;
1166                 }
1167                 break;
1168         }
1169     }
1170
1171     /* calculate new coords for geometrical center */
1172     if (!bSetCenter)
1173     {
1174         calc_box_center(ecenterDEF, box, center);
1175     }
1176
1177     /* center molecule on 'center' */
1178     if (bCenter)
1179     {
1180         center_conf(natom, x, center, gc);
1181     }
1182
1183     /* print some */
1184     if (bCalcGeom)
1185     {
1186         calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1187         printf("new center      :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
1188     }
1189     if (bOrient || bScale || bDist || bSetSize)
1190     {
1191         printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n",
1192                norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
1193         printf("new box angles  :%7.2f%7.2f%7.2f (degrees)\n",
1194                norm2(box[ZZ]) == 0 ? 0 :
1195                RAD2DEG*gmx_angle(box[YY], box[ZZ]),
1196                norm2(box[ZZ]) == 0 ? 0 :
1197                RAD2DEG*gmx_angle(box[XX], box[ZZ]),
1198                norm2(box[YY]) == 0 ? 0 :
1199                RAD2DEG*gmx_angle(box[XX], box[YY]));
1200         printf("new box volume  :%7.2f               (nm^3)\n", det(box));
1201     }
1202
1203     if (check_box(epbcXYZ, box))
1204     {
1205         printf("\nWARNING: %s\n"
1206                "See the GROMACS manual for a description of the requirements that\n"
1207                "must be satisfied by descriptions of simulation cells.\n",
1208                check_box(epbcXYZ, box));
1209     }
1210
1211     if (bDist && btype[0][0] == 't')
1212     {
1213         if (TRICLINIC(box))
1214         {
1215             printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1216                    "distance from the solute to a box surface along the corresponding normal\n"
1217                    "vector might be somewhat smaller than your specified value %f.\n"
1218                    "You can check the actual value with g_mindist -pi\n", dist);
1219         }
1220         else if (!opt2parg_bSet("-bt", NPA, pa))
1221         {
1222             printf("\nWARNING: No boxtype specified - distance condition applied in each dimension.\n"
1223                    "If the molecule rotates the actual distance will be smaller. You might want\n"
1224                    "to use a cubic box instead, or why not try a dodecahedron today?\n");
1225         }
1226     }
1227     if (bCONECT && (outftp == efPDB) && (inftp == efTPR))
1228     {
1229         conect = gmx_conect_generate(top);
1230     }
1231     else
1232     {
1233         conect = NULL;
1234     }
1235
1236     if (bIndex)
1237     {
1238         fprintf(stderr, "\nSelect a group for output:\n");
1239         get_index(&atoms, opt2fn_null("-n", NFILE, fnm),
1240                   1, &isize, &index, &grpname);
1241
1242         if (resnr_start >= 0)
1243         {
1244             renum_resnr(&atoms, isize, index, resnr_start);
1245         }
1246
1247         if (opt2parg_bSet("-label", NPA, pa))
1248         {
1249             for (i = 0; (i < atoms.nr); i++)
1250             {
1251                 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1252             }
1253         }
1254
1255         if (opt2bSet("-bf", NFILE, fnm) || bLegend)
1256         {
1257             gmx_fatal(FARGS, "Sorry, cannot do bfactors with an index group.");
1258         }
1259
1260         if (outftp == efPDB)
1261         {
1262             out = gmx_ffopen(outfile, "w");
1263             write_pdbfile_indexed(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, TRUE);
1264             gmx_ffclose(out);
1265         }
1266         else
1267         {
1268             write_sto_conf_indexed(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : NULL, ePBC, box, isize, index);
1269         }
1270     }
1271     else
1272     {
1273         if (resnr_start >= 0)
1274         {
1275             renum_resnr(&atoms, atoms.nr, NULL, resnr_start);
1276         }
1277
1278         if ((outftp == efPDB) || (outftp == efPQR))
1279         {
1280             out = gmx_ffopen(outfile, "w");
1281             if (bMead)
1282             {
1283                 fprintf(out, "REMARK    "
1284                         "The B-factors in this file hold atomic radii\n");
1285                 fprintf(out, "REMARK    "
1286                         "The occupancy in this file hold atomic charges\n");
1287             }
1288             else if (bGrasp)
1289             {
1290                 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1291                 fprintf(out, "REMARK    "
1292                         "The B-factors in this file hold atomic charges\n");
1293                 fprintf(out, "REMARK    "
1294                         "The occupancy in this file hold atomic radii\n");
1295             }
1296             else if (opt2bSet("-bf", NFILE, fnm))
1297             {
1298                 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1299                 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms,
1300                                   n_bfac, bfac, bfac_nr, peratom);
1301             }
1302             if (opt2parg_bSet("-label", NPA, pa))
1303             {
1304                 for (i = 0; (i < atoms.nr); i++)
1305                 {
1306                     atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1307                 }
1308             }
1309             write_pdbfile(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', -1, conect, TRUE);
1310             if (bLegend)
1311             {
1312                 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1313             }
1314             if (visbox[0] > 0)
1315             {
1316                 visualize_box(out, bLegend ? atoms.nr+12 : atoms.nr,
1317                               bLegend ? atoms.nres = 12 : atoms.nres, box, visbox);
1318             }
1319             gmx_ffclose(out);
1320         }
1321         else
1322         {
1323             write_sto_conf(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : NULL, ePBC, box);
1324         }
1325     }
1326     gmx_atomprop_destroy(aps);
1327
1328     do_view(oenv, outfile, NULL);
1329
1330     return 0;
1331 }