Fix for MSVC linking and release-build OpenMM plugin dir.
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_saltbr.c
blob1bcad58508b547ce4fcb60faf0c39f13c5bd0761
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",
83 *(atoms->resinfo[resnr].name),
84 atoms->resinfo[resnr].nr);
85 cg[ncg].label=strdup(buf);
86 ncg++;
89 *nncg=ncg;
91 for(i=0; (i<ncg); i++) {
92 printf("CG: %10s Q: %6g Atoms:",
93 cg[i].label,cg[i].q);
94 for(j=cgs->index[cg[i].cg]; (j<cgs->index[cg[i].cg+1]); j++)
95 printf(" %4u",j);
96 printf("\n");
99 return cg;
102 static real calc_dist(t_pbc *pbc,rvec x[],t_block *cgs,int icg,int jcg)
104 int i,j;
105 rvec dx;
106 real d2,mindist2=1000;
108 for(i=cgs->index[icg]; (i<cgs->index[icg+1]); i++) {
109 for(j=cgs->index[jcg]; (j<cgs->index[jcg+1]); j++) {
110 pbc_dx(pbc,x[i],x[j],dx);
111 d2 = norm2(dx);
112 if (d2 < mindist2)
113 mindist2 = d2;
116 return sqrt(mindist2);
119 int gmx_saltbr(int argc,char *argv[])
121 const char *desc[] = {
122 "g_saltbr plots the distance between all combination of charged groups",
123 "as a function of time. The groups are combined in different ways."
124 "A minimum distance can be given, (eg. the cut-off), then groups",
125 "that are never closer than that distance will not be plotted.[BR]",
126 "Output will be in a number of fixed filenames, min-min.xvg, plus-min.xvg",
127 "and plus-plus.xvg, or files for every individual ion-pair if selected"
129 static bool bSep=FALSE;
130 static real truncate=1000.0;
131 t_pargs pa[] = {
132 { "-t", FALSE, etREAL, {&truncate},
133 "trunc distance" },
134 { "-sep", FALSE, etBOOL, {&bSep},
135 "Use separate files for each interaction (may be MANY)" }
137 t_filenm fnm[] = {
138 { efTRX, "-f", NULL, ffREAD },
139 { efTPX, NULL, NULL, ffREAD },
141 #define NFILE asize(fnm)
143 FILE *out[3],*fp;
144 static const char *title[3] = {
145 "Distance between positively charged groups",
146 "Distance between negatively charged groups",
147 "Distance between oppositely charged groups"
149 static const char *fn[3] = {
150 "plus-plus.xvg",
151 "min-min.xvg",
152 "plus-min.xvg"
154 int nset[3]={0,0,0};
156 t_topology *top;
157 int ePBC;
158 char *buf;
159 int status,i,j,k,m,nnn,teller,ncg,n1,n2,n3,natoms;
160 real t,*time,qi,qj;
161 t_charge *cg;
162 real ***cgdist;
163 int **nWithin;
165 double t0,dt;
166 char label[234];
167 t_pbc pbc;
168 rvec *x;
169 matrix box;
170 output_env_t oenv;
172 CopyRight(stderr,argv[0]);
174 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_BE_NICE,
175 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
177 top=read_top(ftp2fn(efTPX,NFILE,fnm),&ePBC);
178 cg=mk_charge(&top->atoms,&(top->cgs),&ncg);
179 snew(cgdist,ncg);
180 snew(nWithin,ncg);
181 for(i=0; (i<ncg); i++) {
182 snew(cgdist[i],ncg);
183 snew(nWithin[i],ncg);
186 natoms=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
188 teller=0;
189 time=NULL;
190 do {
191 srenew(time,teller+1);
192 time[teller]=t;
194 set_pbc(&pbc,ePBC,box);
196 for(i=0; (i<ncg); i++) {
197 for(j=i+1; (j<ncg); j++) {
198 srenew(cgdist[i][j],teller+1);
199 cgdist[i][j][teller]=
200 calc_dist(&pbc,x,&(top->cgs),cg[i].cg,cg[j].cg);
201 if (cgdist[i][j][teller] < truncate)
202 nWithin[i][j]=1;
206 teller++;
207 } while (read_next_x(oenv,status,&t,natoms,x,box));
208 fprintf(stderr,"\n");
209 close_trj(status);
211 if (bSep) {
212 snew(buf,256);
213 for(i=0; (i<ncg); i++)
214 for(j=i+1; (j<ncg); j++) {
215 if (nWithin[i][j]) {
216 sprintf(buf,"sb-%s:%s.xvg",cg[i].label,cg[j].label);
217 fp=xvgropen(buf,buf,"Time (ps)","Distance (nm)",oenv);
218 for(k=0; (k<teller); k++)
219 fprintf(fp,"%10g %10g\n",time[k],cgdist[i][j][k]);
220 ffclose(fp);
223 sfree(buf);
225 else {
227 for(m=0; (m<3); m++)
228 out[m]=xvgropen(fn[m],title[m],"Time (ps)","Distance (nm)",oenv);
230 snew(buf,256);
231 for(i=0; (i<ncg); i++) {
232 qi=cg[i].q;
233 for(j=i+1; (j<ncg); j++) {
234 qj=cg[j].q;
235 if (nWithin[i][j]) {
236 sprintf(buf,"%s:%s",cg[i].label,cg[j].label);
237 if (qi*qj < 0)
238 nnn=2;
239 else if (qi+qj > 0)
240 nnn=0;
241 else
242 nnn=1;
244 if (nset[nnn] == 0)
245 xvgr_legend(out[nnn],1,&buf,oenv);
246 else {
247 if (output_env_get_xvg_format(oenv) == exvgXMGR) {
248 fprintf(out[nnn],"@ legend string %d \"%s\"\n",nset[nnn],buf);
249 } else if (output_env_get_xvg_format(oenv) == exvgXMGRACE) {
250 fprintf(out[nnn],"@ s%d legend \"%s\"\n",nset[nnn],buf);
253 nset[nnn]++;
254 nWithin[i][j]=nnn+1;
258 for(k=0; (k<teller); k++) {
259 for(m=0; (m<3); m++)
260 fprintf(out[m],"%10g",time[k]);
262 for(i=0; (i<ncg); i++) {
263 for(j=i+1; (j<ncg); j++) {
264 nnn=nWithin[i][j];
265 if (nnn >0)
266 fprintf(out[nnn-1]," %10g",cgdist[i][j][k]);
269 for(m=0; (m<3); m++)
270 fprintf(out[m],"\n");
272 for(m=0; (m<3); m++) {
273 ffclose(out[m]);
274 if (nset[m] == 0)
275 remove(fn[m]);
278 thanx(stderr);
280 return 0;