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
,
82 std::string(*atoms
->resinfo
[atoms
->atom
[i
].resind
].name
),
83 std::string(*atoms
->atomname
[i
]), &(atoms
->atom
[i
].m
));
85 tmass
+= atoms
->atom
[i
].m
;
91 static real
calc_geom(int isize
, const int *index
, rvec
*x
, rvec geom_center
, rvec minval
,
92 rvec maxval
, gmx_bool bDiam
)
97 clear_rvec(geom_center
);
114 for (j
= 0; j
< DIM
; j
++)
116 minval
[j
] = maxval
[j
] = x
[ii
][j
];
118 for (i
= 0; i
< isize
; i
++)
128 rvec_inc(geom_center
, x
[ii
]);
129 for (j
= 0; j
< DIM
; j
++)
131 if (x
[ii
][j
] < minval
[j
])
133 minval
[j
] = x
[ii
][j
];
135 if (x
[ii
][j
] > maxval
[j
])
137 maxval
[j
] = x
[ii
][j
];
144 for (j
= i
+ 1; j
< isize
; j
++)
146 d
= distance2(x
[ii
], x
[index
[j
]]);
147 diam2
= std::max(d
, diam2
);
152 for (j
= i
+ 1; j
< isize
; j
++)
154 d
= distance2(x
[i
], x
[j
]);
155 diam2
= std::max(d
, diam2
);
160 svmul(1.0 / isize
, geom_center
, geom_center
);
163 return std::sqrt(diam2
);
166 static void center_conf(int natom
, rvec
*x
, rvec center
, rvec geom_cent
)
171 rvec_sub(center
, geom_cent
, shift
);
173 printf(" shift :%7.3f%7.3f%7.3f (nm)\n", shift
[XX
], shift
[YY
],
176 for (i
= 0; (i
< natom
); i
++)
178 rvec_inc(x
[i
], shift
);
182 static void scale_conf(int natom
, rvec x
[], matrix box
, const rvec scale
)
186 for (i
= 0; i
< natom
; i
++)
188 for (j
= 0; j
< DIM
; j
++)
193 for (i
= 0; i
< DIM
; i
++)
195 for (j
= 0; j
< DIM
; j
++)
197 box
[i
][j
] *= scale
[j
];
202 static void read_bfac(const char *fn
, int *n_bfac
, double **bfac_val
, int **bfac_nr
)
207 *n_bfac
= get_lines(fn
, &bfac_lines
);
208 snew(*bfac_val
, *n_bfac
);
209 snew(*bfac_nr
, *n_bfac
);
210 fprintf(stderr
, "Reading %d B-factors from %s\n", *n_bfac
, fn
);
211 for (i
= 0; (i
< *n_bfac
); i
++)
213 sscanf(bfac_lines
[i
], "%d %lf", &(*bfac_nr
)[i
], &(*bfac_val
)[i
]);
218 static 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
;
225 if (n_bfac
> atoms
->nres
)
232 for (i
= 0; (i
< n_bfac
); i
++)
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
)
241 if (bfac
[i
] < bfac_min
)
246 while ((bfac_max
> 99.99) || (bfac_min
< -99.99))
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
++)
258 while ((std::abs(bfac_max
) < 0.5) && (std::abs(bfac_min
) < 0.5))
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
++)
271 for (i
= 0; (i
< natoms
); i
++)
273 atoms
->pdbinfo
[i
].bfac
= 0;
278 fprintf(stderr
, "Will attach %d B-factors to %d residues\n", n_bfac
,
280 for (i
= 0; (i
< n_bfac
); i
++)
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
];
293 gmx_warning("Residue nr %d not found\n", bfac_nr
[i
]);
299 fprintf(stderr
, "Will attach %d B-factors to %d atoms\n", n_bfac
,
301 for (i
= 0; (i
< n_bfac
); i
++)
303 atoms
->pdbinfo
[bfac_nr
[i
] - 1].bfac
= bfac
[i
];
308 static void pdb_legend(FILE *out
, int natoms
, int nres
, t_atoms
*atoms
, rvec x
[])
310 real bfac_min
, bfac_max
, xmin
, ymin
, zmin
;
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
++)
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 static void visualize_images(const char *fn
, int ePBC
, matrix box
)
346 init_t_atoms(&atoms
, nat
, FALSE
);
349 /* FIXME: Constness should not be cast away */
350 c
= const_cast<char*>("C");
351 ala
= const_cast<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
, nullptr, ePBC
, box
);
368 static void visualize_box(FILE *out
, int a0
, int r0
, matrix box
, const rvec gridsize
)
372 int nx
, ny
, nz
, nbox
, nat
;
376 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6,
383 nx
= gmx::roundToInt(gridsize
[XX
]);
384 ny
= gmx::roundToInt(gridsize
[YY
]);
385 nz
= gmx::roundToInt(gridsize
[ZZ
]);
389 nat
= nbox
* NCUCVERT
;
391 calc_compact_unitcell_vertices(ecenterDEF
, box
, vert
);
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
404 for (i
= 0; i
< NCUCVERT
; i
++)
406 rvec_add(vert
[i
], shift
, vert
[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]);
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, "");
446 for (i
= 0; i
< 24; i
+= 2)
448 fprintf(out
, "CONECT%5d%5d\n", a0
+ rectedge
[i
], a0
+ rectedge
[i
+ 1]);
453 static void calc_rotmatrix(rvec principal_axis
, rvec targetvec
, matrix rotmatrix
)
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
]);
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
,
490 int i
, resind_prev
, resind
;
493 for (i
= 0; i
< isize
; i
++)
495 resind
= atoms
->atom
[index
== nullptr ? i
: index
[i
]].resind
;
496 if (resind
!= resind_prev
)
498 atoms
->resinfo
[resind
].nr
= resnr_start
;
501 resind_prev
= resind
;
505 int gmx_editconf(int argc
, char *argv
[])
509 "[THISMODULE] converts generic structure format to [REF].gro[ref], [TT].g96[tt]",
510 "or [REF].pdb[ref].",
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 "The [TT]-center[tt] option can be used to shift the geometric center",
516 "of the system from the default of (x/2, y/2, z/2) implied by [TT]-c[tt]",
517 "to some other value.",
519 "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
520 "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
521 "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
522 "[TT]octahedron[tt] is a truncated octahedron.",
523 "The last two are special cases of a triclinic box.",
524 "The length of the three box vectors of the truncated octahedron is the",
525 "shortest distance between two opposite hexagons.",
526 "Relative to a cubic box with some periodic image distance, the volume of a ",
527 "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
528 "and that of a truncated octahedron is 0.77 times.",
530 "Option [TT]-box[tt] requires only",
531 "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
533 "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, [IT]y[it]-,",
534 "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
535 "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
536 "to the diameter of the system (largest distance between atoms) plus twice",
537 "the specified distance.",
539 "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
540 "a triclinic box and cannot be used with option [TT]-d[tt].",
542 "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
543 "can be selected for calculating the size and the geometric center,",
544 "otherwise the whole system is used.",
546 "[TT]-rotate[tt] rotates the coordinates and velocities.",
548 "[TT]-princ[tt] aligns the principal axes of the system along the",
549 "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
550 "This may allow you to decrease the box volume,",
551 "but beware that molecules can rotate significantly in a nanosecond.",
553 "Scaling is applied before any of the other operations are",
554 "performed. Boxes and coordinates can be scaled to give a certain density (option",
555 "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
556 "file is given as input. A special feature of the scaling option is that when the",
557 "factor -1 is given in one dimension, one obtains a mirror image,",
558 "mirrored in one of the planes. When one uses -1 in three dimensions, ",
559 "a point-mirror image is obtained.[PAR]",
560 "Groups are selected after all operations have been applied.[PAR]",
561 "Periodicity can be removed in a crude manner.",
562 "It is important that the box vectors at the bottom of your input file",
563 "are correct when the periodicity is to be removed.",
565 "When writing [REF].pdb[ref] files, B-factors can be",
566 "added with the [TT]-bf[tt] option. B-factors are read",
567 "from a file with with following format: first line states number of",
568 "entries in the file, next lines state an index",
569 "followed by a B-factor. The B-factors will be attached per residue",
570 "unless the number of B-factors is larger than the number of the residues or unless the",
571 "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
572 "be added instead of B-factors. [TT]-legend[tt] will produce",
573 "a row of CA atoms with B-factors ranging from the minimum to the",
574 "maximum value found, effectively making a legend for viewing.",
576 "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
577 "file for the MEAD electrostatics",
578 "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
579 "is that the input file is a run input file.",
580 "The B-factor field is then filled with the Van der Waals radius",
581 "of the atoms while the occupancy field will hold the charge.",
583 "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
584 "and the radius in the occupancy.",
586 "Option [TT]-align[tt] allows alignment",
587 "of the principal axis of a specified group against the given vector, ",
588 "with an optional center of rotation specified by [TT]-aligncenter[tt].",
590 "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
591 "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
593 "To convert a truncated octrahedron file produced by a package which uses",
594 "a cubic box with the corners cut off (such as GROMOS), use::",
596 " gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
598 "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
602 "For complex molecules, the periodicity removal routine may break down, "
603 "in that case you can use [gmx-trjconv]."
605 static real dist
= 0.0;
606 static gmx_bool bNDEF
= FALSE
, bRMPBC
= FALSE
, bCenter
= FALSE
, bReadVDW
=
607 FALSE
, bCONECT
= FALSE
;
608 static gmx_bool peratom
= FALSE
, bLegend
= FALSE
, bOrient
= FALSE
, bMead
=
609 FALSE
, bGrasp
= FALSE
, bSig56
= FALSE
;
611 { 1, 1, 1 }, newbox
=
612 { 0, 0, 0 }, newang
=
614 static real rho
= 1000.0, rvdw
= 0.12;
616 { 0, 0, 0 }, translation
=
617 { 0, 0, 0 }, rotangles
=
618 { 0, 0, 0 }, aligncenter
=
619 { 0, 0, 0 }, targetvec
=
621 static const char *btype
[] =
622 { nullptr, "triclinic", "cubic", "dodecahedron", "octahedron", nullptr },
626 static int resnr_start
= -1;
630 { "-ndef", FALSE
, etBOOL
,
631 { &bNDEF
}, "Choose output from default index groups" },
632 { "-visbox", FALSE
, etRVEC
,
634 "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
635 { "-bt", FALSE
, etENUM
,
636 { btype
}, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
637 { "-box", FALSE
, etRVEC
,
638 { newbox
}, "Box vector lengths (a,b,c)" },
639 { "-angles", FALSE
, etRVEC
,
640 { newang
}, "Angles between the box vectors (bc,ac,ab)" },
641 { "-d", FALSE
, etREAL
,
642 { &dist
}, "Distance between the solute and the box" },
643 { "-c", FALSE
, etBOOL
,
645 "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
646 { "-center", FALSE
, etRVEC
,
647 { center
}, "Shift the geometrical center to (x,y,z)" },
648 { "-aligncenter", FALSE
, etRVEC
,
649 { aligncenter
}, "Center of rotation for alignment" },
650 { "-align", FALSE
, etRVEC
,
652 "Align to target vector" },
653 { "-translate", FALSE
, etRVEC
,
654 { translation
}, "Translation" },
655 { "-rotate", FALSE
, etRVEC
,
657 "Rotation around the X, Y and Z axes in degrees" },
658 { "-princ", FALSE
, etBOOL
,
660 "Orient molecule(s) along their principal axes" },
661 { "-scale", FALSE
, etRVEC
,
662 { scale
}, "Scaling factor" },
663 { "-density", FALSE
, etREAL
,
665 "Density (g/L) of the output box achieved by scaling" },
666 { "-pbc", FALSE
, etBOOL
,
668 "Remove the periodicity (make molecule whole again)" },
669 { "-resnr", FALSE
, etINT
,
671 " Renumber residues starting from resnr" },
672 { "-grasp", FALSE
, etBOOL
,
674 "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" },
676 "-rvdw", FALSE
, etREAL
,
678 "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"
680 { "-sig56", FALSE
, etBOOL
,
682 "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
684 "-vdwread", FALSE
, etBOOL
,
686 "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing the radii based on the force field"
688 { "-atom", FALSE
, etBOOL
,
689 { &peratom
}, "Force B-factor attachment per atom" },
690 { "-legend", FALSE
, etBOOL
,
691 { &bLegend
}, "Make B-factor legend" },
692 { "-label", FALSE
, etSTR
,
693 { &label
}, "Add chain label for all residues" },
695 "-conect", FALSE
, etBOOL
,
697 "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a topology is present"
700 #define NPA asize(pa)
703 const char *infile
, *outfile
;
704 int outftp
, inftp
, natom
, i
, j
, n_bfac
, itype
, ntype
;
705 double *bfac
= nullptr, c6
, c12
;
706 int *bfac_nr
= nullptr;
707 t_topology
*top
= nullptr;
708 char *grpname
, *sgrpname
, *agrpname
;
709 int isize
, ssize
, numAlignmentAtoms
;
710 int *index
, *sindex
, *aindex
;
711 rvec
*x
, *v
, gc
, rmin
, rmax
, size
;
713 matrix box
, rotmatrix
, trans
;
715 gmx_bool bIndex
, bSetSize
, bSetAng
, bDist
, bSetCenter
, bAlign
;
716 gmx_bool bHaveV
, bScale
, bRho
, bTranslate
, bRotate
, bCalcGeom
, bCalcDiam
;
717 real diam
= 0, mass
= 0, d
, vdw
;
719 gmx_output_env_t
*oenv
;
722 { efSTX
, "-f", nullptr, ffREAD
},
723 { efNDX
, "-n", nullptr, ffOPTRD
},
724 { efSTO
, nullptr, nullptr, ffOPTWR
},
725 { efPQR
, "-mead", "mead", ffOPTWR
},
726 { efDAT
, "-bf", "bfact", ffOPTRD
}
728 #define NFILE asize(fnm)
730 if (!parse_common_args(&argc
, argv
, PCA_CAN_VIEW
, NFILE
, fnm
, NPA
, pa
,
731 asize(desc
), desc
, asize(bugs
), bugs
, &oenv
))
735 fprintf(stdout
, "Note that major changes are planned in future for "
736 "editconf, to improve usability and utility.\n");
738 bIndex
= opt2bSet("-n", NFILE
, fnm
) || bNDEF
;
739 bMead
= opt2bSet("-mead", NFILE
, fnm
);
740 bSetSize
= opt2parg_bSet("-box", NPA
, pa
);
741 bSetAng
= opt2parg_bSet("-angles", NPA
, pa
);
742 bSetCenter
= opt2parg_bSet("-center", NPA
, pa
);
743 bDist
= opt2parg_bSet("-d", NPA
, pa
);
744 bAlign
= opt2parg_bSet("-align", NPA
, pa
);
745 /* Only automatically turn on centering without -noc */
746 if ((bDist
|| bSetSize
|| bSetCenter
) && !opt2parg_bSet("-c", NPA
, pa
))
750 bScale
= opt2parg_bSet("-scale", NPA
, pa
);
751 bRho
= opt2parg_bSet("-density", NPA
, pa
);
752 bTranslate
= opt2parg_bSet("-translate", NPA
, pa
);
753 bRotate
= opt2parg_bSet("-rotate", NPA
, pa
);
756 fprintf(stderr
, "WARNING: setting -density overrides -scale\n");
758 bScale
= bScale
|| bRho
;
759 bCalcGeom
= bCenter
|| bRotate
|| bOrient
|| bScale
;
761 GMX_RELEASE_ASSERT(btype
[0] != nullptr, "Option setting inconsistency; btype[0] is NULL");
763 bCalcDiam
= (btype
[0][0] == 'c' || btype
[0][0] == 'd' || btype
[0][0] == 'o');
765 infile
= ftp2fn(efSTX
, NFILE
, fnm
);
768 outfile
= ftp2fn(efPQR
, NFILE
, fnm
);
772 outfile
= ftp2fn(efSTO
, NFILE
, fnm
);
774 outftp
= fn2ftp(outfile
);
775 inftp
= fn2ftp(infile
);
781 printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
784 if (bGrasp
&& (outftp
!= efPDB
))
786 gmx_fatal(FARGS
, "Output file should be a .pdb file"
787 " when using the -grasp option\n");
789 if ((bMead
|| bGrasp
) && (fn2ftp(infile
) != efTPR
))
791 gmx_fatal(FARGS
, "Input file should be a .tpr file"
792 " when using the -mead option\n");
798 open_symtab(&symtab
);
799 readConfAndAtoms(infile
, &symtab
, &name
, &atoms
, &ePBC
, &x
, &v
, box
);
801 if (atoms
.pdbinfo
== nullptr)
803 snew(atoms
.pdbinfo
, atoms
.nr
);
805 atoms
.havePdbInfo
= TRUE
;
807 if (fn2ftp(infile
) == efPDB
)
809 get_pdb_atomnumber(&atoms
, &aps
);
811 printf("Read %d atoms\n", atoms
.nr
);
813 /* Get the element numbers if available in a pdb file */
814 if (fn2ftp(infile
) == efPDB
)
816 get_pdb_atomnumber(&atoms
, &aps
);
819 if (ePBC
!= epbcNONE
)
822 printf("Volume: %g nm^3, corresponds to roughly %d electrons\n",
823 vol
, 100*(static_cast<int>(vol
*4.5)));
826 if (bMead
|| bGrasp
|| bCONECT
)
828 top
= read_top(infile
, nullptr);
833 if (atoms
.nr
!= top
->atoms
.nr
)
835 gmx_fatal(FARGS
, "Atom numbers don't match (%d vs. %d)", atoms
.nr
, top
->atoms
.nr
);
837 snew(atoms
.pdbinfo
, top
->atoms
.nr
);
838 ntype
= top
->idef
.atnr
;
839 for (i
= 0; (i
< atoms
.nr
); i
++)
841 /* Determine the Van der Waals radius from the force field */
844 if (!aps
.setAtomProperty(epropVDW
,
845 *top
->atoms
.resinfo
[top
->atoms
.atom
[i
].resind
].name
,
846 *top
->atoms
.atomname
[i
], &vdw
))
853 itype
= top
->atoms
.atom
[i
].type
;
854 c12
= top
->idef
.iparams
[itype
*ntype
+itype
].lj
.c12
;
855 c6
= top
->idef
.iparams
[itype
*ntype
+itype
].lj
.c6
;
856 if ((c6
!= 0) && (c12
!= 0))
867 vdw
= 0.5*gmx::sixthroot(sig6
);
874 /* Factor of 10 for nm -> Angstroms */
879 atoms
.pdbinfo
[i
].occup
= top
->atoms
.atom
[i
].q
;
880 atoms
.pdbinfo
[i
].bfac
= vdw
;
884 atoms
.pdbinfo
[i
].occup
= vdw
;
885 atoms
.pdbinfo
[i
].bfac
= top
->atoms
.atom
[i
].q
;
890 for (i
= 0; (i
< natom
) && !bHaveV
; i
++)
892 for (j
= 0; (j
< DIM
) && !bHaveV
; j
++)
894 bHaveV
= bHaveV
|| (v
[i
][j
] != 0);
897 printf("%selocities found\n", bHaveV
? "V" : "No v");
903 gmx_fatal(FARGS
, "Sorry, can not visualize box with index groups");
907 gmx_fatal(FARGS
, "Sorry, can only visualize box with a pdb file");
910 else if (visbox
[0] == -1)
912 visualize_images("images.pdb", ePBC
, box
);
918 rm_gropbc(&atoms
, x
, box
);
925 fprintf(stderr
, "\nSelect a group for determining the system size:\n");
926 get_index(&atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
927 1, &ssize
, &sindex
, &sgrpname
);
934 diam
= calc_geom(ssize
, sindex
, x
, gc
, rmin
, rmax
, bCalcDiam
);
935 rvec_sub(rmax
, rmin
, size
);
936 printf(" system size :%7.3f%7.3f%7.3f (nm)\n",
937 size
[XX
], size
[YY
], size
[ZZ
]);
940 printf(" diameter :%7.3f (nm)\n", diam
);
942 printf(" center :%7.3f%7.3f%7.3f (nm)\n", gc
[XX
], gc
[YY
], gc
[ZZ
]);
943 printf(" box vectors :%7.3f%7.3f%7.3f (nm)\n",
944 norm(box
[XX
]), norm(box
[YY
]), norm(box
[ZZ
]));
945 printf(" box angles :%7.2f%7.2f%7.2f (degrees)\n",
946 norm2(box
[ZZ
]) == 0 ? 0 :
947 RAD2DEG
*gmx_angle(box
[YY
], box
[ZZ
]),
948 norm2(box
[ZZ
]) == 0 ? 0 :
949 RAD2DEG
*gmx_angle(box
[XX
], box
[ZZ
]),
950 norm2(box
[YY
]) == 0 ? 0 :
951 RAD2DEG
*gmx_angle(box
[XX
], box
[YY
]));
952 printf(" box volume :%7.2f (nm^3)\n", det(box
));
955 if (bRho
|| bOrient
|| bAlign
)
957 mass
= calc_mass(&atoms
, !fn2bTPX(infile
), &aps
);
965 /* Get a group for principal component analysis */
966 fprintf(stderr
, "\nSelect group for the determining the orientation\n");
967 get_index(&atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
), 1, &isize
, &index
, &grpnames
);
969 /* Orient the principal axes along the coordinate axes */
970 orient_princ(&atoms
, isize
, index
, natom
, x
, bHaveV
? v
: nullptr, nullptr);
977 /* scale coordinates and box */
980 /* Compute scaling constant */
984 dens
= (mass
*AMU
)/(vol
*NANO
*NANO
*NANO
);
985 fprintf(stderr
, "Volume of input %g (nm^3)\n", vol
);
986 fprintf(stderr
, "Mass of input %g (a.m.u.)\n", mass
);
987 fprintf(stderr
, "Density of input %g (g/l)\n", dens
);
988 if (vol
== 0 || mass
== 0)
990 gmx_fatal(FARGS
, "Cannot scale density with "
991 "zero mass (%g) or volume (%g)\n", mass
, vol
);
994 scale
[XX
] = scale
[YY
] = scale
[ZZ
] = std::cbrt(dens
/rho
);
995 fprintf(stderr
, "Scaling all box vectors by %g\n", scale
[XX
]);
997 scale_conf(atoms
.nr
, x
, box
, scale
);
1004 fprintf(stderr
, "\nSelect a group that you want to align:\n");
1005 get_index(&atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
1006 1, &numAlignmentAtoms
, &aindex
, &agrpname
);
1010 numAlignmentAtoms
= atoms
.nr
;
1011 snew(aindex
, numAlignmentAtoms
);
1012 for (i
= 0; i
< numAlignmentAtoms
; i
++)
1017 printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n", numAlignmentAtoms
, natom
,
1018 targetvec
[XX
], targetvec
[YY
], targetvec
[ZZ
],
1019 aligncenter
[XX
], aligncenter
[YY
], aligncenter
[ZZ
]);
1020 /*subtract out pivot point*/
1021 for (i
= 0; i
< numAlignmentAtoms
; i
++)
1023 rvec_dec(x
[aindex
[i
]], aligncenter
);
1025 /*now determine transform and rotate*/
1027 principal_comp(numAlignmentAtoms
, aindex
, atoms
.atom
, x
, trans
, princd
);
1029 unitv(targetvec
, targetvec
);
1030 printf("Using %g %g %g as principal axis\n", trans
[0][2], trans
[1][2], trans
[2][2]);
1031 tmpvec
[XX
] = trans
[0][2]; tmpvec
[YY
] = trans
[1][2]; tmpvec
[ZZ
] = trans
[2][2];
1032 calc_rotmatrix(tmpvec
, targetvec
, rotmatrix
);
1033 /* rotmatrix finished */
1035 for (i
= 0; i
< numAlignmentAtoms
; ++i
)
1037 mvmul(rotmatrix
, x
[aindex
[i
]], tmpvec
);
1038 copy_rvec(tmpvec
, x
[aindex
[i
]]);
1041 /*add pivot point back*/
1042 for (i
= 0; i
< numAlignmentAtoms
; i
++)
1044 rvec_inc(x
[aindex
[i
]], aligncenter
);
1056 fprintf(stderr
, "\nSelect a group that you want to translate:\n");
1057 get_index(&atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
),
1058 1, &ssize
, &sindex
, &sgrpname
);
1065 printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize
, natom
,
1066 translation
[XX
], translation
[YY
], translation
[ZZ
]);
1069 for (i
= 0; i
< ssize
; i
++)
1071 rvec_inc(x
[sindex
[i
]], translation
);
1076 for (i
= 0; i
< natom
; i
++)
1078 rvec_inc(x
[i
], translation
);
1085 printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n", rotangles
[XX
], rotangles
[YY
], rotangles
[ZZ
]);
1086 for (i
= 0; i
< DIM
; i
++)
1088 rotangles
[i
] *= DEG2RAD
;
1090 rotate_conf(natom
, x
, v
, rotangles
[XX
], rotangles
[YY
], rotangles
[ZZ
]);
1095 /* recalc geometrical center and max and min coordinates and size */
1096 calc_geom(ssize
, sindex
, x
, gc
, rmin
, rmax
, FALSE
);
1097 rvec_sub(rmax
, rmin
, size
);
1098 if (bScale
|| bOrient
|| bRotate
)
1100 printf("new system size : %6.3f %6.3f %6.3f\n",
1101 size
[XX
], size
[YY
], size
[ZZ
]);
1105 if ((btype
[0] != nullptr) && (bSetSize
|| bDist
|| (btype
[0][0] == 't' && bSetAng
)))
1108 if (!(bSetSize
|| bDist
))
1110 for (i
= 0; i
< DIM
; i
++)
1112 newbox
[i
] = norm(box
[i
]);
1116 /* calculate new boxsize */
1117 switch (btype
[0][0])
1122 for (i
= 0; i
< DIM
; i
++)
1124 newbox
[i
] = size
[i
]+2*dist
;
1129 box
[XX
][XX
] = newbox
[XX
];
1130 box
[YY
][YY
] = newbox
[YY
];
1131 box
[ZZ
][ZZ
] = newbox
[ZZ
];
1135 matrix_convert(box
, newbox
, newang
);
1149 if (btype
[0][0] == 'c')
1151 for (i
= 0; i
< DIM
; i
++)
1156 else if (btype
[0][0] == 'd')
1162 box
[ZZ
][ZZ
] = d
*std::sqrt(2.0)/2.0;
1168 box
[YY
][YY
] = d
*std::sqrt(2.0)*2.0/3.0;
1170 box
[ZZ
][YY
] = d
*std::sqrt(2.0)/3.0;
1171 box
[ZZ
][ZZ
] = d
*std::sqrt(6.0)/3.0;
1177 /* calculate new coords for geometrical center */
1180 calc_box_center(ecenterDEF
, box
, center
);
1183 /* center molecule on 'center' */
1186 center_conf(natom
, x
, center
, gc
);
1192 calc_geom(ssize
, sindex
, x
, gc
, rmin
, rmax
, FALSE
);
1193 printf("new center :%7.3f%7.3f%7.3f (nm)\n", gc
[XX
], gc
[YY
], gc
[ZZ
]);
1195 if (bOrient
|| bScale
|| bDist
|| bSetSize
)
1197 printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n",
1198 norm(box
[XX
]), norm(box
[YY
]), norm(box
[ZZ
]));
1199 printf("new box angles :%7.2f%7.2f%7.2f (degrees)\n",
1200 norm2(box
[ZZ
]) == 0 ? 0 :
1201 RAD2DEG
*gmx_angle(box
[YY
], box
[ZZ
]),
1202 norm2(box
[ZZ
]) == 0 ? 0 :
1203 RAD2DEG
*gmx_angle(box
[XX
], box
[ZZ
]),
1204 norm2(box
[YY
]) == 0 ? 0 :
1205 RAD2DEG
*gmx_angle(box
[XX
], box
[YY
]));
1206 printf("new box volume :%7.2f (nm^3)\n", det(box
));
1209 if (check_box(epbcXYZ
, box
))
1211 printf("\nWARNING: %s\n"
1212 "See the GROMACS manual for a description of the requirements that\n"
1213 "must be satisfied by descriptions of simulation cells.\n",
1214 check_box(epbcXYZ
, box
));
1217 if (bDist
&& btype
[0][0] == 't')
1221 printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1222 "distance from the solute to a box surface along the corresponding normal\n"
1223 "vector might be somewhat smaller than your specified value %f.\n"
1224 "You can check the actual value with g_mindist -pi\n", dist
);
1226 else if (!opt2parg_bSet("-bt", NPA
, pa
))
1228 printf("\nWARNING: No boxtype specified - distance condition applied in each dimension.\n"
1229 "If the molecule rotates the actual distance will be smaller. You might want\n"
1230 "to use a cubic box instead, or why not try a dodecahedron today?\n");
1233 if (bCONECT
&& (outftp
== efPDB
) && (inftp
== efTPR
))
1235 conect
= gmx_conect_generate(top
);
1244 fprintf(stderr
, "\nSelect a group for output:\n");
1245 get_index(&atoms
, opt2fn_null("-n", NFILE
, fnm
),
1246 1, &isize
, &index
, &grpname
);
1248 if (resnr_start
>= 0)
1250 renum_resnr(&atoms
, isize
, index
, resnr_start
);
1253 if (opt2parg_bSet("-label", NPA
, pa
))
1255 for (i
= 0; (i
< atoms
.nr
); i
++)
1257 atoms
.resinfo
[atoms
.atom
[i
].resind
].chainid
= label
[0];
1261 if (opt2bSet("-bf", NFILE
, fnm
) || bLegend
)
1263 gmx_fatal(FARGS
, "Sorry, cannot do bfactors with an index group.");
1266 if (outftp
== efPDB
)
1268 out
= gmx_ffopen(outfile
, "w");
1269 write_pdbfile_indexed(out
, name
, &atoms
, x
, ePBC
, box
, ' ', 1, isize
, index
, conect
, FALSE
);
1274 write_sto_conf_indexed(outfile
, name
, &atoms
, x
, bHaveV
? v
: nullptr, ePBC
, box
, isize
, index
);
1281 if (resnr_start
>= 0)
1283 renum_resnr(&atoms
, atoms
.nr
, nullptr, resnr_start
);
1286 if ((outftp
== efPDB
) || (outftp
== efPQR
))
1288 out
= gmx_ffopen(outfile
, "w");
1291 fprintf(out
, "REMARK "
1292 "The B-factors in this file hold atomic radii\n");
1293 fprintf(out
, "REMARK "
1294 "The occupancy in this file hold atomic charges\n");
1298 fprintf(out
, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1299 fprintf(out
, "REMARK "
1300 "The B-factors in this file hold atomic charges\n");
1301 fprintf(out
, "REMARK "
1302 "The occupancy in this file hold atomic radii\n");
1304 else if (opt2bSet("-bf", NFILE
, fnm
))
1306 read_bfac(opt2fn("-bf", NFILE
, fnm
), &n_bfac
, &bfac
, &bfac_nr
);
1307 set_pdb_conf_bfac(atoms
.nr
, atoms
.nres
, &atoms
,
1308 n_bfac
, bfac
, bfac_nr
, peratom
);
1310 if (opt2parg_bSet("-label", NPA
, pa
))
1312 for (i
= 0; (i
< atoms
.nr
); i
++)
1314 atoms
.resinfo
[atoms
.atom
[i
].resind
].chainid
= label
[0];
1317 /* Need to bypass the regular write_pdbfile because I don't want to change
1318 * all instances to include the boolean flag for writing out PQR files.
1321 snew(index
, atoms
.nr
);
1322 for (int i
= 0; i
< atoms
.nr
; i
++)
1326 write_pdbfile_indexed(out
, name
, &atoms
, x
, ePBC
, box
, ' ', -1, atoms
.nr
, index
, conect
,
1331 pdb_legend(out
, atoms
.nr
, atoms
.nres
, &atoms
, x
);
1335 visualize_box(out
, bLegend
? atoms
.nr
+12 : atoms
.nr
,
1336 bLegend
? atoms
.nres
= 12 : atoms
.nres
, box
, visbox
);
1342 write_sto_conf(outfile
, name
, &atoms
, x
, bHaveV
? v
: nullptr, ePBC
, box
);
1346 done_symtab(&symtab
);
1356 do_view(oenv
, outfile
, nullptr);
1357 output_env_done(oenv
);