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