Removed legacyheaders/typedefs.h
[gromacs.git] / src / gromacs / gmxana / gmx_editconf.cpp
blob727b1daae73771ed6cb943802363cd3b734ce58b
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/strdb.h"
49 #include "gromacs/fileio/tpxio.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/gmxana/princ.h"
53 #include "gromacs/gmxlib/conformation-utilities.h"
54 #include "gromacs/legacyheaders/txtdump.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/pbcutil/rmpbc.h"
59 #include "gromacs/topology/atomprop.h"
60 #include "gromacs/topology/index.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/utility/arraysize.h"
63 #include "gromacs/utility/cstringutil.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/futil.h"
66 #include "gromacs/utility/gmxassert.h"
67 #include "gromacs/utility/smalloc.h"
70 real calc_mass(t_atoms *atoms, gmx_bool bGetMass, gmx_atomprop_t aps)
72 real tmass;
73 int i;
75 tmass = 0;
76 for (i = 0; (i < atoms->nr); i++)
78 if (bGetMass)
80 gmx_atomprop_query(aps, epropMass,
81 *atoms->resinfo[atoms->atom[i].resind].name,
82 *atoms->atomname[i], &(atoms->atom[i].m));
84 tmass += atoms->atom[i].m;
87 return tmass;
90 real calc_geom(int isize, atom_id *index, rvec *x, rvec geom_center, rvec minval,
91 rvec maxval, gmx_bool bDiam)
93 real diam2, d;
94 int ii, i, j;
96 clear_rvec(geom_center);
97 diam2 = 0;
98 if (isize == 0)
100 clear_rvec(minval);
101 clear_rvec(maxval);
103 else
105 if (index)
107 ii = index[0];
109 else
111 ii = 0;
113 for (j = 0; j < DIM; j++)
115 minval[j] = maxval[j] = x[ii][j];
117 for (i = 0; i < isize; i++)
119 if (index)
121 ii = index[i];
123 else
125 ii = i;
127 rvec_inc(geom_center, x[ii]);
128 for (j = 0; j < DIM; j++)
130 if (x[ii][j] < minval[j])
132 minval[j] = x[ii][j];
134 if (x[ii][j] > maxval[j])
136 maxval[j] = x[ii][j];
139 if (bDiam)
141 if (index)
143 for (j = i + 1; j < isize; j++)
145 d = distance2(x[ii], x[index[j]]);
146 diam2 = std::max(d, diam2);
149 else
151 for (j = i + 1; j < isize; j++)
153 d = distance2(x[i], x[j]);
154 diam2 = std::max(d, diam2);
159 svmul(1.0 / isize, geom_center, geom_center);
162 return std::sqrt(diam2);
165 void center_conf(int natom, rvec *x, rvec center, rvec geom_cent)
167 int i;
168 rvec shift;
170 rvec_sub(center, geom_cent, shift);
172 printf(" shift :%7.3f%7.3f%7.3f (nm)\n", shift[XX], shift[YY],
173 shift[ZZ]);
175 for (i = 0; (i < natom); i++)
177 rvec_inc(x[i], shift);
181 void scale_conf(int natom, rvec x[], matrix box, rvec scale)
183 int i, j;
185 for (i = 0; i < natom; i++)
187 for (j = 0; j < DIM; j++)
189 x[i][j] *= scale[j];
192 for (i = 0; i < DIM; i++)
194 for (j = 0; j < DIM; j++)
196 box[i][j] *= scale[j];
201 void read_bfac(const char *fn, int *n_bfac, double **bfac_val, int **bfac_nr)
203 int i;
204 char **bfac_lines;
206 *n_bfac = get_lines(fn, &bfac_lines);
207 snew(*bfac_val, *n_bfac);
208 snew(*bfac_nr, *n_bfac);
209 fprintf(stderr, "Reading %d B-factors from %s\n", *n_bfac, fn);
210 for (i = 0; (i < *n_bfac); i++)
212 /*fprintf(stderr, "Line %d: %s",i,bfac_lines[i]);*/
213 sscanf(bfac_lines[i], "%d %lf", &(*bfac_nr)[i], &(*bfac_val)[i]);
214 /*fprintf(stderr," nr %d val %g\n",(*bfac_nr)[i],(*bfac_val)[i]);*/
219 void set_pdb_conf_bfac(int natoms, int nres, t_atoms *atoms, int n_bfac,
220 double *bfac, int *bfac_nr, gmx_bool peratom)
222 real bfac_min, bfac_max;
223 int i, n;
224 gmx_bool found;
226 bfac_max = -1e10;
227 bfac_min = 1e10;
228 for (i = 0; (i < n_bfac); i++)
230 if (bfac_nr[i] - 1 >= atoms->nres)
232 peratom = TRUE;
234 /* if ((bfac_nr[i]-1<0) || (bfac_nr[i]-1>=atoms->nr))
235 gmx_fatal(FARGS,"Index of B-Factor %d is out of range: %d (%g)",
236 i+1,bfac_nr[i],bfac[i]); */
237 if (bfac[i] > bfac_max)
239 bfac_max = bfac[i];
241 if (bfac[i] < bfac_min)
243 bfac_min = bfac[i];
246 while ((bfac_max > 99.99) || (bfac_min < -99.99))
248 fprintf(stderr,
249 "Range of values for B-factors too large (min %g, max %g) "
250 "will scale down a factor 10\n", bfac_min, bfac_max);
251 for (i = 0; (i < n_bfac); i++)
253 bfac[i] /= 10;
255 bfac_max /= 10;
256 bfac_min /= 10;
258 while ((std::abs(bfac_max) < 0.5) && (std::abs(bfac_min) < 0.5))
260 fprintf(stderr,
261 "Range of values for B-factors too small (min %g, max %g) "
262 "will scale up a factor 10\n", bfac_min, bfac_max);
263 for (i = 0; (i < n_bfac); i++)
265 bfac[i] *= 10;
267 bfac_max *= 10;
268 bfac_min *= 10;
271 for (i = 0; (i < natoms); i++)
273 atoms->pdbinfo[i].bfac = 0;
276 if (!peratom)
278 fprintf(stderr, "Will attach %d B-factors to %d residues\n", n_bfac,
279 nres);
280 for (i = 0; (i < n_bfac); i++)
282 found = FALSE;
283 for (n = 0; (n < natoms); n++)
285 if (bfac_nr[i] == atoms->resinfo[atoms->atom[n].resind].nr)
287 atoms->pdbinfo[n].bfac = bfac[i];
288 found = TRUE;
291 if (!found)
293 gmx_warning("Residue nr %d not found\n", bfac_nr[i]);
297 else
299 fprintf(stderr, "Will attach %d B-factors to %d atoms\n", n_bfac,
300 natoms);
301 for (i = 0; (i < n_bfac); i++)
303 atoms->pdbinfo[bfac_nr[i] - 1].bfac = bfac[i];
308 void pdb_legend(FILE *out, int natoms, int nres, t_atoms *atoms, rvec x[])
310 real bfac_min, bfac_max, xmin, ymin, zmin;
311 int i;
312 int space = ' ';
314 bfac_max = -1e10;
315 bfac_min = 1e10;
316 xmin = 1e10;
317 ymin = 1e10;
318 zmin = 1e10;
319 for (i = 0; (i < natoms); i++)
321 xmin = std::min(xmin, x[i][XX]);
322 ymin = std::min(ymin, x[i][YY]);
323 zmin = std::min(zmin, x[i][ZZ]);
324 bfac_min = std::min(bfac_min, atoms->pdbinfo[i].bfac);
325 bfac_max = std::max(bfac_max, atoms->pdbinfo[i].bfac);
327 fprintf(stderr, "B-factors range from %g to %g\n", bfac_min, bfac_max);
328 for (i = 1; (i < 12); i++)
330 fprintf(out,
331 "%-6s%5d %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f\n",
332 "ATOM ", natoms + 1 + i, "CA", "LEG", space, nres + 1, space,
333 (xmin + (i * 0.12)) * 10, ymin * 10, zmin * 10, 1.0, bfac_min
334 + ((i - 1.0) * (bfac_max - bfac_min) / 10));
338 void visualize_images(const char *fn, int ePBC, matrix box)
340 t_atoms atoms;
341 rvec *img;
342 char *c, *ala;
343 int nat, i;
345 nat = NTRICIMG + 1;
346 init_t_atoms(&atoms, nat, FALSE);
347 atoms.nr = nat;
348 snew(img, nat);
349 /* FIXME: Constness should not be cast away */
350 c = (char *) "C";
351 ala = (char *) "ALA";
352 for (i = 0; i < nat; i++)
354 atoms.atomname[i] = &c;
355 atoms.atom[i].resind = i;
356 atoms.resinfo[i].name = &ala;
357 atoms.resinfo[i].nr = i + 1;
358 atoms.resinfo[i].chainid = 'A' + i / NCUCVERT;
360 calc_triclinic_images(box, img + 1);
362 write_sto_conf(fn, "Images", &atoms, img, NULL, ePBC, box);
364 done_atom(&atoms);
365 sfree(img);
368 void visualize_box(FILE *out, int a0, int r0, matrix box, rvec gridsize)
370 int *edge;
371 rvec *vert, shift;
372 int nx, ny, nz, nbox, nat;
373 int i, j, x, y, z;
374 int rectedge[24] =
376 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6,
380 a0++;
381 r0++;
383 nx = static_cast<int>(gridsize[XX] + 0.5);
384 ny = static_cast<int>(gridsize[YY] + 0.5);
385 nz = static_cast<int>(gridsize[ZZ] + 0.5);
386 nbox = nx * ny * nz;
387 if (TRICLINIC(box))
389 nat = nbox * NCUCVERT;
390 snew(vert, nat);
391 calc_compact_unitcell_vertices(ecenterDEF, box, vert);
392 j = 0;
393 for (z = 0; z < nz; z++)
395 for (y = 0; y < ny; y++)
397 for (x = 0; x < nx; x++)
399 for (i = 0; i < DIM; i++)
401 shift[i] = x * box[0][i] + y * box[1][i] + z
402 * box[2][i];
404 for (i = 0; i < NCUCVERT; i++)
406 rvec_add(vert[i], shift, vert[j]);
407 j++;
413 for (i = 0; i < nat; i++)
415 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i / NCUCVERT, r0 + i, ' ',
416 10*vert[i][XX], 10*vert[i][YY], 10*vert[i][ZZ], 1.0, 0.0, "");
419 edge = compact_unitcell_edges();
420 for (j = 0; j < nbox; j++)
422 for (i = 0; i < NCUCEDGE; i++)
424 fprintf(out, "CONECT%5d%5d\n", a0 + j * NCUCVERT + edge[2 * i],
425 a0 + j * NCUCVERT + edge[2 * i + 1]);
429 sfree(vert);
431 else
433 i = 0;
434 for (z = 0; z <= 1; z++)
436 for (y = 0; y <= 1; y++)
438 for (x = 0; x <= 1; x++)
440 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i/8, r0+i, ' ',
441 x * 10 * box[XX][XX], y * 10 * box[YY][YY], z * 10 * box[ZZ][ZZ], 1.0, 0.0, "");
442 i++;
446 for (i = 0; i < 24; i += 2)
448 fprintf(out, "CONECT%5d%5d\n", a0 + rectedge[i], a0 + rectedge[i + 1]);
453 void calc_rotmatrix(rvec principal_axis, rvec targetvec, matrix rotmatrix)
455 rvec rotvec;
456 real ux, uy, uz, costheta, sintheta;
458 costheta = cos_angle(principal_axis, targetvec);
459 sintheta = std::sqrt(1.0-costheta*costheta); /* sign is always positive since 0<theta<pi */
461 /* Determine rotation from cross product with target vector */
462 cprod(principal_axis, targetvec, rotvec);
463 unitv(rotvec, rotvec);
464 printf("Aligning %g %g %g to %g %g %g : xprod %g %g %g\n",
465 principal_axis[XX], principal_axis[YY], principal_axis[ZZ], targetvec[XX], targetvec[YY], targetvec[ZZ],
466 rotvec[XX], rotvec[YY], rotvec[ZZ]);
468 ux = rotvec[XX];
469 uy = rotvec[YY];
470 uz = rotvec[ZZ];
471 rotmatrix[0][0] = ux*ux + (1.0-ux*ux)*costheta;
472 rotmatrix[0][1] = ux*uy*(1-costheta)-uz*sintheta;
473 rotmatrix[0][2] = ux*uz*(1-costheta)+uy*sintheta;
474 rotmatrix[1][0] = ux*uy*(1-costheta)+uz*sintheta;
475 rotmatrix[1][1] = uy*uy + (1.0-uy*uy)*costheta;
476 rotmatrix[1][2] = uy*uz*(1-costheta)-ux*sintheta;
477 rotmatrix[2][0] = ux*uz*(1-costheta)-uy*sintheta;
478 rotmatrix[2][1] = uy*uz*(1-costheta)+ux*sintheta;
479 rotmatrix[2][2] = uz*uz + (1.0-uz*uz)*costheta;
481 printf("Rotation matrix: \n%g %g %g\n%g %g %g\n%g %g %g\n",
482 rotmatrix[0][0], rotmatrix[0][1], rotmatrix[0][2],
483 rotmatrix[1][0], rotmatrix[1][1], rotmatrix[1][2],
484 rotmatrix[2][0], rotmatrix[2][1], rotmatrix[2][2]);
487 static void renum_resnr(t_atoms *atoms, int isize, const int *index,
488 int resnr_start)
490 int i, resind_prev, resind;
492 resind_prev = -1;
493 for (i = 0; i < isize; i++)
495 resind = atoms->atom[index == NULL ? i : index[i]].resind;
496 if (resind != resind_prev)
498 atoms->resinfo[resind].nr = resnr_start;
499 resnr_start++;
501 resind_prev = resind;
505 int gmx_editconf(int argc, char *argv[])
507 const char *desc[] =
509 "[THISMODULE] converts generic structure format to [REF].gro[ref], [TT].g96[tt]",
510 "or [REF].pdb[ref].",
511 "[PAR]",
512 "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and",
513 "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]",
514 "will center the system in the box, unless [TT]-noc[tt] is used.",
515 "[PAR]",
516 "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
517 "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
518 "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
519 "[TT]octahedron[tt] is a truncated octahedron.",
520 "The last two are special cases of a triclinic box.",
521 "The length of the three box vectors of the truncated octahedron is the",
522 "shortest distance between two opposite hexagons.",
523 "Relative to a cubic box with some periodic image distance, the volume of a ",
524 "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
525 "and that of a truncated octahedron is 0.77 times.",
526 "[PAR]",
527 "Option [TT]-box[tt] requires only",
528 "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
529 "[PAR]",
530 "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, [IT]y[it]-,",
531 "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
532 "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
533 "to the diameter of the system (largest distance between atoms) plus twice",
534 "the specified distance.",
535 "[PAR]",
536 "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
537 "a triclinic box and cannot be used with option [TT]-d[tt].",
538 "[PAR]",
539 "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
540 "can be selected for calculating the size and the geometric center,",
541 "otherwise the whole system is used.",
542 "[PAR]",
543 "[TT]-rotate[tt] rotates the coordinates and velocities.",
544 "[PAR]",
545 "[TT]-princ[tt] aligns the principal axes of the system along the",
546 "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
547 "This may allow you to decrease the box volume,",
548 "but beware that molecules can rotate significantly in a nanosecond.",
549 "[PAR]",
550 "Scaling is applied before any of the other operations are",
551 "performed. Boxes and coordinates can be scaled to give a certain density (option",
552 "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
553 "file is given as input. A special feature of the scaling option is that when the",
554 "factor -1 is given in one dimension, one obtains a mirror image,",
555 "mirrored in one of the planes. When one uses -1 in three dimensions, ",
556 "a point-mirror image is obtained.[PAR]",
557 "Groups are selected after all operations have been applied.[PAR]",
558 "Periodicity can be removed in a crude manner.",
559 "It is important that the box vectors at the bottom of your input file",
560 "are correct when the periodicity is to be removed.",
561 "[PAR]",
562 "When writing [REF].pdb[ref] files, B-factors can be",
563 "added with the [TT]-bf[tt] option. B-factors are read",
564 "from a file with with following format: first line states number of",
565 "entries in the file, next lines state an index",
566 "followed by a B-factor. The B-factors will be attached per residue",
567 "unless an index is larger than the number of residues or unless the",
568 "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
569 "be added instead of B-factors. [TT]-legend[tt] will produce",
570 "a row of CA atoms with B-factors ranging from the minimum to the",
571 "maximum value found, effectively making a legend for viewing.",
572 "[PAR]",
573 "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
574 "file for the MEAD electrostatics",
575 "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
576 "is that the input file is a run input file.",
577 "The B-factor field is then filled with the Van der Waals radius",
578 "of the atoms while the occupancy field will hold the charge.",
579 "[PAR]",
580 "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
581 "and the radius in the occupancy.",
582 "[PAR]",
583 "Option [TT]-align[tt] allows alignment",
584 "of the principal axis of a specified group against the given vector, ",
585 "with an optional center of rotation specified by [TT]-aligncenter[tt].",
586 "[PAR]",
587 "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
588 "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
589 "[PAR]",
590 "To convert a truncated octrahedron file produced by a package which uses",
591 "a cubic box with the corners cut off (such as GROMOS), use::",
593 " gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
595 "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
597 const char *bugs[] =
599 "For complex molecules, the periodicity removal routine may break down, "
600 "in that case you can use [gmx-trjconv]."
602 static real dist = 0.0;
603 static gmx_bool bNDEF = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW =
604 FALSE, bCONECT = FALSE;
605 static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead =
606 FALSE, bGrasp = FALSE, bSig56 = FALSE;
607 static rvec scale =
608 { 1, 1, 1 }, newbox =
609 { 0, 0, 0 }, newang =
610 { 90, 90, 90 };
611 static real rho = 1000.0, rvdw = 0.12;
612 static rvec center =
613 { 0, 0, 0 }, translation =
614 { 0, 0, 0 }, rotangles =
615 { 0, 0, 0 }, aligncenter =
616 { 0, 0, 0 }, targetvec =
617 { 0, 0, 0 };
618 static const char *btype[] =
619 { NULL, "triclinic", "cubic", "dodecahedron", "octahedron", NULL },
620 *label = "A";
621 static rvec visbox =
622 { 0, 0, 0 };
623 static int resnr_start = -1;
624 t_pargs
625 pa[] =
627 { "-ndef", FALSE, etBOOL,
628 { &bNDEF }, "Choose output from default index groups" },
629 { "-visbox", FALSE, etRVEC,
630 { visbox },
631 "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
632 { "-bt", FALSE, etENUM,
633 { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
634 { "-box", FALSE, etRVEC,
635 { newbox }, "Box vector lengths (a,b,c)" },
636 { "-angles", FALSE, etRVEC,
637 { newang }, "Angles between the box vectors (bc,ac,ab)" },
638 { "-d", FALSE, etREAL,
639 { &dist }, "Distance between the solute and the box" },
640 { "-c", FALSE, etBOOL,
641 { &bCenter },
642 "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
643 { "-center", FALSE, etRVEC,
644 { center }, "Coordinates of geometrical center" },
645 { "-aligncenter", FALSE, etRVEC,
646 { aligncenter }, "Center of rotation for alignment" },
647 { "-align", FALSE, etRVEC,
648 { targetvec },
649 "Align to target vector" },
650 { "-translate", FALSE, etRVEC,
651 { translation }, "Translation" },
652 { "-rotate", FALSE, etRVEC,
653 { rotangles },
654 "Rotation around the X, Y and Z axes in degrees" },
655 { "-princ", FALSE, etBOOL,
656 { &bOrient },
657 "Orient molecule(s) along their principal axes" },
658 { "-scale", FALSE, etRVEC,
659 { scale }, "Scaling factor" },
660 { "-density", FALSE, etREAL,
661 { &rho },
662 "Density (g/L) of the output box achieved by scaling" },
663 { "-pbc", FALSE, etBOOL,
664 { &bRMPBC },
665 "Remove the periodicity (make molecule whole again)" },
666 { "-resnr", FALSE, etINT,
667 { &resnr_start },
668 " Renumber residues starting from resnr" },
669 { "-grasp", FALSE, etBOOL,
670 { &bGrasp },
671 "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" },
673 "-rvdw", FALSE, etREAL,
674 { &rvdw },
675 "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"
677 { "-sig56", FALSE, etBOOL,
678 { &bSig56 },
679 "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
681 "-vdwread", FALSE, etBOOL,
682 { &bReadVDW },
683 "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing the radii based on the force field"
685 { "-atom", FALSE, etBOOL,
686 { &peratom }, "Force B-factor attachment per atom" },
687 { "-legend", FALSE, etBOOL,
688 { &bLegend }, "Make B-factor legend" },
689 { "-label", FALSE, etSTR,
690 { &label }, "Add chain label for all residues" },
692 "-conect", FALSE, etBOOL,
693 { &bCONECT },
694 "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a topology is present"
697 #define NPA asize(pa)
699 FILE *out;
700 const char *infile, *outfile;
701 int outftp, inftp, natom, i, j, n_bfac, itype, ntype;
702 double *bfac = NULL, c6, c12;
703 int *bfac_nr = NULL;
704 t_topology *top = NULL;
705 char *grpname, *sgrpname, *agrpname;
706 int isize, ssize, asize;
707 atom_id *index, *sindex, *aindex;
708 rvec *x, *v, gc, rmin, rmax, size;
709 int ePBC;
710 matrix box, rotmatrix, trans;
711 rvec princd, tmpvec;
712 gmx_bool bIndex, bSetSize, bSetAng, bDist, bSetCenter, bAlign;
713 gmx_bool bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
714 real diam = 0, mass = 0, d, vdw;
715 gmx_atomprop_t aps;
716 gmx_conect conect;
717 gmx_output_env_t *oenv;
718 t_filenm fnm[] =
720 { efSTX, "-f", NULL, ffREAD },
721 { efNDX, "-n", NULL, ffOPTRD },
722 { efSTO, NULL, NULL, ffOPTWR },
723 { efPQR, "-mead", "mead", ffOPTWR },
724 { efDAT, "-bf", "bfact", ffOPTRD }
726 #define NFILE asize(fnm)
728 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa,
729 asize(desc), desc, asize(bugs), bugs, &oenv))
731 return 0;
734 bIndex = opt2bSet("-n", NFILE, fnm) || bNDEF;
735 bMead = opt2bSet("-mead", NFILE, fnm);
736 bSetSize = opt2parg_bSet("-box", NPA, pa);
737 bSetAng = opt2parg_bSet("-angles", NPA, pa);
738 bSetCenter = opt2parg_bSet("-center", NPA, pa);
739 bDist = opt2parg_bSet("-d", NPA, pa);
740 bAlign = opt2parg_bSet("-align", NPA, pa);
741 /* Only automatically turn on centering without -noc */
742 if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
744 bCenter = TRUE;
746 bScale = opt2parg_bSet("-scale", NPA, pa);
747 bRho = opt2parg_bSet("-density", NPA, pa);
748 bTranslate = opt2parg_bSet("-translate", NPA, pa);
749 bRotate = opt2parg_bSet("-rotate", NPA, pa);
750 if (bScale && bRho)
752 fprintf(stderr, "WARNING: setting -density overrides -scale\n");
754 bScale = bScale || bRho;
755 bCalcGeom = bCenter || bRotate || bOrient || bScale;
757 GMX_RELEASE_ASSERT(btype[0] != NULL, "Option setting inconsistency; btype[0] is NULL");
759 bCalcDiam = (btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o');
761 infile = ftp2fn(efSTX, NFILE, fnm);
762 if (bMead)
764 outfile = ftp2fn(efPQR, NFILE, fnm);
766 else
768 outfile = ftp2fn(efSTO, NFILE, fnm);
770 outftp = fn2ftp(outfile);
771 inftp = fn2ftp(infile);
773 aps = gmx_atomprop_init();
775 if (bMead && bGrasp)
777 printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
778 bGrasp = FALSE;
780 if (bGrasp && (outftp != efPDB))
782 gmx_fatal(FARGS, "Output file should be a .pdb file"
783 " when using the -grasp option\n");
785 if ((bMead || bGrasp) && (fn2ftp(infile) != efTPR))
787 gmx_fatal(FARGS, "Input file should be a .tpr file"
788 " when using the -mead option\n");
791 t_topology *top_tmp;
792 snew(top_tmp, 1);
793 read_tps_conf(infile, top_tmp, &ePBC, &x, &v, box, FALSE);
794 t_atoms &atoms = top_tmp->atoms;
795 natom = atoms.nr;
796 if (atoms.pdbinfo == NULL)
798 snew(atoms.pdbinfo, atoms.nr);
800 if (fn2ftp(infile) == efPDB)
802 get_pdb_atomnumber(&atoms, aps);
804 printf("Read %d atoms\n", atoms.nr);
806 /* Get the element numbers if available in a pdb file */
807 if (fn2ftp(infile) == efPDB)
809 get_pdb_atomnumber(&atoms, aps);
812 if (ePBC != epbcNONE)
814 real vol = det(box);
815 printf("Volume: %g nm^3, corresponds to roughly %d electrons\n",
816 vol, 100*(static_cast<int>(vol*4.5)));
819 if (bMead || bGrasp || bCONECT)
821 top = read_top(infile, NULL);
824 if (bMead || bGrasp)
826 if (atoms.nr != top->atoms.nr)
828 gmx_fatal(FARGS, "Atom numbers don't match (%d vs. %d)", atoms.nr, top->atoms.nr);
830 snew(atoms.pdbinfo, top->atoms.nr);
831 ntype = top->idef.atnr;
832 for (i = 0; (i < atoms.nr); i++)
834 /* Determine the Van der Waals radius from the force field */
835 if (bReadVDW)
837 if (!gmx_atomprop_query(aps, epropVDW,
838 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
839 *top->atoms.atomname[i], &vdw))
841 vdw = rvdw;
844 else
846 itype = top->atoms.atom[i].type;
847 c12 = top->idef.iparams[itype*ntype+itype].lj.c12;
848 c6 = top->idef.iparams[itype*ntype+itype].lj.c6;
849 if ((c6 != 0) && (c12 != 0))
851 real sig6;
852 if (bSig56)
854 sig6 = 2*c12/c6;
856 else
858 sig6 = c12/c6;
860 vdw = 0.5*std::pow(sig6, static_cast<real>(1.0/6.0));
862 else
864 vdw = rvdw;
867 /* Factor of 10 for nm -> Angstroms */
868 vdw *= 10;
870 if (bMead)
872 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
873 atoms.pdbinfo[i].bfac = vdw;
875 else
877 atoms.pdbinfo[i].occup = vdw;
878 atoms.pdbinfo[i].bfac = top->atoms.atom[i].q;
882 bHaveV = FALSE;
883 for (i = 0; (i < natom) && !bHaveV; i++)
885 for (j = 0; (j < DIM) && !bHaveV; j++)
887 bHaveV = bHaveV || (v[i][j] != 0);
890 printf("%selocities found\n", bHaveV ? "V" : "No v");
892 if (visbox[0] > 0)
894 if (bIndex)
896 gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
898 if (outftp != efPDB)
900 gmx_fatal(FARGS, "Sorry, can only visualize box with a pdb file");
903 else if (visbox[0] == -1)
905 visualize_images("images.pdb", ePBC, box);
908 /* remove pbc */
909 if (bRMPBC)
911 rm_gropbc(&atoms, x, box);
914 if (bCalcGeom)
916 if (bIndex)
918 fprintf(stderr, "\nSelect a group for determining the system size:\n");
919 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
920 1, &ssize, &sindex, &sgrpname);
922 else
924 ssize = atoms.nr;
925 sindex = NULL;
927 diam = calc_geom(ssize, sindex, x, gc, rmin, rmax, bCalcDiam);
928 rvec_sub(rmax, rmin, size);
929 printf(" system size :%7.3f%7.3f%7.3f (nm)\n",
930 size[XX], size[YY], size[ZZ]);
931 if (bCalcDiam)
933 printf(" diameter :%7.3f (nm)\n", diam);
935 printf(" center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
936 printf(" box vectors :%7.3f%7.3f%7.3f (nm)\n",
937 norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
938 printf(" box angles :%7.2f%7.2f%7.2f (degrees)\n",
939 norm2(box[ZZ]) == 0 ? 0 :
940 RAD2DEG*std::acos(cos_angle_no_table(box[YY], box[ZZ])),
941 norm2(box[ZZ]) == 0 ? 0 :
942 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[ZZ])),
943 norm2(box[YY]) == 0 ? 0 :
944 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[YY])));
945 printf(" box volume :%7.2f (nm^3)\n", det(box));
948 if (bRho || bOrient || bAlign)
950 mass = calc_mass(&atoms, !fn2bTPX(infile), aps);
953 if (bOrient)
955 atom_id *index;
956 char *grpnames;
958 /* Get a group for principal component analysis */
959 fprintf(stderr, "\nSelect group for the determining the orientation\n");
960 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
962 /* Orient the principal axes along the coordinate axes */
963 orient_princ(&atoms, isize, index, natom, x, bHaveV ? v : NULL, NULL);
964 sfree(index);
965 sfree(grpnames);
968 if (bScale)
970 /* scale coordinates and box */
971 if (bRho)
973 /* Compute scaling constant */
974 real vol, dens;
976 vol = det(box);
977 dens = (mass*AMU)/(vol*NANO*NANO*NANO);
978 fprintf(stderr, "Volume of input %g (nm^3)\n", vol);
979 fprintf(stderr, "Mass of input %g (a.m.u.)\n", mass);
980 fprintf(stderr, "Density of input %g (g/l)\n", dens);
981 if (vol == 0 || mass == 0)
983 gmx_fatal(FARGS, "Cannot scale density with "
984 "zero mass (%g) or volume (%g)\n", mass, vol);
987 scale[XX] = scale[YY] = scale[ZZ] = std::pow(dens/rho, static_cast<real>(1.0/3.0));
988 fprintf(stderr, "Scaling all box vectors by %g\n", scale[XX]);
990 scale_conf(atoms.nr, x, box, scale);
993 if (bAlign)
995 if (bIndex)
997 fprintf(stderr, "\nSelect a group that you want to align:\n");
998 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
999 1, &asize, &aindex, &agrpname);
1001 else
1003 asize = atoms.nr;
1004 snew(aindex, asize);
1005 for (i = 0; i < asize; i++)
1007 aindex[i] = i;
1010 printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n", asize, natom,
1011 targetvec[XX], targetvec[YY], targetvec[ZZ],
1012 aligncenter[XX], aligncenter[YY], aligncenter[ZZ]);
1013 /*subtract out pivot point*/
1014 for (i = 0; i < asize; i++)
1016 rvec_dec(x[aindex[i]], aligncenter);
1018 /*now determine transform and rotate*/
1019 /*will this work?*/
1020 principal_comp(asize, aindex, atoms.atom, x, trans, princd);
1022 unitv(targetvec, targetvec);
1023 printf("Using %g %g %g as principal axis\n", trans[0][2], trans[1][2], trans[2][2]);
1024 tmpvec[XX] = trans[0][2]; tmpvec[YY] = trans[1][2]; tmpvec[ZZ] = trans[2][2];
1025 calc_rotmatrix(tmpvec, targetvec, rotmatrix);
1026 /* rotmatrix finished */
1028 for (i = 0; i < asize; ++i)
1030 mvmul(rotmatrix, x[aindex[i]], tmpvec);
1031 copy_rvec(tmpvec, x[aindex[i]]);
1034 /*add pivot point back*/
1035 for (i = 0; i < asize; i++)
1037 rvec_inc(x[aindex[i]], aligncenter);
1039 if (!bIndex)
1041 sfree(aindex);
1045 if (bTranslate)
1047 if (bIndex)
1049 fprintf(stderr, "\nSelect a group that you want to translate:\n");
1050 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
1051 1, &ssize, &sindex, &sgrpname);
1053 else
1055 ssize = atoms.nr;
1056 sindex = NULL;
1058 printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize, natom,
1059 translation[XX], translation[YY], translation[ZZ]);
1060 if (sindex)
1062 for (i = 0; i < ssize; i++)
1064 rvec_inc(x[sindex[i]], translation);
1067 else
1069 for (i = 0; i < natom; i++)
1071 rvec_inc(x[i], translation);
1075 if (bRotate)
1077 /* Rotate */
1078 printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n", rotangles[XX], rotangles[YY], rotangles[ZZ]);
1079 for (i = 0; i < DIM; i++)
1081 rotangles[i] *= DEG2RAD;
1083 rotate_conf(natom, x, v, rotangles[XX], rotangles[YY], rotangles[ZZ]);
1086 if (bCalcGeom)
1088 /* recalc geometrical center and max and min coordinates and size */
1089 calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1090 rvec_sub(rmax, rmin, size);
1091 if (bScale || bOrient || bRotate)
1093 printf("new system size : %6.3f %6.3f %6.3f\n",
1094 size[XX], size[YY], size[ZZ]);
1098 if ((btype[0] != NULL) && (bSetSize || bDist || (btype[0][0] == 't' && bSetAng)))
1100 ePBC = epbcXYZ;
1101 if (!(bSetSize || bDist))
1103 for (i = 0; i < DIM; i++)
1105 newbox[i] = norm(box[i]);
1108 clear_mat(box);
1109 /* calculate new boxsize */
1110 switch (btype[0][0])
1112 case 't':
1113 if (bDist)
1115 for (i = 0; i < DIM; i++)
1117 newbox[i] = size[i]+2*dist;
1120 if (!bSetAng)
1122 box[XX][XX] = newbox[XX];
1123 box[YY][YY] = newbox[YY];
1124 box[ZZ][ZZ] = newbox[ZZ];
1126 else
1128 matrix_convert(box, newbox, newang);
1130 break;
1131 case 'c':
1132 case 'd':
1133 case 'o':
1134 if (bSetSize)
1136 d = newbox[0];
1138 else
1140 d = diam+2*dist;
1142 if (btype[0][0] == 'c')
1144 for (i = 0; i < DIM; i++)
1146 box[i][i] = d;
1149 else if (btype[0][0] == 'd')
1151 box[XX][XX] = d;
1152 box[YY][YY] = d;
1153 box[ZZ][XX] = d/2;
1154 box[ZZ][YY] = d/2;
1155 box[ZZ][ZZ] = d*std::sqrt(2.0)/2.0;
1157 else
1159 box[XX][XX] = d;
1160 box[YY][XX] = d/3;
1161 box[YY][YY] = d*std::sqrt(2.0)*2.0/3.0;
1162 box[ZZ][XX] = -d/3;
1163 box[ZZ][YY] = d*std::sqrt(2.0)/3.0;
1164 box[ZZ][ZZ] = d*std::sqrt(6.0)/3.0;
1166 break;
1170 /* calculate new coords for geometrical center */
1171 if (!bSetCenter)
1173 calc_box_center(ecenterDEF, box, center);
1176 /* center molecule on 'center' */
1177 if (bCenter)
1179 center_conf(natom, x, center, gc);
1182 /* print some */
1183 if (bCalcGeom)
1185 calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1186 printf("new center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
1188 if (bOrient || bScale || bDist || bSetSize)
1190 printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n",
1191 norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
1192 printf("new box angles :%7.2f%7.2f%7.2f (degrees)\n",
1193 norm2(box[ZZ]) == 0 ? 0 :
1194 RAD2DEG*std::acos(cos_angle_no_table(box[YY], box[ZZ])),
1195 norm2(box[ZZ]) == 0 ? 0 :
1196 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[ZZ])),
1197 norm2(box[YY]) == 0 ? 0 :
1198 RAD2DEG*std::acos(cos_angle_no_table(box[XX], box[YY])));
1199 printf("new box volume :%7.2f (nm^3)\n", det(box));
1202 if (check_box(epbcXYZ, box))
1204 printf("\nWARNING: %s\n"
1205 "See the GROMACS manual for a description of the requirements that\n"
1206 "must be satisfied by descriptions of simulation cells.\n",
1207 check_box(epbcXYZ, box));
1210 if (bDist && btype[0][0] == 't')
1212 if (TRICLINIC(box))
1214 printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1215 "distance from the solute to a box surface along the corresponding normal\n"
1216 "vector might be somewhat smaller than your specified value %f.\n"
1217 "You can check the actual value with g_mindist -pi\n", dist);
1219 else if (!opt2parg_bSet("-bt", NPA, pa))
1221 printf("\nWARNING: No boxtype specified - distance condition applied in each dimension.\n"
1222 "If the molecule rotates the actual distance will be smaller. You might want\n"
1223 "to use a cubic box instead, or why not try a dodecahedron today?\n");
1226 if (bCONECT && (outftp == efPDB) && (inftp == efTPR))
1228 conect = gmx_conect_generate(top);
1230 else
1232 conect = NULL;
1235 if (bIndex)
1237 fprintf(stderr, "\nSelect a group for output:\n");
1238 get_index(&atoms, opt2fn_null("-n", NFILE, fnm),
1239 1, &isize, &index, &grpname);
1241 if (resnr_start >= 0)
1243 renum_resnr(&atoms, isize, index, resnr_start);
1246 if (opt2parg_bSet("-label", NPA, pa))
1248 for (i = 0; (i < atoms.nr); i++)
1250 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1254 if (opt2bSet("-bf", NFILE, fnm) || bLegend)
1256 gmx_fatal(FARGS, "Sorry, cannot do bfactors with an index group.");
1259 if (outftp == efPDB)
1261 out = gmx_ffopen(outfile, "w");
1262 write_pdbfile_indexed(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, TRUE);
1263 gmx_ffclose(out);
1265 else
1267 write_sto_conf_indexed(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : NULL, ePBC, box, isize, index);
1270 else
1272 if (resnr_start >= 0)
1274 renum_resnr(&atoms, atoms.nr, NULL, resnr_start);
1277 if ((outftp == efPDB) || (outftp == efPQR))
1279 out = gmx_ffopen(outfile, "w");
1280 if (bMead)
1282 fprintf(out, "REMARK "
1283 "The B-factors in this file hold atomic radii\n");
1284 fprintf(out, "REMARK "
1285 "The occupancy in this file hold atomic charges\n");
1287 else if (bGrasp)
1289 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1290 fprintf(out, "REMARK "
1291 "The B-factors in this file hold atomic charges\n");
1292 fprintf(out, "REMARK "
1293 "The occupancy in this file hold atomic radii\n");
1295 else if (opt2bSet("-bf", NFILE, fnm))
1297 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1298 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms,
1299 n_bfac, bfac, bfac_nr, peratom);
1301 if (opt2parg_bSet("-label", NPA, pa))
1303 for (i = 0; (i < atoms.nr); i++)
1305 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1308 write_pdbfile(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', -1, conect, TRUE);
1309 if (bLegend)
1311 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1313 if (visbox[0] > 0)
1315 visualize_box(out, bLegend ? atoms.nr+12 : atoms.nr,
1316 bLegend ? atoms.nres = 12 : atoms.nres, box, visbox);
1318 gmx_ffclose(out);
1320 else
1322 write_sto_conf(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : NULL, ePBC, box);
1325 gmx_atomprop_destroy(aps);
1327 do_view(oenv, outfile, NULL);
1329 return 0;