Fixed a bug in the pdb-writing code.
[gromacs.git] / src / tools / eigio.c
blob02466019cd20c9de4fc94b3dd9850119adfd0ed7
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.1
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
29 * And Hey:
30 * Gromacs Runs One Microsecond At Cannonball Speeds
32 static char *SRCID_eigio_c = "$Id$";
33 #include "smalloc.h"
34 #include "vec.h"
35 #include "eigio.h"
36 #include "trnio.h"
37 #include "tpxio.h"
39 void read_eigenvectors(char *file,int *natoms,bool *bFit,
40 rvec **xref,bool *bDMR,
41 rvec **xav,bool *bDMA,
42 int *nvec, int **eignr, rvec ***eigvec)
44 t_trnheader head;
45 int status,i,snew_size;
46 rvec *x;
47 matrix box;
48 bool bOK;
50 *bDMR=FALSE;
52 /* read (reference (t=-1) and) average (t=0) structure */
53 status=open_trn(file,"r");
54 fread_trnheader(status,&head,&bOK);
55 *natoms=head.natoms;
56 snew(*xav,*natoms);
57 fread_htrn(status,&head,box,*xav,NULL,NULL);
58 if ((head.t>=-1.1) && (head.t<=-0.9)) {
59 snew(*xref,*natoms);
60 for(i=0; i<*natoms; i++)
61 copy_rvec((*xav)[i],(*xref)[i]);
62 *bDMR = (head.lambda > 0.5);
63 *bFit = (head.lambda > -0.5);
64 if (*bFit)
65 fprintf(stderr,"Read %smass weighted reference structure with %d atoms from %s\n", *bDMR ? "" : "non ",*natoms,file);
66 else {
67 fprintf(stderr,"Eigenvectors in %s were determined without fitting\n",
68 file);
69 sfree(*xref);
70 *xref=NULL;
72 fread_trnheader(status,&head,&bOK);
73 fread_htrn(status,&head,box,*xav,NULL,NULL);
75 else {
76 *bFit=TRUE;
77 *xref=NULL;
79 *bDMA = (head.lambda > 0.5);
80 if ((head.t<=-0.01) || (head.t>=0.01))
81 fprintf(stderr,"WARNING: %s does not start with t=0, which should be the "
82 "average structure. This might not be a eigenvector file. "
83 "Some things might go wrong.\n",
84 file);
85 else
86 fprintf(stderr,
87 "Read %smass weighted average/minimum structure with %d atoms from %s\n",
88 *bDMA ? "" : "non ",*natoms,file);
90 snew(x,*natoms);
91 snew_size=0;
92 *nvec=0;
93 while (fread_trnheader(status,&head,&bOK)) {
94 fread_htrn(status,&head,box,x,NULL,NULL);
95 if (*nvec >= snew_size) {
96 snew_size+=10;
97 srenew(*eignr,snew_size);
98 srenew(*eigvec,snew_size);
100 i=(int)(head.t+0.01);
101 if ((head.t-i<=-0.01) || (head.t-i>=0.01))
102 fatal_error(0,"%s contains a frame with non-integer time (%f), this "
103 "time should be an eigenvector index. "
104 "This might not be a eigenvector file.",file,head.t);
105 (*eignr)[*nvec]=i-1;
106 snew((*eigvec)[*nvec],*natoms);
107 for(i=0; i<*natoms; i++)
108 copy_rvec(x[i],(*eigvec)[*nvec][i]);
109 (*nvec)++;
111 sfree(x);
112 fprintf(stderr,"Read %d eigenvectors (dim=%d)\n\n",*nvec,*natoms*DIM);
115 void write_eigenvectors(char *trnname,int natoms,real mat[],
116 bool bReverse,int begin,int end,
117 int WriteXref,rvec *xref,bool bDMR,
118 rvec xav[],bool bDMA)
120 int trnout;
121 int ndim,i,j,d,vec;
122 matrix zerobox;
123 rvec *x;
125 ndim = natoms*DIM;
126 clear_mat(zerobox);
127 snew(x,natoms);
129 fprintf (stderr,
130 "\nWriting %saverage structure\nand eigenvectors 1 to %d to %s\n",
131 (WriteXref==eWXR_YES) ? "reference and " : "",
132 end,trnname);
134 trnout = open_tpx(trnname,"w");
135 if (WriteXref == eWXR_YES)
136 /* misuse lambda: 0/1 mass weighted fit no/yes */
137 fwrite_trn(trnout,-1,-1,bDMR ? 1.0 : 0.0,zerobox,natoms,xref,NULL,NULL);
138 else if (WriteXref == eWXR_NOFIT)
139 /* misuse lambda: -1 no fit */
140 fwrite_trn(trnout,-1,-1,-1.0,zerobox,natoms,x,NULL,NULL);
142 /* misuse lambda: 0/1 mass weighted analysis no/yes */
143 fwrite_trn(trnout,0,0,bDMA ? 1.0 : 0.0,zerobox,natoms,xav,NULL,NULL);
145 for(i=begin; i<=end; i++) {
146 if (!bReverse)
147 vec = i-1;
148 else
149 vec = ndim-i;
150 for (j=0; j<natoms; j++)
151 for(d=0; d<DIM; d++)
152 x[j][d]=mat[vec*ndim+DIM*j+d];
153 fwrite_trn(trnout,i,(real)i,0,zerobox,natoms,x,NULL,NULL);
155 close_trn(trnout);
157 sfree(x);