1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Green Red Orange Magenta Azure Cyan Skyblue
82 real calc_mass(t_atoms *atoms, gmx_bool bGetMass, gmx_atomprop_t aps)
88 for (i = 0; (i < atoms->nr); i++)
92 gmx_atomprop_query(aps, epropMass,
93 *atoms->resinfo[atoms->atom[i].resind].name,
94 *atoms->atomname[i], &(atoms->atom[i].m));
96 tmass += atoms->atom[i].m;
102 real calc_geom(int isize, atom_id *index, rvec *x, rvec geom_center, rvec min,
103 rvec max, gmx_bool bDiam)
109 clear_rvec(geom_center);
126 for (j = 0; j < DIM; j++)
128 min[j] = max[j] = x[ii][j];
130 for (i = 0; i < isize; i++)
140 rvec_inc(geom_center, x[ii]);
141 for (j = 0; j < DIM; j++)
143 if (x[ii][j] < min[j])
147 if (x[ii][j] > max[j])
156 for (j = i + 1; j < isize; j++)
158 d = distance2(x[ii], x[index[j]]);
159 diam2 = max(d, diam2);
164 for (j = i + 1; j < isize; j++)
166 d = distance2(x[i], x[j]);
167 diam2 = max(d, diam2);
172 svmul(1.0 / isize, geom_center, geom_center);
178 void center_conf(int natom, rvec *x, rvec center, rvec geom_cent)
183 rvec_sub(center, geom_cent, shift);
185 printf(" shift :%7.3f%7.3f%7.3f (nm)\n", shift[XX], shift[YY],
188 for (i = 0; (i < natom); i++)
190 rvec_inc(x[i], shift);
194 void scale_conf(int natom, rvec x[], matrix box, rvec scale)
198 for (i = 0; i < natom; i++)
200 for (j = 0; j < DIM; j++)
205 for (i = 0; i < DIM; i++)
207 for (j = 0; j < DIM; j++)
209 box[i][j] *= scale[j];
214 void read_bfac(const char *fn, int *n_bfac, double **bfac_val, int **bfac_nr)
219 *n_bfac = get_lines(fn, &bfac_lines);
220 snew(*bfac_val, *n_bfac);
221 snew(*bfac_nr, *n_bfac);
222 fprintf(stderr, "Reading %d B-factors from %s\n", *n_bfac, fn);
223 for (i = 0; (i < *n_bfac); i++)
225 /*fprintf(stderr, "Line %d: %s",i,bfac_lines[i]);*/
226 sscanf(bfac_lines[i], "%d %lf", &(*bfac_nr)[i], &(*bfac_val)[i]);
227 /*fprintf(stderr," nr %d val %g\n",(*bfac_nr)[i],(*bfac_val)[i]);*/
232 void set_pdb_conf_bfac(int natoms, int nres, t_atoms *atoms, int n_bfac,
233 double *bfac, int *bfac_nr, gmx_bool peratom)
236 real bfac_min, bfac_max;
242 for (i = 0; (i < n_bfac); i++)
244 if (bfac_nr[i] - 1 >= atoms->nres)
248 /* if ((bfac_nr[i]-1<0) || (bfac_nr[i]-1>=atoms->nr))
249 gmx_fatal(FARGS,"Index of B-Factor %d is out of range: %d (%g)",
250 i+1,bfac_nr[i],bfac[i]); */
251 if (bfac[i] > bfac_max)
255 if (bfac[i] < bfac_min)
260 while ((bfac_max > 99.99) || (bfac_min < -99.99))
263 "Range of values for B-factors too large (min %g, max %g) "
264 "will scale down a factor 10\n", bfac_min, bfac_max);
265 for (i = 0; (i < n_bfac); i++)
272 while ((fabs(bfac_max) < 0.5) && (fabs(bfac_min) < 0.5))
275 "Range of values for B-factors too small (min %g, max %g) "
276 "will scale up a factor 10\n", bfac_min, bfac_max);
277 for (i = 0; (i < n_bfac); i++)
285 for (i = 0; (i < natoms); i++)
287 atoms->pdbinfo[i].bfac = 0;
292 fprintf(stderr, "Will attach %d B-factors to %d residues\n", n_bfac,
294 for (i = 0; (i < n_bfac); i++)
297 for (n = 0; (n < natoms); n++)
299 if (bfac_nr[i] == atoms->resinfo[atoms->atom[n].resind].nr)
301 atoms->pdbinfo[n].bfac = bfac[i];
307 gmx_warning("Residue nr %d not found\n", bfac_nr[i]);
313 fprintf(stderr, "Will attach %d B-factors to %d atoms\n", n_bfac,
315 for (i = 0; (i < n_bfac); i++)
317 atoms->pdbinfo[bfac_nr[i] - 1].bfac = bfac[i];
322 void pdb_legend(FILE *out, int natoms, int nres, t_atoms *atoms, rvec x[])
324 real bfac_min, bfac_max, xmin, ymin, zmin;
333 for (i = 0; (i < natoms); i++)
335 xmin = min(xmin, x[i][XX]);
336 ymin = min(ymin, x[i][YY]);
337 zmin = min(zmin, x[i][ZZ]);
338 bfac_min = min(bfac_min, atoms->pdbinfo[i].bfac);
339 bfac_max = max(bfac_max, atoms->pdbinfo[i].bfac);
341 fprintf(stderr, "B-factors range from %g to %g\n", bfac_min, bfac_max);
342 for (i = 1; (i < 12); i++)
345 "%-6s%5u %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f\n",
346 "ATOM ", natoms + 1 + i, "CA", "LEG", space, nres + 1, space,
347 (xmin + (i * 0.12)) * 10, ymin * 10, zmin * 10, 1.0, bfac_min
348 + ((i - 1.0) * (bfac_max - bfac_min) / 10));
352 void visualize_images(const char *fn, int ePBC, matrix box)
360 init_t_atoms(&atoms, nat, FALSE);
363 /* FIXME: Constness should not be cast away */
365 ala = (char *) "ALA";
366 for (i = 0; i < nat; i++)
368 atoms.atomname[i] = &c;
369 atoms.atom[i].resind = i;
370 atoms.resinfo[i].name = &ala;
371 atoms.resinfo[i].nr = i + 1;
372 atoms.resinfo[i].chainid = 'A' + i / NCUCVERT;
374 calc_triclinic_images(box, img + 1);
376 write_sto_conf(fn, "Images", &atoms, img, NULL, ePBC, box);
378 free_t_atoms(&atoms, FALSE);
382 void visualize_box(FILE *out, int a0, int r0, matrix box, rvec gridsize)
386 int nx, ny, nz, nbox, nat;
390 0, 1, 1, 3, 3, 2, 0, 2, 0, 4, 1, 5, 3, 7, 2, 6, 4, 5, 5, 7, 7, 6, 6,
397 nx = (int) (gridsize[XX] + 0.5);
398 ny = (int) (gridsize[YY] + 0.5);
399 nz = (int) (gridsize[ZZ] + 0.5);
403 nat = nbox * NCUCVERT;
405 calc_compact_unitcell_vertices(ecenterDEF, box, vert);
407 for (z = 0; z < nz; z++)
409 for (y = 0; y < ny; y++)
411 for (x = 0; x < nx; x++)
413 for (i = 0; i < DIM; i++)
415 shift[i] = x * box[0][i] + y * box[1][i] + z
418 for (i = 0; i < NCUCVERT; i++)
420 rvec_add(vert[i], shift, vert[j]);
427 for (i = 0; i < nat; i++)
429 fprintf(out, get_pdbformat(), "ATOM", a0 + i, "C", "BOX", 'K' + i
430 / NCUCVERT, r0 + i, 10 * vert[i][XX], 10 * vert[i][YY], 10
435 edge = compact_unitcell_edges();
436 for (j = 0; j < nbox; j++)
438 for (i = 0; i < NCUCEDGE; i++)
440 fprintf(out, "CONECT%5d%5d\n", a0 + j * NCUCVERT + edge[2 * i],
441 a0 + j * NCUCVERT + edge[2 * i + 1]);
450 for (z = 0; z <= 1; z++)
452 for (y = 0; y <= 1; y++)
454 for (x = 0; x <= 1; x++)
456 fprintf(out, get_pdbformat(), "ATOM", a0 + i, "C", "BOX", 'K' + i
457 / 8, r0 + i, x * 10 * box[XX][XX],
458 y * 10 * box[YY][YY], z * 10 * box[ZZ][ZZ]);
464 for (i = 0; i < 24; i += 2)
466 fprintf(out, "CONECT%5d%5d\n", a0 + rectedge[i], a0 + rectedge[i
472 void calc_rotmatrix(rvec principal_axis, rvec targetvec, matrix rotmatrix)
475 real ux, uy, uz, costheta, sintheta;
477 costheta = cos_angle(principal_axis, targetvec);
478 sintheta = sqrt(1.0-costheta*costheta); /* sign is always positive since 0<theta<pi */
480 /* Determine rotation from cross product with target vector */
481 cprod(principal_axis, targetvec, rotvec);
482 unitv(rotvec, rotvec);
483 printf("Aligning %g %g %g to %g %g %g : xprod %g %g %g\n",
484 principal_axis[XX], principal_axis[YY], principal_axis[ZZ], targetvec[XX], targetvec[YY], targetvec[ZZ],
485 rotvec[XX], rotvec[YY], rotvec[ZZ]);
490 rotmatrix[0][0] = ux*ux + (1.0-ux*ux)*costheta;
491 rotmatrix[0][1] = ux*uy*(1-costheta)-uz*sintheta;
492 rotmatrix[0][2] = ux*uz*(1-costheta)+uy*sintheta;
493 rotmatrix[1][0] = ux*uy*(1-costheta)+uz*sintheta;
494 rotmatrix[1][1] = uy*uy + (1.0-uy*uy)*costheta;
495 rotmatrix[1][2] = uy*uz*(1-costheta)-ux*sintheta;
496 rotmatrix[2][0] = ux*uz*(1-costheta)-uy*sintheta;
497 rotmatrix[2][1] = uy*uz*(1-costheta)+ux*sintheta;
498 rotmatrix[2][2] = uz*uz + (1.0-uz*uz)*costheta;
500 printf("Rotation matrix: \n%g %g %g\n%g %g %g\n%g %g %g\n",
501 rotmatrix[0][0], rotmatrix[0][1], rotmatrix[0][2],
502 rotmatrix[1][0], rotmatrix[1][1], rotmatrix[1][2],
503 rotmatrix[2][0], rotmatrix[2][1], rotmatrix[2][2]);
506 static void renum_resnr(t_atoms *atoms, int isize, const int *index,
509 int i, resind_prev, resind;
512 for (i = 0; i < isize; i++)
514 resind = atoms->atom[index == NULL ? i : index[i]].resind;
515 if (resind != resind_prev)
517 atoms->resinfo[resind].nr = resnr_start;
520 resind_prev = resind;
524 int gmx_editconf(int argc, char *argv[])
529 "[TT]editconf[tt] converts generic structure format to [TT].gro[tt], [TT].g96[tt]",
532 "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and",
533 "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]",
534 "will center the system in the box, unless [TT]-noc[tt] is used.",
536 "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a",
537 "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal",
538 "[TT]dodecahedron[tt] represents a rhombic dodecahedron and",
539 "[TT]octahedron[tt] is a truncated octahedron.",
540 "The last two are special cases of a triclinic box.",
541 "The length of the three box vectors of the truncated octahedron is the",
542 "shortest distance between two opposite hexagons.",
543 "Relative to a cubic box with some periodic image distance, the volume of a ",
544 "dodecahedron with this same periodic distance is 0.71 times that of the cube, ",
545 "and that of a truncated octahedron is 0.77 times.",
547 "Option [TT]-box[tt] requires only",
548 "one value for a cubic, rhombic dodecahedral, or truncated octahedral box.",
550 "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the [IT]x[it]-, [IT]y[it]-,",
551 "and [IT]z[it]-directions is used. With [TT]-d[tt] and [TT]cubic[tt],",
552 "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set",
553 "to the diameter of the system (largest distance between atoms) plus twice",
554 "the specified distance.",
556 "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and",
557 "a triclinic box and cannot be used with option [TT]-d[tt].",
559 "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group",
560 "can be selected for calculating the size and the geometric center,",
561 "otherwise the whole system is used.",
563 "[TT]-rotate[tt] rotates the coordinates and velocities.",
565 "[TT]-princ[tt] aligns the principal axes of the system along the",
566 "coordinate axes, with the longest axis aligned with the [IT]x[it]-axis. ",
567 "This may allow you to decrease the box volume,",
568 "but beware that molecules can rotate significantly in a nanosecond.",
570 "Scaling is applied before any of the other operations are",
571 "performed. Boxes and coordinates can be scaled to give a certain density (option",
572 "[TT]-density[tt]). Note that this may be inaccurate in case a [TT].gro[tt]",
573 "file is given as input. A special feature of the scaling option is that when the",
574 "factor -1 is given in one dimension, one obtains a mirror image,",
575 "mirrored in one of the planes. When one uses -1 in three dimensions, ",
576 "a point-mirror image is obtained.[PAR]",
577 "Groups are selected after all operations have been applied.[PAR]",
578 "Periodicity can be removed in a crude manner.",
579 "It is important that the box vectors at the bottom of your input file",
580 "are correct when the periodicity is to be removed.",
582 "When writing [TT].pdb[tt] files, B-factors can be",
583 "added with the [TT]-bf[tt] option. B-factors are read",
584 "from a file with with following format: first line states number of",
585 "entries in the file, next lines state an index",
586 "followed by a B-factor. The B-factors will be attached per residue",
587 "unless an index is larger than the number of residues or unless the",
588 "[TT]-atom[tt] option is set. Obviously, any type of numeric data can",
589 "be added instead of B-factors. [TT]-legend[tt] will produce",
590 "a row of CA atoms with B-factors ranging from the minimum to the",
591 "maximum value found, effectively making a legend for viewing.",
593 "With the option [TT]-mead[tt] a special [TT].pdb[tt] ([TT].pqr[tt])",
594 "file for the MEAD electrostatics",
595 "program (Poisson-Boltzmann solver) can be made. A further prerequisite",
596 "is that the input file is a run input file.",
597 "The B-factor field is then filled with the Van der Waals radius",
598 "of the atoms while the occupancy field will hold the charge.",
600 "The option [TT]-grasp[tt] is similar, but it puts the charges in the B-factor",
601 "and the radius in the occupancy.",
603 "Option [TT]-align[tt] allows alignment",
604 "of the principal axis of a specified group against the given vector, ",
605 "with an optional center of rotation specified by [TT]-aligncenter[tt].",
607 "Finally, with option [TT]-label[tt], [TT]editconf[tt] can add a chain identifier",
608 "to a [TT].pdb[tt] file, which can be useful for analysis with e.g. Rasmol.",
610 "To convert a truncated octrahedron file produced by a package which uses",
611 "a cubic box with the corners cut off (such as GROMOS), use:[BR]",
612 "[TT]editconf -f in -rotate 0 45 35.264 -bt o -box veclen -o out[tt][BR]",
613 "where [TT]veclen[tt] is the size of the cubic box times [SQRT]3[sqrt]/2."
617 "For complex molecules, the periodicity removal routine may break down, "
618 "in that case you can use [TT]trjconv[tt]."
620 static real dist = 0.0, rbox = 0.0, to_diam = 0.0;
621 static gmx_bool bNDEF = FALSE, bRMPBC = FALSE, bCenter = FALSE, bReadVDW =
622 FALSE, bCONECT = FALSE;
623 static gmx_bool peratom = FALSE, bLegend = FALSE, bOrient = FALSE, bMead =
624 FALSE, bGrasp = FALSE, bSig56 = FALSE;
626 { 1, 1, 1 }, newbox =
627 { 0, 0, 0 }, newang =
629 static real rho = 1000.0, rvdw = 0.12;
631 { 0, 0, 0 }, translation =
632 { 0, 0, 0 }, rotangles =
633 { 0, 0, 0 }, aligncenter =
634 { 0, 0, 0 }, targetvec =
636 static const char *btype[] =
637 { NULL, "triclinic", "cubic", "dodecahedron", "octahedron", NULL },
641 static int resnr_start = -1;
645 { "-ndef", FALSE, etBOOL,
646 { &bNDEF }, "Choose output from default index groups" },
647 { "-visbox", FALSE, etRVEC,
649 "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" },
650 { "-bt", FALSE, etENUM,
651 { btype }, "Box type for [TT]-box[tt] and [TT]-d[tt]" },
652 { "-box", FALSE, etRVEC,
653 { newbox }, "Box vector lengths (a,b,c)" },
654 { "-angles", FALSE, etRVEC,
655 { newang }, "Angles between the box vectors (bc,ac,ab)" },
656 { "-d", FALSE, etREAL,
657 { &dist }, "Distance between the solute and the box" },
658 { "-c", FALSE, etBOOL,
660 "Center molecule in box (implied by [TT]-box[tt] and [TT]-d[tt])" },
661 { "-center", FALSE, etRVEC,
662 { center }, "Coordinates of geometrical center" },
663 { "-aligncenter", FALSE, etRVEC,
664 { aligncenter }, "Center of rotation for alignment" },
665 { "-align", FALSE, etRVEC,
667 "Align to target vector" },
668 { "-translate", FALSE, etRVEC,
669 { translation }, "Translation" },
670 { "-rotate", FALSE, etRVEC,
672 "Rotation around the X, Y and Z axes in degrees" },
673 { "-princ", FALSE, etBOOL,
675 "Orient molecule(s) along their principal axes" },
676 { "-scale", FALSE, etRVEC,
677 { scale }, "Scaling factor" },
678 { "-density", FALSE, etREAL,
680 "Density (g/L) of the output box achieved by scaling" },
681 { "-pbc", FALSE, etBOOL,
683 "Remove the periodicity (make molecule whole again)" },
684 { "-resnr", FALSE, etINT,
686 " Renumber residues starting from resnr" },
687 { "-grasp", FALSE, etBOOL,
689 "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" },
691 "-rvdw", FALSE, etREAL,
693 "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"
695 { "-sig56", FALSE, etBOOL,
697 "Use rmin/2 (minimum in the Van der Waals potential) rather than [GRK]sigma[grk]/2 " },
699 "-vdwread", FALSE, etBOOL,
701 "Read the Van der Waals radii from the file [TT]vdwradii.dat[tt] rather than computing the radii based on the force field"
703 { "-atom", FALSE, etBOOL,
704 { &peratom }, "Force B-factor attachment per atom" },
705 { "-legend", FALSE, etBOOL,
706 { &bLegend }, "Make B-factor legend" },
707 { "-label", FALSE, etSTR,
708 { &label }, "Add chain label for all residues" },
710 "-conect", FALSE, etBOOL,
712 "Add CONECT records to a [TT].pdb[tt] file when written. Can only be done when a topology is present"
715 #define NPA asize(pa)
718 const char *infile, *outfile;
720 int outftp, inftp, natom, i, j, n_bfac, itype, ntype;
721 double *bfac = NULL, c6, c12;
723 t_topology *top = NULL;
725 char *grpname, *sgrpname, *agrpname;
726 int isize, ssize, tsize, asize;
727 atom_id *index, *sindex, *tindex, *aindex;
728 rvec *x, *v, gc, min, max, size;
730 matrix box, rotmatrix, trans;
732 gmx_bool bIndex, bSetSize, bSetAng, bCubic, bDist, bSetCenter, bAlign;
733 gmx_bool bHaveV, bScale, bRho, bTranslate, bRotate, bCalcGeom, bCalcDiam;
734 real xs, ys, zs, xcent, ycent, zcent, diam = 0, mass = 0, d, vdw;
740 { efSTX, "-f", NULL, ffREAD },
741 { efNDX, "-n", NULL, ffOPTRD },
742 { efSTO, NULL, NULL, ffOPTWR },
743 { efPQR, "-mead", "mead", ffOPTWR },
744 { efDAT, "-bf", "bfact", ffOPTRD }
746 #define NFILE asize(fnm)
748 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, NPA, pa,
749 asize(desc), desc, asize(bugs), bugs, &oenv))
754 bIndex = opt2bSet("-n", NFILE, fnm) || bNDEF;
755 bMead = opt2bSet("-mead", NFILE, fnm);
756 bSetSize = opt2parg_bSet("-box", NPA, pa);
757 bSetAng = opt2parg_bSet("-angles", NPA, pa);
758 bSetCenter = opt2parg_bSet("-center", NPA, pa);
759 bDist = opt2parg_bSet("-d", NPA, pa);
760 bAlign = opt2parg_bSet("-align", NPA, pa);
761 /* Only automatically turn on centering without -noc */
762 if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c", NPA, pa))
766 bScale = opt2parg_bSet("-scale", NPA, pa);
767 bRho = opt2parg_bSet("-density", NPA, pa);
768 bTranslate = opt2parg_bSet("-translate", NPA, pa);
769 bRotate = opt2parg_bSet("-rotate", NPA, pa);
772 fprintf(stderr, "WARNING: setting -density overrides -scale\n");
774 bScale = bScale || bRho;
775 bCalcGeom = bCenter || bRotate || bOrient || bScale;
776 bCalcDiam = btype[0][0] == 'c' || btype[0][0] == 'd' || btype[0][0] == 'o';
778 infile = ftp2fn(efSTX, NFILE, fnm);
781 outfile = ftp2fn(efPQR, NFILE, fnm);
785 outfile = ftp2fn(efSTO, NFILE, fnm);
787 outftp = fn2ftp(outfile);
788 inftp = fn2ftp(infile);
790 aps = gmx_atomprop_init();
794 printf("Incompatible options -mead and -grasp. Turning off -grasp\n");
797 if (bGrasp && (outftp != efPDB))
799 gmx_fatal(FARGS, "Output file should be a .pdb file"
800 " when using the -grasp option\n");
802 if ((bMead || bGrasp) && !((fn2ftp(infile) == efTPR) ||
803 (fn2ftp(infile) == efTPA) ||
804 (fn2ftp(infile) == efTPB)))
806 gmx_fatal(FARGS, "Input file should be a .tp[abr] file"
807 " when using the -mead option\n");
810 get_stx_coordnum(infile, &natom);
811 init_t_atoms(&atoms, natom, TRUE);
814 read_stx_conf(infile, title, &atoms, x, v, &ePBC, box);
815 if (fn2ftp(infile) == efPDB)
817 get_pdb_atomnumber(&atoms, aps);
819 printf("Read %d atoms\n", atoms.nr);
821 /* Get the element numbers if available in a pdb file */
822 if (fn2ftp(infile) == efPDB)
824 get_pdb_atomnumber(&atoms, aps);
827 if (ePBC != epbcNONE)
830 printf("Volume: %g nm^3, corresponds to roughly %d electrons\n",
831 vol, 100*((int)(vol*4.5)));
834 if (bMead || bGrasp || bCONECT)
836 top = read_top(infile, NULL);
841 if (atoms.nr != top->atoms.nr)
843 gmx_fatal(FARGS, "Atom numbers don't match (%d vs. %d)", atoms.nr, top->atoms.nr);
845 snew(atoms.pdbinfo, top->atoms.nr);
846 ntype = top->idef.atnr;
847 for (i = 0; (i < atoms.nr); i++)
849 /* Determine the Van der Waals radius from the force field */
852 if (!gmx_atomprop_query(aps, epropVDW,
853 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
854 *top->atoms.atomname[i], &vdw))
861 itype = top->atoms.atom[i].type;
862 c12 = top->idef.iparams[itype*ntype+itype].lj.c12;
863 c6 = top->idef.iparams[itype*ntype+itype].lj.c6;
864 if ((c6 != 0) && (c12 != 0))
875 vdw = 0.5*pow(sig6, 1.0/6.0);
882 /* Factor of 10 for nm -> Angstroms */
887 atoms.pdbinfo[i].occup = top->atoms.atom[i].q;
888 atoms.pdbinfo[i].bfac = vdw;
892 atoms.pdbinfo[i].occup = vdw;
893 atoms.pdbinfo[i].bfac = top->atoms.atom[i].q;
898 for (i = 0; (i < natom) && !bHaveV; i++)
900 for (j = 0; (j < DIM) && !bHaveV; j++)
902 bHaveV = bHaveV || (v[i][j] != 0);
905 printf("%selocities found\n", bHaveV ? "V" : "No v");
911 gmx_fatal(FARGS, "Sorry, can not visualize box with index groups");
915 gmx_fatal(FARGS, "Sorry, can only visualize box with a pdb file");
918 else if (visbox[0] == -1)
920 visualize_images("images.pdb", ePBC, box);
926 rm_gropbc(&atoms, x, box);
933 fprintf(stderr, "\nSelect a group for determining the system size:\n");
934 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
935 1, &ssize, &sindex, &sgrpname);
942 diam = calc_geom(ssize, sindex, x, gc, min, max, bCalcDiam);
943 rvec_sub(max, min, size);
944 printf(" system size :%7.3f%7.3f%7.3f (nm)\n",
945 size[XX], size[YY], size[ZZ]);
948 printf(" diameter :%7.3f (nm)\n", diam);
950 printf(" center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
951 printf(" box vectors :%7.3f%7.3f%7.3f (nm)\n",
952 norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
953 printf(" box angles :%7.2f%7.2f%7.2f (degrees)\n",
954 norm2(box[ZZ]) == 0 ? 0 :
955 RAD2DEG*acos(cos_angle_no_table(box[YY], box[ZZ])),
956 norm2(box[ZZ]) == 0 ? 0 :
957 RAD2DEG*acos(cos_angle_no_table(box[XX], box[ZZ])),
958 norm2(box[YY]) == 0 ? 0 :
959 RAD2DEG*acos(cos_angle_no_table(box[XX], box[YY])));
960 printf(" box volume :%7.2f (nm^3)\n", det(box));
963 if (bRho || bOrient || bAlign)
965 mass = calc_mass(&atoms, !fn2bTPX(infile), aps);
973 /* Get a group for principal component analysis */
974 fprintf(stderr, "\nSelect group for the determining the orientation\n");
975 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
977 /* Orient the principal axes along the coordinate axes */
978 orient_princ(&atoms, isize, index, natom, x, bHaveV ? v : NULL, NULL);
985 /* scale coordinates and box */
988 /* Compute scaling constant */
992 dens = (mass*AMU)/(vol*NANO*NANO*NANO);
993 fprintf(stderr, "Volume of input %g (nm^3)\n", vol);
994 fprintf(stderr, "Mass of input %g (a.m.u.)\n", mass);
995 fprintf(stderr, "Density of input %g (g/l)\n", dens);
996 if (vol == 0 || mass == 0)
998 gmx_fatal(FARGS, "Cannot scale density with "
999 "zero mass (%g) or volume (%g)\n", mass, vol);
1002 scale[XX] = scale[YY] = scale[ZZ] = pow(dens/rho, 1.0/3.0);
1003 fprintf(stderr, "Scaling all box vectors by %g\n", scale[XX]);
1005 scale_conf(atoms.nr, x, box, scale);
1012 fprintf(stderr, "\nSelect a group that you want to align:\n");
1013 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
1014 1, &asize, &aindex, &agrpname);
1019 snew(aindex, asize);
1020 for (i = 0; i < asize; i++)
1025 printf("Aligning %d atoms (out of %d) to %g %g %g, center of rotation %g %g %g\n", asize, natom,
1026 targetvec[XX], targetvec[YY], targetvec[ZZ],
1027 aligncenter[XX], aligncenter[YY], aligncenter[ZZ]);
1028 /*subtract out pivot point*/
1029 for (i = 0; i < asize; i++)
1031 rvec_dec(x[aindex[i]], aligncenter);
1033 /*now determine transform and rotate*/
1035 principal_comp(asize, aindex, atoms.atom, x, trans, princd);
1037 unitv(targetvec, targetvec);
1038 printf("Using %g %g %g as principal axis\n", trans[0][2], trans[1][2], trans[2][2]);
1039 tmpvec[XX] = trans[0][2]; tmpvec[YY] = trans[1][2]; tmpvec[ZZ] = trans[2][2];
1040 calc_rotmatrix(tmpvec, targetvec, rotmatrix);
1041 /* rotmatrix finished */
1043 for (i = 0; i < asize; ++i)
1045 mvmul(rotmatrix, x[aindex[i]], tmpvec);
1046 copy_rvec(tmpvec, x[aindex[i]]);
1049 /*add pivot point back*/
1050 for (i = 0; i < asize; i++)
1052 rvec_inc(x[aindex[i]], aligncenter);
1064 fprintf(stderr, "\nSelect a group that you want to translate:\n");
1065 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm),
1066 1, &ssize, &sindex, &sgrpname);
1073 printf("Translating %d atoms (out of %d) by %g %g %g nm\n", ssize, natom,
1074 translation[XX], translation[YY], translation[ZZ]);
1077 for (i = 0; i < ssize; i++)
1079 rvec_inc(x[sindex[i]], translation);
1084 for (i = 0; i < natom; i++)
1086 rvec_inc(x[i], translation);
1093 printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n", rotangles[XX], rotangles[YY], rotangles[ZZ]);
1094 for (i = 0; i < DIM; i++)
1096 rotangles[i] *= DEG2RAD;
1098 rotate_conf(natom, x, v, rotangles[XX], rotangles[YY], rotangles[ZZ]);
1103 /* recalc geometrical center and max and min coordinates and size */
1104 calc_geom(ssize, sindex, x, gc, min, max, FALSE);
1105 rvec_sub(max, min, size);
1106 if (bScale || bOrient || bRotate)
1108 printf("new system size : %6.3f %6.3f %6.3f\n",
1109 size[XX], size[YY], size[ZZ]);
1113 if (bSetSize || bDist || (btype[0][0] == 't' && bSetAng))
1116 if (!(bSetSize || bDist))
1118 for (i = 0; i < DIM; i++)
1120 newbox[i] = norm(box[i]);
1124 /* calculate new boxsize */
1125 switch (btype[0][0])
1130 for (i = 0; i < DIM; i++)
1132 newbox[i] = size[i]+2*dist;
1137 box[XX][XX] = newbox[XX];
1138 box[YY][YY] = newbox[YY];
1139 box[ZZ][ZZ] = newbox[ZZ];
1143 matrix_convert(box, newbox, newang);
1157 if (btype[0][0] == 'c')
1159 for (i = 0; i < DIM; i++)
1164 else if (btype[0][0] == 'd')
1170 box[ZZ][ZZ] = d*sqrt(2)/2;
1176 box[YY][YY] = d*sqrt(2)*2/3;
1178 box[ZZ][YY] = d*sqrt(2)/3;
1179 box[ZZ][ZZ] = d*sqrt(6)/3;
1185 /* calculate new coords for geometrical center */
1188 calc_box_center(ecenterDEF, box, center);
1191 /* center molecule on 'center' */
1194 center_conf(natom, x, center, gc);
1200 calc_geom(ssize, sindex, x, gc, min, max, FALSE);
1201 printf("new center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]);
1203 if (bOrient || bScale || bDist || bSetSize)
1205 printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n",
1206 norm(box[XX]), norm(box[YY]), norm(box[ZZ]));
1207 printf("new box angles :%7.2f%7.2f%7.2f (degrees)\n",
1208 norm2(box[ZZ]) == 0 ? 0 :
1209 RAD2DEG*acos(cos_angle_no_table(box[YY], box[ZZ])),
1210 norm2(box[ZZ]) == 0 ? 0 :
1211 RAD2DEG*acos(cos_angle_no_table(box[XX], box[ZZ])),
1212 norm2(box[YY]) == 0 ? 0 :
1213 RAD2DEG*acos(cos_angle_no_table(box[XX], box[YY])));
1214 printf("new box volume :%7.2f (nm^3)\n", det(box));
1217 if (check_box(epbcXYZ, box))
1219 printf("\nWARNING: %s\n"
1220 "See the GROMACS manual for a description of the requirements that\n"
1221 "must be satisfied by descriptions of simulation cells.\n",
1222 check_box(epbcXYZ, box));
1225 if (bDist && btype[0][0] == 't')
1229 printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n"
1230 "distance from the solute to a box surface along the corresponding normal\n"
1231 "vector might be somewhat smaller than your specified value %f.\n"
1232 "You can check the actual value with g_mindist -pi\n", dist);
1234 else if (!opt2parg_bSet("-bt", NPA, pa))
1236 printf("\nWARNING: No boxtype specified - distance condition applied in each dimension.\n"
1237 "If the molecule rotates the actual distance will be smaller. You might want\n"
1238 "to use a cubic box instead, or why not try a dodecahedron today?\n");
1241 if (bCONECT && (outftp == efPDB) && (inftp == efTPR))
1243 conect = gmx_conect_generate(top);
1252 fprintf(stderr, "\nSelect a group for output:\n");
1253 get_index(&atoms, opt2fn_null("-n", NFILE, fnm),
1254 1, &isize, &index, &grpname);
1256 if (resnr_start >= 0)
1258 renum_resnr(&atoms, isize, index, resnr_start);
1261 if (opt2parg_bSet("-label", NPA, pa))
1263 for (i = 0; (i < atoms.nr); i++)
1265 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1269 if (opt2bSet("-bf", NFILE, fnm) || bLegend)
1271 gmx_fatal(FARGS, "Sorry, cannot do bfactors with an index group.");
1274 if (outftp == efPDB)
1276 out = ffopen(outfile, "w");
1277 write_pdbfile_indexed(out, title, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, TRUE);
1282 write_sto_conf_indexed(outfile, title, &atoms, x, bHaveV ? v : NULL, ePBC, box, isize, index);
1287 if (resnr_start >= 0)
1289 renum_resnr(&atoms, atoms.nr, NULL, resnr_start);
1292 if ((outftp == efPDB) || (outftp == efPQR))
1294 out = ffopen(outfile, "w");
1297 set_pdb_wide_format(TRUE);
1298 fprintf(out, "REMARK "
1299 "The B-factors in this file hold atomic radii\n");
1300 fprintf(out, "REMARK "
1301 "The occupancy in this file hold atomic charges\n");
1305 fprintf(out, "GRASP PDB FILE\nFORMAT NUMBER=1\n");
1306 fprintf(out, "REMARK "
1307 "The B-factors in this file hold atomic charges\n");
1308 fprintf(out, "REMARK "
1309 "The occupancy in this file hold atomic radii\n");
1311 else if (opt2bSet("-bf", NFILE, fnm))
1313 read_bfac(opt2fn("-bf", NFILE, fnm), &n_bfac, &bfac, &bfac_nr);
1314 set_pdb_conf_bfac(atoms.nr, atoms.nres, &atoms,
1315 n_bfac, bfac, bfac_nr, peratom);
1317 if (opt2parg_bSet("-label", NPA, pa))
1319 for (i = 0; (i < atoms.nr); i++)
1321 atoms.resinfo[atoms.atom[i].resind].chainid = label[0];
1324 write_pdbfile(out, title, &atoms, x, ePBC, box, ' ', -1, conect, TRUE);
1327 pdb_legend(out, atoms.nr, atoms.nres, &atoms, x);
1331 visualize_box(out, bLegend ? atoms.nr+12 : atoms.nr,
1332 bLegend ? atoms.nres = 12 : atoms.nres, box, visbox);
1338 write_sto_conf(outfile, title, &atoms, x, bHaveV ? v : NULL, ePBC, box);
1341 gmx_atomprop_destroy(aps);
1343 do_view(oenv, outfile, NULL);