Enforced rotation: fixed torque calculation for FLEX potential when using mass-weighting
[gromacs/adressmacs.git] / src / tools / gmx_saltbr.c
blob0d59616405dbe769e9f5a0f40a7ae68447bc3c49
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 truncation/cut-off distance can be supplied with -t, ",
126 "such that only groups interacting within this distance will be plotted in the ",
127 "output.[BR]",
128 "Output will be in a number of fixed filenames, min-min.xvg, plus-min.xvg",
129 "and plus-plus.xvg, or files for every individual ion-pair if the [TT]-sep[tt]",
130 "option is selected. In this case files are named as [TT]sb-ResnameResnr-Atomnr[tt].",
131 "There may be many such files."
133 static gmx_bool bSep=FALSE;
134 static real truncate=1000.0;
135 t_pargs pa[] = {
136 { "-t", FALSE, etREAL, {&truncate},
137 "trunc distance" },
138 { "-sep", FALSE, etBOOL, {&bSep},
139 "Use separate files for each interaction (may be MANY)" }
141 t_filenm fnm[] = {
142 { efTRX, "-f", NULL, ffREAD },
143 { efTPX, NULL, NULL, ffREAD },
145 #define NFILE asize(fnm)
147 FILE *out[3],*fp;
148 static const char *title[3] = {
149 "Distance between positively charged groups",
150 "Distance between negatively charged groups",
151 "Distance between oppositely charged groups"
153 static const char *fn[3] = {
154 "plus-plus.xvg",
155 "min-min.xvg",
156 "plus-min.xvg"
158 int nset[3]={0,0,0};
160 t_topology *top;
161 int ePBC;
162 char *buf;
163 t_trxstatus *status;
164 int i,j,k,m,nnn,teller,ncg,n1,n2,n3,natoms;
165 real t,*time,qi,qj;
166 t_charge *cg;
167 real ***cgdist;
168 int **nWithin;
170 double t0,dt;
171 char label[234];
172 t_pbc pbc;
173 rvec *x;
174 matrix box;
175 output_env_t oenv;
177 CopyRight(stderr,argv[0]);
179 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_BE_NICE,
180 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
182 top=read_top(ftp2fn(efTPX,NFILE,fnm),&ePBC);
183 cg=mk_charge(&top->atoms,&(top->cgs),&ncg);
184 snew(cgdist,ncg);
185 snew(nWithin,ncg);
186 for(i=0; (i<ncg); i++) {
187 snew(cgdist[i],ncg);
188 snew(nWithin[i],ncg);
191 natoms=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
193 teller=0;
194 time=NULL;
195 do {
196 srenew(time,teller+1);
197 time[teller]=t;
199 set_pbc(&pbc,ePBC,box);
201 for(i=0; (i<ncg); i++) {
202 for(j=i+1; (j<ncg); j++) {
203 srenew(cgdist[i][j],teller+1);
204 cgdist[i][j][teller]=
205 calc_dist(&pbc,x,&(top->cgs),cg[i].cg,cg[j].cg);
206 if (cgdist[i][j][teller] < truncate)
207 nWithin[i][j]=1;
211 teller++;
212 } while (read_next_x(oenv,status,&t,natoms,x,box));
213 fprintf(stderr,"\n");
214 close_trj(status);
216 if (bSep) {
217 snew(buf,256);
218 for(i=0; (i<ncg); i++)
219 for(j=i+1; (j<ncg); j++) {
220 if (nWithin[i][j]) {
221 sprintf(buf,"sb-%s:%s.xvg",cg[i].label,cg[j].label);
222 fp=xvgropen(buf,buf,"Time (ps)","Distance (nm)",oenv);
223 for(k=0; (k<teller); k++)
224 fprintf(fp,"%10g %10g\n",time[k],cgdist[i][j][k]);
225 ffclose(fp);
228 sfree(buf);
230 else {
232 for(m=0; (m<3); m++)
233 out[m]=xvgropen(fn[m],title[m],"Time (ps)","Distance (nm)",oenv);
235 snew(buf,256);
236 for(i=0; (i<ncg); i++) {
237 qi=cg[i].q;
238 for(j=i+1; (j<ncg); j++) {
239 qj=cg[j].q;
240 if (nWithin[i][j]) {
241 sprintf(buf,"%s:%s",cg[i].label,cg[j].label);
242 if (qi*qj < 0)
243 nnn=2;
244 else if (qi+qj > 0)
245 nnn=0;
246 else
247 nnn=1;
249 if (nset[nnn] == 0)
250 xvgr_legend(out[nnn],1,(const char**)&buf,oenv);
251 else {
252 if (output_env_get_xvg_format(oenv) == exvgXMGR) {
253 fprintf(out[nnn],"@ legend string %d \"%s\"\n",nset[nnn],buf);
254 } else if (output_env_get_xvg_format(oenv) == exvgXMGRACE) {
255 fprintf(out[nnn],"@ s%d legend \"%s\"\n",nset[nnn],buf);
258 nset[nnn]++;
259 nWithin[i][j]=nnn+1;
263 for(k=0; (k<teller); k++) {
264 for(m=0; (m<3); m++)
265 fprintf(out[m],"%10g",time[k]);
267 for(i=0; (i<ncg); i++) {
268 for(j=i+1; (j<ncg); j++) {
269 nnn=nWithin[i][j];
270 if (nnn >0)
271 fprintf(out[nnn-1]," %10g",cgdist[i][j][k]);
274 for(m=0; (m<3); m++)
275 fprintf(out[m],"\n");
277 for(m=0; (m<3); m++) {
278 ffclose(out[m]);
279 if (nset[m] == 0)
280 remove(fn[m]);
283 thanx(stderr);
285 return 0;