Remove ctest section called simple
[gromacs.git] / src / gromacs / gmxana / gmx_editconf.cpp
blob933b718c65d581cddca883a5068a7226cba65acf
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,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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstring>
42 #include <algorithm>
43 #include <string>
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/pdbio.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/math/functions.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"
68 #include "gromacs/utility/strdb.h"
70 static real calc_mass(t_atoms *atoms, gmx_bool bGetMass, AtomProperties *aps)
72 real tmass;
73 int i;
75 tmass = 0;
76 for (i = 0; (i < atoms->nr); i++)
78 if (bGetMass)
80 aps->setAtomProperty(epropMass,
81 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;
87 return tmass;
90 static real calc_geom(int isize, const int *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 static 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 static void scale_conf(int natom, rvec x[], matrix box, const 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 static 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 static 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 if (n_bfac > atoms->nres)
228 peratom = TRUE;
231 bfac_max = -1e10;
232 bfac_min = 1e10;
233 for (i = 0; (i < n_bfac); i++)
235 /* if ((bfac_nr[i]-1<0) || (bfac_nr[i]-1>=atoms->nr))
236 gmx_fatal(FARGS,"Index of B-Factor %d is out of range: %d (%g)",
237 i+1,bfac_nr[i],bfac[i]); */
238 if (bfac[i] > bfac_max)
240 bfac_max = bfac[i];
242 if (bfac[i] < bfac_min)
244 bfac_min = bfac[i];
247 while ((bfac_max > 99.99) || (bfac_min < -99.99))
249 fprintf(stderr,
250 "Range of values for B-factors too large (min %g, max %g) "
251 "will scale down a factor 10\n", bfac_min, bfac_max);
252 for (i = 0; (i < n_bfac); i++)
254 bfac[i] /= 10;
256 bfac_max /= 10;
257 bfac_min /= 10;
259 while ((std::abs(bfac_max) < 0.5) && (std::abs(bfac_min) < 0.5))
261 fprintf(stderr,
262 "Range of values for B-factors too small (min %g, max %g) "
263 "will scale up a factor 10\n", bfac_min, bfac_max);
264 for (i = 0; (i < n_bfac); i++)
266 bfac[i] *= 10;
268 bfac_max *= 10;
269 bfac_min *= 10;
272 for (i = 0; (i < natoms); i++)
274 atoms->pdbinfo[i].bfac = 0;
277 if (!peratom)
279 fprintf(stderr, "Will attach %d B-factors to %d residues\n", n_bfac,
280 nres);
281 for (i = 0; (i < n_bfac); i++)
283 found = FALSE;
284 for (n = 0; (n < natoms); n++)
286 if (bfac_nr[i] == atoms->resinfo[atoms->atom[n].resind].nr)
288 atoms->pdbinfo[n].bfac = bfac[i];
289 found = TRUE;
292 if (!found)
294 gmx_warning("Residue nr %d not found\n", bfac_nr[i]);
298 else
300 fprintf(stderr, "Will attach %d B-factors to %d atoms\n", n_bfac,
301 natoms);
302 for (i = 0; (i < n_bfac); i++)
304 atoms->pdbinfo[bfac_nr[i] - 1].bfac = bfac[i];
309 static void pdb_legend(FILE *out, int natoms, int nres, t_atoms *atoms, rvec x[])
311 real bfac_min, bfac_max, xmin, ymin, zmin;
312 int i;
313 int space = ' ';
315 bfac_max = -1e10;
316 bfac_min = 1e10;
317 xmin = 1e10;
318 ymin = 1e10;
319 zmin = 1e10;
320 for (i = 0; (i < natoms); i++)
322 xmin = std::min(xmin, x[i][XX]);
323 ymin = std::min(ymin, x[i][YY]);
324 zmin = std::min(zmin, x[i][ZZ]);
325 bfac_min = std::min(bfac_min, atoms->pdbinfo[i].bfac);
326 bfac_max = std::max(bfac_max, atoms->pdbinfo[i].bfac);
328 fprintf(stderr, "B-factors range from %g to %g\n", bfac_min, bfac_max);
329 for (i = 1; (i < 12); i++)
331 fprintf(out,
332 "%-6s%5d %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f\n",
333 "ATOM ", natoms + 1 + i, "CA", "LEG", space, nres + 1, space,
334 (xmin + (i * 0.12)) * 10, ymin * 10, zmin * 10, 1.0, bfac_min
335 + ((i - 1.0) * (bfac_max - bfac_min) / 10));
339 static void visualize_images(const char *fn, int ePBC, matrix box)
341 t_atoms atoms;
342 rvec *img;
343 char *c, *ala;
344 int nat, i;
346 nat = NTRICIMG + 1;
347 init_t_atoms(&atoms, nat, FALSE);
348 atoms.nr = nat;
349 snew(img, nat);
350 /* FIXME: Constness should not be cast away */
351 c = const_cast<char*>("C");
352 ala = const_cast<char*>("ALA");
353 for (i = 0; i < nat; i++)
355 atoms.atomname[i] = &c;
356 atoms.atom[i].resind = i;
357 atoms.resinfo[i].name = &ala;
358 atoms.resinfo[i].nr = i + 1;
359 atoms.resinfo[i].chainid = 'A' + i / NCUCVERT;
361 calc_triclinic_images(box, img + 1);
363 write_sto_conf(fn, "Images", &atoms, img, nullptr, ePBC, box);
365 done_atom(&atoms);
366 sfree(img);
369 static void visualize_box(FILE *out, int a0, int r0, matrix box, const rvec gridsize)
371 int *edge;
372 rvec *vert, shift;
373 int nx, ny, nz, nbox, nat;
374 int i, j, x, y, z;
375 int rectedge[24] =
377 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6,
381 a0++;
382 r0++;
384 nx = gmx::roundToInt(gridsize[XX]);
385 ny = gmx::roundToInt(gridsize[YY]);
386 nz = gmx::roundToInt(gridsize[ZZ]);
387 nbox = nx * ny * nz;
388 if (TRICLINIC(box))
390 nat = nbox * NCUCVERT;
391 snew(vert, nat);
392 calc_compact_unitcell_vertices(ecenterDEF, box, vert);
393 j = 0;
394 for (z = 0; z < nz; z++)
396 for (y = 0; y < ny; y++)
398 for (x = 0; x < nx; x++)
400 for (i = 0; i < DIM; i++)
402 shift[i] = x * box[0][i] + y * box[1][i] + z
403 * box[2][i];
405 for (i = 0; i < NCUCVERT; i++)
407 rvec_add(vert[i], shift, vert[j]);
408 j++;
414 for (i = 0; i < nat; i++)
416 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i / NCUCVERT, r0 + i, ' ',
417 10*vert[i][XX], 10*vert[i][YY], 10*vert[i][ZZ], 1.0, 0.0, "");
420 edge = compact_unitcell_edges();
421 for (j = 0; j < nbox; j++)
423 for (i = 0; i < NCUCEDGE; i++)
425 fprintf(out, "CONECT%5d%5d\n", a0 + j * NCUCVERT + edge[2 * i],
426 a0 + j * NCUCVERT + edge[2 * i + 1]);
430 sfree(vert);
432 else
434 i = 0;
435 for (z = 0; z <= 1; z++)
437 for (y = 0; y <= 1; y++)
439 for (x = 0; x <= 1; x++)
441 gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i/8, r0+i, ' ',
442 x * 10 * box[XX][XX], y * 10 * box[YY][YY], z * 10 * box[ZZ][ZZ], 1.0, 0.0, "");
443 i++;
447 for (i = 0; i < 24; i += 2)
449 fprintf(out, "CONECT%5d%5d\n", a0 + rectedge[i], a0 + rectedge[i + 1]);
454 static void calc_rotmatrix(rvec principal_axis, rvec targetvec, matrix rotmatrix)
456 rvec rotvec;
457 real ux, uy, uz, costheta, sintheta;
459 costheta = cos_angle(principal_axis, targetvec);
460 sintheta = std::sqrt(1.0-costheta*costheta); /* sign is always positive since 0<theta<pi */
462 /* Determine rotation from cross product with target vector */
463 cprod(principal_axis, targetvec, rotvec);
464 unitv(rotvec, rotvec);
465 printf("Aligning %g %g %g to %g %g %g : xprod %g %g %g\n",
466 principal_axis[XX], principal_axis[YY], principal_axis[ZZ], targetvec[XX], targetvec[YY], targetvec[ZZ],
467 rotvec[XX], rotvec[YY], rotvec[ZZ]);
469 ux = rotvec[XX];
470 uy = rotvec[YY];
471 uz = rotvec[ZZ];
472 rotmatrix[0][0] = ux*ux + (1.0-ux*ux)*costheta;
473 rotmatrix[0][1] = ux*uy*(1-costheta)-uz*sintheta;
474 rotmatrix[0][2] = ux*uz*(1-costheta)+uy*sintheta;
475 rotmatrix[1][0] = ux*uy*(1-costheta)+uz*sintheta;
476 rotmatrix[1][1] = uy*uy + (1.0-uy*uy)*costheta;
477 rotmatrix[1][2] = uy*uz*(1-costheta)-ux*sintheta;
478 rotmatrix[2][0] = ux*uz*(1-costheta)-uy*sintheta;
479 rotmatrix[2][1] = uy*uz*(1-costheta)+ux*sintheta;
480 rotmatrix[2][2] = uz*uz + (1.0-uz*uz)*costheta;
482 printf("Rotation matrix: \n%g %g %g\n%g %g %g\n%g %g %g\n",
483 rotmatrix[0][0], rotmatrix[0][1], rotmatrix[0][2],
484 rotmatrix[1][0], rotmatrix[1][1], rotmatrix[1][2],
485 rotmatrix[2][0], rotmatrix[2][1], rotmatrix[2][2]);
488 static void renum_resnr(t_atoms *atoms, int isize, const int *index,
489 int resnr_start)
491 int i, resind_prev, resind;
493 resind_prev = -1;
494 for (i = 0; i < isize; i++)
496 resind = atoms->atom[index == nullptr ? i : index[i]].resind;
497 if (resind != resind_prev)
499 atoms->resinfo[resind].nr = resnr_start;
500 resnr_start++;
502 resind_prev = resind;
506 int gmx_editconf(int argc, char *argv[])
508 const char *desc[] =
510 "[THISMODULE] converts generic structure format to [REF].gro[ref], [TT].g96[tt]",
511 "or [REF].pdb[ref].",
512 "[PAR]",
513 "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and",
514 "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]",
515 "will center the system in the box, unless [TT]-noc[tt] is used.",
516 "The [TT]-center[tt] option can be used to shift the geometric center",
517 "of the system from the default of (x/2, y/2, z/2) implied by [TT]-c[tt]",
518 "to some other value.",
519 "[PAR]",
520 "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
521 "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
522 "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
523 "[TT]octahedron[tt] is a truncated octahedron.",
524 "The last two are special cases of a triclinic box.",
525 "The length of the three box vectors of the truncated octahedron is the",
526 "shortest distance between two opposite hexagons.",
527 "Relative to a cubic box with some periodic image distance, the volume of a ",
528 "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
529 "and that of a truncated octahedron is 0.77 times.",
530 "[PAR]",
531 "Option [TT]-box[tt] requires only",
532 "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
533 "[PAR]",
534 "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, [IT]y[it]-,",
535 "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
536 "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
537 "to the diameter of the system (largest distance between atoms) plus twice",
538 "the specified distance.",
539 "[PAR]",
540 "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
541 "a triclinic box and cannot be used with option [TT]-d[tt].",
542 "[PAR]",
543 "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
544 "can be selected for calculating the size and the geometric center,",
545 "otherwise the whole system is used.",
546 "[PAR]",
547 "[TT]-rotate[tt] rotates the coordinates and velocities.",
548 "[PAR]",
549 "[TT]-princ[tt] aligns the principal axes of the system along the",
550 "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
551 "This may allow you to decrease the box volume,",
552 "but beware that molecules can rotate significantly in a nanosecond.",
553 "[PAR]",
554 "Scaling is applied before any of the other operations are",
555 "performed. Boxes and coordinates can be scaled to give a certain density (option",
556 "[TT]-density[tt]). Note that this may be inaccurate in case a [REF].gro[ref]",
557 "file is given as input. A special feature of the scaling option is that when the",
558 "factor -1 is given in one dimension, one obtains a mirror image,",
559 "mirrored in one of the planes. When one uses -1 in three dimensions, ",
560 "a point-mirror image is obtained.[PAR]",
561 "Groups are selected after all operations have been applied.[PAR]",
562 "Periodicity can be removed in a crude manner.",
563 "It is important that the box vectors at the bottom of your input file",
564 "are correct when the periodicity is to be removed.",
565 "[PAR]",
566 "When writing [REF].pdb[ref] files, B-factors can be",
567 "added with the [TT]-bf[tt] option. B-factors are read",
568 "from a file with with following format: first line states number of",
569 "entries in the file, next lines state an index",
570 "followed by a B-factor. The B-factors will be attached per residue",
571 "unless the number of B-factors is larger than the number of the residues or unless the",
572 "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
573 "be added instead of B-factors. [TT]-legend[tt] will produce",
574 "a row of CA atoms with B-factors ranging from the minimum to the",
575 "maximum value found, effectively making a legend for viewing.",
576 "[PAR]",
577 "With the option [TT]-mead[tt] a special [REF].pdb[ref] ([REF].pqr[ref])",
578 "file for the MEAD electrostatics",
579 "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
580 "is that the input file is a run input file.",
581 "The B-factor field is then filled with the Van der Waals radius",
582 "of the atoms while the occupancy field will hold the charge.",
583 "[PAR]",
584 "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
585 "and the radius in the occupancy.",
586 "[PAR]",
587 "Option [TT]-align[tt] allows alignment",
588 "of the principal axis of a specified group against the given vector, ",
589 "with an optional center of rotation specified by [TT]-aligncenter[tt].",
590 "[PAR]",
591 "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
592 "to a [REF].pdb[ref] file, which can be useful for analysis with e.g. Rasmol.",
593 "[PAR]",
594 "To convert a truncated octrahedron file produced by a package which uses",
595 "a cubic box with the corners cut off (such as GROMOS), use::",
597 " gmx editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out",
599 "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
601 const char *bugs[] =
603 "For complex molecules, the periodicity removal routine may break down, "
604 "in that case you can use [gmx-trjconv]."
606 static real dist = 0.0;
607 static gmx_bool bNDEF = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW =
608 FALSE, bCONECT = FALSE;
609 static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead =
610 FALSE, bGrasp = FALSE, bSig56 = FALSE;
611 static rvec scale =
612 { 1, 1, 1 }, newbox =
613 { 0, 0, 0 }, newang =
614 { 90, 90, 90 };
615 static real rho = 1000.0, rvdw = 0.12;
616 static rvec center =
617 { 0, 0, 0 }, translation =
618 { 0, 0, 0 }, rotangles =
619 { 0, 0, 0 }, aligncenter =
620 { 0, 0, 0 }, targetvec =
621 { 0, 0, 0 };
622 static const char *btype[] =
623 { nullptr, "triclinic", "cubic", "dodecahedron", "octahedron", nullptr },
624 *label = "A";
625 static rvec visbox =
626 { 0, 0, 0 };
627 static int resnr_start = -1;
628 t_pargs
629 pa[] =
631 { "-ndef", FALSE, etBOOL,
632 { &bNDEF }, "Choose output from default index groups" },
633 { "-visbox", FALSE, etRVEC,
634 { visbox },
635 "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
636 { "-bt", FALSE, etENUM,
637 { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
638 { "-box", FALSE, etRVEC,
639 { newbox }, "Box vector lengths (a,b,c)" },
640 { "-angles", FALSE, etRVEC,
641 { newang }, "Angles between the box vectors (bc,ac,ab)" },
642 { "-d", FALSE, etREAL,
643 { &dist }, "Distance between the solute and the box" },
644 { "-c", FALSE, etBOOL,
645 { &bCenter },
646 "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
647 { "-center", FALSE, etRVEC,
648 { center }, "Shift the geometrical center to (x,y,z)" },
649 { "-aligncenter", FALSE, etRVEC,
650 { aligncenter }, "Center of rotation for alignment" },
651 { "-align", FALSE, etRVEC,
652 { targetvec },
653 "Align to target vector" },
654 { "-translate", FALSE, etRVEC,
655 { translation }, "Translation" },
656 { "-rotate", FALSE, etRVEC,
657 { rotangles },
658 "Rotation around the X, Y and Z axes in degrees" },
659 { "-princ", FALSE, etBOOL,
660 { &bOrient },
661 "Orient molecule(s) along their principal axes" },
662 { "-scale", FALSE, etRVEC,
663 { scale }, "Scaling factor" },
664 { "-density", FALSE, etREAL,
665 { &rho },
666 "Density (g/L) of the output box achieved by scaling" },
667 { "-pbc", FALSE, etBOOL,
668 { &bRMPBC },
669 "Remove the periodicity (make molecule whole again)" },
670 { "-resnr", FALSE, etINT,
671 { &resnr_start },
672 " Renumber residues starting from resnr" },
673 { "-grasp", FALSE, etBOOL,
674 { &bGrasp },
675 "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" },
677 "-rvdw", FALSE, etREAL,
678 { &rvdw },
679 "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"
681 { "-sig56", FALSE, etBOOL,
682 { &bSig56 },
683 "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
685 "-vdwread", FALSE, etBOOL,
686 { &bReadVDW },
687 "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing the radii based on the force field"
689 { "-atom", FALSE, etBOOL,
690 { &peratom }, "Force B-factor attachment per atom" },
691 { "-legend", FALSE, etBOOL,
692 { &bLegend }, "Make B-factor legend" },
693 { "-label", FALSE, etSTR,
694 { &label }, "Add chain label for all residues" },
696 "-conect", FALSE, etBOOL,
697 { &bCONECT },
698 "Add CONECT records to a [REF].pdb[ref] file when written. Can only be done when a topology is present"
701 #define NPA asize(pa)
703 FILE *out;
704 const char *infile, *outfile;
705 int outftp, inftp, natom, i, j, n_bfac, itype, ntype;
706 double *bfac = nullptr, c6, c12;
707 int *bfac_nr = nullptr;
708 t_topology *top = nullptr;
709 char *grpname, *sgrpname, *agrpname;
710 int isize, ssize, numAlignmentAtoms;
711 int *index, *sindex, *aindex;
712 rvec *x, *v, gc, rmin, rmax, size;
713 int ePBC;
714 matrix box, rotmatrix, trans;
715 rvec princd, tmpvec;
716 gmx_bool bIndex, bSetSize, bSetAng, bDist, bSetCenter, bAlign;
717 gmx_bool bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
718 real diam = 0, mass = 0, d, vdw;
719 gmx_conect conect;
720 gmx_output_env_t *oenv;
721 t_filenm fnm[] =
723 { efSTX, "-f", nullptr, ffREAD },
724 { efNDX, "-n", nullptr, ffOPTRD },
725 { efSTO, nullptr, nullptr, ffOPTWR },
726 { efPQR, "-mead", "mead", ffOPTWR },
727 { efDAT, "-bf", "bfact", ffOPTRD }
729 #define NFILE asize(fnm)
731 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa,
732 asize(desc), desc, asize(bugs), bugs, &oenv))
734 return 0;
736 fprintf(stdout, "Note that major changes are planned in future for "
737 "editconf, to improve usability and utility.");
739 bIndex = opt2bSet("-n", NFILE, fnm) || bNDEF;
740 bMead = opt2bSet("-mead", NFILE, fnm);
741 bSetSize = opt2parg_bSet("-box", NPA, pa);
742 bSetAng = opt2parg_bSet("-angles", NPA, pa);
743 bSetCenter = opt2parg_bSet("-center", NPA, pa);
744 bDist = opt2parg_bSet("-d", NPA, pa);
745 bAlign = opt2parg_bSet("-align", NPA, pa);
746 /* Only automatically turn on centering without -noc */
747 if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
749 bCenter = TRUE;
751 bScale = opt2parg_bSet("-scale", NPA, pa);
752 bRho = opt2parg_bSet("-density", NPA, pa);
753 bTranslate = opt2parg_bSet("-translate", NPA, pa);
754 bRotate = opt2parg_bSet("-rotate", NPA, pa);
755 if (bScale && bRho)
757 fprintf(stderr, "WARNING: setting -density overrides -scale\n");
759 bScale = bScale || bRho;
760 bCalcGeom = bCenter || bRotate || bOrient || bScale;
762 GMX_RELEASE_ASSERT(btype[0] != nullptr, "Option setting inconsistency; btype[0] is NULL");
764 bCalcDiam = (btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o');
766 infile = ftp2fn(efSTX, NFILE, fnm);
767 if (bMead)
769 outfile = ftp2fn(efPQR, NFILE, fnm);
771 else
773 outfile = ftp2fn(efSTO, NFILE, fnm);
775 outftp = fn2ftp(outfile);
776 inftp = fn2ftp(infile);
778 AtomProperties aps;
780 if (bMead && bGrasp)
782 printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
783 bGrasp = FALSE;
785 if (bGrasp && (outftp != efPDB))
787 gmx_fatal(FARGS, "Output file should be a .pdb file"
788 " when using the -grasp option\n");
790 if ((bMead || bGrasp) && (fn2ftp(infile) != efTPR))
792 gmx_fatal(FARGS, "Input file should be a .tpr file"
793 " when using the -mead option\n");
796 t_topology *top_tmp;
797 snew(top_tmp, 1);
798 read_tps_conf(infile, top_tmp, &ePBC, &x, &v, box, FALSE);
799 t_atoms &atoms = top_tmp->atoms;
800 natom = atoms.nr;
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)
821 real vol = det(box);
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);
831 if (bMead || bGrasp)
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 */
842 if (bReadVDW)
844 if (!aps.setAtomProperty(epropVDW,
845 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
846 *top->atoms.atomname[i], &vdw))
848 vdw = rvdw;
851 else
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))
858 real sig6;
859 if (bSig56)
861 sig6 = 2*c12/c6;
863 else
865 sig6 = c12/c6;
867 vdw = 0.5*gmx::sixthroot(sig6);
869 else
871 vdw = rvdw;
874 /* Factor of 10 for nm -> Angstroms */
875 vdw *= 10;
877 if (bMead)
879 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
880 atoms.pdbinfo[i].bfac = vdw;
882 else
884 atoms.pdbinfo[i].occup = vdw;
885 atoms.pdbinfo[i].bfac = top->atoms.atom[i].q;
889 bHaveV = FALSE;
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");
899 if (visbox[0] > 0)
901 if (bIndex)
903 gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
905 if (outftp != efPDB)
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);
915 /* remove pbc */
916 if (bRMPBC)
918 rm_gropbc(&atoms, x, box);
921 if (bCalcGeom)
923 if (bIndex)
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);
929 else
931 ssize = atoms.nr;
932 sindex = nullptr;
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]);
938 if (bCalcDiam)
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);
960 if (bOrient)
962 int *index;
963 char *grpnames;
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);
971 sfree(index);
972 sfree(grpnames);
975 if (bScale)
977 /* scale coordinates and box */
978 if (bRho)
980 /* Compute scaling constant */
981 real vol, dens;
983 vol = det(box);
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);
1000 if (bAlign)
1002 if (bIndex)
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);
1008 else
1010 numAlignmentAtoms = atoms.nr;
1011 snew(aindex, numAlignmentAtoms);
1012 for (i = 0; i < numAlignmentAtoms; i++)
1014 aindex[i] = 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*/
1026 /*will this work?*/
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);
1046 if (!bIndex)
1048 sfree(aindex);
1052 if (bTranslate)
1054 if (bIndex)
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);
1060 else
1062 ssize = atoms.nr;
1063 sindex = nullptr;
1065 printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize, natom,
1066 translation[XX], translation[YY], translation[ZZ]);
1067 if (sindex)
1069 for (i = 0; i < ssize; i++)
1071 rvec_inc(x[sindex[i]], translation);
1074 else
1076 for (i = 0; i < natom; i++)
1078 rvec_inc(x[i], translation);
1082 if (bRotate)
1084 /* Rotate */
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]);
1093 if (bCalcGeom)
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)))
1107 ePBC = epbcXYZ;
1108 if (!(bSetSize || bDist))
1110 for (i = 0; i < DIM; i++)
1112 newbox[i] = norm(box[i]);
1115 clear_mat(box);
1116 /* calculate new boxsize */
1117 switch (btype[0][0])
1119 case 't':
1120 if (bDist)
1122 for (i = 0; i < DIM; i++)
1124 newbox[i] = size[i]+2*dist;
1127 if (!bSetAng)
1129 box[XX][XX] = newbox[XX];
1130 box[YY][YY] = newbox[YY];
1131 box[ZZ][ZZ] = newbox[ZZ];
1133 else
1135 matrix_convert(box, newbox, newang);
1137 break;
1138 case 'c':
1139 case 'd':
1140 case 'o':
1141 if (bSetSize)
1143 d = newbox[0];
1145 else
1147 d = diam+2*dist;
1149 if (btype[0][0] == 'c')
1151 for (i = 0; i < DIM; i++)
1153 box[i][i] = d;
1156 else if (btype[0][0] == 'd')
1158 box[XX][XX] = d;
1159 box[YY][YY] = d;
1160 box[ZZ][XX] = d/2;
1161 box[ZZ][YY] = d/2;
1162 box[ZZ][ZZ] = d*std::sqrt(2.0)/2.0;
1164 else
1166 box[XX][XX] = d;
1167 box[YY][XX] = d/3;
1168 box[YY][YY] = d*std::sqrt(2.0)*2.0/3.0;
1169 box[ZZ][XX] = -d/3;
1170 box[ZZ][YY] = d*std::sqrt(2.0)/3.0;
1171 box[ZZ][ZZ] = d*std::sqrt(6.0)/3.0;
1173 break;
1177 /* calculate new coords for geometrical center */
1178 if (!bSetCenter)
1180 calc_box_center(ecenterDEF, box, center);
1183 /* center molecule on 'center' */
1184 if (bCenter)
1186 center_conf(natom, x, center, gc);
1189 /* print some */
1190 if (bCalcGeom)
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')
1219 if (TRICLINIC(box))
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);
1237 else
1239 conect = nullptr;
1242 if (bIndex)
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, *top_tmp->name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, TRUE, FALSE);
1270 gmx_ffclose(out);
1272 else
1274 write_sto_conf_indexed(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : nullptr, ePBC, box, isize, index);
1277 else
1279 if (resnr_start >= 0)
1281 renum_resnr(&atoms, atoms.nr, nullptr, resnr_start);
1284 if ((outftp == efPDB) || (outftp == efPQR))
1286 out = gmx_ffopen(outfile, "w");
1287 if (bMead)
1289 fprintf(out, "REMARK "
1290 "The B-factors in this file hold atomic radii\n");
1291 fprintf(out, "REMARK "
1292 "The occupancy in this file hold atomic charges\n");
1294 else if (bGrasp)
1296 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1297 fprintf(out, "REMARK "
1298 "The B-factors in this file hold atomic charges\n");
1299 fprintf(out, "REMARK "
1300 "The occupancy in this file hold atomic radii\n");
1302 else if (opt2bSet("-bf", NFILE, fnm))
1304 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1305 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms,
1306 n_bfac, bfac, bfac_nr, peratom);
1308 if (opt2parg_bSet("-label", NPA, pa))
1310 for (i = 0; (i < atoms.nr); i++)
1312 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1315 /* Need to bypass the regular write_pdbfile because I don't want to change
1316 * all instances to include the boolean flag for writing out PQR files.
1318 int *index;
1319 snew(index, atoms.nr);
1320 for (int i = 0; i < atoms.nr; i++)
1322 index[i] = i;
1324 write_pdbfile_indexed(out, *top_tmp->name, &atoms, x, ePBC, box, ' ', -1, atoms.nr, index, conect,
1325 TRUE, outftp == efPQR);
1326 sfree(index);
1327 if (bLegend)
1329 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1331 if (visbox[0] > 0)
1333 visualize_box(out, bLegend ? atoms.nr+12 : atoms.nr,
1334 bLegend ? atoms.nres = 12 : atoms.nres, box, visbox);
1336 gmx_ffclose(out);
1338 else
1340 write_sto_conf(outfile, *top_tmp->name, &atoms, x, bHaveV ? v : nullptr, ePBC, box);
1343 do_view(oenv, outfile, nullptr);
1345 return 0;