2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/viewit.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
61 static void periodic_dist(int ePBC
,
62 matrix box
, rvec x
[], int n
, atom_id index
[],
63 real
*rmin
, real
*rmax
, int *min_ind
)
66 int nsz
, nshift
, sx
, sy
, sz
, i
, j
, s
;
67 real sqr_box
, r2min
, r2max
, r2
;
68 rvec shift
[NSHIFT_MAX
], d0
, d
;
70 sqr_box
= min(norm2(box
[XX
]), norm2(box
[YY
]));
73 sqr_box
= min(sqr_box
, norm2(box
[ZZ
]));
76 else if (ePBC
== epbcXY
)
82 gmx_fatal(FARGS
, "pbc = %s is not supported by g_mindist",
84 nsz
= 0; /* Keep compilers quiet */
88 for (sz
= -nsz
; sz
<= nsz
; sz
++)
90 for (sy
= -1; sy
<= 1; sy
++)
92 for (sx
= -1; sx
<= 1; sx
++)
94 if (sx
!= 0 || sy
!= 0 || sz
!= 0)
96 for (i
= 0; i
< DIM
; i
++)
99 sx
*box
[XX
][i
] + sy
*box
[YY
][i
] + sz
*box
[ZZ
][i
];
110 for (i
= 0; i
< n
; i
++)
112 for (j
= i
+1; j
< n
; j
++)
114 rvec_sub(x
[index
[i
]], x
[index
[j
]], d0
);
120 for (s
= 0; s
< nshift
; s
++)
122 rvec_add(d0
, shift
[s
], d
);
138 static void periodic_mindist_plot(const char *trxfn
, const char *outfn
,
139 t_topology
*top
, int ePBC
,
140 int n
, atom_id index
[], gmx_bool bSplit
,
141 const output_env_t oenv
)
144 const char *leg
[5] = { "min per.", "max int.", "box1", "box2", "box3" };
149 int natoms
, ind_min
[2] = {0, 0}, ind_mini
= 0, ind_minj
= 0;
150 real r
, rmin
, rmax
, rmint
, tmint
;
152 gmx_rmpbc_t gpbc
= NULL
;
154 natoms
= read_first_x(oenv
, &status
, trxfn
, &t
, &x
, box
);
156 check_index(NULL
, n
, index
, NULL
, natoms
);
158 out
= xvgropen(outfn
, "Minimum distance to periodic image",
159 output_env_get_time_label(oenv
), "Distance (nm)", oenv
);
160 if (output_env_get_print_xvgr_codes(oenv
))
162 fprintf(out
, "@ subtitle \"and maximum internal distance\"\n");
164 xvgr_legend(out
, 5, leg
, oenv
);
171 gpbc
= gmx_rmpbc_init(&top
->idef
, ePBC
, natoms
);
179 gmx_rmpbc(gpbc
, natoms
, box
, x
);
182 periodic_dist(ePBC
, box
, x
, n
, index
, &rmin
, &rmax
, ind_min
);
187 ind_mini
= ind_min
[0];
188 ind_minj
= ind_min
[1];
190 if (bSplit
&& !bFirst
&& fabs(t
/output_env_get_time_factor(oenv
)) < 1e-5)
192 fprintf(out
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
194 fprintf(out
, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
195 output_env_conv_time(oenv
, t
), rmin
, rmax
, norm(box
[0]), norm(box
[1]), norm(box
[2]));
198 while (read_next_x(oenv
, status
, &t
, x
, box
));
202 gmx_rmpbc_done(gpbc
);
208 "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
209 "between atoms %d and %d\n",
210 rmint
, output_env_conv_time(oenv
, tmint
), output_env_get_time_unit(oenv
),
211 index
[ind_mini
]+1, index
[ind_minj
]+1);
214 static void calc_dist(real rcut
, gmx_bool bPBC
, int ePBC
, matrix box
, rvec x
[],
215 int nx1
, int nx2
, atom_id index1
[], atom_id index2
[],
217 real
*rmin
, real
*rmax
, int *nmin
, int *nmax
,
218 int *ixmin
, int *jxmin
, int *ixmax
, int *jxmax
)
220 int i
, j
, i0
= 0, j1
;
224 real r2
, rmin2
, rmax2
, rcut2
;
237 /* Must init pbc every step because of pressure coupling */
240 set_pbc(&pbc
, ePBC
, box
);
257 for (j
= 0; (j
< j1
); j
++)
266 for (i
= i0
; (i
< nx1
); i
++)
273 pbc_dx(&pbc
, x
[ix
], x
[jx
], dx
);
277 rvec_sub(x
[ix
], x
[jx
], dx
);
323 void dist_plot(const char *fn
, const char *afile
, const char *dfile
,
324 const char *nfile
, const char *rfile
, const char *xfile
,
325 real rcut
, gmx_bool bMat
, t_atoms
*atoms
,
326 int ng
, atom_id
*index
[], int gnx
[], char *grpn
[], gmx_bool bSplit
,
327 gmx_bool bMin
, int nres
, atom_id
*residue
, gmx_bool bPBC
, int ePBC
,
328 gmx_bool bGroup
, gmx_bool bEachResEachTime
, gmx_bool bPrintResName
,
329 const output_env_t oenv
)
331 FILE *atm
, *dist
, *num
;
335 real t
, dmin
, dmax
, **mindres
= NULL
, **maxdres
= NULL
;
338 int i
= -1, j
, k
, natoms
;
339 int min1
, min2
, max1
, max2
, min1r
, min2r
, max1r
, max2r
;
345 FILE *respertime
= NULL
;
347 if ((natoms
= read_first_x(oenv
, &status
, fn
, &t
, &x0
, box
)) == 0)
349 gmx_fatal(FARGS
, "Could not read coordinates from statusfile\n");
352 sprintf(buf
, "%simum Distance", bMin
? "Min" : "Max");
353 dist
= xvgropen(dfile
, buf
, output_env_get_time_label(oenv
), "Distance (nm)", oenv
);
354 sprintf(buf
, "Number of Contacts %s %g nm", bMin
? "<" : ">", rcut
);
355 num
= nfile
? xvgropen(nfile
, buf
, output_env_get_time_label(oenv
), "Number", oenv
) : NULL
;
356 atm
= afile
? gmx_ffopen(afile
, "w") : NULL
;
357 trxout
= xfile
? open_trx(xfile
, "w") : NULL
;
364 sprintf(buf
, "Internal in %s", grpn
[0]);
365 leg
[0] = gmx_strdup(buf
);
366 xvgr_legend(dist
, 0, (const char**)leg
, oenv
);
369 xvgr_legend(num
, 0, (const char**)leg
, oenv
);
374 snew(leg
, (ng
*(ng
-1))/2);
375 for (i
= j
= 0; (i
< ng
-1); i
++)
377 for (k
= i
+1; (k
< ng
); k
++, j
++)
379 sprintf(buf
, "%s-%s", grpn
[i
], grpn
[k
]);
380 leg
[j
] = gmx_strdup(buf
);
383 xvgr_legend(dist
, j
, (const char**)leg
, oenv
);
386 xvgr_legend(num
, j
, (const char**)leg
, oenv
);
393 for (i
= 0; (i
< ng
-1); i
++)
395 sprintf(buf
, "%s-%s", grpn
[0], grpn
[i
+1]);
396 leg
[i
] = gmx_strdup(buf
);
398 xvgr_legend(dist
, ng
-1, (const char**)leg
, oenv
);
401 xvgr_legend(num
, ng
-1, (const char**)leg
, oenv
);
405 if (bEachResEachTime
)
407 sprintf(buf
, "%simum Distance", bMin
? "Min" : "Max");
408 respertime
= xvgropen(rfile
, buf
, output_env_get_time_label(oenv
), "Distance (nm)", oenv
);
409 xvgr_legend(respertime
, ng
-1, (const char**)leg
, oenv
);
412 fprintf(respertime
, "# ");
414 for (j
= 0; j
< nres
; j
++)
416 fprintf(respertime
, "%s%d ", *(atoms
->resinfo
[atoms
->atom
[index
[0][residue
[j
]]].resind
].name
), atoms
->atom
[index
[0][residue
[j
]]].resind
);
418 fprintf(respertime
, "\n");
427 for (i
= 1; i
< ng
; i
++)
429 snew(mindres
[i
-1], nres
);
430 snew(maxdres
[i
-1], nres
);
431 for (j
= 0; j
< nres
; j
++)
433 mindres
[i
-1][j
] = 1e6
;
435 /* maxdres[*][*] is already 0 */
441 if (bSplit
&& !bFirst
&& fabs(t
/output_env_get_time_factor(oenv
)) < 1e-5)
443 fprintf(dist
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
446 fprintf(num
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
450 fprintf(atm
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
453 fprintf(dist
, "%12e", output_env_conv_time(oenv
, t
));
456 fprintf(num
, "%12e", output_env_conv_time(oenv
, t
));
463 calc_dist(rcut
, bPBC
, ePBC
, box
, x0
, gnx
[0], gnx
[0], index
[0], index
[0], bGroup
,
464 &dmin
, &dmax
, &nmin
, &nmax
, &min1
, &min2
, &max1
, &max2
);
465 fprintf(dist
, " %12e", bMin
? dmin
: dmax
);
468 fprintf(num
, " %8d", bMin
? nmin
: nmax
);
473 for (i
= 0; (i
< ng
-1); i
++)
475 for (k
= i
+1; (k
< ng
); k
++)
477 calc_dist(rcut
, bPBC
, ePBC
, box
, x0
, gnx
[i
], gnx
[k
], index
[i
], index
[k
],
478 bGroup
, &dmin
, &dmax
, &nmin
, &nmax
, &min1
, &min2
, &max1
, &max2
);
479 fprintf(dist
, " %12e", bMin
? dmin
: dmax
);
482 fprintf(num
, " %8d", bMin
? nmin
: nmax
);
490 for (i
= 1; (i
< ng
); i
++)
492 calc_dist(rcut
, bPBC
, ePBC
, box
, x0
, gnx
[0], gnx
[i
], index
[0], index
[i
], bGroup
,
493 &dmin
, &dmax
, &nmin
, &nmax
, &min1
, &min2
, &max1
, &max2
);
494 fprintf(dist
, " %12e", bMin
? dmin
: dmax
);
497 fprintf(num
, " %8d", bMin
? nmin
: nmax
);
501 for (j
= 0; j
< nres
; j
++)
503 calc_dist(rcut
, bPBC
, ePBC
, box
, x0
, residue
[j
+1]-residue
[j
], gnx
[i
],
504 &(index
[0][residue
[j
]]), index
[i
], bGroup
,
505 &dmin
, &dmax
, &nmin
, &nmax
, &min1r
, &min2r
, &max1r
, &max2r
);
506 mindres
[i
-1][j
] = min(mindres
[i
-1][j
], dmin
);
507 maxdres
[i
-1][j
] = max(maxdres
[i
-1][j
], dmax
);
517 if ( (bMin
? min1
: max1
) != -1)
521 fprintf(atm
, "%12e %12d %12d\n",
522 output_env_conv_time(oenv
, t
), 1+(bMin
? min1
: max1
),
523 1+(bMin
? min2
: max2
));
529 oindex
[0] = bMin
? min1
: max1
;
530 oindex
[1] = bMin
? min2
: max2
;
531 write_trx(trxout
, 2, oindex
, atoms
, i
, t
, box
, x0
, NULL
, NULL
);
534 /*dmin should be minimum distance for residue and group*/
535 if (bEachResEachTime
)
537 fprintf(respertime
, "%12e", t
);
538 for (i
= 1; i
< ng
; i
++)
540 for (j
= 0; j
< nres
; j
++)
542 fprintf(respertime
, " %7g", bMin
? mindres
[i
-1][j
] : maxdres
[i
-1][j
]);
543 /*reset distances for next time point*/
544 mindres
[i
-1][j
] = 1e6
;
548 fprintf(respertime
, "\n");
551 while (read_next_x(oenv
, status
, &t
, x0
, box
));
569 xvgrclose(respertime
);
572 if (nres
&& !bEachResEachTime
)
576 sprintf(buf
, "%simum Distance", bMin
? "Min" : "Max");
577 res
= xvgropen(rfile
, buf
, "Residue (#)", "Distance (nm)", oenv
);
578 xvgr_legend(res
, ng
-1, (const char**)leg
, oenv
);
579 for (j
= 0; j
< nres
; j
++)
581 fprintf(res
, "%4d", j
+1);
582 for (i
= 1; i
< ng
; i
++)
584 fprintf(res
, " %7g", bMin
? mindres
[i
-1][j
] : maxdres
[i
-1][j
]);
594 int find_residues(t_atoms
*atoms
, int n
, atom_id index
[], atom_id
**resindex
)
597 int nres
= 0, resnr
, presnr
;
600 /* build index of first atom numbers for each residue */
602 snew(residx
, atoms
->nres
+1);
603 for (i
= 0; i
< n
; i
++)
605 resnr
= atoms
->atom
[index
[i
]].resind
;
615 printf("Found %d residues out of %d (%d/%d atoms)\n",
616 nres
, atoms
->nres
, atoms
->nr
, n
);
618 srenew(residx
, nres
+1);
619 /* mark end of last residue */
625 void dump_res(FILE *out
, int nres
, atom_id
*resindex
, atom_id index
[])
629 for (i
= 0; i
< nres
-1; i
++)
631 fprintf(out
, "Res %d (%d):", i
, resindex
[i
+1]-resindex
[i
]);
632 for (j
= resindex
[i
]; j
< resindex
[i
+1]; j
++)
634 fprintf(out
, " %d(%d)", j
, index
[j
]);
640 int gmx_mindist(int argc
, char *argv
[])
642 const char *desc
[] = {
643 "[THISMODULE] computes the distance between one group and a number of",
644 "other groups. Both the minimum distance",
645 "(between any pair of atoms from the respective groups)",
646 "and the number of contacts within a given",
647 "distance are written to two separate output files.",
648 "With the [TT]-group[tt] option a contact of an atom in another group",
649 "with multiple atoms in the first group is counted as one contact",
650 "instead of as multiple contacts.",
651 "With [TT]-or[tt], minimum distances to each residue in the first",
652 "group are determined and plotted as a function of residue number.[PAR]",
653 "With option [TT]-pi[tt] the minimum distance of a group to its",
654 "periodic image is plotted. This is useful for checking if a protein",
655 "has seen its periodic image during a simulation. Only one shift in",
656 "each direction is considered, giving a total of 26 shifts.",
657 "It also plots the maximum distance within the group and the lengths",
658 "of the three box vectors.[PAR]",
659 "Also [gmx-distance] and [gmx-pairdist] calculate distances."
661 const char *bugs
[] = {
662 "The [TT]-pi[tt] option is very slow."
665 static gmx_bool bMat
= FALSE
, bPI
= FALSE
, bSplit
= FALSE
, bMax
= FALSE
, bPBC
= TRUE
;
666 static gmx_bool bGroup
= FALSE
;
667 static real rcutoff
= 0.6;
669 static gmx_bool bEachResEachTime
= FALSE
, bPrintResName
= FALSE
;
671 { "-matrix", FALSE
, etBOOL
, {&bMat
},
672 "Calculate half a matrix of group-group distances" },
673 { "-max", FALSE
, etBOOL
, {&bMax
},
674 "Calculate *maximum* distance instead of minimum" },
675 { "-d", FALSE
, etREAL
, {&rcutoff
},
676 "Distance for contacts" },
677 { "-group", FALSE
, etBOOL
, {&bGroup
},
678 "Count contacts with multiple atoms in the first group as one" },
679 { "-pi", FALSE
, etBOOL
, {&bPI
},
680 "Calculate minimum distance with periodic images" },
681 { "-split", FALSE
, etBOOL
, {&bSplit
},
682 "Split graph where time is zero" },
683 { "-ng", FALSE
, etINT
, {&ng
},
684 "Number of secondary groups to compute distance to a central group" },
685 { "-pbc", FALSE
, etBOOL
, {&bPBC
},
686 "Take periodic boundary conditions into account" },
687 { "-respertime", FALSE
, etBOOL
, {&bEachResEachTime
},
688 "When writing per-residue distances, write distance for each time point" },
689 { "-printresname", FALSE
, etBOOL
, {&bPrintResName
},
690 "Write residue names" }
693 t_topology
*top
= NULL
;
699 gmx_bool bTop
= FALSE
;
703 const char *trxfnm
, *tpsfnm
, *ndxfnm
, *distfnm
, *numfnm
, *atmfnm
, *oxfnm
, *resfnm
;
706 atom_id
**index
, *residues
= NULL
;
708 { efTRX
, "-f", NULL
, ffREAD
},
709 { efTPS
, NULL
, NULL
, ffOPTRD
},
710 { efNDX
, NULL
, NULL
, ffOPTRD
},
711 { efXVG
, "-od", "mindist", ffWRITE
},
712 { efXVG
, "-on", "numcont", ffOPTWR
},
713 { efOUT
, "-o", "atm-pair", ffOPTWR
},
714 { efTRO
, "-ox", "mindist", ffOPTWR
},
715 { efXVG
, "-or", "mindistres", ffOPTWR
}
717 #define NFILE asize(fnm)
719 if (!parse_common_args(&argc
, argv
,
720 PCA_CAN_VIEW
| PCA_CAN_TIME
| PCA_TIME_UNIT
,
721 NFILE
, fnm
, asize(pa
), pa
, asize(desc
), desc
, 0, NULL
, &oenv
))
726 trxfnm
= ftp2fn(efTRX
, NFILE
, fnm
);
727 ndxfnm
= ftp2fn_null(efNDX
, NFILE
, fnm
);
728 distfnm
= opt2fn("-od", NFILE
, fnm
);
729 numfnm
= opt2fn_null("-on", NFILE
, fnm
);
730 atmfnm
= ftp2fn_null(efOUT
, NFILE
, fnm
);
731 oxfnm
= opt2fn_null("-ox", NFILE
, fnm
);
732 resfnm
= opt2fn_null("-or", NFILE
, fnm
);
733 if (bPI
|| resfnm
!= NULL
)
735 /* We need a tps file */
736 tpsfnm
= ftp2fn(efTPS
, NFILE
, fnm
);
740 tpsfnm
= ftp2fn_null(efTPS
, NFILE
, fnm
);
743 if (!tpsfnm
&& !ndxfnm
)
745 gmx_fatal(FARGS
, "You have to specify either the index file or a tpr file");
751 fprintf(stderr
, "Choose a group for distance calculation\n");
762 if (tpsfnm
|| resfnm
|| !ndxfnm
)
765 bTop
= read_tps_conf(tpsfnm
, title
, top
, &ePBC
, &x
, NULL
, box
, FALSE
);
768 printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
771 get_index(top
? &(top
->atoms
) : NULL
, ndxfnm
, ng
, gnx
, index
, grpname
);
773 if (bMat
&& (ng
== 1))
776 printf("Special case: making distance matrix between all atoms in group %s\n",
781 for (i
= 1; (i
< ng
); i
++)
784 grpname
[i
] = grpname
[0];
786 index
[i
][0] = index
[0][i
];
793 nres
= find_residues(top
? &(top
->atoms
) : NULL
,
794 gnx
[0], index
[0], &residues
);
797 dump_res(debug
, nres
, residues
, index
[0]);
803 periodic_mindist_plot(trxfnm
, distfnm
, top
, ePBC
, gnx
[0], index
[0], bSplit
, oenv
);
807 dist_plot(trxfnm
, atmfnm
, distfnm
, numfnm
, resfnm
, oxfnm
,
808 rcutoff
, bMat
, top
? &(top
->atoms
) : NULL
,
809 ng
, index
, gnx
, grpname
, bSplit
, !bMax
, nres
, residues
, bPBC
, ePBC
,
810 bGroup
, bEachResEachTime
, bPrintResName
, oenv
);
813 do_view(oenv
, distfnm
, "-nxy");
816 do_view(oenv
, numfnm
, "-nxy");