Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / tools / pp2shift.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43 #include <math.h>
44 #include "typedefs.h"
45 #include "futil.h"
46 #include "macros.h"
47 #include "physics.h"
48 #include "smalloc.h"
49 #include "statutil.h"
50 #include "gstat.h"
51 #include "matio.h"
52 #include "copyrite.h"
53 #include "gmx_fatal.h"
54
55 typedef struct {
56   int  nx,ny;
57   real dx,dy;
58   real **data;
59 } t_shiftdata;
60
61 static real interpolate(real phi,real psi,t_shiftdata *sd)
62 {
63   int  iphi,ipsi,iphi1,ipsi1;
64   real fphi,fpsi,wx0,wx1,wy0,wy1;
65
66   /*phi  += M_PI;
67   if (phi > 2*M_PI) phi -= 2*M_PI;
68   psi  += M_PI;  
69   if (psi > 2*M_PI) psi -= 2*M_PI;
70   */
71   while(phi < 0)
72     phi += 2*M_PI;
73   while(psi < 0)
74     psi += 2*M_PI;
75   phi    = 2*M_PI-phi;
76       
77   fphi  = phi*sd->dx;
78   fpsi  = psi*sd->dy;
79   
80   iphi  = (int)fphi;
81   ipsi  = (int)fpsi;
82   fphi -= iphi; /* Fraction (offset from gridpoint) */
83   fpsi -= ipsi;
84     
85   wx0   = 1.0-fphi;
86   wx1   = fphi;
87   wy0   = 1.0-fpsi;
88   wy1   = fpsi;
89   iphi  = iphi % sd->nx;
90   ipsi  = ipsi % sd->ny;
91   iphi1 = (iphi+1) % sd->nx;
92   ipsi1 = (ipsi+1) % sd->ny;
93   
94   return (sd->data[iphi]  [ipsi]  * wx0*wy0 + 
95           sd->data[iphi1] [ipsi]  * wx1*wy0 +
96           sd->data[iphi]  [ipsi1] * wx0*wy1 +
97           sd->data[iphi1] [ipsi1] * wx1*wy1);
98 }
99
100 static void dump_sd(const char *fn,t_shiftdata *sd)
101 {
102   FILE  *fp;
103   int   i,j;
104   char  buf[256];
105   int   nnx,nny,nfac=4,nlevels=20;
106   real  phi,psi,*x_phi,*y_psi,**newdata;
107   real  lo,hi;
108   t_rgb rlo = { 1, 0, 0 }, rhi = { 0, 0, 1 };
109   
110   nnx = sd->nx*nfac+1;
111   nny = sd->ny*nfac+1;
112   snew(x_phi,nnx);
113   snew(y_psi,nny);
114   snew(newdata,nnx);
115   lo =  100000;
116   hi = -100000;
117   for(i=0; (i<nnx); i++) {
118     snew(newdata[i],nny);
119     phi = i*2*M_PI/(nnx-1);
120     x_phi[i] =phi*RAD2DEG-180;
121     for(j=0; (j<nny); j++) {
122       psi = j*2*M_PI/(nny-1);
123       if (i == 0)
124         y_psi[j] = psi*RAD2DEG-180;
125       /*if (((i % nfac) == 0) && ((j % nfac) == 0))
126         newdata[i][j] = sd->data[i/nfac][j/nfac];
127       else*/
128         newdata[i][j] = interpolate(phi,psi,sd);
129       lo = min(lo,newdata[i][j]);
130       hi = max(hi,newdata[i][j]);
131     }
132   }
133   sprintf(buf,"%s.xpm",fn);
134   fp=ffopen(buf,"w");
135   write_xpm(fp,0,fn,fn,"Phi","Psi",nnx,nny,
136             x_phi,y_psi,newdata,lo,hi,rlo,rhi,&nlevels);
137   for(i=0; (i<nnx); i++)
138     sfree(newdata[i]);
139   sfree(newdata);
140   sfree(x_phi);
141   sfree(y_psi);
142 }
143
144 static t_shiftdata *read_shifts(const char *fn)
145 {
146   FILE        *fp;
147   double      xx;
148   int         i,j,nx,ny;
149   t_shiftdata *sd;
150   
151   snew(sd,1);
152   fp=libopen(fn);
153   if(2 != fscanf(fp,"%d%d",&nx,&ny))
154   {
155     gmx_fatal(FARGS,"Error reading from file %s",fn);
156   }
157  
158   sd->nx = nx;
159   sd->ny = ny;
160   sd->dx = nx/(2*M_PI);
161   sd->dy = ny/(2*M_PI);
162   snew(sd->data,nx+1);
163   for(i=0; (i<=nx); i++) {
164     snew(sd->data[i],ny+1);
165     for(j=0; (j<ny); j++) {
166       if (i == nx)
167         sd->data[i][j] = sd->data[0][j];
168       else {
169         if(1 != fscanf(fp,"%lf",&xx))
170         {
171             gmx_fatal(FARGS,"Error reading from file %s",fn);
172         }
173         sd->data[i][j] = xx;
174       }
175     }
176     sd->data[i][j] = sd->data[i][0];
177   }
178   ffclose(fp);
179   
180   if (bDebugMode()) 
181     dump_sd(fn,sd);
182     
183   return sd;
184 }
185
186
187 static void done_shifts(t_shiftdata *sd)
188 {
189   int i;
190   
191   for(i=0; (i<=sd->nx); i++)
192     sfree(sd->data[i]);
193   sfree(sd->data);
194   sfree(sd);
195 }
196
197 void do_pp2shifts(FILE *fp,int nf,int nlist,t_dlist dlist[],real **dih)
198 {
199   t_shiftdata *ca_sd,*co_sd,*ha_sd,*cb_sd;
200   int         i,j,Phi,Psi;
201   real        phi,psi;
202   real        ca,co,ha,cb;
203     
204   /* Read the shift files */
205   ca_sd = read_shifts("ca-shift.dat");
206   cb_sd = read_shifts("cb-shift.dat");
207   ha_sd = read_shifts("ha-shift.dat");
208   co_sd = read_shifts("co-shift.dat");
209
210   fprintf(fp,"\n *** Chemical shifts from the chemical shift index ***\n");  
211   please_cite(fp,"Wishart98a");
212   fprintf(fp,"%12s  %10s  %10s  %10s  %10s\n",
213           "Residue","delta Ca","delta Ha","delta CO","delta Cb");
214   for(i=0; (i<nlist); i++) {
215     if ((has_dihedral(edPhi,&(dlist[i]))) &&
216         (has_dihedral(edPsi,&(dlist[i])))) {
217       Phi = dlist[i].j0[edPhi];   
218       Psi = dlist[i].j0[edPsi];
219       ca = cb = co = ha = 0;
220       for(j=0; (j<nf); j++) {
221         phi = dih[Phi][j];
222         psi = dih[Psi][j];
223         
224         ca += interpolate(phi,psi,ca_sd);
225         cb += interpolate(phi,psi,cb_sd);
226         co += interpolate(phi,psi,co_sd);
227         ha += interpolate(phi,psi,ha_sd);
228       }
229       fprintf(fp,"%12s  %10g  %10g  %10g  %10g\n",
230               dlist[i].name,ca/nf,ha/nf,co/nf,cb/nf);
231     }
232   }
233   fprintf(fp,"\n");
234   
235   /* Free memory */
236   done_shifts(ca_sd);
237   done_shifts(cb_sd);
238   done_shifts(co_sd);
239   done_shifts(ha_sd);
240 }
241