48eb4e8ffdd236658bcdf0b34af7d000b7609bec
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rmsdist.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstring>
41
42 #include <algorithm>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/matio.h"
48 #include "gromacs/fileio/trxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/strdb.h"
63
64
65 static void calc_dist(int nind, const int index[], const rvec x[], int ePBC, matrix box, real** d)
66 {
67     int   i, j;
68     rvec  dx;
69     real  temp2;
70     t_pbc pbc;
71
72     set_pbc(&pbc, ePBC, box);
73     for (i = 0; (i < nind - 1); i++)
74     {
75         const real* xi = x[index[i]];
76         for (j = i + 1; (j < nind); j++)
77         {
78             pbc_dx(&pbc, xi, x[index[j]], dx);
79             temp2   = norm2(dx);
80             d[i][j] = std::sqrt(temp2);
81         }
82     }
83 }
84
85 static void calc_dist_tot(int       nind,
86                           const int index[],
87                           rvec      x[],
88                           int       ePBC,
89                           matrix    box,
90                           real**    d,
91                           real**    dtot,
92                           real**    dtot2,
93                           gmx_bool  bNMR,
94                           real**    dtot1_3,
95                           real**    dtot1_6)
96 {
97     int   i, j;
98     real* xi;
99     real  temp, temp2, temp1_3;
100     rvec  dx;
101     t_pbc pbc;
102
103     set_pbc(&pbc, ePBC, box);
104     for (i = 0; (i < nind - 1); i++)
105     {
106         xi = x[index[i]];
107         for (j = i + 1; (j < nind); j++)
108         {
109             pbc_dx(&pbc, xi, x[index[j]], dx);
110             temp2   = norm2(dx);
111             temp    = std::sqrt(temp2);
112             d[i][j] = temp;
113             dtot[i][j] += temp;
114             dtot2[i][j] += temp2;
115             if (bNMR)
116             {
117                 temp1_3 = 1.0 / (temp * temp2);
118                 dtot1_3[i][j] += temp1_3;
119                 dtot1_6[i][j] += temp1_3 * temp1_3;
120             }
121         }
122     }
123 }
124
125 static void calc_nmr(int nind, int nframes, real** dtot1_3, real** dtot1_6, real* max1_3, real* max1_6)
126 {
127     int  i, j;
128     real temp1_3, temp1_6;
129
130     for (i = 0; (i < nind - 1); i++)
131     {
132         for (j = i + 1; (j < nind); j++)
133         {
134             temp1_3 = gmx::invcbrt(dtot1_3[i][j] / static_cast<real>(nframes));
135             temp1_6 = gmx::invsixthroot(dtot1_6[i][j] / static_cast<real>(nframes));
136             if (temp1_3 > *max1_3)
137             {
138                 *max1_3 = temp1_3;
139             }
140             if (temp1_6 > *max1_6)
141             {
142                 *max1_6 = temp1_6;
143             }
144             dtot1_3[i][j] = temp1_3;
145             dtot1_6[i][j] = temp1_6;
146             dtot1_3[j][i] = temp1_3;
147             dtot1_6[j][i] = temp1_6;
148         }
149     }
150 }
151
152 static char Hnum[] = "123";
153
154 typedef struct
155 {
156     int  nr;
157     real r_3;
158     real r_6;
159     real i_3;
160     real i_6;
161 } t_noe;
162
163 typedef struct
164 {
165     int   anr;
166     int   ianr;
167     int   rnr;
168     char* aname;
169     char* rname;
170 } t_noe_gr;
171
172 typedef struct
173 {
174     bool  set;
175     int   rnr;
176     char* nname;
177     char* rname;
178     char* aname;
179 } t_equiv;
180
181 static int read_equiv(const char* eq_fn, t_equiv*** equivptr)
182 {
183     FILE*     fp;
184     char      line[STRLEN], resname[10], atomname[10], *lp;
185     int       neq, na, n, resnr;
186     t_equiv** equiv;
187
188     fp    = gmx_ffopen(eq_fn, "r");
189     neq   = 0;
190     equiv = nullptr;
191     while (get_a_line(fp, line, STRLEN))
192     {
193         lp = line;
194         /* this is not efficient, but I'm lazy */
195         srenew(equiv, neq + 1);
196         equiv[neq] = nullptr;
197         na         = 0;
198         if (sscanf(lp, "%s %n", atomname, &n) == 1)
199         {
200             lp += n;
201             snew(equiv[neq], 1);
202             equiv[neq][0].nname = gmx_strdup(atomname);
203             while (sscanf(lp, "%d %s %s %n", &resnr, resname, atomname, &n) == 3)
204             {
205                 /* this is not efficient, but I'm lazy (again) */
206                 srenew(equiv[neq], na + 1);
207                 equiv[neq][na].set   = true;
208                 equiv[neq][na].rnr   = resnr - 1;
209                 equiv[neq][na].rname = gmx_strdup(resname);
210                 equiv[neq][na].aname = gmx_strdup(atomname);
211                 if (na > 0)
212                 {
213                     equiv[neq][na].nname = nullptr;
214                 }
215                 na++;
216                 lp += n;
217             }
218         }
219         /* make empty element as flag for end of array */
220         srenew(equiv[neq], na + 1);
221         equiv[neq][na].set   = false;
222         equiv[neq][na].rnr   = 0;
223         equiv[neq][na].rname = nullptr;
224         equiv[neq][na].aname = nullptr;
225
226         /* next */
227         neq++;
228     }
229     gmx_ffclose(fp);
230
231     *equivptr = equiv;
232
233     return neq;
234 }
235
236 static void dump_equiv(FILE* out, int neq, t_equiv** equiv)
237 {
238     int i, j;
239
240     fprintf(out, "Dumping equivalent list\n");
241     for (i = 0; i < neq; i++)
242     {
243         fprintf(out, "%s", equiv[i][0].nname);
244         for (j = 0; equiv[i][j].set; j++)
245         {
246             fprintf(out, " %d %s %s", equiv[i][j].rnr, equiv[i][j].rname, equiv[i][j].aname);
247         }
248         fprintf(out, "\n");
249     }
250 }
251
252 static gmx_bool
253 is_equiv(int neq, t_equiv** equiv, char** nname, int rnr1, char* rname1, char* aname1, int rnr2, char* rname2, char* aname2)
254 {
255     int      i, j;
256     gmx_bool bFound;
257
258     bFound = FALSE;
259     /* we can terminate each loop when bFound is true! */
260     for (i = 0; i < neq && !bFound; i++)
261     {
262         /* find first atom */
263         for (j = 0; equiv[i][j].set && !bFound; j++)
264         {
265             bFound = (equiv[i][j].rnr == rnr1 && std::strcmp(equiv[i][j].rname, rname1) == 0
266                       && std::strcmp(equiv[i][j].aname, aname1) == 0);
267         }
268         if (bFound)
269         {
270             /* find second atom */
271             bFound = FALSE;
272             for (j = 0; equiv[i][j].set && !bFound; j++)
273             {
274                 bFound = (equiv[i][j].rnr == rnr2 && std::strcmp(equiv[i][j].rname, rname2) == 0
275                           && std::strcmp(equiv[i][j].aname, aname2) == 0);
276             }
277         }
278     }
279     if (bFound)
280     {
281         *nname = gmx_strdup(equiv[i - 1][0].nname);
282     }
283
284     return bFound;
285 }
286
287 static int analyze_noe_equivalent(const char*    eq_fn,
288                                   const t_atoms* atoms,
289                                   int            isize,
290                                   const int*     index,
291                                   gmx_bool       bSumH,
292                                   int*           noe_index,
293                                   t_noe_gr*      noe_gr)
294 {
295     int       i, j, anmil, anmjl, rnri, rnrj, gi, groupnr, neq;
296     char *    anmi, *anmj, **nnm;
297     gmx_bool  bMatch, bEquiv;
298     t_equiv** equiv;
299
300     snew(nnm, isize);
301     if (bSumH)
302     {
303         if (eq_fn)
304         {
305             neq = read_equiv(eq_fn, &equiv);
306             if (debug)
307             {
308                 dump_equiv(debug, neq, equiv);
309             }
310         }
311         else
312         {
313             neq   = 0;
314             equiv = nullptr;
315         }
316
317         groupnr = 0;
318         for (i = 0; i < isize; i++)
319         {
320             if (equiv && i < isize - 1)
321             {
322                 /* check explicit list of equivalent atoms */
323                 do
324                 {
325                     j      = i + 1;
326                     rnri   = atoms->atom[index[i]].resind;
327                     rnrj   = atoms->atom[index[j]].resind;
328                     bEquiv = is_equiv(neq, equiv, &nnm[i], rnri, *atoms->resinfo[rnri].name,
329                                       *atoms->atomname[index[i]], rnrj, *atoms->resinfo[rnrj].name,
330                                       *atoms->atomname[index[j]]);
331                     if (nnm[i] && bEquiv)
332                     {
333                         nnm[j] = gmx_strdup(nnm[i]);
334                     }
335                     if (bEquiv)
336                     {
337                         /* set index for matching atom */
338                         noe_index[i] = groupnr;
339                         /* skip matching atom */
340                         i = j;
341                     }
342                 } while (bEquiv && i < isize - 1);
343             }
344             else
345             {
346                 bEquiv = FALSE;
347             }
348             if (!bEquiv)
349             {
350                 /* look for triplets of consecutive atoms with name XX?,
351                    X are any number of letters or digits and ? goes from 1 to 3
352                    This is supposed to cover all CH3 groups and the like */
353                 anmi   = *atoms->atomname[index[i]];
354                 anmil  = std::strlen(anmi);
355                 bMatch = i <= isize - 3 && anmi[anmil - 1] == '1';
356                 if (bMatch)
357                 {
358                     for (j = 1; j < 3; j++)
359                     {
360                         anmj   = *atoms->atomname[index[i + j]];
361                         anmjl  = std::strlen(anmj);
362                         bMatch = bMatch
363                                  && (anmil == anmjl && anmj[anmjl - 1] == Hnum[j]
364                                      && std::strncmp(anmi, anmj, anmil - 1) == 0);
365                     }
366                 }
367                 /* set index for this atom */
368                 noe_index[i] = groupnr;
369                 if (bMatch)
370                 {
371                     /* set index for next two matching atoms */
372                     for (j = 1; j < 3; j++)
373                     {
374                         noe_index[i + j] = groupnr;
375                     }
376                     /* skip matching atoms */
377                     i += 2;
378                 }
379             }
380             groupnr++;
381         }
382     }
383     else
384     {
385         /* make index without looking for equivalent atoms */
386         for (i = 0; i < isize; i++)
387         {
388             noe_index[i] = i;
389         }
390         groupnr = isize;
391     }
392     noe_index[isize] = groupnr;
393
394     if (debug)
395     {
396         /* dump new names */
397         for (i = 0; i < isize; i++)
398         {
399             rnri = atoms->atom[index[i]].resind;
400             fprintf(debug, "%s %s %d -> %s\n", *atoms->atomname[index[i]],
401                     *atoms->resinfo[rnri].name, rnri, nnm[i] ? nnm[i] : "");
402         }
403     }
404
405     for (i = 0; i < isize; i++)
406     {
407         gi = noe_index[i];
408         if (!noe_gr[gi].aname)
409         {
410             noe_gr[gi].ianr = i;
411             noe_gr[gi].anr  = index[i];
412             if (nnm[i])
413             {
414                 noe_gr[gi].aname = gmx_strdup(nnm[i]);
415             }
416             else
417             {
418                 noe_gr[gi].aname = gmx_strdup(*atoms->atomname[index[i]]);
419                 if (noe_index[i] == noe_index[i + 1])
420                 {
421                     noe_gr[gi].aname[std::strlen(noe_gr[gi].aname) - 1] = '*';
422                 }
423             }
424             noe_gr[gi].rnr   = atoms->atom[index[i]].resind;
425             noe_gr[gi].rname = gmx_strdup(*atoms->resinfo[noe_gr[gi].rnr].name);
426             /* dump group definitions */
427             if (debug)
428             {
429                 fprintf(debug, "%d %d %d %d %s %s %d\n", i, gi, noe_gr[gi].ianr, noe_gr[gi].anr,
430                         noe_gr[gi].aname, noe_gr[gi].rname, noe_gr[gi].rnr);
431             }
432         }
433     }
434     for (i = 0; i < isize; i++)
435     {
436         sfree(nnm[i]);
437     }
438     sfree(nnm);
439
440     return groupnr;
441 }
442
443 /* #define NSCALE 3 */
444 /*  static char *noe_scale[NSCALE+1] = { */
445 /*    "strong", "medium", "weak", "none" */
446 /*  }; */
447 #define NSCALE 6
448
449 static char* noe2scale(real r3, real r6, real rmax)
450 {
451     int         i, s3, s6;
452     static char buf[NSCALE + 1];
453
454     /* r goes from 0 to rmax
455        NSCALE*r/rmax goes from 0 to NSCALE
456        NSCALE - NSCALE*r/rmax goes from NSCALE to 0 */
457     s3 = NSCALE - std::min(NSCALE, static_cast<int>(NSCALE * r3 / rmax));
458     s6 = NSCALE - std::min(NSCALE, static_cast<int>(NSCALE * r6 / rmax));
459
460     for (i = 0; i < s3; i++)
461     {
462         buf[i] = '=';
463     }
464     for (; i < s6; i++)
465     {
466         buf[i] = '-';
467     }
468     buf[i] = '\0';
469
470     return buf;
471 }
472
473 static void calc_noe(int isize, const int* noe_index, real** dtot1_3, real** dtot1_6, int gnr, t_noe** noe)
474 {
475     int i, j, gi, gj;
476
477     /* make half matrix of noe-group distances from atom distances */
478     for (i = 0; i < isize; i++)
479     {
480         gi = noe_index[i];
481         for (j = i; j < isize; j++)
482         {
483             gj = noe_index[j];
484             noe[gi][gj].nr++;
485             noe[gi][gj].i_3 += 1.0 / gmx::power3(dtot1_3[i][j]);
486             noe[gi][gj].i_6 += 1.0 / gmx::power6(dtot1_6[i][j]);
487         }
488     }
489
490     /* make averages */
491     for (i = 0; i < gnr; i++)
492     {
493         for (j = i + 1; j < gnr; j++)
494         {
495             noe[i][j].r_3 = gmx::invcbrt(noe[i][j].i_3 / static_cast<real>(noe[i][j].nr));
496             noe[i][j].r_6 = gmx::invsixthroot(noe[i][j].i_6 / static_cast<real>(noe[i][j].nr));
497             noe[j][i]     = noe[i][j];
498         }
499     }
500 }
501
502 static void write_noe(FILE* fp, int gnr, t_noe** noe, t_noe_gr* noe_gr, real rmax)
503 {
504     int      i, j;
505     real     r3, r6, min3, min6;
506     char     buf[10], b3[10], b6[10];
507     t_noe_gr gri, grj;
508
509     min3 = min6 = 1e6;
510     fprintf(fp, ";%4s %3s %4s %4s%3s %4s %4s %4s %4s%3s %5s %5s %8s %2s %2s %s\n", "ianr", "anr",
511             "anm", "rnm", "rnr", "ianr", "anr", "anm", "rnm", "rnr", "1/r^3", "1/r^6", "intnsty",
512             "Dr", "Da", "scale");
513     for (i = 0; i < gnr; i++)
514     {
515         gri = noe_gr[i];
516         for (j = i + 1; j < gnr; j++)
517         {
518             grj  = noe_gr[j];
519             r3   = noe[i][j].r_3;
520             r6   = noe[i][j].r_6;
521             min3 = std::min(r3, min3);
522             min6 = std::min(r6, min6);
523             if (r3 < rmax || r6 < rmax)
524             {
525                 if (grj.rnr == gri.rnr)
526                 {
527                     sprintf(buf, "%2d", grj.anr - gri.anr);
528                 }
529                 else
530                 {
531                     buf[0] = '\0';
532                 }
533                 if (r3 < rmax)
534                 {
535                     sprintf(b3, "%-5.3f", r3);
536                 }
537                 else
538                 {
539                     std::strcpy(b3, "-");
540                 }
541                 if (r6 < rmax)
542                 {
543                     sprintf(b6, "%-5.3f", r6);
544                 }
545                 else
546                 {
547                     std::strcpy(b6, "-");
548                 }
549                 fprintf(fp, "%4d %4d %4s %4s%3d %4d %4d %4s %4s%3d %5s %5s %8d %2d %2s %s\n",
550                         gri.ianr + 1, gri.anr + 1, gri.aname, gri.rname, gri.rnr + 1, grj.ianr + 1,
551                         grj.anr + 1, grj.aname, grj.rname, grj.rnr + 1, b3, b6,
552                         gmx::roundToInt(noe[i][j].i_6), grj.rnr - gri.rnr, buf, noe2scale(r3, r6, rmax));
553             }
554         }
555     }
556 #define MINI ((i == 3) ? min3 : min6)
557     for (i = 3; i <= 6; i += 3)
558     {
559         if (MINI > rmax)
560         {
561             fprintf(stdout,
562                     "NOTE: no 1/r^%d averaged distances found below %g, "
563                     "smallest was %g\n",
564                     i, rmax, MINI);
565         }
566         else
567         {
568             fprintf(stdout, "Smallest 1/r^%d averaged distance was %g\n", i, MINI);
569         }
570     }
571 #undef MINI
572 }
573
574 static void calc_rms(int    nind,
575                      int    nframes,
576                      real** dtot,
577                      real** dtot2,
578                      real** rmsmat,
579                      real*  rmsmax,
580                      real** rmscmat,
581                      real*  rmscmax,
582                      real** meanmat,
583                      real*  meanmax)
584 {
585     int  i, j;
586     real mean, mean2, rms, rmsc;
587     /* N.B. dtot and dtot2 contain the total distance and the total squared
588      * distance respectively, BUT they return RMS and the scaled RMS resp.
589      */
590
591     *rmsmax  = -1000;
592     *rmscmax = -1000;
593     *meanmax = -1000;
594
595     for (i = 0; (i < nind - 1); i++)
596     {
597         for (j = i + 1; (j < nind); j++)
598         {
599             mean  = dtot[i][j] / static_cast<real>(nframes);
600             mean2 = dtot2[i][j] / static_cast<real>(nframes);
601             rms   = std::sqrt(std::max(0.0_real, mean2 - mean * mean));
602             rmsc  = rms / mean;
603             if (mean > *meanmax)
604             {
605                 *meanmax = mean;
606             }
607             if (rms > *rmsmax)
608             {
609                 *rmsmax = rms;
610             }
611             if (rmsc > *rmscmax)
612             {
613                 *rmscmax = rmsc;
614             }
615             meanmat[i][j] = meanmat[j][i] = mean;
616             rmsmat[i][j] = rmsmat[j][i] = rms;
617             rmscmat[i][j] = rmscmat[j][i] = rmsc;
618         }
619     }
620 }
621
622 static real rms_diff(int natom, real** d, real** d_r)
623 {
624     int  i, j;
625     real r, r2;
626
627     r2 = 0.0;
628     for (i = 0; (i < natom - 1); i++)
629     {
630         for (j = i + 1; (j < natom); j++)
631         {
632             r = d[i][j] - d_r[i][j];
633             r2 += r * r;
634         }
635     }
636     r2 /= gmx::exactDiv(natom * (natom - 1), 2);
637
638     return std::sqrt(r2);
639 }
640
641 int gmx_rmsdist(int argc, char* argv[])
642 {
643     const char* desc[] = {
644         "[THISMODULE] computes the root mean square deviation of atom distances,",
645         "which has the advantage that no fit is needed like in standard RMS",
646         "deviation as computed by [gmx-rms].",
647         "The reference structure is taken from the structure file.",
648         "The RMSD at time t is calculated as the RMS",
649         "of the differences in distance between atom-pairs in the reference",
650         "structure and the structure at time t.[PAR]",
651         "[THISMODULE] can also produce matrices of the rms distances, rms distances",
652         "scaled with the mean distance and the mean distances and matrices with",
653         "NMR averaged distances (1/r^3 and 1/r^6 averaging). Finally, lists",
654         "of atom pairs with 1/r^3 and 1/r^6 averaged distance below the",
655         "maximum distance ([TT]-max[tt], which will default to 0.6 in this case)",
656         "can be generated, by default averaging over equivalent hydrogens",
657         "(all triplets of hydrogens named \\*[123]). Additionally a list of",
658         "equivalent atoms can be supplied ([TT]-equiv[tt]), each line containing",
659         "a set of equivalent atoms specified as residue number and name and",
660         "atom name; e.g.:[PAR]",
661         "[TT]HB* 3 SER  HB1 3 SER  HB2[tt][PAR]",
662         "Residue and atom names must exactly match those in the structure",
663         "file, including case. Specifying non-sequential atoms is undefined."
664
665     };
666
667     int  i, teller;
668     real t;
669
670     t_topology top;
671     int        ePBC;
672     t_atoms*   atoms;
673     matrix     box;
674     rvec*      x;
675     FILE*      fp;
676
677     t_trxstatus* status;
678     int          isize, gnr = 0;
679     int *        index, *noe_index;
680     char*        grpname;
681     real **      d_r, **d, **dtot, **dtot2, **mean, **rms, **rmsc, *resnr;
682     real **      dtot1_3 = nullptr, **dtot1_6 = nullptr;
683     real         rmsnow, meanmax, rmsmax, rmscmax;
684     real         max1_3, max1_6;
685     t_noe_gr*    noe_gr = nullptr;
686     t_noe**      noe    = nullptr;
687     t_rgb        rlo, rhi;
688     gmx_bool     bRMS, bScale, bMean, bNOE, bNMR3, bNMR6, bNMR;
689
690     static int        nlevels  = 40;
691     static real       scalemax = -1.0;
692     static gmx_bool   bSumH    = TRUE;
693     static gmx_bool   bPBC     = TRUE;
694     gmx_output_env_t* oenv;
695
696     t_pargs pa[] = {
697         { "-nlevels", FALSE, etINT, { &nlevels }, "Discretize RMS in this number of levels" },
698         { "-max", FALSE, etREAL, { &scalemax }, "Maximum level in matrices" },
699         { "-sumh", FALSE, etBOOL, { &bSumH }, "Average distance over equivalent hydrogens" },
700         { "-pbc",
701           FALSE,
702           etBOOL,
703           { &bPBC },
704           "Use periodic boundary conditions when computing distances" }
705     };
706     t_filenm fnm[] = {
707         { efTRX, "-f", nullptr, ffREAD },        { efTPS, nullptr, nullptr, ffREAD },
708         { efNDX, nullptr, nullptr, ffOPTRD },    { efDAT, "-equiv", "equiv", ffOPTRD },
709         { efXVG, nullptr, "distrmsd", ffWRITE }, { efXPM, "-rms", "rmsdist", ffOPTWR },
710         { efXPM, "-scl", "rmsscale", ffOPTWR },  { efXPM, "-mean", "rmsmean", ffOPTWR },
711         { efXPM, "-nmr3", "nmr3", ffOPTWR },     { efXPM, "-nmr6", "nmr6", ffOPTWR },
712         { efDAT, "-noe", "noe", ffOPTWR },
713     };
714 #define NFILE asize(fnm)
715
716     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, asize(pa), pa,
717                            asize(desc), desc, 0, nullptr, &oenv))
718     {
719         return 0;
720     }
721
722     bRMS   = opt2bSet("-rms", NFILE, fnm);
723     bScale = opt2bSet("-scl", NFILE, fnm);
724     bMean  = opt2bSet("-mean", NFILE, fnm);
725     bNOE   = opt2bSet("-noe", NFILE, fnm);
726     bNMR3  = opt2bSet("-nmr3", NFILE, fnm);
727     bNMR6  = opt2bSet("-nmr6", NFILE, fnm);
728     bNMR   = bNMR3 || bNMR6 || bNOE;
729
730     max1_3 = 0;
731     max1_6 = 0;
732
733     /* check input */
734     if (bNOE && scalemax < 0)
735     {
736         scalemax = 0.6;
737         fprintf(stderr,
738                 "WARNING: using -noe without -max "
739                 "makes no sense, setting -max to %g\n\n",
740                 scalemax);
741     }
742
743     /* get topology and index */
744     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &x, nullptr, box, FALSE);
745
746     if (!bPBC)
747     {
748         ePBC = epbcNONE;
749     }
750     atoms = &(top.atoms);
751
752     get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
753
754     /* initialize arrays */
755     snew(d, isize);
756     snew(dtot, isize);
757     snew(dtot2, isize);
758     if (bNMR)
759     {
760         snew(dtot1_3, isize);
761         snew(dtot1_6, isize);
762     }
763     snew(mean, isize);
764     snew(rms, isize);
765     snew(rmsc, isize);
766     snew(d_r, isize);
767     snew(resnr, isize);
768     for (i = 0; (i < isize); i++)
769     {
770         snew(d[i], isize);
771         snew(dtot[i], isize);
772         snew(dtot2[i], isize);
773         if (bNMR)
774         {
775             snew(dtot1_3[i], isize);
776             snew(dtot1_6[i], isize);
777         }
778         snew(mean[i], isize);
779         snew(rms[i], isize);
780         snew(rmsc[i], isize);
781         snew(d_r[i], isize);
782         resnr[i] = i + 1;
783     }
784
785     /*set box type*/
786     calc_dist(isize, index, x, ePBC, box, d_r);
787     sfree(x);
788
789     /*open output files*/
790     fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS Deviation", "Time (ps)", "RMSD (nm)", oenv);
791     if (output_env_get_print_xvgr_codes(oenv))
792     {
793         fprintf(fp, "@ subtitle \"of distances between %s atoms\"\n", grpname);
794     }
795
796     /*do a first step*/
797     read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
798     teller = 0;
799
800     do
801     {
802         calc_dist_tot(isize, index, x, ePBC, box, d, dtot, dtot2, bNMR, dtot1_3, dtot1_6);
803
804         rmsnow = rms_diff(isize, d, d_r);
805         fprintf(fp, "%g  %g\n", t, rmsnow);
806         teller++;
807     } while (read_next_x(oenv, status, &t, x, box));
808     fprintf(stderr, "\n");
809
810     xvgrclose(fp);
811
812     close_trx(status);
813
814     calc_rms(isize, teller, dtot, dtot2, rms, &rmsmax, rmsc, &rmscmax, mean, &meanmax);
815     fprintf(stderr, "rmsmax = %g, rmscmax = %g\n", rmsmax, rmscmax);
816
817     if (bNMR)
818     {
819         calc_nmr(isize, teller, dtot1_3, dtot1_6, &max1_3, &max1_6);
820     }
821
822     if (scalemax > -1.0)
823     {
824         rmsmax  = scalemax;
825         rmscmax = scalemax;
826         meanmax = scalemax;
827         max1_3  = scalemax;
828         max1_6  = scalemax;
829     }
830
831     if (bNOE)
832     {
833         /* make list of noe atom groups */
834         snew(noe_index, isize + 1);
835         snew(noe_gr, isize);
836         gnr = analyze_noe_equivalent(opt2fn_null("-equiv", NFILE, fnm), atoms, isize, index, bSumH,
837                                      noe_index, noe_gr);
838         fprintf(stdout, "Found %d non-equivalent atom-groups in %d atoms\n", gnr, isize);
839         /* make half matrix of of noe-group distances from atom distances */
840         snew(noe, gnr);
841         for (i = 0; i < gnr; i++)
842         {
843             snew(noe[i], gnr);
844         }
845         calc_noe(isize, noe_index, dtot1_3, dtot1_6, gnr, noe);
846     }
847
848     rlo.r = 1.0;
849     rlo.g = 1.0;
850     rlo.b = 1.0;
851     rhi.r = 0.0;
852     rhi.g = 0.0;
853     rhi.b = 0.0;
854
855     if (bRMS)
856     {
857         write_xpm(opt2FILE("-rms", NFILE, fnm, "w"), 0, "RMS of distance", "RMS (nm)", "Atom Index",
858                   "Atom Index", isize, isize, resnr, resnr, rms, 0.0, rmsmax, rlo, rhi, &nlevels);
859     }
860
861     if (bScale)
862     {
863         write_xpm(opt2FILE("-scl", NFILE, fnm, "w"), 0, "Relative RMS", "RMS", "Atom Index",
864                   "Atom Index", isize, isize, resnr, resnr, rmsc, 0.0, rmscmax, rlo, rhi, &nlevels);
865     }
866
867     if (bMean)
868     {
869         write_xpm(opt2FILE("-mean", NFILE, fnm, "w"), 0, "Mean Distance", "Distance (nm)",
870                   "Atom Index", "Atom Index", isize, isize, resnr, resnr, mean, 0.0, meanmax, rlo,
871                   rhi, &nlevels);
872     }
873
874     if (bNMR3)
875     {
876         write_xpm(opt2FILE("-nmr3", NFILE, fnm, "w"), 0, "1/r^3 averaged distances",
877                   "Distance (nm)", "Atom Index", "Atom Index", isize, isize, resnr, resnr, dtot1_3,
878                   0.0, max1_3, rlo, rhi, &nlevels);
879     }
880     if (bNMR6)
881     {
882         write_xpm(opt2FILE("-nmr6", NFILE, fnm, "w"), 0, "1/r^6 averaged distances",
883                   "Distance (nm)", "Atom Index", "Atom Index", isize, isize, resnr, resnr, dtot1_6,
884                   0.0, max1_6, rlo, rhi, &nlevels);
885     }
886
887     if (bNOE)
888     {
889         write_noe(opt2FILE("-noe", NFILE, fnm, "w"), gnr, noe, noe_gr, scalemax);
890     }
891
892     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), nullptr);
893
894     return 0;
895 }