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,2016,2017,2018,2019, 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.
47 #include "gromacs/commandline/pargs.h"
48 #include "gromacs/commandline/viewit.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/fileio/pdbio.h"
51 #include "gromacs/fileio/tpxio.h"
52 #include "gromacs/fileio/trxio.h"
53 #include "gromacs/gmxana/princ.h"
54 #include "gromacs/gmxlib/conformation_utilities.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/pbcutil/rmpbc.h"
60 #include "gromacs/topology/atomprop.h"
61 #include "gromacs/topology/index.h"
62 #include "gromacs/topology/topology.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/futil.h"
67 #include "gromacs/utility/gmxassert.h"
68 #include "gromacs/utility/smalloc.h"
69 #include "gromacs/utility/strdb.h"
71 static real calc_mass(t_atoms* atoms, gmx_bool bGetMass, AtomProperties* aps)
77 for (i = 0; (i < atoms->nr); i++)
81 aps->setAtomProperty(epropMass, std::string(*atoms->resinfo[atoms->atom[i].resind].name),
82 std::string(*atoms->atomname[i]), &(atoms->atom[i].m));
84 tmass += atoms->atom[i].m;
90 static real calc_geom(int isize, const int* index, rvec* x, rvec geom_center, rvec minval, rvec maxval, gmx_bool bDiam)
95 clear_rvec(geom_center);
112 for (j = 0; j < DIM; j++)
114 minval[j] = maxval[j] = x[ii][j];
116 for (i = 0; i < isize; 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];
142 for (j = i + 1; j < isize; j++)
144 d = distance2(x[ii], x[index[j]]);
145 diam2 = std::max(d, diam2);
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 static void center_conf(int natom, rvec* x, rvec center, rvec geom_cent)
169 rvec_sub(center, geom_cent, shift);
171 printf(" shift :%7.3f%7.3f%7.3f (nm)\n", shift[XX], shift[YY], shift[ZZ]);
173 for (i = 0; (i < natom); i++)
175 rvec_inc(x[i], shift);
179 static void scale_conf(int natom, rvec x[], matrix box, const rvec scale)
183 for (i = 0; i < natom; i++)
185 for (j = 0; j < DIM; j++)
190 for (i = 0; i < DIM; i++)
192 for (j = 0; j < DIM; j++)
194 box[i][j] *= scale[j];
199 static void read_bfac(const char* fn, int* n_bfac, double** bfac_val, int** bfac_nr)
204 *n_bfac = get_lines(fn, &bfac_lines);
205 snew(*bfac_val, *n_bfac);
206 snew(*bfac_nr, *n_bfac);
207 fprintf(stderr, "Reading %d B-factors from %s\n", *n_bfac, fn);
208 for (i = 0; (i < *n_bfac); i++)
210 sscanf(bfac_lines[i], "%d %lf", &(*bfac_nr)[i], &(*bfac_val)[i]);
214 static void set_pdb_conf_bfac(int natoms, int nres, t_atoms* atoms, int n_bfac, double* bfac, int* bfac_nr, gmx_bool peratom)
216 real bfac_min, bfac_max;
220 if (n_bfac > atoms->nres)
227 for (i = 0; (i < n_bfac); i++)
229 /* if ((bfac_nr[i]-1<0) || (bfac_nr[i]-1>=atoms->nr))
230 gmx_fatal(FARGS,"Index of B-Factor %d is out of range: %d (%g)",
231 i+1,bfac_nr[i],bfac[i]); */
232 if (bfac[i] > bfac_max)
236 if (bfac[i] < bfac_min)
241 while ((bfac_max > 99.99) || (bfac_min < -99.99))
244 "Range of values for B-factors too large (min %g, max %g) "
245 "will scale down a factor 10\n",
247 for (i = 0; (i < n_bfac); i++)
254 while ((std::abs(bfac_max) < 0.5) && (std::abs(bfac_min) < 0.5))
257 "Range of values for B-factors too small (min %g, max %g) "
258 "will scale up a factor 10\n",
260 for (i = 0; (i < n_bfac); i++)
268 for (i = 0; (i < natoms); i++)
270 atoms->pdbinfo[i].bfac = 0;
275 fprintf(stderr, "Will attach %d B-factors to %d residues\n", n_bfac, nres);
276 for (i = 0; (i < n_bfac); i++)
279 for (n = 0; (n < natoms); n++)
281 if (bfac_nr[i] == atoms->resinfo[atoms->atom[n].resind].nr)
283 atoms->pdbinfo[n].bfac = bfac[i];
289 gmx_warning("Residue nr %d not found\n", bfac_nr[i]);
295 fprintf(stderr, "Will attach %d B-factors to %d atoms\n", n_bfac, natoms);
296 for (i = 0; (i < n_bfac); i++)
298 atoms->pdbinfo[bfac_nr[i] - 1].bfac = bfac[i];
303 static void pdb_legend(FILE* out, int natoms, int nres, t_atoms* atoms, rvec x[])
305 real bfac_min, bfac_max, xmin, ymin, zmin;
314 for (i = 0; (i < natoms); i++)
316 xmin = std::min(xmin, x[i][XX]);
317 ymin = std::min(ymin, x[i][YY]);
318 zmin = std::min(zmin, x[i][ZZ]);
319 bfac_min = std::min(bfac_min, atoms->pdbinfo[i].bfac);
320 bfac_max = std::max(bfac_max, atoms->pdbinfo[i].bfac);
322 fprintf(stderr, "B-factors range from %g to %g\n", bfac_min, bfac_max);
323 for (i = 1; (i < 12); i++)
325 fprintf(out, "%-6s%5d %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f\n", "ATOM ",
326 natoms + 1 + i, "CA", "LEG", space, nres + 1, space, (xmin + (i * 0.12)) * 10,
327 ymin * 10, zmin * 10, 1.0, bfac_min + ((i - 1.0) * (bfac_max - bfac_min) / 10));
331 static void visualize_images(const char* fn, int ePBC, matrix box)
339 init_t_atoms(&atoms, nat, FALSE);
342 /* FIXME: Constness should not be cast away */
343 c = const_cast<char*>("C");
344 ala = const_cast<char*>("ALA");
345 for (i = 0; i < nat; i++)
347 atoms.atomname[i] = &c;
348 atoms.atom[i].resind = i;
349 atoms.resinfo[i].name = &ala;
350 atoms.resinfo[i].nr = i + 1;
351 atoms.resinfo[i].chainid = 'A' + i / NCUCVERT;
353 calc_triclinic_images(box, img + 1);
355 write_sto_conf(fn, "Images", &atoms, img, nullptr, ePBC, box);
361 static void visualize_box(FILE* out, int a0, int r0, matrix box, const rvec gridsize)
365 int nx, ny, nz, nbox, nat;
367 int rectedge[24] = { 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6, 4 };
372 nx = gmx::roundToInt(gridsize[XX]);
373 ny = gmx::roundToInt(gridsize[YY]);
374 nz = gmx::roundToInt(gridsize[ZZ]);
378 nat = nbox * NCUCVERT;
380 calc_compact_unitcell_vertices(ecenterDEF, box, vert);
382 for (z = 0; z < nz; z++)
384 for (y = 0; y < ny; y++)
386 for (x = 0; x < nx; x++)
388 for (i = 0; i < DIM; i++)
390 shift[i] = x * box[0][i] + y * box[1][i] + z * box[2][i];
392 for (i = 0; i < NCUCVERT; i++)
394 rvec_add(vert[i], shift, vert[j]);
401 for (i = 0; i < nat; i++)
403 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i / NCUCVERT,
404 r0 + i, ' ', 10 * vert[i][XX], 10 * vert[i][YY],
405 10 * vert[i][ZZ], 1.0, 0.0, "");
408 edge = compact_unitcell_edges();
409 for (j = 0; j < nbox; j++)
411 for (i = 0; i < NCUCEDGE; i++)
413 fprintf(out, "CONECT%5d%5d\n", a0 + j * NCUCVERT + edge[2 * i],
414 a0 + j * NCUCVERT + edge[2 * i + 1]);
423 for (z = 0; z <= 1; z++)
425 for (y = 0; y <= 1; y++)
427 for (x = 0; x <= 1; x++)
429 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i / 8,
430 r0 + i, ' ', x * 10 * box[XX][XX], y * 10 * box[YY][YY],
431 z * 10 * box[ZZ][ZZ], 1.0, 0.0, "");
436 for (i = 0; i < 24; i += 2)
438 fprintf(out, "CONECT%5d%5d\n", a0 + rectedge[i], a0 + rectedge[i + 1]);
443 static void calc_rotmatrix(rvec principal_axis, rvec targetvec, matrix rotmatrix)
446 real ux, uy, uz, costheta, sintheta;
448 costheta = cos_angle(principal_axis, targetvec);
449 sintheta = std::sqrt(1.0 - costheta * costheta); /* sign is always positive since 0<theta<pi */
451 /* Determine rotation from cross product with target vector */
452 cprod(principal_axis, targetvec, rotvec);
453 unitv(rotvec, rotvec);
454 printf("Aligning %g %g %g to %g %g %g : xprod %g %g %g\n", principal_axis[XX],
455 principal_axis[YY], principal_axis[ZZ], targetvec[XX], targetvec[YY], targetvec[ZZ],
456 rotvec[XX], rotvec[YY], rotvec[ZZ]);
461 rotmatrix[0][0] = ux * ux + (1.0 - ux * ux) * costheta;
462 rotmatrix[0][1] = ux * uy * (1 - costheta) - uz * sintheta;
463 rotmatrix[0][2] = ux * uz * (1 - costheta) + uy * sintheta;
464 rotmatrix[1][0] = ux * uy * (1 - costheta) + uz * sintheta;
465 rotmatrix[1][1] = uy * uy + (1.0 - uy * uy) * costheta;
466 rotmatrix[1][2] = uy * uz * (1 - costheta) - ux * sintheta;
467 rotmatrix[2][0] = ux * uz * (1 - costheta) - uy * sintheta;
468 rotmatrix[2][1] = uy * uz * (1 - costheta) + ux * sintheta;
469 rotmatrix[2][2] = uz * uz + (1.0 - uz * uz) * costheta;
471 printf("Rotation matrix: \n%g %g %g\n%g %g %g\n%g %g %g\n", rotmatrix[0][0], rotmatrix[0][1],
472 rotmatrix[0][2], rotmatrix[1][0], rotmatrix[1][1], rotmatrix[1][2], rotmatrix[2][0],
473 rotmatrix[2][1], rotmatrix[2][2]);
476 static void renum_resnr(t_atoms* atoms, int isize, const int* index, int resnr_start)
478 int i, resind_prev, resind;
481 for (i = 0; i < isize; i++)
483 resind = atoms->atom[index == nullptr ? i : index[i]].resind;
484 if (resind != resind_prev)
486 atoms->resinfo[resind].nr = resnr_start;
489 resind_prev = resind;
493 int gmx_editconf(int argc, char* argv[])
495 const char* desc[] = {
496 "[THISMODULE] converts generic structure format to [REF].gro[ref], [TT].g96[tt]",
497 "or [REF].pdb[ref].",
499 "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and",
500 "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]",
501 "will center the system in the box, unless [TT]-noc[tt] is used.",
502 "The [TT]-center[tt] option can be used to shift the geometric center",
503 "of the system from the default of (x/2, y/2, z/2) implied by [TT]-c[tt]",
504 "to some other value.",
506 "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
507 "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
508 "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
509 "[TT]octahedron[tt] is a truncated octahedron.",
510 "The last two are special cases of a triclinic box.",
511 "The length of the three box vectors of the truncated octahedron is the",
512 "shortest distance between two opposite hexagons.",
513 "Relative to a cubic box with some periodic image distance, the volume of a ",
514 "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
515 "and that of a truncated octahedron is 0.77 times.",
517 "Option [TT]-box[tt] requires only",
518 "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
520 "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, ",
522 "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
523 "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
524 "to the diameter of the system (largest distance between atoms) plus twice",
525 "the specified distance.",
527 "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
528 "a triclinic box and cannot be used with option [TT]-d[tt].",
530 "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
531 "can be selected for calculating the size and the geometric center,",
532 "otherwise the whole system is used.",
534 "[TT]-rotate[tt] rotates the coordinates and velocities.",
536 "[TT]-princ[tt] aligns the principal axes of the system along the",
537 "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
538 "This may allow you to decrease the box volume,",
539 "but beware that molecules can rotate significantly in a nanosecond.",
541 "Scaling is applied before any of the other operations are",
542 "performed. Boxes and coordinates can be scaled to give a certain density (option",
543 "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
544 "file is given as input. A special feature of the scaling option is that when the",
545 "factor -1 is given in one dimension, one obtains a mirror image,",
546 "mirrored in one of the planes. When one uses -1 in three dimensions, ",
547 "a point-mirror image is obtained.[PAR]",
548 "Groups are selected after all operations have been applied.[PAR]",
549 "Periodicity can be removed in a crude manner.",
550 "It is important that the box vectors at the bottom of your input file",
551 "are correct when the periodicity is to be removed.",
553 "When writing [REF].pdb[ref] files, B-factors can be",
554 "added with the [TT]-bf[tt] option. B-factors are read",
555 "from a file with with following format: first line states number of",
556 "entries in the file, next lines state an index",
557 "followed by a B-factor. The B-factors will be attached per residue",
558 "unless the number of B-factors is larger than the number of the residues or unless the",
559 "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
560 "be added instead of B-factors. [TT]-legend[tt] will produce",
561 "a row of CA atoms with B-factors ranging from the minimum to the",
562 "maximum value found, effectively making a legend for viewing.",
564 "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
565 "file for the MEAD electrostatics",
566 "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
567 "is that the input file is a run input file.",
568 "The B-factor field is then filled with the Van der Waals radius",
569 "of the atoms while the occupancy field will hold the charge.",
571 "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
572 "and the radius in the occupancy.",
574 "Option [TT]-align[tt] allows alignment",
575 "of the principal axis of a specified group against the given vector, ",
576 "with an optional center of rotation specified by [TT]-aligncenter[tt].",
578 "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
579 "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
581 "To convert a truncated octrahedron file produced by a package which uses",
582 "a cubic box with the corners cut off (such as GROMOS), use::",
584 " gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
586 "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
588 const char* bugs[] = {
589 "For complex molecules, the periodicity removal routine may break down, ",
590 "in that case you can use [gmx-trjconv]."
592 static real dist = 0.0;
593 static gmx_bool bNDEF = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW = FALSE, bCONECT = FALSE;
594 static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead = FALSE,
595 bGrasp = FALSE, bSig56 = FALSE;
596 static rvec scale = { 1, 1, 1 }, newbox = { 0, 0, 0 }, newang = { 90, 90, 90 };
597 static real rho = 1000.0, rvdw = 0.12;
598 static rvec center = { 0, 0, 0 }, translation = { 0, 0, 0 }, rotangles = { 0, 0, 0 },
599 aligncenter = { 0, 0, 0 }, targetvec = { 0, 0, 0 };
600 static const char *btype[] = { nullptr, "triclinic", "cubic",
601 "dodecahedron", "octahedron", nullptr },
603 static rvec visbox = { 0, 0, 0 };
604 static int resnr_start = -1;
606 { "-ndef", FALSE, etBOOL, { &bNDEF }, "Choose output from default index groups" },
611 "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
612 { "-bt", FALSE, etENUM, { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
613 { "-box", FALSE, etRVEC, { newbox }, "Box vector lengths (a,b,c)" },
614 { "-angles", FALSE, etRVEC, { newang }, "Angles between the box vectors (bc,ac,ab)" },
615 { "-d", FALSE, etREAL, { &dist }, "Distance between the solute and the box" },
620 "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
621 { "-center", FALSE, etRVEC, { center }, "Shift the geometrical center to (x,y,z)" },
622 { "-aligncenter", FALSE, etRVEC, { aligncenter }, "Center of rotation for alignment" },
623 { "-align", FALSE, etRVEC, { targetvec }, "Align to target vector" },
624 { "-translate", FALSE, etRVEC, { translation }, "Translation" },
629 "Rotation around the X, Y and Z axes in degrees" },
630 { "-princ", FALSE, etBOOL, { &bOrient }, "Orient molecule(s) along their principal axes" },
631 { "-scale", FALSE, etRVEC, { scale }, "Scaling factor" },
636 "Density (g/L) of the output box achieved by scaling" },
637 { "-pbc", FALSE, etBOOL, { &bRMPBC }, "Remove the periodicity (make molecule whole again)" },
638 { "-resnr", FALSE, etINT, { &resnr_start }, " Renumber residues starting from resnr" },
643 "Store the charge of the atom in the B-factor field and the radius of the atom in the "
649 "Default Van der Waals radius (in nm) if one can not be found in the database or if no "
650 "parameters are present in the topology file" },
655 "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
660 "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing "
661 "the radii based on the force field" },
662 { "-atom", FALSE, etBOOL, { &peratom }, "Force B-factor attachment per atom" },
663 { "-legend", FALSE, etBOOL, { &bLegend }, "Make B-factor legend" },
664 { "-label", FALSE, etSTR, { &label }, "Add chain label for all residues" },
669 "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a "
670 "topology is present" }
672 #define NPA asize(pa)
675 const char * infile, *outfile;
676 int outftp, inftp, natom, i, j, n_bfac, itype, ntype;
677 double * bfac = nullptr, c6, c12;
678 int* bfac_nr = nullptr;
679 t_topology* top = nullptr;
680 char * grpname, *sgrpname, *agrpname;
681 int isize, ssize, numAlignmentAtoms;
682 int * index, *sindex, *aindex;
683 rvec * x, *v, gc, rmin, rmax, size;
685 matrix box, rotmatrix, trans;
687 gmx_bool bIndex, bSetSize, bSetAng, bDist, bSetCenter, bAlign;
688 gmx_bool bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
689 real diam = 0, mass = 0, d, vdw;
691 gmx_output_env_t* oenv;
692 t_filenm fnm[] = { { efSTX, "-f", nullptr, ffREAD },
693 { efNDX, "-n", nullptr, ffOPTRD },
694 { efSTO, nullptr, nullptr, ffOPTWR },
695 { efPQR, "-mead", "mead", ffOPTWR },
696 { efDAT, "-bf", "bfact", ffOPTRD } };
697 #define NFILE asize(fnm)
699 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa, asize(desc), desc,
700 asize(bugs), bugs, &oenv))
705 "Note that major changes are planned in future for "
706 "editconf, to improve usability and utility.\n");
708 bIndex = opt2bSet("-n", NFILE, fnm) || bNDEF;
709 bMead = opt2bSet("-mead", NFILE, fnm);
710 bSetSize = opt2parg_bSet("-box", NPA, pa);
711 bSetAng = opt2parg_bSet("-angles", NPA, pa);
712 bSetCenter = opt2parg_bSet("-center", NPA, pa);
713 bDist = opt2parg_bSet("-d", NPA, pa);
714 bAlign = opt2parg_bSet("-align", NPA, pa);
715 /* Only automatically turn on centering without -noc */
716 if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
720 bScale = opt2parg_bSet("-scale", NPA, pa);
721 bRho = opt2parg_bSet("-density", NPA, pa);
722 bTranslate = opt2parg_bSet("-translate", NPA, pa);
723 bRotate = opt2parg_bSet("-rotate", NPA, pa);
726 fprintf(stderr, "WARNING: setting -density overrides -scale\n");
728 bScale = bScale || bRho;
729 bCalcGeom = bCenter || bRotate || bOrient || bScale;
731 GMX_RELEASE_ASSERT(btype[0] != nullptr, "Option setting inconsistency; btype[0] is NULL");
733 bCalcDiam = (btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o');
735 infile = ftp2fn(efSTX, NFILE, fnm);
738 outfile = ftp2fn(efPQR, NFILE, fnm);
742 outfile = ftp2fn(efSTO, NFILE, fnm);
744 outftp = fn2ftp(outfile);
745 inftp = fn2ftp(infile);
751 printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
754 if (bGrasp && (outftp != efPDB))
757 "Output file should be a .pdb file"
758 " when using the -grasp option\n");
760 if ((bMead || bGrasp) && (fn2ftp(infile) != efTPR))
763 "Input file should be a .tpr file"
764 " when using the -mead option\n");
770 open_symtab(&symtab);
771 readConfAndAtoms(infile, &symtab, &name, &atoms, &ePBC, &x, &v, box);
773 if (atoms.pdbinfo == nullptr)
775 snew(atoms.pdbinfo, atoms.nr);
777 atoms.havePdbInfo = TRUE;
779 if (fn2ftp(infile) == efPDB)
781 get_pdb_atomnumber(&atoms, &aps);
783 printf("Read %d atoms\n", atoms.nr);
785 /* Get the element numbers if available in a pdb file */
786 if (fn2ftp(infile) == efPDB)
788 get_pdb_atomnumber(&atoms, &aps);
791 if (ePBC != epbcNONE)
794 printf("Volume: %g nm^3, corresponds to roughly %d electrons\n", vol,
795 100 * (static_cast<int>(vol * 4.5)));
798 if (bMead || bGrasp || bCONECT)
800 top = read_top(infile, nullptr);
805 if (atoms.nr != top->atoms.nr)
807 gmx_fatal(FARGS, "Atom numbers don't match (%d vs. %d)", atoms.nr, top->atoms.nr);
809 snew(atoms.pdbinfo, top->atoms.nr);
810 ntype = top->idef.atnr;
811 for (i = 0; (i < atoms.nr); i++)
813 /* Determine the Van der Waals radius from the force field */
816 if (!aps.setAtomProperty(epropVDW, *top->atoms.resinfo[top->atoms.atom[i].resind].name,
817 *top->atoms.atomname[i], &vdw))
824 itype = top->atoms.atom[i].type;
825 c12 = top->idef.iparams[itype * ntype + itype].lj.c12;
826 c6 = top->idef.iparams[itype * ntype + itype].lj.c6;
827 if ((c6 != 0) && (c12 != 0))
838 vdw = 0.5 * gmx::sixthroot(sig6);
845 /* Factor of 10 for nm -> Angstroms */
850 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
851 atoms.pdbinfo[i].bfac = vdw;
855 atoms.pdbinfo[i].occup = vdw;
856 atoms.pdbinfo[i].bfac = top->atoms.atom[i].q;
861 for (i = 0; (i < natom) && !bHaveV; i++)
863 for (j = 0; (j < DIM) && !bHaveV; j++)
865 bHaveV = bHaveV || (v[i][j] != 0);
868 printf("%selocities found\n", bHaveV ? "V" : "No v");
874 gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
878 gmx_fatal(FARGS, "Sorry, can only visualize box with a pdb file");
881 else if (visbox[0] == -1)
883 visualize_images("images.pdb", ePBC, box);
889 rm_gropbc(&atoms, x, box);
896 fprintf(stderr, "\nSelect a group for determining the system size:\n");
897 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ssize, &sindex, &sgrpname);
904 diam = calc_geom(ssize, sindex, x, gc, rmin, rmax, bCalcDiam);
905 rvec_sub(rmax, rmin, size);
906 printf(" system size :%7.3f%7.3f%7.3f (nm)\n", size[XX], size[YY], size[ZZ]);
909 printf(" diameter :%7.3f (nm)\n", diam);
911 printf(" center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
912 printf(" box vectors :%7.3f%7.3f%7.3f (nm)\n", norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
913 printf(" box angles :%7.2f%7.2f%7.2f (degrees)\n",
914 norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[YY], box[ZZ]),
915 norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[ZZ]),
916 norm2(box[YY]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[YY]));
917 printf(" box volume :%7.2f (nm^3)\n", det(box));
920 if (bRho || bOrient || bAlign)
922 mass = calc_mass(&atoms, !fn2bTPX(infile), &aps);
930 /* Get a group for principal component analysis */
931 fprintf(stderr, "\nSelect group for the determining the orientation\n");
932 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
934 /* Orient the principal axes along the coordinate axes */
935 orient_princ(&atoms, isize, index, natom, x, bHaveV ? v : nullptr, nullptr);
942 /* scale coordinates and box */
945 /* Compute scaling constant */
949 dens = (mass * AMU) / (vol * NANO * NANO * NANO);
950 fprintf(stderr, "Volume of input %g (nm^3)\n", vol);
951 fprintf(stderr, "Mass of input %g (a.m.u.)\n", mass);
952 fprintf(stderr, "Density of input %g (g/l)\n", dens);
953 if (vol == 0 || mass == 0)
956 "Cannot scale density with "
957 "zero mass (%g) or volume (%g)\n",
961 scale[XX] = scale[YY] = scale[ZZ] = std::cbrt(dens / rho);
962 fprintf(stderr, "Scaling all box vectors by %g\n", scale[XX]);
964 scale_conf(atoms.nr, x, box, scale);
971 fprintf(stderr, "\nSelect a group that you want to align:\n");
972 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &numAlignmentAtoms, &aindex, &agrpname);
976 numAlignmentAtoms = atoms.nr;
977 snew(aindex, numAlignmentAtoms);
978 for (i = 0; i < numAlignmentAtoms; i++)
983 printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n",
984 numAlignmentAtoms, natom, targetvec[XX], targetvec[YY], targetvec[ZZ],
985 aligncenter[XX], aligncenter[YY], aligncenter[ZZ]);
986 /*subtract out pivot point*/
987 for (i = 0; i < numAlignmentAtoms; i++)
989 rvec_dec(x[aindex[i]], aligncenter);
991 /*now determine transform and rotate*/
993 principal_comp(numAlignmentAtoms, aindex, atoms.atom, x, trans, princd);
995 unitv(targetvec, targetvec);
996 printf("Using %g %g %g as principal axis\n", trans[0][2], trans[1][2], trans[2][2]);
997 tmpvec[XX] = trans[0][2];
998 tmpvec[YY] = trans[1][2];
999 tmpvec[ZZ] = trans[2][2];
1000 calc_rotmatrix(tmpvec, targetvec, rotmatrix);
1001 /* rotmatrix finished */
1003 for (i = 0; i < numAlignmentAtoms; ++i)
1005 mvmul(rotmatrix, x[aindex[i]], tmpvec);
1006 copy_rvec(tmpvec, x[aindex[i]]);
1009 /*add pivot point back*/
1010 for (i = 0; i < numAlignmentAtoms; i++)
1012 rvec_inc(x[aindex[i]], aligncenter);
1024 fprintf(stderr, "\nSelect a group that you want to translate:\n");
1025 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ssize, &sindex, &sgrpname);
1032 printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize, natom, translation[XX],
1033 translation[YY], translation[ZZ]);
1036 for (i = 0; i < ssize; i++)
1038 rvec_inc(x[sindex[i]], translation);
1043 for (i = 0; i < natom; i++)
1045 rvec_inc(x[i], translation);
1052 printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n",
1053 rotangles[XX], rotangles[YY], rotangles[ZZ]);
1054 for (i = 0; i < DIM; i++)
1056 rotangles[i] *= DEG2RAD;
1058 rotate_conf(natom, x, v, rotangles[XX], rotangles[YY], rotangles[ZZ]);
1063 /* recalc geometrical center and max and min coordinates and size */
1064 calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1065 rvec_sub(rmax, rmin, size);
1066 if (bScale || bOrient || bRotate)
1068 printf("new system size : %6.3f %6.3f %6.3f\n", size[XX], size[YY], size[ZZ]);
1072 if ((btype[0] != nullptr) && (bSetSize || bDist || (btype[0][0] == 't' && bSetAng)))
1075 if (!(bSetSize || bDist))
1077 for (i = 0; i < DIM; i++)
1079 newbox[i] = norm(box[i]);
1083 /* calculate new boxsize */
1084 switch (btype[0][0])
1089 for (i = 0; i < DIM; i++)
1091 newbox[i] = size[i] + 2 * dist;
1096 box[XX][XX] = newbox[XX];
1097 box[YY][YY] = newbox[YY];
1098 box[ZZ][ZZ] = newbox[ZZ];
1102 matrix_convert(box, newbox, newang);
1114 d = diam + 2 * dist;
1116 if (btype[0][0] == 'c')
1118 for (i = 0; i < DIM; i++)
1123 else if (btype[0][0] == 'd')
1127 box[ZZ][XX] = d / 2;
1128 box[ZZ][YY] = d / 2;
1129 box[ZZ][ZZ] = d * std::sqrt(2.0) / 2.0;
1134 box[YY][XX] = d / 3;
1135 box[YY][YY] = d * std::sqrt(2.0) * 2.0 / 3.0;
1136 box[ZZ][XX] = -d / 3;
1137 box[ZZ][YY] = d * std::sqrt(2.0) / 3.0;
1138 box[ZZ][ZZ] = d * std::sqrt(6.0) / 3.0;
1144 /* calculate new coords for geometrical center */
1147 calc_box_center(ecenterDEF, box, center);
1150 /* center molecule on 'center' */
1153 center_conf(natom, x, center, gc);
1159 calc_geom(ssize, sindex, x, gc, rmin, rmax, FALSE);
1160 printf("new center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
1162 if (bOrient || bScale || bDist || bSetSize)
1164 printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n", norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
1165 printf("new box angles :%7.2f%7.2f%7.2f (degrees)\n",
1166 norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[YY], box[ZZ]),
1167 norm2(box[ZZ]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[ZZ]),
1168 norm2(box[YY]) == 0 ? 0 : RAD2DEG * gmx_angle(box[XX], box[YY]));
1169 printf("new box volume :%7.2f (nm^3)\n", det(box));
1172 if (check_box(epbcXYZ, box))
1174 printf("\nWARNING: %s\n"
1175 "See the GROMACS manual for a description of the requirements that\n"
1176 "must be satisfied by descriptions of simulation cells.\n",
1177 check_box(epbcXYZ, box));
1180 if (bDist && btype[0][0] == 't')
1184 printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1185 "distance from the solute to a box surface along the corresponding normal\n"
1186 "vector might be somewhat smaller than your specified value %f.\n"
1187 "You can check the actual value with g_mindist -pi\n",
1190 else if (!opt2parg_bSet("-bt", NPA, pa))
1192 printf("\nWARNING: No boxtype specified - distance condition applied in each "
1194 "If the molecule rotates the actual distance will be smaller. You might want\n"
1195 "to use a cubic box instead, or why not try a dodecahedron today?\n");
1198 if (bCONECT && (outftp == efPDB) && (inftp == efTPR))
1200 conect = gmx_conect_generate(top);
1209 fprintf(stderr, "\nSelect a group for output:\n");
1210 get_index(&atoms, opt2fn_null("-n", NFILE, fnm), 1, &isize, &index, &grpname);
1212 if (resnr_start >= 0)
1214 renum_resnr(&atoms, isize, index, resnr_start);
1217 if (opt2parg_bSet("-label", NPA, pa))
1219 for (i = 0; (i < atoms.nr); i++)
1221 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1225 if (opt2bSet("-bf", NFILE, fnm) || bLegend)
1227 gmx_fatal(FARGS, "Sorry, cannot do bfactors with an index group.");
1230 if (outftp == efPDB)
1232 out = gmx_ffopen(outfile, "w");
1233 write_pdbfile_indexed(out, name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, FALSE);
1238 write_sto_conf_indexed(outfile, name, &atoms, x, bHaveV ? v : nullptr, ePBC, box, isize, index);
1245 if (resnr_start >= 0)
1247 renum_resnr(&atoms, atoms.nr, nullptr, resnr_start);
1250 if ((outftp == efPDB) || (outftp == efPQR))
1252 out = gmx_ffopen(outfile, "w");
1257 "The B-factors in this file hold atomic radii\n");
1260 "The occupancy in this file hold atomic charges\n");
1264 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1267 "The B-factors in this file hold atomic charges\n");
1270 "The occupancy in this file hold atomic radii\n");
1272 else if (opt2bSet("-bf", NFILE, fnm))
1274 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1275 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms, n_bfac, bfac, bfac_nr, peratom);
1277 if (opt2parg_bSet("-label", NPA, pa))
1279 for (i = 0; (i < atoms.nr); i++)
1281 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1284 /* Need to bypass the regular write_pdbfile because I don't want to change
1285 * all instances to include the boolean flag for writing out PQR files.
1288 snew(index, atoms.nr);
1289 for (int i = 0; i < atoms.nr; i++)
1293 write_pdbfile_indexed(out, name, &atoms, x, ePBC, box, ' ', -1, atoms.nr, index, conect,
1298 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1302 visualize_box(out, bLegend ? atoms.nr + 12 : atoms.nr,
1303 bLegend ? atoms.nres = 12 : atoms.nres, box, visbox);
1309 write_sto_conf(outfile, name, &atoms, x, bHaveV ? v : nullptr, ePBC, box);
1313 done_symtab(&symtab);
1323 do_view(oenv, outfile, nullptr);
1324 output_env_done(oenv);