Upped the version to 3.2.0
[gromacs.git] / src / tools / eigio.c
bloba89ede886c0cea3da9be882d7386620525911ab2
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Green Red Orange Magenta Azure Cyan Skyblue
36 #include "smalloc.h"
37 #include "vec.h"
38 #include "eigio.h"
39 #include "trnio.h"
40 #include "tpxio.h"
42 void read_eigenvectors(char *file,int *natoms,bool *bFit,
43 rvec **xref,bool *bDMR,
44 rvec **xav,bool *bDMA,
45 int *nvec, int **eignr, rvec ***eigvec)
47 t_trnheader head;
48 int status,i,snew_size;
49 rvec *x;
50 matrix box;
51 bool bOK;
53 *bDMR=FALSE;
55 /* read (reference (t=-1) and) average (t=0) structure */
56 status=open_trn(file,"r");
57 fread_trnheader(status,&head,&bOK);
58 *natoms=head.natoms;
59 snew(*xav,*natoms);
60 fread_htrn(status,&head,box,*xav,NULL,NULL);
61 if ((head.t>=-1.1) && (head.t<=-0.9)) {
62 snew(*xref,*natoms);
63 for(i=0; i<*natoms; i++)
64 copy_rvec((*xav)[i],(*xref)[i]);
65 *bDMR = (head.lambda > 0.5);
66 *bFit = (head.lambda > -0.5);
67 if (*bFit)
68 fprintf(stderr,"Read %smass weighted reference structure with %d atoms from %s\n", *bDMR ? "" : "non ",*natoms,file);
69 else {
70 fprintf(stderr,"Eigenvectors in %s were determined without fitting\n",
71 file);
72 sfree(*xref);
73 *xref=NULL;
75 fread_trnheader(status,&head,&bOK);
76 fread_htrn(status,&head,box,*xav,NULL,NULL);
78 else {
79 *bFit=TRUE;
80 *xref=NULL;
82 *bDMA = (head.lambda > 0.5);
83 if ((head.t<=-0.01) || (head.t>=0.01))
84 fprintf(stderr,"WARNING: %s does not start with t=0, which should be the "
85 "average structure. This might not be a eigenvector file. "
86 "Some things might go wrong.\n",
87 file);
88 else
89 fprintf(stderr,
90 "Read %smass weighted average/minimum structure with %d atoms from %s\n",
91 *bDMA ? "" : "non ",*natoms,file);
93 snew(x,*natoms);
94 snew_size=0;
95 *nvec=0;
96 while (fread_trnheader(status,&head,&bOK)) {
97 fread_htrn(status,&head,box,x,NULL,NULL);
98 if (*nvec >= snew_size) {
99 snew_size+=10;
100 srenew(*eignr,snew_size);
101 srenew(*eigvec,snew_size);
103 i=(int)(head.t+0.01);
104 if ((head.t-i<=-0.01) || (head.t-i>=0.01))
105 fatal_error(0,"%s contains a frame with non-integer time (%f), this "
106 "time should be an eigenvector index. "
107 "This might not be a eigenvector file.",file,head.t);
108 (*eignr)[*nvec]=i-1;
109 snew((*eigvec)[*nvec],*natoms);
110 for(i=0; i<*natoms; i++)
111 copy_rvec(x[i],(*eigvec)[*nvec][i]);
112 (*nvec)++;
114 sfree(x);
115 fprintf(stderr,"Read %d eigenvectors (dim=%d)\n\n",*nvec,*natoms*DIM);
118 void write_eigenvectors(char *trnname,int natoms,real mat[],
119 bool bReverse,int begin,int end,
120 int WriteXref,rvec *xref,bool bDMR,
121 rvec xav[],bool bDMA)
123 int trnout;
124 int ndim,i,j,d,vec;
125 matrix zerobox;
126 rvec *x;
128 ndim = natoms*DIM;
129 clear_mat(zerobox);
130 snew(x,natoms);
132 fprintf (stderr,
133 "\nWriting %saverage structure\nand eigenvectors 1 to %d to %s\n",
134 (WriteXref==eWXR_YES) ? "reference and " : "",
135 end,trnname);
137 trnout = open_tpx(trnname,"w");
138 if (WriteXref == eWXR_YES)
139 /* misuse lambda: 0/1 mass weighted fit no/yes */
140 fwrite_trn(trnout,-1,-1,bDMR ? 1.0 : 0.0,zerobox,natoms,xref,NULL,NULL);
141 else if (WriteXref == eWXR_NOFIT)
142 /* misuse lambda: -1 no fit */
143 fwrite_trn(trnout,-1,-1,-1.0,zerobox,natoms,x,NULL,NULL);
145 /* misuse lambda: 0/1 mass weighted analysis no/yes */
146 fwrite_trn(trnout,0,0,bDMA ? 1.0 : 0.0,zerobox,natoms,xav,NULL,NULL);
148 for(i=begin; i<=end; i++) {
149 if (!bReverse)
150 vec = i-1;
151 else
152 vec = ndim-i;
153 for (j=0; j<natoms; j++)
154 for(d=0; d<DIM; d++)
155 x[j][d]=mat[vec*ndim+DIM*j+d];
156 fwrite_trn(trnout,i,(real)i,0,zerobox,natoms,x,NULL,NULL);
158 close_trn(trnout);
160 sfree(x);