3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Green Red Orange Magenta Azure Cyan Skyblue
50 #include "gmx_fatal.h"
58 static real interpolate(real phi,real psi,t_shiftdata *sd)
60 int iphi,ipsi,iphi1,ipsi1;
61 real fphi,fpsi,wx0,wx1,wy0,wy1;
64 if (phi > 2*M_PI) phi -= 2*M_PI;
66 if (psi > 2*M_PI) psi -= 2*M_PI;
79 fphi -= iphi; /* Fraction (offset from gridpoint) */
88 iphi1 = (iphi+1) % sd->nx;
89 ipsi1 = (ipsi+1) % sd->ny;
91 return (sd->data[iphi] [ipsi] * wx0*wy0 +
92 sd->data[iphi1] [ipsi] * wx1*wy0 +
93 sd->data[iphi] [ipsi1] * wx0*wy1 +
94 sd->data[iphi1] [ipsi1] * wx1*wy1);
97 static void dump_sd(const char *fn,t_shiftdata *sd)
102 int nnx,nny,nfac=4,nlevels=20;
103 real phi,psi,*x_phi,*y_psi,**newdata;
105 t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
114 for(i=0; (i<nnx); i++) {
115 snew(newdata[i],nny);
116 phi = i*2*M_PI/(nnx-1);
117 x_phi[i] =phi*RAD2DEG-180;
118 for(j=0; (j<nny); j++) {
119 psi = j*2*M_PI/(nny-1);
121 y_psi[j] = psi*RAD2DEG-180;
122 /*if (((i % nfac) == 0) && ((j % nfac) == 0))
123 newdata[i][j] = sd->data[i/nfac][j/nfac];
125 newdata[i][j] = interpolate(phi,psi,sd);
126 lo = min(lo,newdata[i][j]);
127 hi = max(hi,newdata[i][j]);
130 sprintf(buf,"%s.xpm",fn);
132 write_xpm(fp,0,fn,fn,"Phi","Psi",nnx,nny,
133 x_phi,y_psi,newdata,lo,hi,rlo,rhi,&nlevels);
134 for(i=0; (i<nnx); i++)
141 static t_shiftdata *read_shifts(const char *fn)
150 if(2 != fscanf(fp,"%d%d",&nx,&ny))
152 gmx_fatal(FARGS,"Error reading from file %s",fn);
157 sd->dx = nx/(2*M_PI);
158 sd->dy = ny/(2*M_PI);
160 for(i=0; (i<=nx); i++) {
161 snew(sd->data[i],ny+1);
162 for(j=0; (j<ny); j++) {
164 sd->data[i][j] = sd->data[0][j];
166 if(1 != fscanf(fp,"%lf",&xx))
168 gmx_fatal(FARGS,"Error reading from file %s",fn);
173 sd->data[i][j] = sd->data[i][0];
184 static void done_shifts(t_shiftdata *sd)
188 for(i=0; (i<=sd->nx); i++)
194 void do_pp2shifts(FILE *fp,int nf,int nlist,t_dlist dlist[],real **dih)
196 t_shiftdata *ca_sd,*co_sd,*ha_sd,*cb_sd;
201 /* Read the shift files */
202 ca_sd = read_shifts("ca-shift.dat");
203 cb_sd = read_shifts("cb-shift.dat");
204 ha_sd = read_shifts("ha-shift.dat");
205 co_sd = read_shifts("co-shift.dat");
207 fprintf(fp,"\n *** Chemical shifts from the chemical shift index ***\n");
208 please_cite(fp,"Wishart98a");
209 fprintf(fp,"%12s %10s %10s %10s %10s\n",
210 "Residue","delta Ca","delta Ha","delta CO","delta Cb");
211 for(i=0; (i<nlist); i++) {
212 if ((has_dihedral(edPhi,&(dlist[i]))) &&
213 (has_dihedral(edPsi,&(dlist[i])))) {
214 Phi = dlist[i].j0[edPhi];
215 Psi = dlist[i].j0[edPsi];
216 ca = cb = co = ha = 0;
217 for(j=0; (j<nf); j++) {
221 ca += interpolate(phi,psi,ca_sd);
222 cb += interpolate(phi,psi,cb_sd);
223 co += interpolate(phi,psi,co_sd);
224 ha += interpolate(phi,psi,ha_sd);
226 fprintf(fp,"%12s %10g %10g %10g %10g\n",
227 dlist[i].name,ca/nf,ha/nf,co/nf,cb/nf);