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