Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / gmx_mindist.cpp
blobed51f851d9a415b1f73a8feec5181a1d4a6896ee
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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include <cmath>
41 #include <cstdlib>
42 #include <cstring>
44 #include <algorithm>
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/commandline/viewit.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/pbcutil/rmpbc.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/smalloc.h"
67 static void
68 periodic_dist(PbcType pbcType, matrix box, rvec x[], int n, const int index[], real* rmin, real* rmax, int* min_ind)
70 #define NSHIFT_MAX 26
71 int nsz, nshift, sx, sy, sz, i, j, s;
72 real sqr_box, r2min, r2max, r2;
73 rvec shift[NSHIFT_MAX], d0, d;
75 sqr_box = std::min(norm2(box[XX]), norm2(box[YY]));
76 if (pbcType == PbcType::Xyz)
78 sqr_box = std::min(sqr_box, norm2(box[ZZ]));
79 nsz = 1;
81 else if (pbcType == PbcType::XY)
83 nsz = 0;
85 else
87 gmx_fatal(FARGS, "pbc = %s is not supported by g_mindist", c_pbcTypeNames[pbcType].c_str());
90 nshift = 0;
91 for (sz = -nsz; sz <= nsz; sz++)
93 for (sy = -1; sy <= 1; sy++)
95 for (sx = -1; sx <= 1; sx++)
97 if (sx != 0 || sy != 0 || sz != 0)
99 for (i = 0; i < DIM; i++)
101 shift[nshift][i] = sx * box[XX][i] + sy * box[YY][i] + sz * box[ZZ][i];
103 nshift++;
109 r2min = sqr_box;
110 r2max = 0;
112 for (i = 0; i < n; i++)
114 for (j = i + 1; j < n; j++)
116 rvec_sub(x[index[i]], x[index[j]], d0);
117 r2 = norm2(d0);
118 if (r2 > r2max)
120 r2max = r2;
122 for (s = 0; s < nshift; s++)
124 rvec_add(d0, shift[s], d);
125 r2 = norm2(d);
126 if (r2 < r2min)
128 r2min = r2;
129 min_ind[0] = i;
130 min_ind[1] = j;
136 *rmin = std::sqrt(r2min);
137 *rmax = std::sqrt(r2max);
140 static void periodic_mindist_plot(const char* trxfn,
141 const char* outfn,
142 const t_topology* top,
143 PbcType pbcType,
144 int n,
145 int index[],
146 gmx_bool bSplit,
147 const gmx_output_env_t* oenv)
149 FILE* out;
150 const char* leg[5] = { "min per.", "max int.", "box1", "box2", "box3" };
151 t_trxstatus* status;
152 real t;
153 rvec* x;
154 matrix box;
155 int natoms, ind_min[2] = { 0, 0 }, ind_mini = 0, ind_minj = 0;
156 real rmin, rmax, rmint, tmint;
157 gmx_bool bFirst;
158 gmx_rmpbc_t gpbc = nullptr;
160 natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
162 check_index(nullptr, n, index, nullptr, natoms);
164 out = xvgropen(outfn, "Minimum distance to periodic image", output_env_get_time_label(oenv),
165 "Distance (nm)", oenv);
166 if (output_env_get_print_xvgr_codes(oenv))
168 fprintf(out, "@ subtitle \"and maximum internal distance\"\n");
170 xvgr_legend(out, 5, leg, oenv);
172 rmint = box[XX][XX];
173 tmint = 0;
175 if (nullptr != top)
177 gpbc = gmx_rmpbc_init(&top->idef, pbcType, natoms);
180 bFirst = TRUE;
183 if (nullptr != top)
185 gmx_rmpbc(gpbc, natoms, box, x);
188 periodic_dist(pbcType, box, x, n, index, &rmin, &rmax, ind_min);
189 if (rmin < rmint)
191 rmint = rmin;
192 tmint = t;
193 ind_mini = ind_min[0];
194 ind_minj = ind_min[1];
196 if (bSplit && !bFirst && std::abs(t / output_env_get_time_factor(oenv)) < 1e-5)
198 fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
200 fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n", output_env_conv_time(oenv, t), rmin,
201 rmax, norm(box[0]), norm(box[1]), norm(box[2]));
202 bFirst = FALSE;
203 } while (read_next_x(oenv, status, &t, x, box));
205 if (nullptr != top)
207 gmx_rmpbc_done(gpbc);
210 xvgrclose(out);
212 fprintf(stdout,
213 "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
214 "between atoms %d and %d\n",
215 rmint, output_env_conv_time(oenv, tmint), output_env_get_time_unit(oenv).c_str(),
216 index[ind_mini] + 1, index[ind_minj] + 1);
219 static void calc_dist(real rcut,
220 gmx_bool bPBC,
221 PbcType pbcType,
222 matrix box,
223 rvec x[],
224 int nx1,
225 int nx2,
226 int index1[],
227 int index2[],
228 gmx_bool bGroup,
229 real* rmin,
230 real* rmax,
231 int* nmin,
232 int* nmax,
233 int* ixmin,
234 int* jxmin,
235 int* ixmax,
236 int* jxmax)
238 int i, j, i0 = 0, j1;
239 int ix, jx;
240 int* index3;
241 rvec dx;
242 real r2, rmin2, rmax2, rcut2;
243 t_pbc pbc;
244 int nmin_j, nmax_j;
246 *ixmin = -1;
247 *jxmin = -1;
248 *ixmax = -1;
249 *jxmax = -1;
250 *nmin = 0;
251 *nmax = 0;
253 rcut2 = gmx::square(rcut);
255 /* Must init pbc every step because of pressure coupling */
256 if (bPBC)
258 set_pbc(&pbc, pbcType, box);
260 if (index2)
262 i0 = 0;
263 j1 = nx2;
264 index3 = index2;
266 else
268 j1 = nx1;
269 index3 = index1;
271 GMX_RELEASE_ASSERT(index1 != nullptr, "Need a valid index for plotting distances");
273 rmin2 = 1e12;
274 rmax2 = -1e12;
276 for (j = 0; (j < j1); j++)
278 jx = index3[j];
279 if (index2 == nullptr)
281 i0 = j + 1;
283 nmin_j = 0;
284 nmax_j = 0;
285 for (i = i0; (i < nx1); i++)
287 ix = index1[i];
288 if (ix != jx)
290 if (bPBC)
292 pbc_dx(&pbc, x[ix], x[jx], dx);
294 else
296 rvec_sub(x[ix], x[jx], dx);
298 r2 = iprod(dx, dx);
299 if (r2 < rmin2)
301 rmin2 = r2;
302 *ixmin = ix;
303 *jxmin = jx;
305 if (r2 > rmax2)
307 rmax2 = r2;
308 *ixmax = ix;
309 *jxmax = jx;
311 if (r2 <= rcut2)
313 nmin_j++;
315 else
317 nmax_j++;
321 if (bGroup)
323 if (nmin_j > 0)
325 (*nmin)++;
327 if (nmax_j > 0)
329 (*nmax)++;
332 else
334 *nmin += nmin_j;
335 *nmax += nmax_j;
338 *rmin = std::sqrt(rmin2);
339 *rmax = std::sqrt(rmax2);
342 static void dist_plot(const char* fn,
343 const char* afile,
344 const char* dfile,
345 const char* nfile,
346 const char* rfile,
347 const char* xfile,
348 real rcut,
349 gmx_bool bMat,
350 const t_atoms* atoms,
351 int ng,
352 int* index[],
353 int gnx[],
354 char* grpn[],
355 gmx_bool bSplit,
356 gmx_bool bMin,
357 int nres,
358 int* residue,
359 gmx_bool bPBC,
360 PbcType pbcType,
361 gmx_bool bGroup,
362 gmx_bool bEachResEachTime,
363 gmx_bool bPrintResName,
364 const gmx_output_env_t* oenv)
366 FILE * atm, *dist, *num;
367 t_trxstatus* trxout;
368 char buf[256];
369 char** leg;
370 real t, dmin, dmax, **mindres = nullptr, **maxdres = nullptr;
371 int nmin, nmax;
372 t_trxstatus* status;
373 int i = -1, j, k;
374 int min2, max2, min1r, min2r, max1r, max2r;
375 int min1 = 0;
376 int max1 = 0;
377 int oindex[2];
378 rvec* x0;
379 matrix box;
380 gmx_bool bFirst;
381 FILE* respertime = nullptr;
383 if (read_first_x(oenv, &status, fn, &t, &x0, box) == 0)
385 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
388 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
389 dist = xvgropen(dfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
390 sprintf(buf, "Number of Contacts %s %g nm", bMin ? "<" : ">", rcut);
391 num = nfile ? xvgropen(nfile, buf, output_env_get_time_label(oenv), "Number", oenv) : nullptr;
392 atm = afile ? gmx_ffopen(afile, "w") : nullptr;
393 trxout = xfile ? open_trx(xfile, "w") : nullptr;
395 if (bMat)
397 if (ng == 1)
399 snew(leg, 1);
400 sprintf(buf, "Internal in %s", grpn[0]);
401 leg[0] = gmx_strdup(buf);
402 xvgr_legend(dist, 0, leg, oenv);
403 if (num)
405 xvgr_legend(num, 0, leg, oenv);
408 else
410 GMX_RELEASE_ASSERT(ng > 1, "Must have more than one group with bMat");
411 snew(leg, (ng * (ng - 1)) / 2);
412 for (i = j = 0; (i < ng - 1); i++)
414 for (k = i + 1; (k < ng); k++, j++)
416 sprintf(buf, "%s-%s", grpn[i], grpn[k]);
417 leg[j] = gmx_strdup(buf);
420 xvgr_legend(dist, j, leg, oenv);
421 if (num)
423 xvgr_legend(num, j, leg, oenv);
427 else
429 snew(leg, ng - 1);
430 for (i = 0; (i < ng - 1); i++)
432 sprintf(buf, "%s-%s", grpn[0], grpn[i + 1]);
433 leg[i] = gmx_strdup(buf);
435 xvgr_legend(dist, ng - 1, leg, oenv);
436 if (num)
438 xvgr_legend(num, ng - 1, leg, oenv);
442 if (bEachResEachTime)
444 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
445 respertime = xvgropen(rfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
446 xvgr_legend(respertime, ng - 1, leg, oenv);
447 if (bPrintResName && output_env_get_print_xvgr_codes(oenv))
449 fprintf(respertime, "# ");
451 for (j = 0; j < nres; j++)
453 fprintf(respertime, "%s%d ",
454 *(atoms->resinfo[atoms->atom[index[0][residue[j]]].resind].name),
455 atoms->atom[index[0][residue[j]]].resind);
457 fprintf(respertime, "\n");
461 if (nres)
463 snew(mindres, ng - 1);
464 snew(maxdres, ng - 1);
465 for (i = 1; i < ng; i++)
467 snew(mindres[i - 1], nres);
468 snew(maxdres[i - 1], nres);
469 for (j = 0; j < nres; j++)
471 mindres[i - 1][j] = 1e6;
473 /* maxdres[*][*] is already 0 */
476 bFirst = TRUE;
479 if (bSplit && !bFirst && std::abs(t / output_env_get_time_factor(oenv)) < 1e-5)
481 fprintf(dist, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
482 if (num)
484 fprintf(num, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
486 if (atm)
488 fprintf(atm, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
491 fprintf(dist, "%12e", output_env_conv_time(oenv, t));
492 if (num)
494 fprintf(num, "%12e", output_env_conv_time(oenv, t));
497 if (bMat)
499 if (ng == 1)
501 calc_dist(rcut, bPBC, pbcType, box, x0, gnx[0], gnx[0], index[0], index[0], bGroup,
502 &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
503 fprintf(dist, " %12e", bMin ? dmin : dmax);
504 if (num)
506 fprintf(num, " %8d", bMin ? nmin : nmax);
509 else
511 for (i = 0; (i < ng - 1); i++)
513 for (k = i + 1; (k < ng); k++)
515 calc_dist(rcut, bPBC, pbcType, box, x0, gnx[i], gnx[k], index[i], index[k],
516 bGroup, &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
517 fprintf(dist, " %12e", bMin ? dmin : dmax);
518 if (num)
520 fprintf(num, " %8d", bMin ? nmin : nmax);
526 else
528 GMX_RELEASE_ASSERT(ng > 1, "Must have more than one group when not using -matrix");
529 for (i = 1; (i < ng); i++)
531 calc_dist(rcut, bPBC, pbcType, box, x0, gnx[0], gnx[i], index[0], index[i], bGroup,
532 &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
533 fprintf(dist, " %12e", bMin ? dmin : dmax);
534 if (num)
536 fprintf(num, " %8d", bMin ? nmin : nmax);
538 if (nres)
540 for (j = 0; j < nres; j++)
542 calc_dist(rcut, bPBC, pbcType, box, x0, residue[j + 1] - residue[j], gnx[i],
543 &(index[0][residue[j]]), index[i], bGroup, &dmin, &dmax, &nmin,
544 &nmax, &min1r, &min2r, &max1r, &max2r);
545 mindres[i - 1][j] = std::min(mindres[i - 1][j], dmin);
546 maxdres[i - 1][j] = std::max(maxdres[i - 1][j], dmax);
551 fprintf(dist, "\n");
552 if (num)
554 fprintf(num, "\n");
556 if ((bMin ? min1 : max1) != -1)
558 if (atm)
560 fprintf(atm, "%12e %12d %12d\n", output_env_conv_time(oenv, t),
561 1 + (bMin ? min1 : max1), 1 + (bMin ? min2 : max2));
565 if (trxout)
567 oindex[0] = bMin ? min1 : max1;
568 oindex[1] = bMin ? min2 : max2;
569 write_trx(trxout, 2, oindex, atoms, i, t, box, x0, nullptr, nullptr);
571 bFirst = FALSE;
572 /*dmin should be minimum distance for residue and group*/
573 if (bEachResEachTime)
575 fprintf(respertime, "%12e", t);
576 for (i = 1; i < ng; i++)
578 for (j = 0; j < nres; j++)
580 fprintf(respertime, " %7g", bMin ? mindres[i - 1][j] : maxdres[i - 1][j]);
581 /*reset distances for next time point*/
582 mindres[i - 1][j] = 1e6;
583 maxdres[i - 1][j] = 0;
586 fprintf(respertime, "\n");
588 } while (read_next_x(oenv, status, &t, x0, box));
590 close_trx(status);
591 xvgrclose(dist);
592 if (num)
594 xvgrclose(num);
596 if (atm)
598 gmx_ffclose(atm);
600 if (trxout)
602 close_trx(trxout);
604 if (respertime)
606 xvgrclose(respertime);
609 if (nres && !bEachResEachTime)
611 FILE* res;
613 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
614 res = xvgropen(rfile, buf, "Residue (#)", "Distance (nm)", oenv);
615 xvgr_legend(res, ng - 1, leg, oenv);
616 for (j = 0; j < nres; j++)
618 fprintf(res, "%4d", j + 1);
619 for (i = 1; i < ng; i++)
621 fprintf(res, " %7g", bMin ? mindres[i - 1][j] : maxdres[i - 1][j]);
623 fprintf(res, "\n");
625 xvgrclose(res);
628 if (x0)
630 sfree(x0);
633 int freeLeg = bMat ? (ng == 1 ? 1 : (ng * (ng - 1)) / 2) : ng - 1;
634 for (int i = 0; i < freeLeg; i++)
636 sfree(leg[i]);
638 sfree(leg);
641 static int find_residues(const t_atoms* atoms, int n, const int index[], int** resindex)
643 int i;
644 int nres = 0, resnr, presnr = 0;
645 bool presFound = false;
646 int* residx;
648 /* build index of first atom numbers for each residue */
649 snew(residx, atoms->nres + 1);
650 for (i = 0; i < n; i++)
652 resnr = atoms->atom[index[i]].resind;
653 if (!presFound || resnr != presnr)
655 residx[nres] = i;
656 nres++;
657 presnr = resnr;
658 presFound = true;
661 if (debug)
663 printf("Found %d residues out of %d (%d/%d atoms)\n", nres, atoms->nres, atoms->nr, n);
665 srenew(residx, nres + 1);
666 /* mark end of last residue */
667 residx[nres] = n;
668 *resindex = residx;
669 return nres;
672 static void dump_res(FILE* out, int nres, int* resindex, int index[])
674 int i, j;
676 for (i = 0; i < nres - 1; i++)
678 fprintf(out, "Res %d (%d):", i, resindex[i + 1] - resindex[i]);
679 for (j = resindex[i]; j < resindex[i + 1]; j++)
681 fprintf(out, " %d(%d)", j, index[j]);
683 fprintf(out, "\n");
687 int gmx_mindist(int argc, char* argv[])
689 const char* desc[] = {
690 "[THISMODULE] computes the distance between one group and a number of",
691 "other groups. Both the minimum distance",
692 "(between any pair of atoms from the respective groups)",
693 "and the number of contacts within a given",
694 "distance are written to two separate output files.",
695 "With the [TT]-group[tt] option a contact of an atom in another group",
696 "with multiple atoms in the first group is counted as one contact",
697 "instead of as multiple contacts.",
698 "With [TT]-or[tt], minimum distances to each residue in the first",
699 "group are determined and plotted as a function of residue number.[PAR]",
700 "With option [TT]-pi[tt] the minimum distance of a group to its",
701 "periodic image is plotted. This is useful for checking if a protein",
702 "has seen its periodic image during a simulation. Only one shift in",
703 "each direction is considered, giving a total of 26 shifts. Note",
704 "that periodicity information is required from the file supplied with",
705 "with [TT]-s[tt], either as a .tpr file or a .pdb file with CRYST1 fields.",
706 "It also plots the maximum distance within the group and the lengths",
707 "of the three box vectors.[PAR]",
708 "Also [gmx-distance] and [gmx-pairdist] calculate distances."
711 gmx_bool bMat = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
712 gmx_bool bGroup = FALSE;
713 real rcutoff = 0.6;
714 int ng = 1;
715 gmx_bool bEachResEachTime = FALSE, bPrintResName = FALSE;
716 t_pargs pa[] = {
717 { "-matrix", FALSE, etBOOL, { &bMat }, "Calculate half a matrix of group-group distances" },
718 { "-max", FALSE, etBOOL, { &bMax }, "Calculate *maximum* distance instead of minimum" },
719 { "-d", FALSE, etREAL, { &rcutoff }, "Distance for contacts" },
720 { "-group",
721 FALSE,
722 etBOOL,
723 { &bGroup },
724 "Count contacts with multiple atoms in the first group as one" },
725 { "-pi", FALSE, etBOOL, { &bPI }, "Calculate minimum distance with periodic images" },
726 { "-split", FALSE, etBOOL, { &bSplit }, "Split graph where time is zero" },
727 { "-ng",
728 FALSE,
729 etINT,
730 { &ng },
731 "Number of secondary groups to compute distance to a central group" },
732 { "-pbc", FALSE, etBOOL, { &bPBC }, "Take periodic boundary conditions into account" },
733 { "-respertime",
734 FALSE,
735 etBOOL,
736 { &bEachResEachTime },
737 "When writing per-residue distances, write distance for each time point" },
738 { "-printresname", FALSE, etBOOL, { &bPrintResName }, "Write residue names" }
740 gmx_output_env_t* oenv;
741 t_topology* top = nullptr;
742 PbcType pbcType = PbcType::Unset;
743 rvec* x = nullptr;
744 matrix box;
745 gmx_bool bTop = FALSE;
747 int i, nres = 0;
748 const char *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
749 char** grpname;
750 int* gnx;
751 int ** index, *residues = nullptr;
752 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffOPTRD },
753 { efNDX, nullptr, nullptr, ffOPTRD }, { efXVG, "-od", "mindist", ffWRITE },
754 { efXVG, "-on", "numcont", ffOPTWR }, { efOUT, "-o", "atm-pair", ffOPTWR },
755 { efTRO, "-ox", "mindist", ffOPTWR }, { efXVG, "-or", "mindistres", ffOPTWR } };
756 #define NFILE asize(fnm)
758 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT, NFILE, fnm,
759 asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
761 return 0;
764 trxfnm = ftp2fn(efTRX, NFILE, fnm);
765 ndxfnm = ftp2fn_null(efNDX, NFILE, fnm);
766 distfnm = opt2fn("-od", NFILE, fnm);
767 numfnm = opt2fn_null("-on", NFILE, fnm);
768 atmfnm = ftp2fn_null(efOUT, NFILE, fnm);
769 oxfnm = opt2fn_null("-ox", NFILE, fnm);
770 resfnm = opt2fn_null("-or", NFILE, fnm);
771 if (bPI || resfnm != nullptr)
773 /* We need a tps file */
774 tpsfnm = ftp2fn(efTPS, NFILE, fnm);
776 else
778 tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
781 if (!tpsfnm && !ndxfnm)
783 gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
786 if (bPI)
788 ng = 1;
789 fprintf(stderr, "Choose a group for distance calculation\n");
791 else if (!bMat)
793 ng++;
796 snew(gnx, ng);
797 snew(index, ng);
798 snew(grpname, ng);
800 if (tpsfnm || resfnm || !ndxfnm)
802 snew(top, 1);
803 bTop = read_tps_conf(tpsfnm, top, &pbcType, &x, nullptr, box, FALSE);
804 if (bPI && !bTop)
806 printf("\nWARNING: Without a run input file a trajectory with broken molecules will "
807 "not give the correct periodic image distance\n\n");
810 get_index(top ? &(top->atoms) : nullptr, ndxfnm, ng, gnx, index, grpname);
812 if (bMat && (ng == 1))
814 ng = gnx[0];
815 printf("Special case: making distance matrix between all atoms in group %s\n", grpname[0]);
816 srenew(gnx, ng);
817 srenew(index, ng);
818 srenew(grpname, ng);
819 for (i = 1; (i < ng); i++)
821 gnx[i] = 1;
822 grpname[i] = grpname[0];
823 snew(index[i], 1);
824 index[i][0] = index[0][i];
826 gnx[0] = 1;
828 GMX_RELEASE_ASSERT(!bMat || ng > 1, "Must have more than one group with bMat");
830 if (resfnm)
832 GMX_RELEASE_ASSERT(top != nullptr, "top pointer cannot be NULL when finding residues");
833 nres = find_residues(&(top->atoms), gnx[0], index[0], &residues);
835 if (debug)
837 dump_res(debug, nres, residues, index[0]);
840 else if (bEachResEachTime || bPrintResName)
842 gmx_fatal(FARGS, "Option -or needs to be set to print residues");
845 if (bPI)
847 periodic_mindist_plot(trxfnm, distfnm, top, pbcType, gnx[0], index[0], bSplit, oenv);
849 else
851 dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm, rcutoff, bMat,
852 top ? &(top->atoms) : nullptr, ng, index, gnx, grpname, bSplit, !bMax, nres,
853 residues, bPBC, pbcType, bGroup, bEachResEachTime, bPrintResName, oenv);
856 do_view(oenv, distfnm, "-nxy");
857 if (!bPI)
859 do_view(oenv, numfnm, "-nxy");
862 output_env_done(oenv);
863 done_top(top);
864 for (int i = 0; i < ng; i++)
866 sfree(index[i]);
868 sfree(index);
869 sfree(gnx);
870 sfree(x);
871 sfree(grpname);
872 sfree(top);
874 return 0;