Fix for MSVC linking and release-build OpenMM plugin dir.
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_rmsdist.c
blob02a52640923c4ab462d711e153ccdd03b45eccf7
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <math.h>
40 #include <ctype.h>
42 #include "macros.h"
43 #include "smalloc.h"
44 #include "typedefs.h"
45 #include "names.h"
46 #include "copyrite.h"
47 #include "statutil.h"
48 #include "tpxio.h"
49 #include "string2.h"
50 #include "strdb.h"
51 #include "vec.h"
52 #include "macros.h"
53 #include "index.h"
54 #include "pbc.h"
55 #include "xvgr.h"
56 #include "futil.h"
57 #include "matio.h"
58 #include "gmx_ana.h"
61 static void calc_dist(int nind,atom_id index[],rvec x[],int ePBC,matrix box,
62 real **d)
64 int i,j;
65 real *xi;
66 rvec dx;
67 t_pbc pbc;
69 set_pbc(&pbc,ePBC,box);
70 for(i=0; (i<nind-1); i++) {
71 xi=x[index[i]];
72 for(j=i+1; (j<nind); j++) {
73 pbc_dx(&pbc,xi,x[index[j]],dx);
74 d[i][j]=norm(dx);
79 static void calc_dist_tot(int nind,atom_id index[],rvec x[],
80 int ePBC,matrix box,
81 real **d, real **dtot, real **dtot2,
82 bool bNMR, real **dtot1_3, real **dtot1_6)
84 int i,j;
85 real *xi;
86 real temp, temp2, temp1_3;
87 rvec dx;
88 t_pbc pbc;
90 set_pbc(&pbc,ePBC,box);
91 for(i=0; (i<nind-1); i++) {
92 xi=x[index[i]];
93 for(j=i+1; (j<nind); j++) {
94 pbc_dx(&pbc,xi,x[index[j]],dx);
95 temp2=dx[XX]*dx[XX]+dx[YY]*dx[YY]+dx[ZZ]*dx[ZZ];
96 temp =sqrt(temp2);
97 d[i][j]=temp;
98 dtot[i][j]+=temp;
99 dtot2[i][j]+=temp2;
100 if (bNMR) {
101 temp1_3 = 1.0/(temp*temp2);
102 dtot1_3[i][j] += temp1_3;
103 dtot1_6[i][j] += temp1_3*temp1_3;
109 static void calc_nmr(int nind, int nframes, real **dtot1_3, real **dtot1_6,
110 real *max1_3, real *max1_6)
112 int i,j;
113 real temp1_3,temp1_6;
115 for(i=0; (i<nind-1); i++) {
116 for(j=i+1; (j<nind); j++) {
117 temp1_3 = pow(dtot1_3[i][j]/nframes,-1.0/3.0);
118 temp1_6 = pow(dtot1_6[i][j]/nframes,-1.0/6.0);
119 if (temp1_3 > *max1_3) *max1_3 = temp1_3;
120 if (temp1_6 > *max1_6) *max1_6 = temp1_6;
121 dtot1_3[i][j] = temp1_3;
122 dtot1_6[i][j] = temp1_6;
123 dtot1_3[j][i] = temp1_3;
124 dtot1_6[j][i] = temp1_6;
129 static char Hnum[] = "123";
131 typedef struct {
132 int nr;
133 real r_3;
134 real r_6;
135 real i_3;
136 real i_6;
137 } t_noe;
139 typedef struct {
140 int anr;
141 int ianr;
142 int rnr;
143 char *aname;
144 char *rname;
145 } t_noe_gr;
147 typedef struct {
148 int rnr;
149 char *nname;
150 char *rname;
151 char *aname;
152 } t_equiv;
154 static int read_equiv(const char *eq_fn, t_equiv ***equivptr)
156 FILE *fp;
157 char line[STRLEN],resname[10],atomname[10],*lp;
158 int neq, na, n, resnr;
159 t_equiv **equiv;
161 fp = ffopen(eq_fn,"r");
162 neq=0;
163 equiv=NULL;
164 while(get_a_line(fp,line,STRLEN)) {
165 lp=line;
166 /* this is not efficient, but I'm lazy */
167 srenew(equiv,neq+1);
168 equiv[neq]=NULL;
169 na=0;
170 if (sscanf(lp,"%s %n",atomname,&n)==1) {
171 lp+=n;
172 snew(equiv[neq], 1);
173 equiv[neq][0].nname=strdup(atomname);
174 while (sscanf(lp,"%d %s %s %n",&resnr,resname,atomname,&n)==3) {
175 /* this is not efficient, but I'm lazy (again) */
176 srenew(equiv[neq], na+1);
177 equiv[neq][na].rnr=resnr-1;
178 equiv[neq][na].rname=strdup(resname);
179 equiv[neq][na].aname=strdup(atomname);
180 if (na>0) equiv[neq][na].nname=NULL;
181 na++;
182 lp+=n;
185 /* make empty element as flag for end of array */
186 srenew(equiv[neq], na+1);
187 equiv[neq][na].rnr=NOTSET;
188 equiv[neq][na].rname=NULL;
189 equiv[neq][na].aname=NULL;
191 /* next */
192 neq++;
194 ffclose(fp);
196 *equivptr=equiv;
198 return neq;
201 static void dump_equiv(FILE *out, int neq, t_equiv **equiv)
203 int i,j;
205 fprintf(out,"Dumping equivalent list\n");
206 for (i=0; i<neq; i++) {
207 fprintf(out,"%s",equiv[i][0].nname);
208 for(j=0; equiv[i][j].rnr!=NOTSET; j++)
209 fprintf(out," %d %s %s",
210 equiv[i][j].rnr,equiv[i][j].rname,equiv[i][j].aname);
211 fprintf(out,"\n");
215 static bool is_equiv(int neq, t_equiv **equiv, char **nname,
216 int rnr1, char *rname1, char *aname1,
217 int rnr2, char *rname2, char *aname2)
219 int i,j;
220 bool bFound;
222 bFound=FALSE;
223 /* we can terminate each loop when bFound is true! */
224 for (i=0; i<neq && !bFound; i++) {
225 /* find first atom */
226 for(j=0; equiv[i][j].rnr!=NOTSET && !bFound; j++)
227 bFound = ( equiv[i][j].rnr==rnr1 &&
228 strcmp(equiv[i][j].rname, rname1)==0 &&
229 strcmp(equiv[i][j].aname, aname1)==0 );
230 if (bFound) {
231 /* find second atom */
232 bFound=FALSE;
233 for(j=0; equiv[i][j].rnr!=NOTSET && !bFound; j++)
234 bFound = ( equiv[i][j].rnr==rnr2 &&
235 strcmp(equiv[i][j].rname, rname2)==0 &&
236 strcmp(equiv[i][j].aname, aname2)==0 );
239 if (bFound)
240 *nname = strdup(equiv[i-1][0].nname);
242 return bFound;
245 static int analyze_noe_equivalent(const char *eq_fn,
246 t_atoms *atoms, int isize, atom_id *index,
247 bool bSumH,
248 atom_id *noe_index, t_noe_gr *noe_gr)
250 FILE *fp;
251 int i, j, anmil, anmjl, rnri, rnrj, gi, groupnr, neq;
252 char *anmi, *anmj, **nnm;
253 bool bMatch,bEquiv;
254 t_equiv **equiv;
256 snew(nnm,isize);
257 if (bSumH) {
258 if (eq_fn) {
259 neq=read_equiv(eq_fn,&equiv);
260 if (debug) dump_equiv(debug,neq,equiv);
261 } else {
262 neq=0;
263 equiv=NULL;
266 groupnr=0;
267 for(i=0; i<isize; i++) {
268 if (equiv && i<isize-1) {
269 /* check explicit list of equivalent atoms */
270 do {
271 j=i+1;
272 rnri=atoms->atom[index[i]].resind;
273 rnrj=atoms->atom[index[j]].resind;
274 bEquiv =
275 is_equiv(neq, equiv, &nnm[i],
276 rnri,*atoms->resinfo[rnri].name,*atoms->atomname[index[i]],
277 rnrj,*atoms->resinfo[rnrj].name,*atoms->atomname[index[j]]);
278 if(nnm[i] && bEquiv)
279 nnm[j]=strdup(nnm[i]);
280 if (bEquiv) {
281 /* set index for matching atom */
282 noe_index[j]=groupnr;
283 /* skip matching atom */
284 i=j;
286 } while ( bEquiv && i<isize-1 );
287 } else
288 bEquiv=FALSE;
289 if (!bEquiv) {
290 /* look for triplets of consecutive atoms with name XX?,
291 X are any number of letters or digits and ? goes from 1 to 3
292 This is supposed to cover all CH3 groups and the like */
293 anmi = *atoms->atomname[index[i]];
294 anmil = strlen(anmi);
295 bMatch = i<isize-3 && anmi[anmil-1]=='1';
296 if (bMatch)
297 for(j=1; j<3; j++) {
298 anmj = *atoms->atomname[index[i+j]];
299 anmjl = strlen(anmj);
300 bMatch = bMatch && ( anmil==anmjl && anmj[anmjl-1]==Hnum[j] &&
301 strncmp(anmi, anmj, anmil-1)==0 );
303 /* set index for this atom */
304 noe_index[i]=groupnr;
305 if (bMatch) {
306 /* set index for next two matching atoms */
307 for(j=1; j<3; j++)
308 noe_index[i+j]=groupnr;
309 /* skip matching atoms */
310 i+=2;
313 groupnr++;
315 } else {
316 /* make index without looking for equivalent atoms */
317 for(i=0; i<isize; i++)
318 noe_index[i]=i;
319 groupnr=isize;
321 noe_index[isize]=groupnr;
323 if (debug)
324 /* dump new names */
325 for(i=0; i<isize; i++) {
326 rnri=atoms->atom[index[i]].resind;
327 fprintf(debug,"%s %s %d -> %s\n",*atoms->atomname[index[i]],
328 *atoms->resinfo[rnri].name,rnri,nnm[i]?nnm[i]:"");
331 for(i=0; i<isize; i++) {
332 gi=noe_index[i];
333 if (!noe_gr[gi].aname) {
334 noe_gr[gi].ianr=i;
335 noe_gr[gi].anr=index[i];
336 if (nnm[i])
337 noe_gr[gi].aname=strdup(nnm[i]);
338 else {
339 noe_gr[gi].aname=strdup(*atoms->atomname[index[i]]);
340 if ( noe_index[i]==noe_index[i+1] )
341 noe_gr[gi].aname[strlen(noe_gr[gi].aname)-1]='*';
343 noe_gr[gi].rnr=atoms->atom[index[i]].resind;
344 noe_gr[gi].rname=strdup(*atoms->resinfo[noe_gr[gi].rnr].name);
345 /* dump group definitions */
346 if (debug) fprintf(debug,"%d %d %d %d %s %s %d\n",i,gi,
347 noe_gr[gi].ianr,noe_gr[gi].anr,noe_gr[gi].aname,
348 noe_gr[gi].rname,noe_gr[gi].rnr);
351 for(i=0; i<isize; i++)
352 sfree(nnm[i]);
353 sfree(nnm);
355 return groupnr;
358 /* #define NSCALE 3 */
359 /* static char *noe_scale[NSCALE+1] = { */
360 /* "strong", "medium", "weak", "none" */
361 /* }; */
362 #define NSCALE 6
364 static char *noe2scale(real r3, real r6, real rmax)
366 int i,s3, s6;
367 static char buf[NSCALE+1];
369 /* r goes from 0 to rmax
370 NSCALE*r/rmax goes from 0 to NSCALE
371 NSCALE - NSCALE*r/rmax goes from NSCALE to 0 */
372 s3 = NSCALE - min(NSCALE, (int)(NSCALE*r3/rmax));
373 s6 = NSCALE - min(NSCALE, (int)(NSCALE*r6/rmax));
375 for(i=0; i<s3; i++)
376 buf[i]='=';
377 for( ; i<s6; i++)
378 buf[i]='-';
379 buf[i]='\0';
381 return buf;
384 static void calc_noe(int isize, atom_id *noe_index,
385 real **dtot1_3, real **dtot1_6, int gnr, t_noe **noe)
387 int i, j, gi, gj;
389 /* make half matrix of noe-group distances from atom distances */
390 for(i=0; i<isize; i++) {
391 gi=noe_index[i];
392 for(j=i; j<isize; j++) {
393 gj=noe_index[j];
394 noe[gi][gj].nr++;
395 noe[gi][gj].i_3+=pow(dtot1_3[i][j],-3);
396 noe[gi][gj].i_6+=pow(dtot1_6[i][j],-6);
400 /* make averages */
401 for(i=0; i<gnr; i++)
402 for(j=i+1; j<gnr; j++) {
403 noe[i][j].r_3 = pow(noe[i][j].i_3/noe[i][j].nr,-1.0/3.0);
404 noe[i][j].r_6 = pow(noe[i][j].i_6/noe[i][j].nr,-1.0/6.0);
405 noe[j][i] = noe[i][j];
409 static void write_noe(FILE *fp,int gnr,t_noe **noe,t_noe_gr *noe_gr,real rmax)
411 int i,j;
412 real r3, r6, min3, min6;
413 char buf[10],b3[10],b6[10];
414 t_noe_gr gri, grj;
416 min3 = min6 = 1e6;
417 fprintf(fp,
418 ";%4s %3s %4s %4s%3s %4s %4s %4s %4s%3s %5s %5s %8s %2s %2s %s\n",
419 "ianr","anr","anm","rnm","rnr","ianr","anr","anm","rnm","rnr",
420 "1/r^3","1/r^6","intnsty","Dr","Da","scale");
421 for(i=0; i<gnr; i++) {
422 gri=noe_gr[i];
423 for(j=i+1; j<gnr; j++) {
424 grj=noe_gr[j];
425 r3 = noe[i][j].r_3;
426 r6 = noe[i][j].r_6;
427 min3 = min(r3,min3);
428 min6 = min(r6,min6);
429 if ( r3 < rmax || r6 < rmax ) {
430 if (grj.rnr == gri.rnr)
431 sprintf(buf,"%2d", grj.anr-gri.anr);
432 else
433 buf[0]='\0';
434 if ( r3 < rmax )
435 sprintf(b3,"%-5.3f",r3);
436 else
437 strcpy(b3,"-");
438 if ( r6 < rmax )
439 sprintf(b6,"%-5.3f",r6);
440 else
441 strcpy(b6,"-");
442 fprintf(fp,
443 "%4d %4d %4s %4s%3d %4d %4d %4s %4s%3d %5s %5s %8d %2d %2s %s\n",
444 gri.ianr+1, gri.anr+1, gri.aname, gri.rname, gri.rnr+1,
445 grj.ianr+1, grj.anr+1, grj.aname, grj.rname, grj.rnr+1,
446 b3, b6, (int)(noe[i][j].i_6+0.5), grj.rnr-gri.rnr, buf,
447 noe2scale(r3, r6, rmax));
451 #define MINI ((i==3)?min3:min6)
452 for(i=3; i<=6; i+=3)
453 if ( MINI > rmax )
454 fprintf(stdout,"NOTE: no 1/r^%d averaged distances found below %g, "
455 "smallest was %g\n", i, rmax, MINI );
456 else
457 fprintf(stdout,"Smallest 1/r^%d averaged distance was %g\n", i, MINI );
458 #undef MINI
461 static void calc_rms(int nind, int nframes,
462 real **dtot, real **dtot2,
463 real **rmsmat, real *rmsmax,
464 real **rmscmat, real *rmscmax,
465 real **meanmat, real *meanmax)
467 int i,j;
468 real mean, mean2, rms, rmsc;
469 /* N.B. dtot and dtot2 contain the total distance and the total squared
470 * distance respectively, BUT they return RMS and the scaled RMS resp.
473 *rmsmax=-1000;
474 *rmscmax=-1000;
475 *meanmax=-1000;
477 for(i=0; (i<nind-1); i++) {
478 for(j=i+1; (j<nind); j++) {
479 mean =dtot[i][j]/nframes;
480 mean2=dtot2[i][j]/nframes;
481 rms=sqrt(max(0,mean2-mean*mean));
482 rmsc=rms/mean;
483 if (mean > *meanmax) *meanmax=mean;
484 if (rms > *rmsmax ) *rmsmax =rms;
485 if (rmsc > *rmscmax) *rmscmax=rmsc;
486 meanmat[i][j]=meanmat[j][i]=mean;
487 rmsmat[i][j] =rmsmat[j][i] =rms;
488 rmscmat[i][j]=rmscmat[j][i]=rmsc;
493 real rms_diff(int natom,real **d,real **d_r)
495 int i,j;
496 real r,r2;
498 r2=0.0;
499 for(i=0; (i<natom-1); i++)
500 for(j=i+1; (j<natom); j++) {
501 r=d[i][j]-d_r[i][j];
502 r2+=r*r;
504 r2/=(natom*(natom-1))/2;
506 return sqrt(r2);
509 int gmx_rmsdist (int argc,char *argv[])
511 const char *desc[] = {
512 "g_rmsdist computes the root mean square deviation of atom distances,",
513 "which has the advantage that no fit is needed like in standard RMS",
514 "deviation as computed by g_rms.",
515 "The reference structure is taken from the structure file.",
516 "The rmsd at time t is calculated as the rms",
517 "of the differences in distance between atom-pairs in the reference",
518 "structure and the structure at time t.[PAR]",
519 "g_rmsdist can also produce matrices of the rms distances, rms distances",
520 "scaled with the mean distance and the mean distances and matrices with",
521 "NMR averaged distances (1/r^3 and 1/r^6 averaging). Finally, lists",
522 "of atom pairs with 1/r^3 and 1/r^6 averaged distance below the",
523 "maximum distance ([TT]-max[tt], which will default to 0.6 in this case)",
524 "can be generated, by default averaging over equivalent hydrogens",
525 "(all triplets of hydrogens named *[123]). Additionally a list of",
526 "equivalent atoms can be supplied ([TT]-equiv[tt]), each line containing",
527 "a set of equivalent atoms specified as residue number and name and",
528 "atom name; e.g.:[PAR]",
529 "[TT]3 SER HB1 3 SER HB2[tt][PAR]",
530 "Residue and atom names must exactly match those in the structure",
531 "file, including case. Specifying non-sequential atoms is undefined."
535 int natom,i,j,teller,gi,gj;
536 real t;
538 t_topology top;
539 int ePBC;
540 t_atoms *atoms;
541 matrix box;
542 rvec *x;
543 FILE *fp;
545 int status,isize,gnr=0;
546 atom_id *index, *noe_index;
547 char *grpname;
548 real **d_r,**d,**dtot,**dtot2,**mean,**rms,**rmsc,*resnr;
549 real **dtot1_3=NULL,**dtot1_6=NULL;
550 real rmsnow,meanmax,rmsmax,rmscmax;
551 real max1_3,max1_6;
552 t_noe_gr *noe_gr=NULL;
553 t_noe **noe=NULL;
554 t_rgb rlo,rhi;
555 char buf[255];
556 bool bRMS, bScale, bMean, bNOE, bNMR3, bNMR6, bNMR;
558 static int nlevels=40;
559 static real scalemax=-1.0;
560 static bool bSumH=TRUE;
562 output_env_t oenv;
564 t_pargs pa[] = {
565 { "-nlevels", FALSE, etINT, {&nlevels},
566 "Discretize rms in # levels" },
567 { "-max", FALSE, etREAL, {&scalemax},
568 "Maximum level in matrices" },
569 { "-sumh", FALSE, etBOOL, {&bSumH},
570 "average distance over equivalent hydrogens" }
572 t_filenm fnm[] = {
573 { efTRX, "-f", NULL, ffREAD },
574 { efTPS, NULL, NULL, ffREAD },
575 { efNDX, NULL, NULL, ffOPTRD },
576 { efDAT, "-equiv","equiv", ffOPTRD },
577 { efXVG, NULL, "distrmsd", ffWRITE },
578 { efXPM, "-rms", "rmsdist", ffOPTWR },
579 { efXPM, "-scl", "rmsscale", ffOPTWR },
580 { efXPM, "-mean","rmsmean", ffOPTWR },
581 { efXPM, "-nmr3","nmr3", ffOPTWR },
582 { efXPM, "-nmr6","nmr6", ffOPTWR },
583 { efDAT, "-noe", "noe", ffOPTWR },
585 #define NFILE asize(fnm)
587 CopyRight(stderr,argv[0]);
588 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE ,
589 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
591 bRMS = opt2bSet("-rms", NFILE,fnm);
592 bScale= opt2bSet("-scl", NFILE,fnm);
593 bMean = opt2bSet("-mean",NFILE,fnm);
594 bNOE = opt2bSet("-noe", NFILE,fnm);
595 bNMR3 = opt2bSet("-nmr3",NFILE,fnm);
596 bNMR6 = opt2bSet("-nmr6",NFILE,fnm);
597 bNMR = bNMR3 || bNMR6 || bNOE;
599 max1_3 = 0;
600 max1_6 = 0;
602 /* check input */
603 if (bNOE && scalemax < 0) {
604 scalemax=0.6;
605 fprintf(stderr,"WARNING: using -noe without -max "
606 "makes no sense, setting -max to %g\n\n",scalemax);
609 /* get topology and index */
610 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),buf,&top,&ePBC,&x,NULL,box,FALSE);
611 atoms=&(top.atoms);
613 get_index(atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
615 /* initialize arrays */
616 snew(d,isize);
617 snew(dtot,isize);
618 snew(dtot2,isize);
619 if (bNMR) {
620 snew(dtot1_3,isize);
621 snew(dtot1_6,isize);
623 snew(mean,isize);
624 snew(rms,isize);
625 snew(rmsc,isize);
626 snew(d_r,isize);
627 snew(resnr,isize);
628 for(i=0; (i<isize); i++) {
629 snew(d[i],isize);
630 snew(dtot[i],isize);
631 snew(dtot2[i],isize);
632 if (bNMR) {
633 snew(dtot1_3[i],isize);
634 snew(dtot1_6[i],isize);
636 snew(mean[i],isize);
637 snew(rms[i],isize);
638 snew(rmsc[i],isize);
639 snew(d_r[i],isize);
640 resnr[i]=i+1;
643 /*set box type*/
644 calc_dist(isize,index,x,ePBC,box,d_r);
645 sfree(x);
647 /*open output files*/
648 fp=xvgropen(ftp2fn(efXVG,NFILE,fnm),"RMS Deviation","Time (ps)","RMSD (nm)",
649 oenv);
650 if (output_env_get_print_xvgr_codes(oenv))
651 fprintf(fp,"@ subtitle \"of distances between %s atoms\"\n",grpname);
653 /*do a first step*/
654 natom=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
656 do {
657 calc_dist_tot(isize,index,x,ePBC,box,d,dtot,dtot2,bNMR,dtot1_3,dtot1_6);
659 rmsnow=rms_diff(isize,d,d_r);
660 fprintf(fp,"%g %g\n",t,rmsnow);
661 } while (read_next_x(oenv,status,&t,natom,x,box));
662 fprintf(stderr, "\n");
664 ffclose(fp);
665 close_trj(status);
667 teller = nframes_read();
668 calc_rms(isize,teller,dtot,dtot2,mean,&meanmax,rms,&rmsmax,rmsc,&rmscmax);
669 fprintf(stderr,"rmsmax = %g, rmscmax = %g\n",rmsmax,rmscmax);
671 if (bNMR)
673 calc_nmr(isize,teller,dtot1_3,dtot1_6,&max1_3,&max1_6);
676 if (scalemax > -1.0) {
677 rmsmax=scalemax;
678 rmscmax=scalemax;
679 meanmax=scalemax;
680 max1_3=scalemax;
681 max1_6=scalemax;
684 if (bNOE) {
685 /* make list of noe atom groups */
686 snew(noe_index, isize+1);
687 snew(noe_gr, isize);
688 gnr=analyze_noe_equivalent(opt2fn_null("-equiv",NFILE,fnm),
689 atoms, isize, index, bSumH, noe_index, noe_gr);
690 fprintf(stdout, "Found %d non-equivalent atom-groups in %d atoms\n",
691 gnr, isize);
692 /* make half matrix of of noe-group distances from atom distances */
693 snew(noe,gnr);
694 for(i=0; i<gnr; i++)
695 snew(noe[i],gnr);
696 calc_noe(isize, noe_index, dtot1_3, dtot1_6, gnr, noe);
699 rlo.r=1.0, rlo.g=1.0, rlo.b=1.0;
700 rhi.r=0.0, rhi.g=0.0, rhi.b=0.0;
702 if ( bRMS )
703 write_xpm(opt2FILE("-rms",NFILE,fnm,"w"),0,
704 "RMS of distance","RMS (nm)","Atom Index","Atom Index",
705 isize,isize,resnr,resnr,rms,0.0,rmsmax,rlo,rhi,&nlevels);
707 if ( bScale )
708 write_xpm(opt2FILE("-scl",NFILE,fnm,"w"),0,
709 "Relative RMS","RMS","Atom Index","Atom Index",
710 isize,isize,resnr,resnr,rmsc,0.0,rmscmax,rlo,rhi,&nlevels);
712 if ( bMean )
713 write_xpm(opt2FILE("-mean",NFILE,fnm,"w"),0,
714 "Mean Distance","Distance (nm)","Atom Index","Atom Index",
715 isize,isize,resnr,resnr,mean,0.0,meanmax,rlo,rhi,&nlevels);
717 if (bNMR3)
718 write_xpm(opt2FILE("-nmr3",NFILE,fnm,"w"),0,"1/r^3 averaged distances",
719 "Distance (nm)","Atom Index","Atom Index",
720 isize,isize,resnr,resnr,dtot1_3,0.0,max1_3,rlo,rhi,&nlevels);
721 if (bNMR6)
722 write_xpm(opt2FILE("-nmr6",NFILE,fnm,"w"),0,"1/r^6 averaged distances",
723 "Distance (nm)","Atom Index","Atom Index",
724 isize,isize,resnr,resnr,dtot1_6,0.0,max1_6,rlo,rhi,&nlevels);
726 if (bNOE)
727 write_noe(opt2FILE("-noe",NFILE,fnm,"w"), gnr, noe, noe_gr, scalemax);
729 do_view(oenv,ftp2fn(efXVG,NFILE,fnm),NULL);
731 thanx(stderr);
732 return 0;