Split txtdump.*, part 1
[gromacs.git] / src / gromacs / gmxana / gmx_editconf.cpp
blob4ce459f3679f62a2367a132cc7fb9fe6c9652c09
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstring>
42 #include <algorithm>
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/pdbio.h"
48 #include "gromacs/fileio/tpxio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/gmxana/princ.h"
52 #include "gromacs/gmxlib/conformation-utilities.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/pbcutil/rmpbc.h"
57 #include "gromacs/topology/atomprop.h"
58 #include "gromacs/topology/index.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/smalloc.h"
66 #include "gromacs/utility/strdb.h"
69 real calc_mass(t_atoms *atoms, gmx_bool bGetMass, gmx_atomprop_t aps)
71 real tmass;
72 int i;
74 tmass = 0;
75 for (i = 0; (i < atoms->nr); i++)
77 if (bGetMass)
79 gmx_atomprop_query(aps, epropMass,
80 *atoms->resinfo[atoms->atom[i].resind].name,
81 *atoms->atomname[i], &(atoms->atom[i].m));
83 tmass += atoms->atom[i].m;
86 return tmass;
89 real calc_geom(int isize, int *index, rvec *x, rvec geom_center, rvec minval,
90 rvec maxval, gmx_bool bDiam)
92 real diam2, d;
93 int ii, i, j;
95 clear_rvec(geom_center);
96 diam2 = 0;
97 if (isize == 0)
99 clear_rvec(minval);
100 clear_rvec(maxval);
102 else
104 if (index)
106 ii = index[0];
108 else
110 ii = 0;
112 for (j = 0; j < DIM; j++)
114 minval[j] = maxval[j] = x[ii][j];
116 for (i = 0; i < isize; i++)
118 if (index)
120 ii = index[i];
122 else
124 ii = i;
126 rvec_inc(geom_center, x[ii]);
127 for (j = 0; j < DIM; j++)
129 if (x[ii][j] < minval[j])
131 minval[j] = x[ii][j];
133 if (x[ii][j] > maxval[j])
135 maxval[j] = x[ii][j];
138 if (bDiam)
140 if (index)
142 for (j = i + 1; j < isize; j++)
144 d = distance2(x[ii], x[index[j]]);
145 diam2 = std::max(d, diam2);
148 else
150 for (j = i + 1; j < isize; j++)
152 d = distance2(x[i], x[j]);
153 diam2 = std::max(d, diam2);
158 svmul(1.0 / isize, geom_center, geom_center);
161 return std::sqrt(diam2);
164 void center_conf(int natom, rvec *x, rvec center, rvec geom_cent)
166 int i;
167 rvec shift;
169 rvec_sub(center, geom_cent, shift);
171 printf(" shift :%7.3f%7.3f%7.3f (nm)\n", shift[XX], shift[YY],
172 shift[ZZ]);
174 for (i = 0; (i < natom); i++)
176 rvec_inc(x[i], shift);
180 void scale_conf(int natom, rvec x[], matrix box, rvec scale)
182 int i, j;
184 for (i = 0; i < natom; i++)
186 for (j = 0; j < DIM; j++)
188 x[i][j] *= scale[j];
191 for (i = 0; i < DIM; i++)
193 for (j = 0; j < DIM; j++)
195 box[i][j] *= scale[j];
200 void read_bfac(const char *fn, int *n_bfac, double **bfac_val, int **bfac_nr)
202 int i;
203 char **bfac_lines;
205 *n_bfac = get_lines(fn, &bfac_lines);
206 snew(*bfac_val, *n_bfac);
207 snew(*bfac_nr, *n_bfac);
208 fprintf(stderr, "Reading %d B-factors from %s\n", *n_bfac, fn);
209 for (i = 0; (i < *n_bfac); i++)
211 /*fprintf(stderr, "Line %d: %s",i,bfac_lines[i]);*/
212 sscanf(bfac_lines[i], "%d %lf", &(*bfac_nr)[i], &(*bfac_val)[i]);
213 /*fprintf(stderr," nr %d val %g\n",(*bfac_nr)[i],(*bfac_val)[i]);*/
218 void set_pdb_conf_bfac(int natoms, int nres, t_atoms *atoms, int n_bfac,
219 double *bfac, int *bfac_nr, gmx_bool peratom)
221 real bfac_min, bfac_max;
222 int i, n;
223 gmx_bool found;
225 bfac_max = -1e10;
226 bfac_min = 1e10;
227 for (i = 0; (i < n_bfac); i++)
229 if (bfac_nr[i] - 1 >= atoms->nres)
231 peratom = TRUE;
233 /* if ((bfac_nr[i]-1<0) || (bfac_nr[i]-1>=atoms->nr))
234 gmx_fatal(FARGS,"Index of B-Factor %d is out of range: %d (%g)",
235 i+1,bfac_nr[i],bfac[i]); */
236 if (bfac[i] > bfac_max)
238 bfac_max = bfac[i];
240 if (bfac[i] < bfac_min)
242 bfac_min = bfac[i];
245 while ((bfac_max > 99.99) || (bfac_min < -99.99))
247 fprintf(stderr,
248 "Range of values for B-factors too large (min %g, max %g) "
249 "will scale down a factor 10\n", bfac_min, bfac_max);
250 for (i = 0; (i < n_bfac); i++)
252 bfac[i] /= 10;
254 bfac_max /= 10;
255 bfac_min /= 10;
257 while ((std::abs(bfac_max) < 0.5) && (std::abs(bfac_min) < 0.5))
259 fprintf(stderr,
260 "Range of values for B-factors too small (min %g, max %g) "
261 "will scale up a factor 10\n", bfac_min, bfac_max);
262 for (i = 0; (i < n_bfac); i++)
264 bfac[i] *= 10;
266 bfac_max *= 10;
267 bfac_min *= 10;
270 for (i = 0; (i < natoms); i++)
272 atoms->pdbinfo[i].bfac = 0;
275 if (!peratom)
277 fprintf(stderr, "Will attach %d B-factors to %d residues\n", n_bfac,
278 nres);
279 for (i = 0; (i < n_bfac); i++)
281 found = FALSE;
282 for (n = 0; (n < natoms); n++)
284 if (bfac_nr[i] == atoms->resinfo[atoms->atom[n].resind].nr)
286 atoms->pdbinfo[n].bfac = bfac[i];
287 found = TRUE;
290 if (!found)
292 gmx_warning("Residue nr %d not found\n", bfac_nr[i]);
296 else
298 fprintf(stderr, "Will attach %d B-factors to %d atoms\n", n_bfac,
299 natoms);
300 for (i = 0; (i < n_bfac); i++)
302 atoms->pdbinfo[bfac_nr[i] - 1].bfac = bfac[i];
307 void pdb_legend(FILE *out, int natoms, int nres, t_atoms *atoms, rvec x[])
309 real bfac_min, bfac_max, xmin, ymin, zmin;
310 int i;
311 int space = ' ';
313 bfac_max = -1e10;
314 bfac_min = 1e10;
315 xmin = 1e10;
316 ymin = 1e10;
317 zmin = 1e10;
318 for (i = 0; (i < natoms); i++)
320 xmin = std::min(xmin, x[i][XX]);
321 ymin = std::min(ymin, x[i][YY]);
322 zmin = std::min(zmin, x[i][ZZ]);
323 bfac_min = std::min(bfac_min, atoms->pdbinfo[i].bfac);
324 bfac_max = std::max(bfac_max, atoms->pdbinfo[i].bfac);
326 fprintf(stderr, "B-factors range from %g to %g\n", bfac_min, bfac_max);
327 for (i = 1; (i < 12); i++)
329 fprintf(out,
330 "%-6s%5d %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f\n",
331 "ATOM ", natoms + 1 + i, "CA", "LEG", space, nres + 1, space,
332 (xmin + (i * 0.12)) * 10, ymin * 10, zmin * 10, 1.0, bfac_min
333 + ((i - 1.0) * (bfac_max - bfac_min) / 10));
337 void visualize_images(const char *fn, int ePBC, matrix box)
339 t_atoms atoms;
340 rvec *img;
341 char *c, *ala;
342 int nat, i;
344 nat = NTRICIMG + 1;
345 init_t_atoms(&atoms, nat, FALSE);
346 atoms.nr = nat;
347 snew(img, nat);
348 /* FIXME: Constness should not be cast away */
349 c = (char *) "C";
350 ala = (char *) "ALA";
351 for (i = 0; i < nat; i++)
353 atoms.atomname[i] = &c;
354 atoms.atom[i].resind = i;
355 atoms.resinfo[i].name = &ala;
356 atoms.resinfo[i].nr = i + 1;
357 atoms.resinfo[i].chainid = 'A' + i / NCUCVERT;
359 calc_triclinic_images(box, img + 1);
361 write_sto_conf(fn, "Images", &atoms, img, NULL, ePBC, box);
363 done_atom(&atoms);
364 sfree(img);
367 void visualize_box(FILE *out, int a0, int r0, matrix box, rvec gridsize)
369 int *edge;
370 rvec *vert, shift;
371 int nx, ny, nz, nbox, nat;
372 int i, j, x, y, z;
373 int rectedge[24] =
375 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6,
379 a0++;
380 r0++;
382 nx = static_cast<int>(gridsize[XX] + 0.5);
383 ny = static_cast<int>(gridsize[YY] + 0.5);
384 nz = static_cast<int>(gridsize[ZZ] + 0.5);
385 nbox = nx * ny * nz;
386 if (TRICLINIC(box))
388 nat = nbox * NCUCVERT;
389 snew(vert, nat);
390 calc_compact_unitcell_vertices(ecenterDEF, box, vert);
391 j = 0;
392 for (z = 0; z < nz; z++)
394 for (y = 0; y < ny; y++)
396 for (x = 0; x < nx; x++)
398 for (i = 0; i < DIM; i++)
400 shift[i] = x * box[0][i] + y * box[1][i] + z
401 * box[2][i];
403 for (i = 0; i < NCUCVERT; i++)
405 rvec_add(vert[i], shift, vert[j]);
406 j++;
412 for (i = 0; i < nat; i++)
414 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i / NCUCVERT, r0 + i, ' ',
415 10*vert[i][XX], 10*vert[i][YY], 10*vert[i][ZZ], 1.0, 0.0, "");
418 edge = compact_unitcell_edges();
419 for (j = 0; j < nbox; j++)
421 for (i = 0; i < NCUCEDGE; i++)
423 fprintf(out, "CONECT%5d%5d\n", a0 + j * NCUCVERT + edge[2 * i],
424 a0 + j * NCUCVERT + edge[2 * i + 1]);
428 sfree(vert);
430 else
432 i = 0;
433 for (z = 0; z <= 1; z++)
435 for (y = 0; y <= 1; y++)
437 for (x = 0; x <= 1; x++)
439 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i/8, r0+i, ' ',
440 x * 10 * box[XX][XX], y * 10 * box[YY][YY], z * 10 * box[ZZ][ZZ], 1.0, 0.0, "");
441 i++;
445 for (i = 0; i < 24; i += 2)
447 fprintf(out, "CONECT%5d%5d\n", a0 + rectedge[i], a0 + rectedge[i + 1]);
452 void calc_rotmatrix(rvec principal_axis, rvec targetvec, matrix rotmatrix)
454 rvec rotvec;
455 real ux, uy, uz, costheta, sintheta;
457 costheta = cos_angle(principal_axis, targetvec);
458 sintheta = std::sqrt(1.0-costheta*costheta); /* sign is always positive since 0<theta<pi */
460 /* Determine rotation from cross product with target vector */
461 cprod(principal_axis, targetvec, rotvec);
462 unitv(rotvec, rotvec);
463 printf("Aligning %g %g %g to %g %g %g : xprod %g %g %g\n",
464 principal_axis[XX], principal_axis[YY], principal_axis[ZZ], targetvec[XX], targetvec[YY], targetvec[ZZ],
465 rotvec[XX], rotvec[YY], rotvec[ZZ]);
467 ux = rotvec[XX];
468 uy = rotvec[YY];
469 uz = rotvec[ZZ];
470 rotmatrix[0][0] = ux*ux + (1.0-ux*ux)*costheta;
471 rotmatrix[0][1] = ux*uy*(1-costheta)-uz*sintheta;
472 rotmatrix[0][2] = ux*uz*(1-costheta)+uy*sintheta;
473 rotmatrix[1][0] = ux*uy*(1-costheta)+uz*sintheta;
474 rotmatrix[1][1] = uy*uy + (1.0-uy*uy)*costheta;
475 rotmatrix[1][2] = uy*uz*(1-costheta)-ux*sintheta;
476 rotmatrix[2][0] = ux*uz*(1-costheta)-uy*sintheta;
477 rotmatrix[2][1] = uy*uz*(1-costheta)+ux*sintheta;
478 rotmatrix[2][2] = uz*uz + (1.0-uz*uz)*costheta;
480 printf("Rotation matrix: \n%g %g %g\n%g %g %g\n%g %g %g\n",
481 rotmatrix[0][0], rotmatrix[0][1], rotmatrix[0][2],
482 rotmatrix[1][0], rotmatrix[1][1], rotmatrix[1][2],
483 rotmatrix[2][0], rotmatrix[2][1], rotmatrix[2][2]);
486 static void renum_resnr(t_atoms *atoms, int isize, const int *index,
487 int resnr_start)
489 int i, resind_prev, resind;
491 resind_prev = -1;
492 for (i = 0; i < isize; i++)
494 resind = atoms->atom[index == NULL ? i : index[i]].resind;
495 if (resind != resind_prev)
497 atoms->resinfo[resind].nr = resnr_start;
498 resnr_start++;
500 resind_prev = resind;
504 int gmx_editconf(int argc, char *argv[])
506 const char *desc[] =
508 "[THISMODULE] converts generic structure format to [REF].gro[ref], [TT].g96[tt]",
509 "or [REF].pdb[ref].",
510 "[PAR]",
511 "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and",
512 "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]",
513 "will center the system in the box, unless [TT]-noc[tt] is used.",
514 "[PAR]",
515 "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
516 "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
517 "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
518 "[TT]octahedron[tt] is a truncated octahedron.",
519 "The last two are special cases of a triclinic box.",
520 "The length of the three box vectors of the truncated octahedron is the",
521 "shortest distance between two opposite hexagons.",
522 "Relative to a cubic box with some periodic image distance, the volume of a ",
523 "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
524 "and that of a truncated octahedron is 0.77 times.",
525 "[PAR]",
526 "Option [TT]-box[tt] requires only",
527 "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
528 "[PAR]",
529 "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, [IT]y[it]-,",
530 "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
531 "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
532 "to the diameter of the system (largest distance between atoms) plus twice",
533 "the specified distance.",
534 "[PAR]",
535 "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
536 "a triclinic box and cannot be used with option [TT]-d[tt].",
537 "[PAR]",
538 "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
539 "can be selected for calculating the size and the geometric center,",
540 "otherwise the whole system is used.",
541 "[PAR]",
542 "[TT]-rotate[tt] rotates the coordinates and velocities.",
543 "[PAR]",
544 "[TT]-princ[tt] aligns the principal axes of the system along the",
545 "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
546 "This may allow you to decrease the box volume,",
547 "but beware that molecules can rotate significantly in a nanosecond.",
548 "[PAR]",
549 "Scaling is applied before any of the other operations are",
550 "performed. Boxes and coordinates can be scaled to give a certain density (option",
551 "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
552 "file is given as input. A special feature of the scaling option is that when the",
553 "factor -1 is given in one dimension, one obtains a mirror image,",
554 "mirrored in one of the planes. When one uses -1 in three dimensions, ",
555 "a point-mirror image is obtained.[PAR]",
556 "Groups are selected after all operations have been applied.[PAR]",
557 "Periodicity can be removed in a crude manner.",
558 "It is important that the box vectors at the bottom of your input file",
559 "are correct when the periodicity is to be removed.",
560 "[PAR]",
561 "When writing [REF].pdb[ref] files, B-factors can be",
562 "added with the [TT]-bf[tt] option. B-factors are read",
563 "from a file with with following format: first line states number of",
564 "entries in the file, next lines state an index",
565 "followed by a B-factor. The B-factors will be attached per residue",
566 "unless an index is larger than the number of residues or unless the",
567 "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
568 "be added instead of B-factors. [TT]-legend[tt] will produce",
569 "a row of CA atoms with B-factors ranging from the minimum to the",
570 "maximum value found, effectively making a legend for viewing.",
571 "[PAR]",
572 "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
573 "file for the MEAD electrostatics",
574 "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
575 "is that the input file is a run input file.",
576 "The B-factor field is then filled with the Van der Waals radius",
577 "of the atoms while the occupancy field will hold the charge.",
578 "[PAR]",
579 "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
580 "and the radius in the occupancy.",
581 "[PAR]",
582 "Option [TT]-align[tt] allows alignment",
583 "of the principal axis of a specified group against the given vector, ",
584 "with an optional center of rotation specified by [TT]-aligncenter[tt].",
585 "[PAR]",
586 "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
587 "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
588 "[PAR]",
589 "To convert a truncated octrahedron file produced by a package which uses",
590 "a cubic box with the corners cut off (such as GROMOS), use::",
592 " gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
594 "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
596 const char *bugs[] =
598 "For complex molecules, the periodicity removal routine may break down, "
599 "in that case you can use [gmx-trjconv]."
601 static real dist = 0.0;
602 static gmx_bool bNDEF = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW =
603 FALSE, bCONECT = FALSE;
604 static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead =
605 FALSE, bGrasp = FALSE, bSig56 = FALSE;
606 static rvec scale =
607 { 1, 1, 1 }, newbox =
608 { 0, 0, 0 }, newang =
609 { 90, 90, 90 };
610 static real rho = 1000.0, rvdw = 0.12;
611 static rvec center =
612 { 0, 0, 0 }, translation =
613 { 0, 0, 0 }, rotangles =
614 { 0, 0, 0 }, aligncenter =
615 { 0, 0, 0 }, targetvec =
616 { 0, 0, 0 };
617 static const char *btype[] =
618 { NULL, "triclinic", "cubic", "dodecahedron", "octahedron", NULL },
619 *label = "A";
620 static rvec visbox =
621 { 0, 0, 0 };
622 static int resnr_start = -1;
623 t_pargs
624 pa[] =
626 { "-ndef", FALSE, etBOOL,
627 { &bNDEF }, "Choose output from default index groups" },
628 { "-visbox", FALSE, etRVEC,
629 { visbox },
630 "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
631 { "-bt", FALSE, etENUM,
632 { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
633 { "-box", FALSE, etRVEC,
634 { newbox }, "Box vector lengths (a,b,c)" },
635 { "-angles", FALSE, etRVEC,
636 { newang }, "Angles between the box vectors (bc,ac,ab)" },
637 { "-d", FALSE, etREAL,
638 { &dist }, "Distance between the solute and the box" },
639 { "-c", FALSE, etBOOL,
640 { &bCenter },
641 "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
642 { "-center", FALSE, etRVEC,
643 { center }, "Coordinates of geometrical center" },
644 { "-aligncenter", FALSE, etRVEC,
645 { aligncenter }, "Center of rotation for alignment" },
646 { "-align", FALSE, etRVEC,
647 { targetvec },
648 "Align to target vector" },
649 { "-translate", FALSE, etRVEC,
650 { translation }, "Translation" },
651 { "-rotate", FALSE, etRVEC,
652 { rotangles },
653 "Rotation around the X, Y and Z axes in degrees" },
654 { "-princ", FALSE, etBOOL,
655 { &bOrient },
656 "Orient molecule(s) along their principal axes" },
657 { "-scale", FALSE, etRVEC,
658 { scale }, "Scaling factor" },
659 { "-density", FALSE, etREAL,
660 { &rho },
661 "Density (g/L) of the output box achieved by scaling" },
662 { "-pbc", FALSE, etBOOL,
663 { &bRMPBC },
664 "Remove the periodicity (make molecule whole again)" },
665 { "-resnr", FALSE, etINT,
666 { &resnr_start },
667 " Renumber residues starting from resnr" },
668 { "-grasp", FALSE, etBOOL,
669 { &bGrasp },
670 "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" },
672 "-rvdw", FALSE, etREAL,
673 { &rvdw },
674 "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"
676 { "-sig56", FALSE, etBOOL,
677 { &bSig56 },
678 "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
680 "-vdwread", FALSE, etBOOL,
681 { &bReadVDW },
682 "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing the radii based on the force field"
684 { "-atom", FALSE, etBOOL,
685 { &peratom }, "Force B-factor attachment per atom" },
686 { "-legend", FALSE, etBOOL,
687 { &bLegend }, "Make B-factor legend" },
688 { "-label", FALSE, etSTR,
689 { &label }, "Add chain label for all residues" },
691 "-conect", FALSE, etBOOL,
692 { &bCONECT },
693 "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a topology is present"
696 #define NPA asize(pa)
698 FILE *out;
699 const char *infile, *outfile;
700 int outftp, inftp, natom, i, j, n_bfac, itype, ntype;
701 double *bfac = NULL, c6, c12;
702 int *bfac_nr = NULL;
703 t_topology *top = NULL;
704 char *grpname, *sgrpname, *agrpname;
705 int isize, ssize, asize;
706 int *index, *sindex, *aindex;
707 rvec *x, *v, gc, rmin, rmax, size;
708 int ePBC;
709 matrix box, rotmatrix, trans;
710 rvec princd, tmpvec;
711 gmx_bool bIndex, bSetSize, bSetAng, bDist, bSetCenter, bAlign;
712 gmx_bool bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
713 real diam = 0, mass = 0, d, vdw;
714 gmx_atomprop_t aps;
715 gmx_conect conect;
716 gmx_output_env_t *oenv;
717 t_filenm fnm[] =
719 { efSTX, "-f", NULL, ffREAD },
720 { efNDX, "-n", NULL, ffOPTRD },
721 { efSTO, NULL, NULL, ffOPTWR },
722 { efPQR, "-mead", "mead", ffOPTWR },
723 { efDAT, "-bf", "bfact", ffOPTRD }
725 #define NFILE asize(fnm)
727 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa,
728 asize(desc), desc, asize(bugs), bugs, &oenv))
730 return 0;
733 bIndex = opt2bSet("-n", NFILE, fnm) || bNDEF;
734 bMead = opt2bSet("-mead", NFILE, fnm);
735 bSetSize = opt2parg_bSet("-box", NPA, pa);
736 bSetAng = opt2parg_bSet("-angles", NPA, pa);
737 bSetCenter = opt2parg_bSet("-center", NPA, pa);
738 bDist = opt2parg_bSet("-d", NPA, pa);
739 bAlign = opt2parg_bSet("-align", NPA, pa);
740 /* Only automatically turn on centering without -noc */
741 if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
743 bCenter = TRUE;
745 bScale = opt2parg_bSet("-scale", NPA, pa);
746 bRho = opt2parg_bSet("-density", NPA, pa);
747 bTranslate = opt2parg_bSet("-translate", NPA, pa);
748 bRotate = opt2parg_bSet("-rotate", NPA, pa);
749 if (bScale && bRho)
751 fprintf(stderr, "WARNING: setting -density overrides -scale\n");
753 bScale = bScale || bRho;
754 bCalcGeom = bCenter || bRotate || bOrient || bScale;
756 GMX_RELEASE_ASSERT(btype[0] != NULL, "Option setting inconsistency; btype[0] is NULL");
758 bCalcDiam = (btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o');
760 infile = ftp2fn(efSTX, NFILE, fnm);
761 if (bMead)
763 outfile = ftp2fn(efPQR, NFILE, fnm);
765 else
767 outfile = ftp2fn(efSTO, NFILE, fnm);
769 outftp = fn2ftp(outfile);
770 inftp = fn2ftp(infile);
772 aps = gmx_atomprop_init();
774 if (bMead && bGrasp)
776 printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
777 bGrasp = FALSE;
779 if (bGrasp && (outftp != efPDB))
781 gmx_fatal(FARGS, "Output file should be a .pdb file"
782 " when using the -grasp option\n");
784 if ((bMead || bGrasp) && (fn2ftp(infile) != efTPR))
786 gmx_fatal(FARGS, "Input file should be a .tpr file"
787 " when using the -mead option\n");
790 t_topology *top_tmp;
791 snew(top_tmp, 1);
792 read_tps_conf(infile, top_tmp, &ePBC, &x, &v, box, FALSE);
793 t_atoms &atoms = top_tmp->atoms;
794 natom = atoms.nr;
795 if (atoms.pdbinfo == NULL)
797 snew(atoms.pdbinfo, atoms.nr);
799 if (fn2ftp(infile) == efPDB)
801 get_pdb_atomnumber(&atoms, aps);
803 printf("Read %d atoms\n", atoms.nr);
805 /* Get the element numbers if available in a pdb file */
806 if (fn2ftp(infile) == efPDB)
808 get_pdb_atomnumber(&atoms, aps);
811 if (ePBC != epbcNONE)
813 real vol = det(box);
814 printf("Volume: %g nm^3, corresponds to roughly %d electrons\n",
815 vol, 100*(static_cast<int>(vol*4.5)));
818 if (bMead || bGrasp || bCONECT)
820 top = read_top(infile, NULL);
823 if (bMead || bGrasp)
825 if (atoms.nr != top->atoms.nr)
827 gmx_fatal(FARGS, "Atom numbers don't match (%d vs. %d)", atoms.nr, top->atoms.nr);
829 snew(atoms.pdbinfo, top->atoms.nr);
830 ntype = top->idef.atnr;
831 for (i = 0; (i < atoms.nr); i++)
833 /* Determine the Van der Waals radius from the force field */
834 if (bReadVDW)
836 if (!gmx_atomprop_query(aps, epropVDW,
837 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
838 *top->atoms.atomname[i], &vdw))
840 vdw = rvdw;
843 else
845 itype = top->atoms.atom[i].type;
846 c12 = top->idef.iparams[itype*ntype+itype].lj.c12;
847 c6 = top->idef.iparams[itype*ntype+itype].lj.c6;
848 if ((c6 != 0) && (c12 != 0))
850 real sig6;
851 if (bSig56)
853 sig6 = 2*c12/c6;
855 else
857 sig6 = c12/c6;
859 vdw = 0.5*std::pow(sig6, static_cast<real>(1.0/6.0));
861 else
863 vdw = rvdw;
866 /* Factor of 10 for nm -> Angstroms */
867 vdw *= 10;
869 if (bMead)
871 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
872 atoms.pdbinfo[i].bfac = vdw;
874 else
876 atoms.pdbinfo[i].occup = vdw;
877 atoms.pdbinfo[i].bfac = top->atoms.atom[i].q;
881 bHaveV = FALSE;
882 for (i = 0; (i < natom) && !bHaveV; i++)
884 for (j = 0; (j < DIM) && !bHaveV; j++)
886 bHaveV = bHaveV || (v[i][j] != 0);
889 printf("%selocities found\n", bHaveV ? "V" : "No v");
891 if (visbox[0] > 0)
893 if (bIndex)
895 gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
897 if (outftp != efPDB)
899 gmx_fatal(FARGS, "Sorry, can only visualize box with a pdb file");
902 else if (visbox[0] == -1)
904 visualize_images("images.pdb", ePBC, box);
907 /* remove pbc */
908 if (bRMPBC)
910 rm_gropbc(&atoms, x, box);
913 if (bCalcGeom)
915 if (bIndex)
917 fprintf(stderr, "\nSelect a group for determining the system size:\n");
918 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
919 1, &ssize, &sindex, &sgrpname);
921 else
923 ssize = atoms.nr;
924 sindex = NULL;
926 diam = calc_geom(ssize, sindex, x, gc, rmin, rmax, bCalcDiam);
927 rvec_sub(rmax, rmin, size);
928 printf(" system size :%7.3f%7.3f%7.3f (nm)\n",
929 size[XX], size[YY], size[ZZ]);
930 if (bCalcDiam)
932 printf(" diameter :%7.3f (nm)\n", diam);
934 printf(" center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
935 printf(" box vectors :%7.3f%7.3f%7.3f (nm)\n",
936 norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
937 printf(" box angles :%7.2f%7.2f%7.2f (degrees)\n",
938 norm2(box[ZZ]) == 0 ? 0 :
939 RAD2DEG*std::acos(cos_angle_no_table(box[YY], box[ZZ])),
940 norm2(box[ZZ]) == 0 ? 0 :
941 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[ZZ])),
942 norm2(box[YY]) == 0 ? 0 :
943 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[YY])));
944 printf(" box volume :%7.2f (nm^3)\n", det(box));
947 if (bRho || bOrient || bAlign)
949 mass = calc_mass(&atoms, !fn2bTPX(infile), aps);
952 if (bOrient)
954 int *index;
955 char *grpnames;
957 /* Get a group for principal component analysis */
958 fprintf(stderr, "\nSelect group for the determining the orientation\n");
959 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
961 /* Orient the principal axes along the coordinate axes */
962 orient_princ(&atoms, isize, index, natom, x, bHaveV ? v : NULL, NULL);
963 sfree(index);
964 sfree(grpnames);
967 if (bScale)
969 /* scale coordinates and box */
970 if (bRho)
972 /* Compute scaling constant */
973 real vol, dens;
975 vol = det(box);
976 dens = (mass*AMU)/(vol*NANO*NANO*NANO);
977 fprintf(stderr, "Volume of input %g (nm^3)\n", vol);
978 fprintf(stderr, "Mass of input %g (a.m.u.)\n", mass);
979 fprintf(stderr, "Density of input %g (g/l)\n", dens);
980 if (vol == 0 || mass == 0)
982 gmx_fatal(FARGS, "Cannot scale density with "
983 "zero mass (%g) or volume (%g)\n", mass, vol);
986 scale[XX] = scale[YY] = scale[ZZ] = std::pow(dens/rho, static_cast<real>(1.0/3.0));
987 fprintf(stderr, "Scaling all box vectors by %g\n", scale[XX]);
989 scale_conf(atoms.nr, x, box, scale);
992 if (bAlign)
994 if (bIndex)
996 fprintf(stderr, "\nSelect a group that you want to align:\n");
997 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
998 1, &asize, &aindex, &agrpname);
1000 else
1002 asize = atoms.nr;
1003 snew(aindex, asize);
1004 for (i = 0; i < asize; i++)
1006 aindex[i] = i;
1009 printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n", asize, natom,
1010 targetvec[XX], targetvec[YY], targetvec[ZZ],
1011 aligncenter[XX], aligncenter[YY], aligncenter[ZZ]);
1012 /*subtract out pivot point*/
1013 for (i = 0; i < asize; i++)
1015 rvec_dec(x[aindex[i]], aligncenter);
1017 /*now determine transform and rotate*/
1018 /*will this work?*/
1019 principal_comp(asize, aindex, atoms.atom, x, trans, princd);
1021 unitv(targetvec, targetvec);
1022 printf("Using %g %g %g as principal axis\n", trans[0][2], trans[1][2], trans[2][2]);
1023 tmpvec[XX] = trans[0][2]; tmpvec[YY] = trans[1][2]; tmpvec[ZZ] = trans[2][2];
1024 calc_rotmatrix(tmpvec, targetvec, rotmatrix);
1025 /* rotmatrix finished */
1027 for (i = 0; i < asize; ++i)
1029 mvmul(rotmatrix, x[aindex[i]], tmpvec);
1030 copy_rvec(tmpvec, x[aindex[i]]);
1033 /*add pivot point back*/
1034 for (i = 0; i < asize; i++)
1036 rvec_inc(x[aindex[i]], aligncenter);
1038 if (!bIndex)
1040 sfree(aindex);
1044 if (bTranslate)
1046 if (bIndex)
1048 fprintf(stderr, "\nSelect a group that you want to translate:\n");
1049 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
1050 1, &ssize, &sindex, &sgrpname);
1052 else
1054 ssize = atoms.nr;
1055 sindex = NULL;
1057 printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize, natom,
1058 translation[XX], translation[YY], translation[ZZ]);
1059 if (sindex)
1061 for (i = 0; i < ssize; i++)
1063 rvec_inc(x[sindex[i]], translation);
1066 else
1068 for (i = 0; i < natom; i++)
1070 rvec_inc(x[i], translation);
1074 if (bRotate)
1076 /* Rotate */
1077 printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n", rotangles[XX], rotangles[YY], rotangles[ZZ]);
1078 for (i = 0; i < DIM; i++)
1080 rotangles[i] *= DEG2RAD;
1082 rotate_conf(natom, x, v, rotangles[XX], rotangles[YY], rotangles[ZZ]);
1085 if (bCalcGeom)
1087 /* recalc geometrical center and max and min coordinates and size */
1088 calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1089 rvec_sub(rmax, rmin, size);
1090 if (bScale || bOrient || bRotate)
1092 printf("new system size : %6.3f %6.3f %6.3f\n",
1093 size[XX], size[YY], size[ZZ]);
1097 if ((btype[0] != NULL) && (bSetSize || bDist || (btype[0][0] == 't' && bSetAng)))
1099 ePBC = epbcXYZ;
1100 if (!(bSetSize || bDist))
1102 for (i = 0; i < DIM; i++)
1104 newbox[i] = norm(box[i]);
1107 clear_mat(box);
1108 /* calculate new boxsize */
1109 switch (btype[0][0])
1111 case 't':
1112 if (bDist)
1114 for (i = 0; i < DIM; i++)
1116 newbox[i] = size[i]+2*dist;
1119 if (!bSetAng)
1121 box[XX][XX] = newbox[XX];
1122 box[YY][YY] = newbox[YY];
1123 box[ZZ][ZZ] = newbox[ZZ];
1125 else
1127 matrix_convert(box, newbox, newang);
1129 break;
1130 case 'c':
1131 case 'd':
1132 case 'o':
1133 if (bSetSize)
1135 d = newbox[0];
1137 else
1139 d = diam+2*dist;
1141 if (btype[0][0] == 'c')
1143 for (i = 0; i < DIM; i++)
1145 box[i][i] = d;
1148 else if (btype[0][0] == 'd')
1150 box[XX][XX] = d;
1151 box[YY][YY] = d;
1152 box[ZZ][XX] = d/2;
1153 box[ZZ][YY] = d/2;
1154 box[ZZ][ZZ] = d*std::sqrt(2.0)/2.0;
1156 else
1158 box[XX][XX] = d;
1159 box[YY][XX] = d/3;
1160 box[YY][YY] = d*std::sqrt(2.0)*2.0/3.0;
1161 box[ZZ][XX] = -d/3;
1162 box[ZZ][YY] = d*std::sqrt(2.0)/3.0;
1163 box[ZZ][ZZ] = d*std::sqrt(6.0)/3.0;
1165 break;
1169 /* calculate new coords for geometrical center */
1170 if (!bSetCenter)
1172 calc_box_center(ecenterDEF, box, center);
1175 /* center molecule on 'center' */
1176 if (bCenter)
1178 center_conf(natom, x, center, gc);
1181 /* print some */
1182 if (bCalcGeom)
1184 calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1185 printf("new center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
1187 if (bOrient || bScale || bDist || bSetSize)
1189 printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n",
1190 norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
1191 printf("new box angles :%7.2f%7.2f%7.2f (degrees)\n",
1192 norm2(box[ZZ]) == 0 ? 0 :
1193 RAD2DEG*std::acos(cos_angle_no_table(box[YY], box[ZZ])),
1194 norm2(box[ZZ]) == 0 ? 0 :
1195 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[ZZ])),
1196 norm2(box[YY]) == 0 ? 0 :
1197 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[YY])));
1198 printf("new box volume :%7.2f (nm^3)\n", det(box));
1201 if (check_box(epbcXYZ, box))
1203 printf("\nWARNING: %s\n"
1204 "See the GROMACS manual for a description of the requirements that\n"
1205 "must be satisfied by descriptions of simulation cells.\n",
1206 check_box(epbcXYZ, box));
1209 if (bDist && btype[0][0] == 't')
1211 if (TRICLINIC(box))
1213 printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1214 "distance from the solute to a box surface along the corresponding normal\n"
1215 "vector might be somewhat smaller than your specified value %f.\n"
1216 "You can check the actual value with g_mindist -pi\n", dist);
1218 else if (!opt2parg_bSet("-bt", NPA, pa))
1220 printf("\nWARNING: No boxtype specified - distance condition applied in each dimension.\n"
1221 "If the molecule rotates the actual distance will be smaller. You might want\n"
1222 "to use a cubic box instead, or why not try a dodecahedron today?\n");
1225 if (bCONECT && (outftp == efPDB) && (inftp == efTPR))
1227 conect = gmx_conect_generate(top);
1229 else
1231 conect = NULL;
1234 if (bIndex)
1236 fprintf(stderr, "\nSelect a group for output:\n");
1237 get_index(&atoms, opt2fn_null("-n", NFILE, fnm),
1238 1, &isize, &index, &grpname);
1240 if (resnr_start >= 0)
1242 renum_resnr(&atoms, isize, index, resnr_start);
1245 if (opt2parg_bSet("-label", NPA, pa))
1247 for (i = 0; (i < atoms.nr); i++)
1249 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1253 if (opt2bSet("-bf", NFILE, fnm) || bLegend)
1255 gmx_fatal(FARGS, "Sorry, cannot do bfactors with an index group.");
1258 if (outftp == efPDB)
1260 out = gmx_ffopen(outfile, "w");
1261 write_pdbfile_indexed(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, TRUE);
1262 gmx_ffclose(out);
1264 else
1266 write_sto_conf_indexed(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : NULL, ePBC, box, isize, index);
1269 else
1271 if (resnr_start >= 0)
1273 renum_resnr(&atoms, atoms.nr, NULL, resnr_start);
1276 if ((outftp == efPDB) || (outftp == efPQR))
1278 out = gmx_ffopen(outfile, "w");
1279 if (bMead)
1281 fprintf(out, "REMARK "
1282 "The B-factors in this file hold atomic radii\n");
1283 fprintf(out, "REMARK "
1284 "The occupancy in this file hold atomic charges\n");
1286 else if (bGrasp)
1288 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1289 fprintf(out, "REMARK "
1290 "The B-factors in this file hold atomic charges\n");
1291 fprintf(out, "REMARK "
1292 "The occupancy in this file hold atomic radii\n");
1294 else if (opt2bSet("-bf", NFILE, fnm))
1296 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1297 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms,
1298 n_bfac, bfac, bfac_nr, peratom);
1300 if (opt2parg_bSet("-label", NPA, pa))
1302 for (i = 0; (i < atoms.nr); i++)
1304 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1307 write_pdbfile(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', -1, conect, TRUE);
1308 if (bLegend)
1310 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1312 if (visbox[0] > 0)
1314 visualize_box(out, bLegend ? atoms.nr+12 : atoms.nr,
1315 bLegend ? atoms.nres = 12 : atoms.nres, box, visbox);
1317 gmx_ffclose(out);
1319 else
1321 write_sto_conf(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : NULL, ePBC, box);
1324 gmx_atomprop_destroy(aps);
1326 do_view(oenv, outfile, NULL);
1328 return 0;