Parallel vs. sequentiual code: I get binary identical result (apart from comment...
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_covar.c
blob57ee9ca47132522b6b19ed185f1c9dfb2ff7f835
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
38 #include <math.h>
39 #include <string.h>
40 #include <time.h>
42 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
43 #include <direct.h>
44 #include <io.h>
45 #endif
47 #include "statutil.h"
48 #include "sysstuff.h"
49 #include "typedefs.h"
50 #include "smalloc.h"
51 #include "macros.h"
52 #include "vec.h"
53 #include "pbc.h"
54 #include "copyrite.h"
55 #include "futil.h"
56 #include "statutil.h"
57 #include "index.h"
58 #include "confio.h"
59 #include "trnio.h"
60 #include "mshift.h"
61 #include "xvgr.h"
62 #include "do_fit.h"
63 #include "rmpbc.h"
64 #include "txtdump.h"
65 #include "matio.h"
66 #include "eigio.h"
67 #include "eigensolver.h"
68 #include "physics.h"
69 #include "gmx_ana.h"
72 int gmx_covar(int argc,char *argv[])
74 const char *desc[] = {
75 "[TT]g_covar[tt] calculates and diagonalizes the (mass-weighted)",
76 "covariance matrix.",
77 "All structures are fitted to the structure in the structure file.",
78 "When this is not a run input file periodicity will not be taken into",
79 "account. When the fit and analysis groups are identical and the analysis",
80 "is non mass-weighted, the fit will also be non mass-weighted.",
81 "[PAR]",
82 "The eigenvectors are written to a trajectory file ([TT]-v[tt]).",
83 "When the same atoms are used for the fit and the covariance analysis,",
84 "the reference structure for the fit is written first with t=-1.",
85 "The average (or reference when [TT]-ref[tt] is used) structure is",
86 "written with t=0, the eigenvectors",
87 "are written as frames with the eigenvector number as timestamp.",
88 "[PAR]",
89 "The eigenvectors can be analyzed with [TT]g_anaeig[tt].",
90 "[PAR]",
91 "Option [TT]-ascii[tt] writes the whole covariance matrix to",
92 "an ASCII file. The order of the elements is: x1x1, x1y1, x1z1, x1x2, ...",
93 "[PAR]",
94 "Option [TT]-xpm[tt] writes the whole covariance matrix to an xpm file.",
95 "[PAR]",
96 "Option [TT]-xpma[tt] writes the atomic covariance matrix to an xpm file,",
97 "i.e. for each atom pair the sum of the xx, yy and zz covariances is",
98 "written."
100 static bool bFit=TRUE,bRef=FALSE,bM=FALSE,bPBC=TRUE;
101 static int end=-1;
102 t_pargs pa[] = {
103 { "-fit", FALSE, etBOOL, {&bFit},
104 "Fit to a reference structure"},
105 { "-ref", FALSE, etBOOL, {&bRef},
106 "Use the deviation from the conformation in the structure file instead of from the average" },
107 { "-mwa", FALSE, etBOOL, {&bM},
108 "Mass-weighted covariance analysis"},
109 { "-last", FALSE, etINT, {&end},
110 "Last eigenvector to write away (-1 is till the last)" },
111 { "-pbc", FALSE, etBOOL, {&bPBC},
112 "Apply corrections for periodic boundary conditions" }
114 FILE *out;
115 int status,trjout;
116 t_topology top;
117 int ePBC;
118 t_atoms *atoms;
119 rvec *x,*xread,*xref,*xav,*xproj;
120 matrix box,zerobox;
121 real *sqrtm,*mat,*eigval,sum,trace,inv_nframes;
122 real t,tstart,tend,**mat2;
123 real xj,*w_rls=NULL;
124 real min,max,*axis;
125 int ntopatoms,step;
126 int natoms,nat,count,nframes0,nframes,nlevels;
127 gmx_large_int_t ndim,i,j,k,l;
128 int WriteXref;
129 const char *fitfile,*trxfile,*ndxfile;
130 const char *eigvalfile,*eigvecfile,*averfile,*logfile;
131 const char *asciifile,*xpmfile,*xpmafile;
132 char str[STRLEN],*fitname,*ananame,*pcwd;
133 int d,dj,nfit;
134 atom_id *index,*ifit;
135 bool bDiffMass1,bDiffMass2;
136 time_t now;
137 t_rgb rlo,rmi,rhi;
138 real *tmp;
139 output_env_t oenv;
141 t_filenm fnm[] = {
142 { efTRX, "-f", NULL, ffREAD },
143 { efTPS, NULL, NULL, ffREAD },
144 { efNDX, NULL, NULL, ffOPTRD },
145 { efXVG, NULL, "eigenval", ffWRITE },
146 { efTRN, "-v", "eigenvec", ffWRITE },
147 { efSTO, "-av", "average.pdb", ffWRITE },
148 { efLOG, NULL, "covar", ffWRITE },
149 { efDAT, "-ascii","covar", ffOPTWR },
150 { efXPM, "-xpm","covar", ffOPTWR },
151 { efXPM, "-xpma","covara", ffOPTWR }
153 #define NFILE asize(fnm)
155 CopyRight(stderr,argv[0]);
156 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
157 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
159 clear_mat(zerobox);
161 fitfile = ftp2fn(efTPS,NFILE,fnm);
162 trxfile = ftp2fn(efTRX,NFILE,fnm);
163 ndxfile = ftp2fn_null(efNDX,NFILE,fnm);
164 eigvalfile = ftp2fn(efXVG,NFILE,fnm);
165 eigvecfile = ftp2fn(efTRN,NFILE,fnm);
166 averfile = ftp2fn(efSTO,NFILE,fnm);
167 logfile = ftp2fn(efLOG,NFILE,fnm);
168 asciifile = opt2fn_null("-ascii",NFILE,fnm);
169 xpmfile = opt2fn_null("-xpm",NFILE,fnm);
170 xpmafile = opt2fn_null("-xpma",NFILE,fnm);
172 read_tps_conf(fitfile,str,&top,&ePBC,&xref,NULL,box,TRUE);
173 atoms=&top.atoms;
175 if (bFit) {
176 printf("\nChoose a group for the least squares fit\n");
177 get_index(atoms,ndxfile,1,&nfit,&ifit,&fitname);
178 if (nfit < 3)
179 gmx_fatal(FARGS,"Need >= 3 points to fit!\n");
180 } else
181 nfit=0;
182 printf("\nChoose a group for the covariance analysis\n");
183 get_index(atoms,ndxfile,1,&natoms,&index,&ananame);
185 bDiffMass1=FALSE;
186 if (bFit) {
187 snew(w_rls,atoms->nr);
188 for(i=0; (i<nfit); i++) {
189 w_rls[ifit[i]]=atoms->atom[ifit[i]].m;
190 if (i)
191 bDiffMass1 = bDiffMass1 || (w_rls[ifit[i]]!=w_rls[ifit[i-1]]);
194 bDiffMass2=FALSE;
195 snew(sqrtm,natoms);
196 for(i=0; (i<natoms); i++)
197 if (bM) {
198 sqrtm[i]=sqrt(atoms->atom[index[i]].m);
199 if (i)
200 bDiffMass2 = bDiffMass2 || (sqrtm[i]!=sqrtm[i-1]);
202 else
203 sqrtm[i]=1.0;
205 if (bFit && bDiffMass1 && !bDiffMass2) {
206 bDiffMass1 = natoms != nfit;
207 i=0;
208 for (i=0; (i<natoms) && !bDiffMass1; i++)
209 bDiffMass1 = index[i] != ifit[i];
210 if (!bDiffMass1) {
211 fprintf(stderr,"\n"
212 "Note: the fit and analysis group are identical,\n"
213 " while the fit is mass weighted and the analysis is not.\n"
214 " Making the fit non mass weighted.\n\n");
215 for(i=0; (i<nfit); i++)
216 w_rls[ifit[i]]=1.0;
220 /* Prepare reference frame */
221 if (bPBC)
222 rm_pbc(&(top.idef),ePBC,atoms->nr,box,xref,xref);
223 if (bFit)
224 reset_x(nfit,ifit,atoms->nr,NULL,xref,w_rls);
226 snew(x,natoms);
227 snew(xav,natoms);
228 ndim=natoms*DIM;
229 if (sqrt(LARGE_INT_MAX)<ndim) {
230 gmx_fatal(FARGS,"Number of degrees of freedoms to large for matrix.\n");
232 snew(mat,ndim*ndim);
234 fprintf(stderr,"Calculating the average structure ...\n");
235 nframes0 = 0;
236 nat=read_first_x(oenv,&status,trxfile,&t,&xread,box);
237 if (nat != atoms->nr)
238 fprintf(stderr,"\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n",natoms,nat);
239 do {
240 nframes0++;
241 /* calculate x: a fitted struture of the selected atoms */
242 if (bPBC)
243 rm_pbc(&(top.idef),ePBC,nat,box,xread,xread);
244 if (bFit) {
245 reset_x(nfit,ifit,nat,NULL,xread,w_rls);
246 do_fit(nat,w_rls,xref,xread);
248 for (i=0; i<natoms; i++)
249 rvec_inc(xav[i],xread[index[i]]);
250 } while (read_next_x(oenv,status,&t,nat,xread,box));
251 close_trj(status);
253 inv_nframes = 1.0/nframes0;
254 for(i=0; i<natoms; i++)
255 for(d=0; d<DIM; d++) {
256 xav[i][d] *= inv_nframes;
257 xread[index[i]][d] = xav[i][d];
259 write_sto_conf_indexed(opt2fn("-av",NFILE,fnm),"Average structure",
260 atoms,xread,NULL,epbcNONE,zerobox,natoms,index);
261 sfree(xread);
263 fprintf(stderr,"Constructing covariance matrix (%dx%d) ...\n",(int)ndim,(int)ndim);
264 nframes=0;
265 nat=read_first_x(oenv,&status,trxfile,&t,&xread,box);
266 tstart = t;
267 do {
268 nframes++;
269 tend = t;
270 /* calculate x: a (fitted) structure of the selected atoms */
271 if (bPBC)
272 rm_pbc(&(top.idef),ePBC,nat,box,xread,xread);
273 if (bFit) {
274 reset_x(nfit,ifit,nat,NULL,xread,w_rls);
275 do_fit(nat,w_rls,xref,xread);
277 if (bRef)
278 for (i=0; i<natoms; i++)
279 rvec_sub(xread[index[i]],xref[index[i]],x[i]);
280 else
281 for (i=0; i<natoms; i++)
282 rvec_sub(xread[index[i]],xav[i],x[i]);
284 for (j=0; j<natoms; j++) {
285 for (dj=0; dj<DIM; dj++) {
286 k=ndim*(DIM*j+dj);
287 xj=x[j][dj];
288 for (i=j; i<natoms; i++) {
289 l=k+DIM*i;
290 for(d=0; d<DIM; d++)
291 mat[l+d] += x[i][d]*xj;
295 } while (read_next_x(oenv,status,&t,nat,xread,box) &&
296 (bRef || nframes < nframes0));
297 close_trj(status);
299 fprintf(stderr,"Read %d frames\n",nframes);
301 if (bRef) {
302 /* copy the reference structure to the ouput array x */
303 snew(xproj,natoms);
304 for (i=0; i<natoms; i++)
305 copy_rvec(xref[index[i]],xproj[i]);
306 } else {
307 xproj = xav;
310 /* correct the covariance matrix for the mass */
311 inv_nframes = 1.0/nframes;
312 for (j=0; j<natoms; j++)
313 for (dj=0; dj<DIM; dj++)
314 for (i=j; i<natoms; i++) {
315 k = ndim*(DIM*j+dj)+DIM*i;
316 for (d=0; d<DIM; d++)
317 mat[k+d] = mat[k+d]*inv_nframes*sqrtm[i]*sqrtm[j];
320 /* symmetrize the matrix */
321 for (j=0; j<ndim; j++)
322 for (i=j; i<ndim; i++)
323 mat[ndim*i+j]=mat[ndim*j+i];
325 trace=0;
326 for(i=0; i<ndim; i++)
327 trace+=mat[i*ndim+i];
328 fprintf(stderr,"\nTrace of the covariance matrix: %g (%snm^2)\n",
329 trace,bM ? "u " : "");
331 if (asciifile) {
332 out = ffopen(asciifile,"w");
333 for (j=0; j<ndim; j++) {
334 for (i=0; i<ndim; i+=3)
335 fprintf(out,"%g %g %g\n",
336 mat[ndim*j+i],mat[ndim*j+i+1],mat[ndim*j+i+2]);
338 ffclose(out);
341 if (xpmfile) {
342 min = 0;
343 max = 0;
344 snew(mat2,ndim);
345 for (j=0; j<ndim; j++) {
346 mat2[j] = &(mat[ndim*j]);
347 for (i=0; i<=j; i++) {
348 if (mat2[j][i] < min)
349 min = mat2[j][i];
350 if (mat2[j][j] > max)
351 max = mat2[j][i];
354 snew(axis,ndim);
355 for(i=0; i<ndim; i++)
356 axis[i] = i+1;
357 rlo.r = 0; rlo.g = 0; rlo.b = 1;
358 rmi.r = 1; rmi.g = 1; rmi.b = 1;
359 rhi.r = 1; rhi.g = 0; rhi.b = 0;
360 out = ffopen(xpmfile,"w");
361 nlevels = 80;
362 write_xpm3(out,0,"Covariance",bM ? "u nm^2" : "nm^2",
363 "dim","dim",ndim,ndim,axis,axis,
364 mat2,min,0.0,max,rlo,rmi,rhi,&nlevels);
365 ffclose(out);
366 sfree(axis);
367 sfree(mat2);
370 if (xpmafile) {
371 min = 0;
372 max = 0;
373 snew(mat2,ndim/DIM);
374 for (i=0; i<ndim/DIM; i++)
375 snew(mat2[i],ndim/DIM);
376 for (j=0; j<ndim/DIM; j++) {
377 for (i=0; i<=j; i++) {
378 mat2[j][i] = 0;
379 for(d=0; d<DIM; d++)
380 mat2[j][i] += mat[ndim*(DIM*j+d)+DIM*i+d];
381 if (mat2[j][i] < min)
382 min = mat2[j][i];
383 if (mat2[j][j] > max)
384 max = mat2[j][i];
385 mat2[i][j] = mat2[j][i];
388 snew(axis,ndim/DIM);
389 for(i=0; i<ndim/DIM; i++)
390 axis[i] = i+1;
391 rlo.r = 0; rlo.g = 0; rlo.b = 1;
392 rmi.r = 1; rmi.g = 1; rmi.b = 1;
393 rhi.r = 1; rhi.g = 0; rhi.b = 0;
394 out = ffopen(xpmafile,"w");
395 nlevels = 80;
396 write_xpm3(out,0,"Covariance",bM ? "u nm^2" : "nm^2",
397 "atom","atom",ndim/DIM,ndim/DIM,axis,axis,
398 mat2,min,0.0,max,rlo,rmi,rhi,&nlevels);
399 ffclose(out);
400 sfree(axis);
401 for (i=0; i<ndim/DIM; i++)
402 sfree(mat2[i]);
403 sfree(mat2);
407 /* call diagonalization routine */
409 fprintf(stderr,"\nDiagonalizing ...\n");
410 fflush(stderr);
412 snew(eigval,ndim);
413 snew(tmp,ndim*ndim);
414 memcpy(tmp,mat,ndim*ndim*sizeof(real));
415 eigensolver(tmp,ndim,0,ndim,eigval,mat);
416 sfree(tmp);
418 /* now write the output */
420 sum=0;
421 for(i=0; i<ndim; i++)
422 sum+=eigval[i];
423 fprintf(stderr,"\nSum of the eigenvalues: %g (%snm^2)\n",
424 sum,bM ? "u " : "");
425 if (fabs(trace-sum)>0.01*trace)
426 fprintf(stderr,"\nWARNING: eigenvalue sum deviates from the trace of the covariance matrix\n");
428 fprintf(stderr,"\nWriting eigenvalues to %s\n",eigvalfile);
430 sprintf(str,"(%snm\\S2\\N)",bM ? "u " : "");
431 out=xvgropen(eigvalfile,
432 "Eigenvalues of the covariance matrix",
433 "Eigenvector index",str,oenv);
434 for (i=0; (i<ndim); i++)
435 fprintf (out,"%10d %g\n",(int)i+1,eigval[ndim-1-i]);
436 ffclose(out);
438 if (end==-1) {
439 if (nframes-1 < ndim)
440 end=nframes-1;
441 else
442 end=ndim;
444 if (bFit) {
445 /* misuse lambda: 0/1 mass weighted analysis no/yes */
446 if (nfit==natoms) {
447 WriteXref = eWXR_YES;
448 for(i=0; i<nfit; i++)
449 copy_rvec(xref[ifit[i]],x[i]);
450 } else
451 WriteXref = eWXR_NO;
452 } else {
453 /* misuse lambda: -1 for no fit */
454 WriteXref = eWXR_NOFIT;
457 write_eigenvectors(eigvecfile,natoms,mat,TRUE,1,end,
458 WriteXref,x,bDiffMass1,xproj,bM,eigval);
460 out = ffopen(logfile,"w");
462 now = time(NULL);
463 fprintf(out,"Covariance analysis log, written %s\n",
464 ctime(&now));
465 fprintf(out,"Program: %s\n",argv[0]);
466 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
467 pcwd=_getcwd(str,STRLEN);
468 #else
469 pcwd=getcwd(str,STRLEN);
470 #endif
471 if(NULL==pcwd)
473 gmx_fatal(FARGS,"Current working directory is undefined");
476 fprintf(out,"Working directory: %s\n\n",str);
478 fprintf(out,"Read %d frames from %s (time %g to %g %s)\n",nframes,trxfile,
479 output_env_conv_time(oenv,tstart),output_env_conv_time(oenv,tend),output_env_get_time_unit(oenv));
480 if (bFit)
481 fprintf(out,"Read reference structure for fit from %s\n",fitfile);
482 if (ndxfile)
483 fprintf(out,"Read index groups from %s\n",ndxfile);
484 fprintf(out,"\n");
486 fprintf(out,"Analysis group is '%s' (%d atoms)\n",ananame,natoms);
487 if (bFit)
488 fprintf(out,"Fit group is '%s' (%d atoms)\n",fitname,nfit);
489 else
490 fprintf(out,"No fit was used\n");
491 fprintf(out,"Analysis is %smass weighted\n", bDiffMass2 ? "":"non-");
492 if (bFit)
493 fprintf(out,"Fit is %smass weighted\n", bDiffMass1 ? "":"non-");
494 fprintf(out,"Diagonalized the %dx%d covariance matrix\n",(int)ndim,(int)ndim);
495 fprintf(out,"Trace of the covariance matrix before diagonalizing: %g\n",
496 trace);
497 fprintf(out,"Trace of the covariance matrix after diagonalizing: %g\n\n",
498 sum);
500 fprintf(out,"Wrote %d eigenvalues to %s\n",(int)ndim,eigvalfile);
501 if (WriteXref == eWXR_YES)
502 fprintf(out,"Wrote reference structure to %s\n",eigvecfile);
503 fprintf(out,"Wrote average structure to %s and %s\n",averfile,eigvecfile);
504 fprintf(out,"Wrote eigenvectors %d to %d to %s\n",1,end,eigvecfile);
506 ffclose(out);
508 fprintf(stderr,"Wrote the log to %s\n",logfile);
510 thanx(stderr);
512 return 0;