Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_saltbr.c
blob057b58241d27ff467c8acb40eb438944c88b678c
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.2.0
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
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38 #include <math.h>
39 #include <string.h>
41 #include "macros.h"
42 #include "vec.h"
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "filenm.h"
46 #include "statutil.h"
47 #include "copyrite.h"
48 #include "futil.h"
49 #include "gmx_fatal.h"
50 #include "smalloc.h"
51 #include "pbc.h"
52 #include "xvgr.h"
53 #include "gmx_ana.h"
56 typedef struct {
57 char *label;
58 int cg;
59 real q;
60 } t_charge;
62 static t_charge *mk_charge(t_atoms *atoms,t_block *cgs,int *nncg)
64 t_charge *cg=NULL;
65 char buf[32];
66 int i,j,ncg,resnr,anr;
67 real qq;
69 /* Find the charged groups */
70 ncg=0;
71 for(i=0; (i<cgs->nr); i++) {
72 qq=0.0;
73 for(j=cgs->index[i]; (j<cgs->index[i+1]); j++) {
74 qq+=atoms->atom[j].q;
76 if (fabs(qq) > 1.0e-5) {
77 srenew(cg,ncg+1);
78 cg[ncg].q=qq;
79 cg[ncg].cg=i;
80 anr=cgs->index[i];
81 resnr=atoms->atom[anr].resind;
82 sprintf(buf,"%s%d-%d",
83 *(atoms->resinfo[resnr].name),
84 atoms->resinfo[resnr].nr,
85 anr+1);
86 cg[ncg].label=strdup(buf);
87 ncg++;
90 *nncg=ncg;
92 for(i=0; (i<ncg); i++) {
93 printf("CG: %10s Q: %6g Atoms:",
94 cg[i].label,cg[i].q);
95 for(j=cgs->index[cg[i].cg]; (j<cgs->index[cg[i].cg+1]); j++)
96 printf(" %4u",j);
97 printf("\n");
100 return cg;
103 static real calc_dist(t_pbc *pbc,rvec x[],t_block *cgs,int icg,int jcg)
105 int i,j;
106 rvec dx;
107 real d2,mindist2=1000;
109 for(i=cgs->index[icg]; (i<cgs->index[icg+1]); i++) {
110 for(j=cgs->index[jcg]; (j<cgs->index[jcg+1]); j++) {
111 pbc_dx(pbc,x[i],x[j],dx);
112 d2 = norm2(dx);
113 if (d2 < mindist2)
114 mindist2 = d2;
117 return sqrt(mindist2);
120 int gmx_saltbr(int argc,char *argv[])
122 const char *desc[] = {
123 "g_saltbr plots the distance between all combination of charged groups",
124 "as a function of time. The groups are combined in different ways.",
125 "A minimum distance can be given, (eg. the cut-off), then groups",
126 "that are never closer than that distance will not be plotted.[BR]",
127 "Output will be in a number of fixed filenames, min-min.xvg, plus-min.xvg",
128 "and plus-plus.xvg, or files for every individual ion-pair if the [TT]-sep[tt]",
129 "option is selected. In this case files are named as [TT]sb-ResnameResnr-Atomnr[tt].",
130 "There may be many such files."
132 static bool bSep=FALSE;
133 static real truncate=1000.0;
134 t_pargs pa[] = {
135 { "-t", FALSE, etREAL, {&truncate},
136 "trunc distance" },
137 { "-sep", FALSE, etBOOL, {&bSep},
138 "Use separate files for each interaction (may be MANY)" }
140 t_filenm fnm[] = {
141 { efTRX, "-f", NULL, ffREAD },
142 { efTPX, NULL, NULL, ffREAD },
144 #define NFILE asize(fnm)
146 FILE *out[3],*fp;
147 static const char *title[3] = {
148 "Distance between positively charged groups",
149 "Distance between negatively charged groups",
150 "Distance between oppositely charged groups"
152 static const char *fn[3] = {
153 "plus-plus.xvg",
154 "min-min.xvg",
155 "plus-min.xvg"
157 int nset[3]={0,0,0};
159 t_topology *top;
160 int ePBC;
161 char *buf;
162 int status,i,j,k,m,nnn,teller,ncg,n1,n2,n3,natoms;
163 real t,*time,qi,qj;
164 t_charge *cg;
165 real ***cgdist;
166 int **nWithin;
168 double t0,dt;
169 char label[234];
170 t_pbc pbc;
171 rvec *x;
172 matrix box;
173 output_env_t oenv;
175 CopyRight(stderr,argv[0]);
177 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_BE_NICE,
178 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
180 top=read_top(ftp2fn(efTPX,NFILE,fnm),&ePBC);
181 cg=mk_charge(&top->atoms,&(top->cgs),&ncg);
182 snew(cgdist,ncg);
183 snew(nWithin,ncg);
184 for(i=0; (i<ncg); i++) {
185 snew(cgdist[i],ncg);
186 snew(nWithin[i],ncg);
189 natoms=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
191 teller=0;
192 time=NULL;
193 do {
194 srenew(time,teller+1);
195 time[teller]=t;
197 set_pbc(&pbc,ePBC,box);
199 for(i=0; (i<ncg); i++) {
200 for(j=i+1; (j<ncg); j++) {
201 srenew(cgdist[i][j],teller+1);
202 cgdist[i][j][teller]=
203 calc_dist(&pbc,x,&(top->cgs),cg[i].cg,cg[j].cg);
204 if (cgdist[i][j][teller] < truncate)
205 nWithin[i][j]=1;
209 teller++;
210 } while (read_next_x(oenv,status,&t,natoms,x,box));
211 fprintf(stderr,"\n");
212 close_trj(status);
214 if (bSep) {
215 snew(buf,256);
216 for(i=0; (i<ncg); i++)
217 for(j=i+1; (j<ncg); j++) {
218 if (nWithin[i][j]) {
219 sprintf(buf,"sb-%s:%s.xvg",cg[i].label,cg[j].label);
220 fp=xvgropen(buf,buf,"Time (ps)","Distance (nm)",oenv);
221 for(k=0; (k<teller); k++)
222 fprintf(fp,"%10g %10g\n",time[k],cgdist[i][j][k]);
223 ffclose(fp);
226 sfree(buf);
228 else {
230 for(m=0; (m<3); m++)
231 out[m]=xvgropen(fn[m],title[m],"Time (ps)","Distance (nm)",oenv);
233 snew(buf,256);
234 for(i=0; (i<ncg); i++) {
235 qi=cg[i].q;
236 for(j=i+1; (j<ncg); j++) {
237 qj=cg[j].q;
238 if (nWithin[i][j]) {
239 sprintf(buf,"%s:%s",cg[i].label,cg[j].label);
240 if (qi*qj < 0)
241 nnn=2;
242 else if (qi+qj > 0)
243 nnn=0;
244 else
245 nnn=1;
247 if (nset[nnn] == 0)
248 xvgr_legend(out[nnn],1,&buf,oenv);
249 else {
250 if (output_env_get_xvg_format(oenv) == exvgXMGR) {
251 fprintf(out[nnn],"@ legend string %d \"%s\"\n",nset[nnn],buf);
252 } else if (output_env_get_xvg_format(oenv) == exvgXMGRACE) {
253 fprintf(out[nnn],"@ s%d legend \"%s\"\n",nset[nnn],buf);
256 nset[nnn]++;
257 nWithin[i][j]=nnn+1;
261 for(k=0; (k<teller); k++) {
262 for(m=0; (m<3); m++)
263 fprintf(out[m],"%10g",time[k]);
265 for(i=0; (i<ncg); i++) {
266 for(j=i+1; (j<ncg); j++) {
267 nnn=nWithin[i][j];
268 if (nnn >0)
269 fprintf(out[nnn-1]," %10g",cgdist[i][j][k]);
272 for(m=0; (m<3); m++)
273 fprintf(out[m],"\n");
275 for(m=0; (m<3); m++) {
276 ffclose(out[m]);
277 if (nset[m] == 0)
278 remove(fn[m]);
281 thanx(stderr);
283 return 0;