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