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