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