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, 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.
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/matio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxana/cmat.h"
52 #include "gromacs/gmxana/gmx_ana.h"
53 #include "gromacs/linearalgebra/eigensolver.h"
54 #include "gromacs/math/do_fit.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/pbcutil/rmpbc.h"
58 #include "gromacs/random/threefry.h"
59 #include "gromacs/random/uniformintdistribution.h"
60 #include "gromacs/random/uniformrealdistribution.h"
61 #include "gromacs/topology/index.h"
62 #include "gromacs/topology/topology.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/futil.h"
67 #include "gromacs/utility/smalloc.h"
68 #include "gromacs/utility/stringutil.h"
70 /* print to two file pointers at once (i.e. stderr and log) */
72 void lo_ffprintf(FILE *fp1, FILE *fp2, const char *buf)
74 fprintf(fp1, "%s", buf);
75 fprintf(fp2, "%s", buf);
78 /* just print a prepared buffer to fp1 and fp2 */
80 void ffprintf(FILE *fp1, FILE *fp2, const char *buf)
82 lo_ffprintf(fp1, fp2, buf);
85 /* prepare buffer with one argument, then print to fp1 and fp2 */
87 void ffprintf_d(FILE *fp1, FILE *fp2, char *buf, const char *fmt, int arg)
89 sprintf(buf, fmt, arg);
90 lo_ffprintf(fp1, fp2, buf);
93 /* prepare buffer with one argument, then print to fp1 and fp2 */
95 void ffprintf_g(FILE *fp1, FILE *fp2, char *buf, const char *fmt, real arg)
97 sprintf(buf, fmt, arg);
98 lo_ffprintf(fp1, fp2, buf);
101 /* prepare buffer with one argument, then print to fp1 and fp2 */
103 void ffprintf_s(FILE *fp1, FILE *fp2, char *buf, const char *fmt, const char *arg)
105 sprintf(buf, fmt, arg);
106 lo_ffprintf(fp1, fp2, buf);
109 /* prepare buffer with two arguments, then print to fp1 and fp2 */
111 void ffprintf_dd(FILE *fp1, FILE *fp2, char *buf, const char *fmt, int arg1, int arg2)
113 sprintf(buf, fmt, arg1, arg2);
114 lo_ffprintf(fp1, fp2, buf);
117 /* prepare buffer with two arguments, then print to fp1 and fp2 */
119 void ffprintf_gg(FILE *fp1, FILE *fp2, char *buf, const char *fmt, real arg1, real arg2)
121 sprintf(buf, fmt, arg1, arg2);
122 lo_ffprintf(fp1, fp2, buf);
125 /* prepare buffer with two arguments, then print to fp1 and fp2 */
127 void ffprintf_ss(FILE *fp1, FILE *fp2, char *buf, const char *fmt, const char *arg1, const char *arg2)
129 sprintf(buf, fmt, arg1, arg2);
130 lo_ffprintf(fp1, fp2, buf);
143 static void mc_optimize(FILE *log, t_mat *m, real *time,
144 int maxiter, int nrandom,
146 const char *conv, gmx_output_env_t *oenv)
149 real ecur, enext, emin, prob, enorm;
150 int i, j, iswap, jswap, nn, nuphill = 0;
155 seed = static_cast<int>(gmx::makeRandomSeed());
157 gmx::DefaultRandomEngine rng(seed);
161 fprintf(stderr, "Can not do Monte Carlo optimization with a non-square matrix.\n");
164 printf("\nDoing Monte Carlo optimization to find the smoothest trajectory\n");
165 printf("by reordering the frames to minimize the path between the two structures\n");
166 printf("that have the largest pairwise RMSD.\n");
167 printf("Using random seed %d.\n", seed);
170 enorm = m->mat[0][0];
171 for (i = 0; (i < m->n1); i++)
173 for (j = 0; (j < m->nn); j++)
175 if (m->mat[i][j] > enorm)
177 enorm = m->mat[i][j];
183 if ((iswap == -1) || (jswap == -1))
185 fprintf(stderr, "Matrix contains identical values in all fields\n");
188 swap_rows(m, 0, iswap);
189 swap_rows(m, m->n1-1, jswap);
190 emin = ecur = mat_energy(m);
191 printf("Largest distance %g between %d and %d. Energy: %g.\n",
192 enorm, iswap, jswap, emin);
196 /* Initiate and store global minimum */
197 minimum = init_mat(nn, m->b1D);
199 copy_t_mat(minimum, m);
203 fp = xvgropen(conv, "Convergence of the MC optimization",
204 "Energy", "Step", oenv);
207 gmx::UniformIntDistribution<int> intDistNN(1, nn-2); // [1,nn-2]
208 gmx::UniformRealDistribution<real> realDistOne; // [0,1)
210 for (i = 0; (i < maxiter); i++)
212 /* Generate new swapping candidates */
215 iswap = intDistNN(rng);
216 jswap = intDistNN(rng);
218 while ((iswap == jswap) || (iswap >= nn-1) || (jswap >= nn-1));
220 /* Apply swap and compute energy */
221 swap_rows(m, iswap, jswap);
222 enext = mat_energy(m);
224 /* Compute probability */
226 if ((enext < ecur) || (i < nrandom))
231 /* Store global minimum */
232 copy_t_mat(minimum, m);
238 /* Try Monte Carlo step */
239 prob = std::exp(-(enext-ecur)/(enorm*kT));
242 if (prob == 1 || realDistOne(rng) < prob)
249 fprintf(log, "Iter: %d Swapped %4d and %4d (energy: %g prob: %g)\n",
250 i, iswap, jswap, enext, prob);
253 fprintf(fp, "%6d %10g\n", i, enext);
259 swap_rows(m, jswap, iswap);
262 fprintf(log, "%d uphill steps were taken during optimization\n",
265 /* Now swap the matrix to get it into global minimum mode */
266 copy_t_mat(m, minimum);
268 fprintf(log, "Global minimum energy %g\n", mat_energy(minimum));
269 fprintf(log, "Global minimum energy %g\n", mat_energy(m));
270 fprintf(log, "Swapped time and frame indices and RMSD to next neighbor:\n");
271 for (i = 0; (i < m->nn); i++)
273 fprintf(log, "%10g %5d %10g\n",
276 (i < m->nn-1) ? m->mat[m->m_ind[i]][m->m_ind[i+1]] : 0);
285 static void calc_dist(int nind, rvec x[], real **d)
291 for (i = 0; (i < nind-1); i++)
294 for (j = i+1; (j < nind); j++)
296 /* Should use pbc_dx when analysing multiple molecueles,
297 * but the box is not stored for every frame.
299 rvec_sub(xi, x[j], dx);
305 static real rms_dist(int isize, real **d, real **d_r)
311 for (i = 0; (i < isize-1); i++)
313 for (j = i+1; (j < isize); j++)
315 r = d[i][j]-d_r[i][j];
319 r2 /= gmx::exactDiv(isize*(isize-1), 2);
321 return std::sqrt(r2);
324 static bool rms_dist_comp(const t_dist &a, const t_dist &b)
326 return a.dist < b.dist;
329 static bool clust_id_comp(const t_clustid &a, const t_clustid &b)
331 return a.clust < b.clust;
334 static bool nrnb_comp(const t_nnb &a, const t_nnb &b)
336 /* return b<a, we want highest first */
340 static void gather(t_mat *m, real cutoff, t_clusters *clust)
344 int i, j, k, nn, cid, n1, diff;
347 /* First we sort the entries in the RMSD matrix */
351 for (i = k = 0; (i < n1); i++)
353 for (j = i+1; (j < n1); j++, k++)
357 d[k].dist = m->mat[i][j];
362 gmx_incons("gather algortihm");
364 std::sort(d, d+nn, rms_dist_comp);
366 /* Now we make a cluster index for all of the conformations */
369 /* Now we check the closest structures, and equalize their cluster numbers */
370 fprintf(stderr, "Linking structures ");
373 fprintf(stderr, "*");
375 for (k = 0; (k < nn) && (d[k].dist < cutoff); k++)
377 diff = c[d[k].j].clust - c[d[k].i].clust;
383 c[d[k].j].clust = c[d[k].i].clust;
387 c[d[k].i].clust = c[d[k].j].clust;
393 fprintf(stderr, "\nSorting and renumbering clusters\n");
394 /* Sort on cluster number */
395 std::sort(c, c+n1, clust_id_comp);
397 /* Renumber clusters */
399 for (k = 1; k < n1; k++)
401 if (c[k].clust != c[k-1].clust)
414 for (k = 0; (k < n1); k++)
416 fprintf(debug, "Cluster index for conformation %d: %d\n",
417 c[k].conf, c[k].clust);
421 for (k = 0; k < n1; k++)
423 clust->cl[c[k].conf] = c[k].clust;
430 static gmx_bool jp_same(int **nnb, int i, int j, int P)
436 for (k = 0; nnb[i][k] >= 0; k++)
438 bIn = bIn || (nnb[i][k] == j);
446 for (k = 0; nnb[j][k] >= 0; k++)
448 bIn = bIn || (nnb[j][k] == i);
456 for (ii = 0; nnb[i][ii] >= 0; ii++)
458 for (jj = 0; nnb[j][jj] >= 0; jj++)
460 if ((nnb[i][ii] == nnb[j][jj]) && (nnb[i][ii] != -1))
470 static void jarvis_patrick(int n1, real **mat, int M, int P,
471 real rmsdcut, t_clusters *clust)
476 int i, j, k, cid, diff, maxval;
478 real **mcpy = nullptr;
485 /* First we sort the entries in the RMSD matrix row by row.
486 * This gives us the nearest neighbor list.
490 for (i = 0; (i < n1); i++)
492 for (j = 0; (j < n1); j++)
495 row[j].dist = mat[i][j];
497 std::sort(row, row+n1, rms_dist_comp);
500 /* Put the M nearest neighbors in the list */
502 for (j = k = 0; (k < M) && (j < n1) && (mat[i][row[j].j] < rmsdcut); j++)
506 nnb[i][k] = row[j].j;
514 /* Put all neighbors nearer than rmsdcut in the list */
517 for (j = 0; (j < n1) && (mat[i][row[j].j] < rmsdcut); j++)
524 srenew(nnb[i], maxval);
526 nnb[i][k] = row[j].j;
532 srenew(nnb[i], maxval+1);
540 fprintf(debug, "Nearest neighborlist. M = %d, P = %d\n", M, P);
541 for (i = 0; (i < n1); i++)
543 fprintf(debug, "i:%5d nbs:", i);
544 for (j = 0; nnb[i][j] >= 0; j++)
546 fprintf(debug, "%5d[%5.3f]", nnb[i][j], mat[i][nnb[i][j]]);
548 fprintf(debug, "\n");
553 fprintf(stderr, "Linking structures ");
554 /* Use mcpy for temporary storage of booleans */
555 mcpy = mk_matrix(n1, n1, FALSE);
556 for (i = 0; i < n1; i++)
558 for (j = i+1; j < n1; j++)
560 mcpy[i][j] = jp_same(nnb, i, j, P);
565 fprintf(stderr, "*");
567 for (i = 0; i < n1; i++)
569 for (j = i+1; j < n1; j++)
573 diff = c[j].clust - c[i].clust;
579 c[j].clust = c[i].clust;
583 c[i].clust = c[j].clust;
592 fprintf(stderr, "\nSorting and renumbering clusters\n");
593 /* Sort on cluster number */
594 std::sort(c, c+n1, clust_id_comp);
596 /* Renumber clusters */
598 for (k = 1; k < n1; k++)
600 if (c[k].clust != c[k-1].clust)
612 for (k = 0; k < n1; k++)
614 clust->cl[c[k].conf] = c[k].clust;
618 for (k = 0; (k < n1); k++)
620 fprintf(debug, "Cluster index for conformation %d: %d\n",
621 c[k].conf, c[k].clust);
625 done_matrix(n1, &mcpy);
628 for (i = 0; (i < n1); i++)
635 static void dump_nnb (FILE *fp, const char *title, int n1, t_nnb *nnb)
639 /* dump neighbor list */
640 fprintf(fp, "%s", title);
641 for (i = 0; (i < n1); i++)
643 fprintf(fp, "i:%5d #:%5d nbs:", i, nnb[i].nr);
644 for (j = 0; j < nnb[i].nr; j++)
646 fprintf(fp, "%5d", nnb[i].nb[j]);
652 static void gromos(int n1, real **mat, real rmsdcut, t_clusters *clust)
656 int i, j, k, j1, maxval;
658 /* Put all neighbors nearer than rmsdcut in the list */
659 fprintf(stderr, "Making list of neighbors within cutoff ");
662 for (i = 0; (i < n1); i++)
666 /* put all neighbors within cut-off in list */
667 for (j = 0; j < n1; j++)
669 if (mat[i][j] < rmsdcut)
674 srenew(nnb[i].nb, maxval);
680 /* store nr of neighbors, we'll need that */
682 if (i%(1+n1/100) == 0)
684 fprintf(stderr, "%3d%%\b\b\b\b", (i*100+1)/n1);
687 fprintf(stderr, "%3d%%\n", 100);
690 /* sort neighbor list on number of neighbors, largest first */
691 std::sort(nnb, nnb+n1, nrnb_comp);
695 dump_nnb(debug, "Nearest neighborlist after sort.\n", n1, nnb);
698 /* turn first structure with all its neighbors (largest) into cluster
699 remove them from pool of structures and repeat for all remaining */
700 fprintf(stderr, "Finding clusters %4d", 0);
701 /* cluster id's start at 1: */
705 /* set cluster id (k) for first item in neighborlist */
706 for (j = 0; j < nnb[0].nr; j++)
708 clust->cl[nnb[0].nb[j]] = k;
714 /* adjust number of neighbors for others, taking removals into account: */
715 for (i = 1; i < n1 && nnb[i].nr; i++)
718 for (j = 0; j < nnb[i].nr; j++)
720 /* if this neighbor wasn't removed */
721 if (clust->cl[nnb[i].nb[j]] == 0)
723 /* shift the rest (j1<=j) */
724 nnb[i].nb[j1] = nnb[i].nb[j];
729 /* now j1 is the new number of neighbors */
732 /* sort again on nnb[].nr, because we have new # neighbors: */
733 /* but we only need to sort upto i, i.e. when nnb[].nr>0 */
734 std::sort(nnb, nnb+i, nrnb_comp);
736 fprintf(stderr, "\b\b\b\b%4d", k);
740 fprintf(stderr, "\n");
744 fprintf(debug, "Clusters (%d):\n", k);
745 for (i = 0; i < n1; i++)
747 fprintf(debug, " %3d", clust->cl[i]);
749 fprintf(debug, "\n");
755 static rvec **read_whole_trj(const char *fn, int isize, const int index[], int skip,
756 int *nframe, real **time, const gmx_output_env_t *oenv, gmx_bool bPBC, gmx_rmpbc_t gpbc)
761 int i, i0, j, max_nf;
769 natom = read_first_x(oenv, &status, fn, &t, &x, box);
776 gmx_rmpbc(gpbc, natom, box, x);
782 srenew(*time, max_nf);
787 /* Store only the interesting atoms */
788 for (j = 0; (j < isize); j++)
790 copy_rvec(x[index[j]], xx[i0][j]);
797 while (read_next_x(oenv, status, &t, x, box));
798 fprintf(stderr, "Allocated %zu bytes for frames\n",
799 (max_nf*isize*sizeof(**xx)));
800 fprintf(stderr, "Read %d frames from trajectory %s\n", i0, fn);
807 static int plot_clusters(int nf, real **mat, t_clusters *clust,
810 int i, j, ncluster, ci;
811 int *cl_id, *nstruct, *strind;
816 for (i = 0; i < nf; i++)
819 cl_id[i] = clust->cl[i];
823 for (i = 0; i < nf; i++)
825 if (nstruct[i] >= minstruct)
828 for (j = 0; (j < nf); j++)
832 strind[j] = ncluster;
838 fprintf(stderr, "There are %d clusters with at least %d conformations\n",
839 ncluster, minstruct);
841 for (i = 0; (i < nf); i++)
844 for (j = 0; j < i; j++)
846 if ((ci == cl_id[j]) && (nstruct[ci] >= minstruct))
848 /* color different clusters with different colors, as long as
849 we don't run out of colors */
850 mat[i][j] = strind[i];
865 static void mark_clusters(int nf, real **mat, real val, t_clusters *clust)
869 for (i = 0; i < nf; i++)
871 for (j = 0; j < i; j++)
873 if (clust->cl[i] == clust->cl[j])
885 static char *parse_filename(const char *fn, int maxnr)
892 if (std::strchr(fn, '%'))
894 gmx_fatal(FARGS, "will not number filename %s containing '%c'", fn, '%');
896 /* number of digits needed in numbering */
897 i = static_cast<int>((std::log(static_cast<real>(maxnr))/std::log(10.0)) + 1);
898 /* split fn and ext */
899 ext = std::strrchr(fn, '.');
902 gmx_fatal(FARGS, "cannot separate extension in filename %s", fn);
905 /* insert e.g. '%03d' between fn and ext */
906 sprintf(buf, "%s%%0%dd.%s", fn, i, ext);
907 snew(fnout, std::strlen(buf)+1);
908 std::strcpy(fnout, buf);
913 static void ana_trans(t_clusters *clust, int nf,
914 const char *transfn, const char *ntransfn, FILE *log,
915 t_rgb rlo, t_rgb rhi, const gmx_output_env_t *oenv)
920 int i, ntranst, maxtrans;
923 snew(ntrans, clust->ncl);
924 snew(trans, clust->ncl);
925 snew(axis, clust->ncl);
926 for (i = 0; i < clust->ncl; i++)
929 snew(trans[i], clust->ncl);
933 for (i = 1; i < nf; i++)
935 if (clust->cl[i] != clust->cl[i-1])
938 ntrans[clust->cl[i-1]-1]++;
939 ntrans[clust->cl[i]-1]++;
940 trans[clust->cl[i-1]-1][clust->cl[i]-1]++;
941 maxtrans = static_cast<int>(std::max(static_cast<real>(maxtrans), trans[clust->cl[i]-1][clust->cl[i-1]-1]));
944 ffprintf_dd(stderr, log, buf, "Counted %d transitions in total, "
945 "max %d between two specific clusters\n", ntranst, maxtrans);
948 fp = gmx_ffopen(transfn, "w");
949 i = std::min(maxtrans+1, 80);
950 write_xpm(fp, 0, "Cluster Transitions", "# transitions",
951 "from cluster", "to cluster",
952 clust->ncl, clust->ncl, axis, axis, trans,
953 0, maxtrans, rlo, rhi, &i);
958 fp = xvgropen(ntransfn, "Cluster Transitions", "Cluster #", "# transitions",
960 for (i = 0; i < clust->ncl; i++)
962 fprintf(fp, "%5d %5d\n", i+1, ntrans[i]);
967 for (i = 0; i < clust->ncl; i++)
975 static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
976 int natom, t_atoms *atoms, rvec *xtps,
977 real *mass, rvec **xx, real *time,
978 int ifsize, int *fitidx,
979 int iosize, int *outidx,
980 const char *trxfn, const char *sizefn,
981 const char *transfn, const char *ntransfn,
982 const char *clustidfn, gmx_bool bAverage,
983 int write_ncl, int write_nst, real rmsmin,
984 gmx_bool bFit, FILE *log, t_rgb rlo, t_rgb rhi,
985 const gmx_output_env_t *oenv)
987 FILE *size_fp = nullptr;
988 char buf[STRLEN], buf1[40], buf2[40], buf3[40], *trxsfn;
989 t_trxstatus *trxout = nullptr;
990 t_trxstatus *trxsout = nullptr;
991 int i, i1, cl, nstr, *structure, first = 0, midstr;
992 gmx_bool *bWrite = nullptr;
993 real r, clrmsd, midrmsd;
999 ffprintf_d(stderr, log, buf, "\nFound %d clusters\n\n", clust->ncl);
1003 /* do we write all structures? */
1006 trxsfn = parse_filename(trxfn, std::max(write_ncl, clust->ncl));
1009 ffprintf_ss(stderr, log, buf, "Writing %s structure for each cluster to %s\n",
1010 bAverage ? "average" : "middle", trxfn);
1013 /* find out what we want to tell the user:
1014 Writing [all structures|structures with rmsd > %g] for
1015 {all|first %d} clusters {with more than %d structures} to %s */
1018 sprintf(buf1, "structures with rmsd > %g", rmsmin);
1022 sprintf(buf1, "all structures");
1024 buf2[0] = buf3[0] = '\0';
1025 if (write_ncl >= clust->ncl)
1029 sprintf(buf2, "all ");
1034 sprintf(buf2, "the first %d ", write_ncl);
1038 sprintf(buf3, " with more than %d structures", write_nst);
1040 sprintf(buf, "Writing %s for %sclusters%s to %s\n", buf1, buf2, buf3, trxsfn);
1041 ffprintf(stderr, log, buf);
1044 /* Prepare a reference structure for the orientation of the clusters */
1047 reset_x(ifsize, fitidx, natom, nullptr, xtps, mass);
1049 trxout = open_trx(trxfn, "w");
1050 /* Calculate the average structure in each cluster, *
1051 * all structures are fitted to the first struture of the cluster */
1053 GMX_ASSERT(xav, "");
1056 if (transfn || ntransfn)
1058 ana_trans(clust, nf, transfn, ntransfn, log, rlo, rhi, oenv);
1063 FILE *fp = xvgropen(clustidfn, "Clusters", output_env_get_xvgr_tlabel(oenv), "Cluster #", oenv);
1064 if (output_env_get_print_xvgr_codes(oenv))
1066 fprintf(fp, "@ s0 symbol 2\n");
1067 fprintf(fp, "@ s0 symbol size 0.2\n");
1068 fprintf(fp, "@ s0 linestyle 0\n");
1070 for (i = 0; i < nf; i++)
1072 fprintf(fp, "%8g %8d\n", time[i], clust->cl[i]);
1078 size_fp = xvgropen(sizefn, "Cluster Sizes", "Cluster #", "# Structures", oenv);
1079 if (output_env_get_print_xvgr_codes(oenv))
1081 fprintf(size_fp, "@g%d type %s\n", 0, "bar");
1084 snew(structure, nf);
1085 fprintf(log, "\n%3s | %3s %4s | %6s %4s | cluster members\n",
1086 "cl.", "#st", "rmsd", "middle", "rmsd");
1087 for (cl = 1; cl <= clust->ncl; cl++)
1089 /* prepare structures (fit, middle, average) */
1092 for (i = 0; i < natom; i++)
1098 for (i1 = 0; i1 < nf; i1++)
1100 if (clust->cl[i1] == cl)
1102 structure[nstr] = i1;
1104 if (trxfn && (bAverage || write_ncl) )
1108 reset_x(ifsize, fitidx, natom, nullptr, xx[i1], mass);
1116 do_fit(natom, mass, xx[first], xx[i1]);
1120 for (i = 0; i < natom; i++)
1122 rvec_inc(xav[i], xx[i1][i]);
1130 fprintf(size_fp, "%8d %8d\n", cl, nstr);
1135 for (i1 = 0; i1 < nstr; i1++)
1140 for (i = 0; i < nstr; i++)
1144 r += rmsd[structure[i]][structure[i1]];
1148 r += rmsd[structure[i1]][structure[i]];
1155 midstr = structure[i1];
1162 /* dump cluster info to logfile */
1165 sprintf(buf1, "%6.3f", clrmsd);
1170 sprintf(buf2, "%5.3f", midrmsd);
1178 sprintf(buf1, "%5s", "");
1179 sprintf(buf2, "%5s", "");
1181 fprintf(log, "%3d | %3d %s | %6g%s |", cl, nstr, buf1, time[midstr], buf2);
1182 for (i = 0; i < nstr; i++)
1184 if ((i % 7 == 0) && i)
1186 sprintf(buf, "\n%3s | %3s %4s | %6s %4s |", "", "", "", "", "");
1193 fprintf(log, "%s %6g", buf, time[i1]);
1197 /* write structures to trajectory file(s) */
1202 for (i = 0; i < nstr; i++)
1207 if (cl < write_ncl+1 && nstr > write_nst)
1209 /* Dump all structures for this cluster */
1210 /* generate numbered filename (there is a %d in trxfn!) */
1211 sprintf(buf, trxsfn, cl);
1212 trxsout = open_trx(buf, "w");
1213 for (i = 0; i < nstr; i++)
1218 for (i1 = 0; i1 < i && bWrite[i]; i1++)
1222 bWrite[i] = rmsd[structure[i1]][structure[i]] > rmsmin;
1228 write_trx(trxsout, iosize, outidx, atoms, i, time[structure[i]], zerobox,
1229 xx[structure[i]], nullptr, nullptr);
1234 /* Dump the average structure for this cluster */
1237 for (i = 0; i < natom; i++)
1239 svmul(1.0/nstr, xav[i], xav[i]);
1244 for (i = 0; i < natom; i++)
1246 copy_rvec(xx[midstr][i], xav[i]);
1250 reset_x(ifsize, fitidx, natom, nullptr, xav, mass);
1255 do_fit(natom, mass, xtps, xav);
1257 write_trx(trxout, iosize, outidx, atoms, cl, time[midstr], zerobox, xav, nullptr, nullptr);
1282 static void convert_mat(t_matrix *mat, t_mat *rms)
1287 matrix2real(mat, rms->mat);
1288 /* free input xpm matrix data */
1289 for (i = 0; i < mat->nx; i++)
1291 sfree(mat->matrix[i]);
1295 for (i = 0; i < mat->nx; i++)
1297 for (j = i; j < mat->nx; j++)
1299 rms->sumrms += rms->mat[i][j];
1300 rms->maxrms = std::max(rms->maxrms, rms->mat[i][j]);
1303 rms->minrms = std::min(rms->minrms, rms->mat[i][j]);
1310 int gmx_cluster(int argc, char *argv[])
1312 const char *desc[] = {
1313 "[THISMODULE] can cluster structures using several different methods.",
1314 "Distances between structures can be determined from a trajectory",
1315 "or read from an [REF].xpm[ref] matrix file with the [TT]-dm[tt] option.",
1316 "RMS deviation after fitting or RMS deviation of atom-pair distances",
1317 "can be used to define the distance between structures.[PAR]",
1319 "single linkage: add a structure to a cluster when its distance to any",
1320 "element of the cluster is less than [TT]cutoff[tt].[PAR]",
1322 "Jarvis Patrick: add a structure to a cluster when this structure",
1323 "and a structure in the cluster have each other as neighbors and",
1324 "they have a least [TT]P[tt] neighbors in common. The neighbors",
1325 "of a structure are the M closest structures or all structures within",
1326 "[TT]cutoff[tt].[PAR]",
1328 "Monte Carlo: reorder the RMSD matrix using Monte Carlo such that",
1329 "the order of the frames is using the smallest possible increments.",
1330 "With this it is possible to make a smooth animation going from one",
1331 "structure to another with the largest possible (e.g.) RMSD between",
1332 "them, however the intermediate steps should be as small as possible.",
1333 "Applications could be to visualize a potential of mean force",
1334 "ensemble of simulations or a pulling simulation. Obviously the user",
1335 "has to prepare the trajectory well (e.g. by not superimposing frames).",
1336 "The final result can be inspect visually by looking at the matrix",
1337 "[REF].xpm[ref] file, which should vary smoothly from bottom to top.[PAR]",
1339 "diagonalization: diagonalize the RMSD matrix.[PAR]",
1341 "gromos: use algorithm as described in Daura [IT]et al.[it]",
1342 "([IT]Angew. Chem. Int. Ed.[it] [BB]1999[bb], [IT]38[it], pp 236-240).",
1343 "Count number of neighbors using cut-off, take structure with",
1344 "largest number of neighbors with all its neighbors as cluster",
1345 "and eliminate it from the pool of clusters. Repeat for remaining",
1346 "structures in pool.[PAR]",
1348 "When the clustering algorithm assigns each structure to exactly one",
1349 "cluster (single linkage, Jarvis Patrick and gromos) and a trajectory",
1350 "file is supplied, the structure with",
1351 "the smallest average distance to the others or the average structure",
1352 "or all structures for each cluster will be written to a trajectory",
1353 "file. When writing all structures, separate numbered files are made",
1354 "for each cluster.[PAR]",
1356 "Two output files are always written:",
1358 " * [TT]-o[tt] writes the RMSD values in the upper left half of the matrix",
1359 " and a graphical depiction of the clusters in the lower right half",
1360 " When [TT]-minstruct[tt] = 1 the graphical depiction is black",
1361 " when two structures are in the same cluster.",
1362 " When [TT]-minstruct[tt] > 1 different colors will be used for each",
1364 " * [TT]-g[tt] writes information on the options used and a detailed list",
1365 " of all clusters and their members.",
1368 "Additionally, a number of optional output files can be written:",
1370 " * [TT]-dist[tt] writes the RMSD distribution.",
1371 " * [TT]-ev[tt] writes the eigenvectors of the RMSD matrix",
1372 " diagonalization.",
1373 " * [TT]-sz[tt] writes the cluster sizes.",
1374 " * [TT]-tr[tt] writes a matrix of the number transitions between",
1376 " * [TT]-ntr[tt] writes the total number of transitions to or from",
1378 " * [TT]-clid[tt] writes the cluster number as a function of time.",
1379 " * [TT]-cl[tt] writes average (with option [TT]-av[tt]) or central",
1380 " structure of each cluster or writes numbered files with cluster members",
1381 " for a selected set of clusters (with option [TT]-wcl[tt], depends on",
1382 " [TT]-nst[tt] and [TT]-rmsmin[tt]). The center of a cluster is the",
1383 " structure with the smallest average RMSD from all other structures",
1388 int nf = 0, i, i1, i2, j;
1392 rvec *xtps, *usextps, *x1, **xx = nullptr;
1393 const char *fn, *trx_out_fn;
1395 t_mat *rms, *orig = nullptr;
1400 t_matrix *readmat = nullptr;
1403 int isize = 0, ifsize = 0, iosize = 0;
1404 int *index = nullptr, *fitidx = nullptr, *outidx = nullptr;
1406 real rmsd, **d1, **d2, *time = nullptr, time_invfac, *mass = nullptr;
1407 char buf[STRLEN], buf1[80];
1408 gmx_bool bAnalyze, bUseRmsdCut, bJP_RMSD = FALSE, bReadMat, bReadTraj, bPBC = TRUE;
1410 int method, ncluster = 0;
1411 static const char *methodname[] = {
1412 nullptr, "linkage", "jarvis-patrick", "monte-carlo",
1413 "diagonalization", "gromos", nullptr
1416 m_null, m_linkage, m_jarvis_patrick,
1417 m_monte_carlo, m_diagonalize, m_gromos, m_nr
1419 /* Set colors for plotting: white = zero RMS, black = maximum */
1420 static t_rgb rlo_top = { 1.0, 1.0, 1.0 };
1421 static t_rgb rhi_top = { 0.0, 0.0, 0.0 };
1422 static t_rgb rlo_bot = { 1.0, 1.0, 1.0 };
1423 static t_rgb rhi_bot = { 0.0, 0.0, 1.0 };
1424 static int nlevels = 40, skip = 1;
1425 static real scalemax = -1.0, rmsdcut = 0.1, rmsmin = 0.0;
1426 gmx_bool bRMSdist = FALSE, bBinary = FALSE, bAverage = FALSE, bFit = TRUE;
1427 static int niter = 10000, nrandom = 0, seed = 0, write_ncl = 0, write_nst = 1, minstruct = 1;
1428 static real kT = 1e-3;
1429 static int M = 10, P = 3;
1430 gmx_output_env_t *oenv;
1431 gmx_rmpbc_t gpbc = nullptr;
1434 { "-dista", FALSE, etBOOL, {&bRMSdist},
1435 "Use RMSD of distances instead of RMS deviation" },
1436 { "-nlevels", FALSE, etINT, {&nlevels},
1437 "Discretize RMSD matrix in this number of levels" },
1438 { "-cutoff", FALSE, etREAL, {&rmsdcut},
1439 "RMSD cut-off (nm) for two structures to be neighbor" },
1440 { "-fit", FALSE, etBOOL, {&bFit},
1441 "Use least squares fitting before RMSD calculation" },
1442 { "-max", FALSE, etREAL, {&scalemax},
1443 "Maximum level in RMSD matrix" },
1444 { "-skip", FALSE, etINT, {&skip},
1445 "Only analyze every nr-th frame" },
1446 { "-av", FALSE, etBOOL, {&bAverage},
1447 "Write average instead of middle structure for each cluster" },
1448 { "-wcl", FALSE, etINT, {&write_ncl},
1449 "Write the structures for this number of clusters to numbered files" },
1450 { "-nst", FALSE, etINT, {&write_nst},
1451 "Only write all structures if more than this number of structures per cluster" },
1452 { "-rmsmin", FALSE, etREAL, {&rmsmin},
1453 "minimum rms difference with rest of cluster for writing structures" },
1454 { "-method", FALSE, etENUM, {methodname},
1455 "Method for cluster determination" },
1456 { "-minstruct", FALSE, etINT, {&minstruct},
1457 "Minimum number of structures in cluster for coloring in the [REF].xpm[ref] file" },
1458 { "-binary", FALSE, etBOOL, {&bBinary},
1459 "Treat the RMSD matrix as consisting of 0 and 1, where the cut-off "
1460 "is given by [TT]-cutoff[tt]" },
1461 { "-M", FALSE, etINT, {&M},
1462 "Number of nearest neighbors considered for Jarvis-Patrick algorithm, "
1463 "0 is use cutoff" },
1464 { "-P", FALSE, etINT, {&P},
1465 "Number of identical nearest neighbors required to form a cluster" },
1466 { "-seed", FALSE, etINT, {&seed},
1467 "Random number seed for Monte Carlo clustering algorithm (0 means generate)" },
1468 { "-niter", FALSE, etINT, {&niter},
1469 "Number of iterations for MC" },
1470 { "-nrandom", FALSE, etINT, {&nrandom},
1471 "The first iterations for MC may be done complete random, to shuffle the frames" },
1472 { "-kT", FALSE, etREAL, {&kT},
1473 "Boltzmann weighting factor for Monte Carlo optimization "
1474 "(zero turns off uphill steps)" },
1475 { "-pbc", FALSE, etBOOL,
1476 { &bPBC }, "PBC check" }
1479 { efTRX, "-f", nullptr, ffOPTRD },
1480 { efTPS, "-s", nullptr, ffREAD },
1481 { efNDX, nullptr, nullptr, ffOPTRD },
1482 { efXPM, "-dm", "rmsd", ffOPTRD },
1483 { efXPM, "-om", "rmsd-raw", ffWRITE },
1484 { efXPM, "-o", "rmsd-clust", ffWRITE },
1485 { efLOG, "-g", "cluster", ffWRITE },
1486 { efXVG, "-dist", "rmsd-dist", ffOPTWR },
1487 { efXVG, "-ev", "rmsd-eig", ffOPTWR },
1488 { efXVG, "-conv", "mc-conv", ffOPTWR },
1489 { efXVG, "-sz", "clust-size", ffOPTWR},
1490 { efXPM, "-tr", "clust-trans", ffOPTWR},
1491 { efXVG, "-ntr", "clust-trans", ffOPTWR},
1492 { efXVG, "-clid", "clust-id", ffOPTWR},
1493 { efTRX, "-cl", "clusters.pdb", ffOPTWR }
1495 #define NFILE asize(fnm)
1497 if (!parse_common_args(&argc, argv,
1498 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
1499 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr,
1506 bReadMat = opt2bSet("-dm", NFILE, fnm);
1507 bReadTraj = opt2bSet("-f", NFILE, fnm) || !bReadMat;
1508 if (opt2parg_bSet("-av", asize(pa), pa) ||
1509 opt2parg_bSet("-wcl", asize(pa), pa) ||
1510 opt2parg_bSet("-nst", asize(pa), pa) ||
1511 opt2parg_bSet("-rmsmin", asize(pa), pa) ||
1512 opt2bSet("-cl", NFILE, fnm) )
1514 trx_out_fn = opt2fn("-cl", NFILE, fnm);
1518 trx_out_fn = nullptr;
1520 if (bReadMat && output_env_get_time_factor(oenv) != 1)
1523 "\nWarning: assuming the time unit in %s is %s\n",
1524 opt2fn("-dm", NFILE, fnm), output_env_get_time_unit(oenv).c_str());
1526 if (trx_out_fn && !bReadTraj)
1528 fprintf(stderr, "\nWarning: "
1529 "cannot write cluster structures without reading trajectory\n"
1530 " ignoring option -cl %s\n", trx_out_fn);
1534 while (method < m_nr && gmx_strcasecmp(methodname[0], methodname[method]) != 0)
1540 gmx_fatal(FARGS, "Invalid method");
1543 bAnalyze = (method == m_linkage || method == m_jarvis_patrick ||
1544 method == m_gromos );
1547 log = ftp2FILE(efLOG, NFILE, fnm, "w");
1549 fprintf(stderr, "Using %s method for clustering\n", methodname[0]);
1550 fprintf(log, "Using %s method for clustering\n", methodname[0]);
1552 /* check input and write parameters to log file */
1553 bUseRmsdCut = FALSE;
1554 if (method == m_jarvis_patrick)
1556 bJP_RMSD = (M == 0) || opt2parg_bSet("-cutoff", asize(pa), pa);
1557 if ((M < 0) || (M == 1))
1559 gmx_fatal(FARGS, "M (%d) must be 0 or larger than 1", M);
1563 sprintf(buf1, "Will use P=%d and RMSD cutoff (%g)", P, rmsdcut);
1570 gmx_fatal(FARGS, "Number of neighbors required (P) must be less than M");
1574 sprintf(buf1, "Will use P=%d, M=%d and RMSD cutoff (%g)", P, M, rmsdcut);
1579 sprintf(buf1, "Will use P=%d, M=%d", P, M);
1582 ffprintf_s(stderr, log, buf, "%s for determining the neighbors\n\n", buf1);
1584 else /* method != m_jarvis */
1586 bUseRmsdCut = ( bBinary || method == m_linkage || method == m_gromos );
1588 if (bUseRmsdCut && method != m_jarvis_patrick)
1590 fprintf(log, "Using RMSD cutoff %g nm\n", rmsdcut);
1592 if (method == m_monte_carlo)
1594 fprintf(log, "Using %d iterations\n", niter);
1599 gmx_fatal(FARGS, "skip (%d) should be >= 1", skip);
1605 /* don't read mass-database as masses (and top) are not used */
1606 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtps, nullptr, box,
1610 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
1613 fprintf(stderr, "\nSelect group for least squares fit%s:\n",
1614 bReadMat ? "" : " and RMSD calculation");
1615 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm),
1616 1, &ifsize, &fitidx, &grpname);
1619 fprintf(stderr, "\nSelect group for output:\n");
1620 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm),
1621 1, &iosize, &outidx, &grpname);
1622 /* merge and convert both index groups: */
1623 /* first copy outidx to index. let outidx refer to elements in index */
1624 snew(index, iosize);
1626 for (i = 0; i < iosize; i++)
1628 index[i] = outidx[i];
1631 /* now lookup elements from fitidx in index, add them if necessary
1632 and also let fitidx refer to elements in index */
1633 for (i = 0; i < ifsize; i++)
1636 while (j < isize && index[j] != fitidx[i])
1642 /* slow this way, but doesn't matter much */
1644 srenew(index, isize);
1646 index[j] = fitidx[i];
1650 else /* !trx_out_fn */
1654 for (i = 0; i < ifsize; i++)
1656 index[i] = fitidx[i];
1664 /* Loop over first coordinate file */
1665 fn = opt2fn("-f", NFILE, fnm);
1667 xx = read_whole_trj(fn, isize, index, skip, &nf, &time, oenv, bPBC, gpbc);
1668 output_env_conv_times(oenv, nf, time);
1669 if (!bRMSdist || bAnalyze)
1671 /* Center all frames on zero */
1673 for (i = 0; i < ifsize; i++)
1675 mass[fitidx[i]] = top.atoms.atom[index[fitidx[i]]].m;
1679 for (i = 0; i < nf; i++)
1681 reset_x(ifsize, fitidx, isize, nullptr, xx[i], mass);
1687 gmx_rmpbc_done(gpbc);
1693 fprintf(stderr, "Reading rms distance matrix ");
1694 read_xpm_matrix(opt2fn("-dm", NFILE, fnm), &readmat);
1695 fprintf(stderr, "\n");
1696 if (readmat[0].nx != readmat[0].ny)
1698 gmx_fatal(FARGS, "Matrix (%dx%d) is not square",
1699 readmat[0].nx, readmat[0].ny);
1701 if (bReadTraj && bAnalyze && (readmat[0].nx != nf))
1703 gmx_fatal(FARGS, "Matrix size (%dx%d) does not match the number of "
1704 "frames (%d)", readmat[0].nx, readmat[0].ny, nf);
1709 time = readmat[0].axis_x;
1710 time_invfac = output_env_get_time_invfactor(oenv);
1711 for (i = 0; i < nf; i++)
1713 time[i] *= time_invfac;
1716 rms = init_mat(readmat[0].nx, method == m_diagonalize);
1717 convert_mat(&(readmat[0]), rms);
1719 nlevels = readmat[0].nmap;
1721 else /* !bReadMat */
1723 rms = init_mat(nf, method == m_diagonalize);
1724 nrms = (static_cast<int64_t>(nf)*static_cast<int64_t>(nf-1))/2;
1727 fprintf(stderr, "Computing %dx%d RMS deviation matrix\n", nf, nf);
1728 /* Initialize work array */
1730 for (i1 = 0; i1 < nf; i1++)
1732 for (i2 = i1+1; i2 < nf; i2++)
1734 for (i = 0; i < isize; i++)
1736 copy_rvec(xx[i1][i], x1[i]);
1740 do_fit(isize, mass, xx[i2], x1);
1742 rmsd = rmsdev(isize, mass, xx[i2], x1);
1743 set_mat_entry(rms, i1, i2, rmsd);
1746 fprintf(stderr, "\r# RMSD calculations left: " "%" PRId64 " ", nrms);
1753 fprintf(stderr, "Computing %dx%d RMS distance deviation matrix\n", nf, nf);
1755 /* Initiate work arrays */
1758 for (i = 0; (i < isize); i++)
1763 for (i1 = 0; i1 < nf; i1++)
1765 calc_dist(isize, xx[i1], d1);
1766 for (i2 = i1+1; (i2 < nf); i2++)
1768 calc_dist(isize, xx[i2], d2);
1769 set_mat_entry(rms, i1, i2, rms_dist(isize, d1, d2));
1772 fprintf(stderr, "\r# RMSD calculations left: " "%" PRId64 " ", nrms);
1775 /* Clean up work arrays */
1776 for (i = 0; (i < isize); i++)
1784 fprintf(stderr, "\n\n");
1786 ffprintf_gg(stderr, log, buf, "The RMSD ranges from %g to %g nm\n",
1787 rms->minrms, rms->maxrms);
1788 ffprintf_g(stderr, log, buf, "Average RMSD is %g\n", 2*rms->sumrms/(nf*(nf-1)));
1789 ffprintf_d(stderr, log, buf, "Number of structures for matrix %d\n", nf);
1790 ffprintf_g(stderr, log, buf, "Energy of the matrix is %g.\n", mat_energy(rms));
1791 if (bUseRmsdCut && (rmsdcut < rms->minrms || rmsdcut > rms->maxrms) )
1793 fprintf(stderr, "WARNING: rmsd cutoff %g is outside range of rmsd values "
1794 "%g to %g\n", rmsdcut, rms->minrms, rms->maxrms);
1796 if (bAnalyze && (rmsmin < rms->minrms) )
1798 fprintf(stderr, "WARNING: rmsd minimum %g is below lowest rmsd value %g\n",
1799 rmsmin, rms->minrms);
1801 if (bAnalyze && (rmsmin > rmsdcut) )
1803 fprintf(stderr, "WARNING: rmsd minimum %g is above rmsd cutoff %g\n",
1807 /* Plot the rmsd distribution */
1808 rmsd_distribution(opt2fn("-dist", NFILE, fnm), rms, oenv);
1812 for (i1 = 0; (i1 < nf); i1++)
1814 for (i2 = 0; (i2 < nf); i2++)
1816 if (rms->mat[i1][i2] < rmsdcut)
1818 rms->mat[i1][i2] = 0;
1822 rms->mat[i1][i2] = 1;
1832 /* Now sort the matrix and write it out again */
1833 gather(rms, rmsdcut, &clust);
1836 /* Do a diagonalization */
1837 snew(eigenvalues, nf);
1838 snew(eigenvectors, nf*nf);
1839 std::memcpy(eigenvectors, rms->mat[0], nf*nf*sizeof(real));
1840 eigensolver(eigenvectors, nf, 0, nf, eigenvalues, rms->mat[0]);
1841 sfree(eigenvectors);
1843 fp = xvgropen(opt2fn("-ev", NFILE, fnm), "RMSD matrix Eigenvalues",
1844 "Eigenvector index", "Eigenvalues (nm\\S2\\N)", oenv);
1845 for (i = 0; (i < nf); i++)
1847 fprintf(fp, "%10d %10g\n", i, eigenvalues[i]);
1852 orig = init_mat(rms->nn, FALSE);
1854 copy_t_mat(orig, rms);
1855 mc_optimize(log, rms, time, niter, nrandom, seed, kT,
1856 opt2fn_null("-conv", NFILE, fnm), oenv);
1858 case m_jarvis_patrick:
1859 jarvis_patrick(rms->nn, rms->mat, M, P, bJP_RMSD ? rmsdcut : -1, &clust);
1862 gromos(rms->nn, rms->mat, rmsdcut, &clust);
1865 gmx_fatal(FARGS, "DEATH HORROR unknown method \"%s\"", methodname[0]);
1868 if (method == m_monte_carlo || method == m_diagonalize)
1870 fprintf(stderr, "Energy of the matrix after clustering is %g.\n",
1878 ncluster = plot_clusters(nf, rms->mat, &clust, minstruct);
1882 mark_clusters(nf, rms->mat, rms->maxrms, &clust);
1884 init_t_atoms(&useatoms, isize, FALSE);
1885 snew(usextps, isize);
1886 useatoms.resinfo = top.atoms.resinfo;
1887 for (i = 0; i < isize; i++)
1889 useatoms.atomname[i] = top.atoms.atomname[index[i]];
1890 useatoms.atom[i].resind = top.atoms.atom[index[i]].resind;
1891 useatoms.nres = std::max(useatoms.nres, useatoms.atom[i].resind+1);
1892 copy_rvec(xtps[index[i]], usextps[i]);
1894 useatoms.nr = isize;
1895 analyze_clusters(nf, &clust, rms->mat, isize, &useatoms, usextps, mass, xx, time,
1896 ifsize, fitidx, iosize, outidx,
1897 bReadTraj ? trx_out_fn : nullptr,
1898 opt2fn_null("-sz", NFILE, fnm),
1899 opt2fn_null("-tr", NFILE, fnm),
1900 opt2fn_null("-ntr", NFILE, fnm),
1901 opt2fn_null("-clid", NFILE, fnm),
1902 bAverage, write_ncl, write_nst, rmsmin, bFit, log,
1903 rlo_bot, rhi_bot, oenv);
1907 if (bBinary && !bAnalyze)
1909 /* Make the clustering visible */
1910 for (i2 = 0; (i2 < nf); i2++)
1912 for (i1 = i2+1; (i1 < nf); i1++)
1914 if (rms->mat[i1][i2])
1916 rms->mat[i1][i2] = rms->maxrms;
1922 fp = opt2FILE("-o", NFILE, fnm, "w");
1923 fprintf(stderr, "Writing rms distance/clustering matrix ");
1926 write_xpm(fp, 0, readmat[0].title, readmat[0].legend, readmat[0].label_x,
1927 readmat[0].label_y, nf, nf, readmat[0].axis_x, readmat[0].axis_y,
1928 rms->mat, 0.0, rms->maxrms, rlo_top, rhi_top, &nlevels);
1932 auto timeLabel = output_env_get_time_label(oenv);
1933 auto title = gmx::formatString("RMS%sDeviation / Cluster Index",
1934 bRMSdist ? " Distance " : " ");
1937 write_xpm_split(fp, 0, title, "RMSD (nm)", timeLabel, timeLabel,
1938 nf, nf, time, time, rms->mat, 0.0, rms->maxrms, &nlevels,
1939 rlo_top, rhi_top, 0.0, ncluster,
1940 &ncluster, TRUE, rlo_bot, rhi_bot);
1944 write_xpm(fp, 0, title, "RMSD (nm)", timeLabel, timeLabel,
1945 nf, nf, time, time, rms->mat, 0.0, rms->maxrms,
1946 rlo_top, rhi_top, &nlevels);
1949 fprintf(stderr, "\n");
1951 if (nullptr != orig)
1953 fp = opt2FILE("-om", NFILE, fnm, "w");
1954 auto timeLabel = output_env_get_time_label(oenv);
1955 auto title = gmx::formatString("RMS%sDeviation", bRMSdist ? " Distance " : " ");
1956 write_xpm(fp, 0, title, "RMSD (nm)", timeLabel, timeLabel,
1957 nf, nf, time, time, orig->mat, 0.0, orig->maxrms,
1958 rlo_top, rhi_top, &nlevels);
1963 /* now show what we've done */
1964 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
1965 do_view(oenv, opt2fn_null("-sz", NFILE, fnm), "-nxy");
1966 if (method == m_diagonalize)
1968 do_view(oenv, opt2fn_null("-ev", NFILE, fnm), "-nxy");
1970 do_view(oenv, opt2fn("-dist", NFILE, fnm), "-nxy");
1973 do_view(oenv, opt2fn_null("-tr", NFILE, fnm), "-nxy");
1974 do_view(oenv, opt2fn_null("-ntr", NFILE, fnm), "-nxy");
1975 do_view(oenv, opt2fn_null("-clid", NFILE, fnm), "-nxy");
1977 do_view(oenv, opt2fn_null("-conv", NFILE, fnm), nullptr);