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