Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_spatial.c
blob1ec61418939ab3c7998ee2b1b94653669e1f3e97
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.0
11 * Copyright (c) 1991-2001
12 * BIOSON Research Institute, Dept. of Biophysical Chemistry
13 * University of Groningen, The Netherlands
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 * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
32 * And Hey:
33 * Gyas ROwers Mature At Cryogenic Speed
37 #include "statutil.h"
38 #include "typedefs.h"
39 #include "smalloc.h"
40 #include "vec.h"
41 #include "copyrite.h"
42 #include "statutil.h"
43 #include "tpxio.h"
44 #include "math.h"
45 #include "index.h"
46 #include "pbc.h"
47 #include "rmpbc.h"
48 #include "gmx_ana.h"
51 static const double bohr=0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
53 static void mequit(void){
54 printf("Memory allocation error\n");
55 exit(1);
58 int gmx_spatial(int argc,char *argv[])
60 const char *desc[] = {
61 "g_spatial calculates the spatial distribution function and ",
62 "outputs it in a form that can be read by VMD as Gaussian98 cube format. ",
63 "This was developed from template.c (gromacs-3.3). ",
64 "For a system of 32K atoms and a 50ns trajectory, the SDF can be generated ",
65 "in about 30 minutes, with most of the time dedicated to the two runs through ",
66 "trjconv that are required to center everything properly. ",
67 "This also takes a whole bunch of space (3 copies of the xtc file). ",
68 "Still, the pictures are pretty and very informative when the fitted selection is properly made. ",
69 "3-4 atoms in a widely mobile group like a free amino acid in solution works ",
70 "well, or select the protein backbone in a stable folded structure to get the SDF ",
71 "of solvent and look at the time-averaged solvation shell. ",
72 "It is also possible using this program to generate the SDF based on some arbitrarty ",
73 "Cartesian coordinate. To do that, simply omit the preliminary trjconv steps. \n",
74 "USAGE: \n",
75 "1. Use make_ndx to create a group containing the atoms around which you want the SDF \n",
76 "2. trjconv -s a.tpr -f a.xtc -o b.xtc -center tric -ur compact -pbc none \n",
77 "3. trjconv -s a.tpr -f b.xtc -o c.xtc -fit rot+trans \n",
78 "4. run g_spatial on the xtc output of step #3. \n",
79 "5. Load grid.cube into VMD and view as an isosurface. \n",
80 "*** Systems such as micelles will require trjconv -pbc cluster between steps 1 and 2\n",
81 "WARNINGS: \n",
82 "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy. ",
83 "However, the preparatory -fit rot+trans option to trjconv implies that your system will be rotating ",
84 "and translating in space (in order that the selected group does not). Therefore the values that are ",
85 "returned will only be valid for some region around your central group/coordinate that has full overlap ",
86 "with system volume throughout the entire translated/rotated system over the course of the trajectory. ",
87 "It is up to the user to ensure that this is the case. \n",
88 "BUGS: \n",
89 "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected ",
90 "and the program is halted prior to the fault while displaying a warning message suggesting the use of the -nab ",
91 "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again ",
92 "with an increased -nab value. \n",
93 "RISKY OPTIONS: \n",
94 "To reduce the amount of space and time required, you can output only the coords ",
95 "that are going to be used in the first and subsequent run through trjconv. ",
96 "However, be sure to set the -nab option to a sufficiently high value since ",
97 "memory is allocated for cube bins based on the initial coords and the -nab ",
98 "(Number of Additional Bins) option value. \n"
101 static bool bPBC=FALSE;
102 static bool bSHIFT=FALSE;
103 static int iIGNOREOUTER=-1; /*Positive values may help if the surface is spikey */
104 static bool bCUTDOWN=TRUE;
105 static real rBINWIDTH=0.05; /* nm */
106 static bool bCALCDIV=TRUE;
107 static int iNAB=4;
109 t_pargs pa[] = {
110 { "-pbc", FALSE, etBOOL, {&bPBC},
111 "Use periodic boundary conditions for computing distances" },
112 { "-div", FALSE, etBOOL, {&bCALCDIV},
113 "Calculate and apply the divisor for bin occupancies based on atoms/minimal cube size. Set as TRUE for visualization and as FALSE (-nodiv) to get accurate counts per frame" },
114 { "-ign", FALSE, etINT, {&iIGNOREOUTER},
115 "Do not display this number of outer cubes (positive values may reduce boundary speckles; -1 ensures outer surface is visible)" },
116 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
117 /* "Display a total cube that is of minimal size" }, */
118 { "-bin", FALSE, etREAL, {&rBINWIDTH},
119 "Width of the bins in nm" },
120 { "-nab", FALSE, etINT, {&iNAB},
121 "Number of additional bins to ensure proper memory allocation" }
124 double MINBIN[3];
125 double MAXBIN[3];
126 t_topology top;
127 int ePBC;
128 char title[STRLEN];
129 t_trxframe fr;
130 rvec *xtop,*shx[26];
131 matrix box,box_pbc;
132 int status;
133 int flags = TRX_READ_X;
134 t_pbc pbc;
135 t_atoms *atoms;
136 int natoms;
137 char *grpnm,*grpnmp;
138 atom_id *index,*indexp;
139 int i,nidx,nidxp;
140 int v;
141 int j,k;
142 long ***bin=(long ***)NULL;
143 long nbin[3];
144 FILE *flp;
145 long x,y,z,minx,miny,minz,maxx,maxy,maxz;
146 long numfr, numcu;
147 long tot,max,min;
148 double norm;
149 output_env_t oenv;
151 t_filenm fnm[] = {
152 { efTPS, NULL, NULL, ffREAD }, /* this is for the topology */
153 { efTRX, "-f", NULL, ffREAD }, /* and this for the trajectory */
154 { efNDX, NULL, NULL, ffOPTRD }
157 #define NFILE asize(fnm)
159 CopyRight(stderr,argv[0]);
161 /* This is the routine responsible for adding default options,
162 * calling the X/motif interface, etc. */
163 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW,
164 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
166 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&xtop,NULL,box,TRUE);
167 sfree(xtop);
169 atoms=&(top.atoms);
170 printf("Select group to generate SDF:\n");
171 get_index(atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&nidx,&index,&grpnm);
172 printf("Select group to output coords (e.g. solute):\n");
173 get_index(atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&nidxp,&indexp,&grpnmp);
175 /* The first time we read data is a little special */
176 natoms=read_first_frame(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&fr,flags);
178 /* Memory Allocation */
179 MINBIN[XX]=MAXBIN[XX]=fr.x[0][XX];
180 MINBIN[YY]=MAXBIN[YY]=fr.x[0][YY];
181 MINBIN[ZZ]=MAXBIN[ZZ]=fr.x[0][ZZ];
182 for(i=1; i<top.atoms.nr; ++i) {
183 if(fr.x[i][XX]<MINBIN[XX])MINBIN[XX]=fr.x[i][XX];
184 if(fr.x[i][XX]>MAXBIN[XX])MAXBIN[XX]=fr.x[i][XX];
185 if(fr.x[i][YY]<MINBIN[YY])MINBIN[YY]=fr.x[i][YY];
186 if(fr.x[i][YY]>MAXBIN[YY])MAXBIN[YY]=fr.x[i][YY];
187 if(fr.x[i][ZZ]<MINBIN[ZZ])MINBIN[ZZ]=fr.x[i][ZZ];
188 if(fr.x[i][ZZ]>MAXBIN[ZZ])MAXBIN[ZZ]=fr.x[i][ZZ];
190 for (i=ZZ; i>=XX; --i){
191 MAXBIN[i]=(ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+(double)iNAB)*rBINWIDTH+MINBIN[i];
192 MINBIN[i]-=(double)iNAB*rBINWIDTH;
193 nbin[i]=(long)ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH);
195 bin=(long ***)malloc(nbin[XX]*sizeof(long **));
196 if(!bin)mequit();
197 for(i=0; i<nbin[XX]; ++i){
198 bin[i]=(long **)malloc(nbin[YY]*sizeof(long *));
199 if(!bin[i])mequit();
200 for(j=0; j<nbin[YY]; ++j){
201 bin[i][j]=(long *)calloc(nbin[ZZ],sizeof(long));
202 if(!bin[i][j])mequit();
205 copy_mat(box,box_pbc);
206 numfr=0;
207 minx=miny=minz=999;
208 maxx=maxy=maxz=0;
209 /* This is the main loop over frames */
210 do {
211 /* Must init pbc every step because of pressure coupling */
213 copy_mat(box,box_pbc);
214 if (bPBC) {
215 rm_pbc(&top.idef,ePBC,natoms,box,fr.x,fr.x);
216 set_pbc(&pbc,ePBC,box_pbc);
219 for(i=0; i<nidx; i++) {
220 if(fr.x[index[i]][XX]<MINBIN[XX]||fr.x[index[i]][XX]>MAXBIN[XX]||
221 fr.x[index[i]][YY]<MINBIN[YY]||fr.x[index[i]][YY]>MAXBIN[YY]||
222 fr.x[index[i]][ZZ]<MINBIN[ZZ]||fr.x[index[i]][ZZ]>MAXBIN[ZZ])
224 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
225 printf("Memory was allocated for [%f,%f,%f]\tto\t[%f,%f,%f]\n",MINBIN[XX],MINBIN[YY],MINBIN[ZZ],MAXBIN[XX],MAXBIN[YY],MAXBIN[ZZ]);
226 printf("Memory was required for [%f,%f,%f]\n",fr.x[index[i]][XX],fr.x[index[i]][YY],fr.x[index[i]][ZZ]);
227 exit(1);
229 x=(long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
230 y=(long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
231 z=(long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
232 ++bin[x][y][z];
233 if(x<minx)minx=x;
234 if(x>maxx)maxx=x;
235 if(y<miny)miny=y;
236 if(y>maxy)maxy=y;
237 if(z<minz)minz=z;
238 if(z>maxz)maxz=z;
240 numfr++;
241 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
243 } while(read_next_frame(oenv,status,&fr));
245 if(!bCUTDOWN){
246 minx=miny=minz=0;
247 maxx=nbin[XX];
248 maxy=nbin[YY];
249 maxz=nbin[ZZ];
252 /* OUTPUT */
253 flp=ffopen("grid.cube","w");
254 fprintf(flp,"Spatial Distribution Function\n");
255 fprintf(flp,"test\n");
256 fprintf(flp,"%5d%12.6f%12.6f%12.6f\n",nidxp,(MINBIN[XX]+(minx+iIGNOREOUTER)*rBINWIDTH)*10./bohr,(MINBIN[YY]+(miny+iIGNOREOUTER)*rBINWIDTH)*10./bohr,(MINBIN[ZZ]+(minz+iIGNOREOUTER)*rBINWIDTH)*10./bohr);
257 fprintf(flp,"%5ld%12.6f%12.6f%12.6f\n",maxx-minx+1-(2*iIGNOREOUTER),rBINWIDTH*10./bohr,0.,0.);
258 fprintf(flp,"%5ld%12.6f%12.6f%12.6f\n",maxy-miny+1-(2*iIGNOREOUTER),0.,rBINWIDTH*10./bohr,0.);
259 fprintf(flp,"%5ld%12.6f%12.6f%12.6f\n",maxz-minz+1-(2*iIGNOREOUTER),0.,0.,rBINWIDTH*10./bohr);
260 for(i=0; i<nidxp; i++){
261 v=2;
262 if(*(top.atoms.atomname[indexp[i]][0])=='C')v=6;
263 if(*(top.atoms.atomname[indexp[i]][0])=='N')v=7;
264 if(*(top.atoms.atomname[indexp[i]][0])=='O')v=8;
265 if(*(top.atoms.atomname[indexp[i]][0])=='H')v=1;
266 if(*(top.atoms.atomname[indexp[i]][0])=='S')v=16;
267 fprintf(flp,"%5d%12.6f%12.6f%12.6f%12.6f\n",v,0.,(double)fr.x[indexp[i]][XX]*10./bohr,(double)fr.x[indexp[i]][YY]*10./bohr,(double)fr.x[indexp[i]][ZZ]*10./bohr);
270 tot=0;
271 for(k=0;k<nbin[XX];k++) {
272 if(!(k<minx||k>maxx))continue;
273 for(j=0;j<nbin[YY];j++) {
274 if(!(j<miny||j>maxy))continue;
275 for(i=0;i<nbin[ZZ];i++) {
276 if(!(i<minz||i>maxz))continue;
277 if(bin[k][j][i]!=0){
278 printf("A bin was not empty when it should have been empty. Programming error.\n");
279 printf("bin[%d][%d][%d] was = %ld\n",k,j,i,bin[k][j][i]);
280 exit(1);
286 min=999;
287 max=0;
288 for(k=0;k<nbin[XX];k++) {
289 if(k<minx+iIGNOREOUTER||k>maxx-iIGNOREOUTER)continue;
290 for(j=0;j<nbin[YY];j++) {
291 if(j<miny+iIGNOREOUTER||j>maxy-iIGNOREOUTER)continue;
292 for(i=0;i<nbin[ZZ];i++) {
293 if(i<minz+iIGNOREOUTER||i>maxz-iIGNOREOUTER)continue;
294 tot+=bin[k][j][i];
295 if(bin[k][j][i]>max)max=bin[k][j][i];
296 if(bin[k][j][i]<min)min=bin[k][j][i];
301 numcu=(maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
302 if(bCALCDIV){
303 norm=((double)numcu*(double)numfr) / (double)tot;
304 }else{
305 norm=1.0;
308 for(k=0;k<nbin[XX];k++) {
309 if(k<minx+iIGNOREOUTER||k>maxx-iIGNOREOUTER)continue;
310 for(j=0;j<nbin[YY];j++) {
311 if(j<miny+iIGNOREOUTER||j>maxy-iIGNOREOUTER)continue;
312 for(i=0;i<nbin[ZZ];i++) {
313 if(i<minz+iIGNOREOUTER||i>maxz-iIGNOREOUTER)continue;
314 fprintf(flp,"%12.6f ",norm*(double)bin[k][j][i]/(double)numfr);
316 fprintf(flp,"\n");
318 fprintf(flp,"\n");
320 ffclose(flp);
322 /* printf("x=%d to %d\n",minx,maxx); */
323 /* printf("y=%d to %d\n",miny,maxy); */
324 /* printf("z=%d to %d\n",minz,maxz); */
326 if(bCALCDIV){
327 printf("Counts per frame in all %ld cubes divided by %le\n",numcu,1.0/norm);
328 printf("Normalized data: average %le, min %le, max %le\n",1.0,norm*(double)min/(double)numfr,norm*(double)max/(double)numfr);
329 }else{
330 printf("grid.cube contains counts per frame in all %ld cubes\n",numcu);
331 printf("Raw data: average %le, min %le, max %le\n",1.0/norm,(double)min/(double)numfr,(double)max/(double)numfr);
334 thanx(stderr);
336 return 0;