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