Fix for MSVC linking and release-build OpenMM plugin dir.
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_rdf.c
blobdddac035f107d56421fff5679526611143803ae8
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
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <math.h>
41 #include <ctype.h>
42 #include "string2.h"
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "macros.h"
46 #include "vec.h"
47 #include "pbc.h"
48 #include "xvgr.h"
49 #include "copyrite.h"
50 #include "futil.h"
51 #include "statutil.h"
52 #include "tpxio.h"
53 #include "physics.h"
54 #include "index.h"
55 #include "smalloc.h"
56 #include "fftgrid.h"
57 #include "calcgrid.h"
58 #include "nrnb.h"
59 #include "coulomb.h"
60 #include "gstat.h"
61 #include "matio.h"
62 #include "gmx_ana.h"
63 #include "names.h"
64 #include "sfactor.h"
67 static void check_box_c(matrix box)
69 if (fabs(box[ZZ][XX]) > GMX_REAL_EPS*box[ZZ][ZZ] ||
70 fabs(box[ZZ][YY]) > GMX_REAL_EPS*box[ZZ][ZZ])
71 gmx_fatal(FARGS,
72 "The last box vector is not parallel to the z-axis: %f %f %f",
73 box[ZZ][XX],box[ZZ][YY],box[ZZ][ZZ]);
76 static void calc_comg(int is,int *coi,int *index,bool bMass,t_atom *atom,
77 rvec *x,rvec *x_comg)
79 int c,i,d;
80 rvec xc;
81 real mtot,m;
83 if (bMass && atom==NULL)
84 gmx_fatal(FARGS,"No masses available while mass weighting was requested");
86 for(c=0; c<is; c++) {
87 clear_rvec(xc);
88 mtot = 0;
89 for(i=coi[c]; i<coi[c+1]; i++) {
90 if (bMass) {
91 m = atom[index[i]].m;
92 for(d=0; d<DIM; d++)
93 xc[d] += m*x[index[i]][d];
94 mtot += m;
95 } else {
96 rvec_inc(xc,x[index[i]]);
97 mtot += 1.0;
100 svmul(1/mtot,xc,x_comg[c]);
104 static void split_group(int isize,int *index,char *grpname,
105 t_topology *top,char type,
106 int *is_out,int **coi_out)
108 t_block *mols=NULL;
109 t_atom *atom=NULL;
110 int is,*coi;
111 int cur,mol,res,i,a,i1;
113 /* Split up the group in molecules or residues */
114 switch (type) {
115 case 'm':
116 mols = &top->mols;
117 break;
118 case 'r':
119 atom = top->atoms.atom;
120 break;
121 default:
122 gmx_fatal(FARGS,"Unknown rdf option '%s'",type);
124 snew(coi,isize+1);
125 is = 0;
126 cur = -1;
127 mol = 0;
128 for(i=0; i<isize; i++) {
129 a = index[i];
130 if (type == 'm') {
131 /* Check if the molecule number has changed */
132 i1 = mols->index[mol+1];
133 while(a >= i1) {
134 mol++;
135 i1 = mols->index[mol+1];
137 if (mol != cur) {
138 coi[is++] = i;
139 cur = mol;
141 } else if (type == 'r') {
142 /* Check if the residue index has changed */
143 res = atom[a].resind;
144 if (res != cur) {
145 coi[is++] = i;
146 cur = res;
150 coi[is] = i;
151 srenew(coi,is+1);
152 printf("Group '%s' of %d atoms consists of %d %s\n",
153 grpname,isize,is,
154 (type=='m' ? "molecules" : "residues"));
156 *is_out = is;
157 *coi_out = coi;
160 static void do_rdf(const char *fnNDX,const char *fnTPS,const char *fnTRX,
161 const char *fnRDF,const char *fnCNRDF, const char *fnHQ,
162 bool bCM,const char *close,
163 const char **rdft,bool bXY,bool bPBC,bool bNormalize,
164 real cutoff,real binwidth,real fade,int ng,
165 const output_env_t oenv)
167 FILE *fp;
168 int status;
169 char outf1[STRLEN],outf2[STRLEN];
170 char title[STRLEN],gtitle[STRLEN],refgt[30];
171 int g,natoms,i,ii,j,k,nbin,j0,j1,n,nframes;
172 int **count;
173 char **grpname;
174 int *isize,isize_cm=0,nrdf=0,max_i,isize0,isize_g;
175 atom_id **index,*index_cm=NULL;
176 #if (defined SIZEOF_LONG_LONG_INT) && (SIZEOF_LONG_LONG_INT >= 8)
177 long long int *sum;
178 #else
179 double *sum;
180 #endif
181 real t,rmax2,cut2,r,r2,r2ii,invhbinw,normfac;
182 real segvol,spherevol,prev_spherevol,**rdf;
183 rvec *x,dx,*x0=NULL,*x_i1,xi;
184 real *inv_segvol,invvol,invvol_sum,rho;
185 bool bClose,*bExcl,bTop,bNonSelfExcl;
186 matrix box,box_pbc;
187 int **npairs;
188 atom_id ix,jx,***pairs;
189 t_topology *top=NULL;
190 int ePBC=-1,ePBCrdf=-1;
191 t_block *mols=NULL;
192 t_blocka *excl;
193 t_atom *atom=NULL;
194 t_pbc pbc;
196 int *is=NULL,**coi=NULL,cur,mol,i1,res,a;
198 excl=NULL;
200 bClose = (close[0] != 'n');
202 if (fnTPS) {
203 snew(top,1);
204 bTop=read_tps_conf(fnTPS,title,top,&ePBC,&x,NULL,box,TRUE);
205 if (bTop && !(bCM || bClose))
206 /* get exclusions from topology */
207 excl = &(top->excls);
209 snew(grpname,ng+1);
210 snew(isize,ng+1);
211 snew(index,ng+1);
212 fprintf(stderr,"\nSelect a reference group and %d group%s\n",
213 ng,ng==1?"":"s");
214 if (fnTPS) {
215 get_index(&(top->atoms),fnNDX,ng+1,isize,index,grpname);
216 atom = top->atoms.atom;
217 } else {
218 rd_index(fnNDX,ng+1,isize,index,grpname);
221 if (bCM || bClose) {
222 snew(is,1);
223 snew(coi,1);
224 if (bClose) {
225 split_group(isize[0],index[0],grpname[0],top,close[0],&is[0],&coi[0]);
228 if (rdft[0][0] != 'a') {
229 /* Split up all the groups in molecules or residues */
230 srenew(is,ng+1);
231 srenew(coi,ng+1);
232 for(g=((bCM || bClose) ? 1 : 0); g<ng+1; g++) {
233 split_group(isize[g],index[g],grpname[g],top,rdft[0][0],&is[g],&coi[g]);
237 if (bCM) {
238 is[0] = 1;
239 snew(coi[0],is[0]+1);
240 coi[0][0] = 0;
241 coi[0][1] = isize[0];
242 isize0 = is[0];
243 snew(x0,isize0);
244 } else if (bClose || rdft[0][0] != 'a') {
245 isize0 = is[0];
246 snew(x0,isize0);
247 } else {
248 isize0 = isize[0];
251 natoms=read_first_x(oenv,&status,fnTRX,&t,&x,box);
252 if ( !natoms )
253 gmx_fatal(FARGS,"Could not read coordinates from statusfile\n");
254 if (fnTPS)
255 /* check with topology */
256 if ( natoms > top->atoms.nr )
257 gmx_fatal(FARGS,"Trajectory (%d atoms) does not match topology (%d atoms)",
258 natoms,top->atoms.nr);
259 /* check with index groups */
260 for (i=0; i<ng+1; i++)
261 for (j=0; j<isize[i]; j++)
262 if ( index[i][j] >= natoms )
263 gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) larger "
264 "than number of atoms in trajectory (%d atoms)",
265 index[i][j],grpname[i],isize[i],natoms);
267 /* initialize some handy things */
268 if (ePBC == -1) {
269 ePBC = guess_ePBC(box);
271 copy_mat(box,box_pbc);
272 if (bXY) {
273 check_box_c(box);
274 switch (ePBC) {
275 case epbcXYZ:
276 case epbcXY: ePBCrdf = epbcXY; break;
277 case epbcNONE: ePBCrdf = epbcNONE; break;
278 default:
279 gmx_fatal(FARGS,"xy-rdf's are not supported for pbc type'%s'",
280 EPBC(ePBC));
281 break;
283 /* Make sure the z-height does not influence the cut-off */
284 box_pbc[ZZ][ZZ] = 2*max(box[XX][XX],box[YY][YY]);
285 } else {
286 ePBCrdf = ePBC;
288 if (bPBC)
289 rmax2 = 0.99*0.99*max_cutoff2(bXY ? epbcXY : epbcXYZ,box_pbc);
290 else
291 rmax2 = sqr(3*max(box[XX][XX],max(box[YY][YY],box[ZZ][ZZ])));
292 if (debug)
293 fprintf(debug,"rmax2 = %g\n",rmax2);
295 /* We use the double amount of bins, so we can correctly
296 * write the rdf and rdf_cn output at i*binwidth values.
298 nbin = (int)(sqrt(rmax2) * 2 / binwidth);
299 invhbinw = 2.0 / binwidth;
300 cut2 = sqr(cutoff);
302 snew(count,ng);
303 snew(pairs,ng);
304 snew(npairs,ng);
306 snew(bExcl,natoms);
307 max_i = 0;
308 for(g=0; g<ng; g++) {
309 if (isize[g+1] > max_i)
310 max_i = isize[g+1];
312 /* this is THE array */
313 snew(count[g],nbin+1);
315 /* make pairlist array for groups and exclusions */
316 snew(pairs[g],isize[0]);
317 snew(npairs[g],isize[0]);
318 for(i=0; i<isize[0]; i++) {
319 /* We can only have exclusions with atomic rdfs */
320 if (!(bCM || bClose || rdft[0][0] != 'a')) {
321 ix = index[0][i];
322 for(j=0; j < natoms; j++)
323 bExcl[j] = FALSE;
324 /* exclusions? */
325 if (excl)
326 for( j = excl->index[ix]; j < excl->index[ix+1]; j++)
327 bExcl[excl->a[j]]=TRUE;
328 k = 0;
329 snew(pairs[g][i], isize[g+1]);
330 bNonSelfExcl = FALSE;
331 for(j=0; j<isize[g+1]; j++) {
332 jx = index[g+1][j];
333 if (!bExcl[jx])
334 pairs[g][i][k++]=jx;
335 else if (ix != jx)
336 /* Check if we have exclusions other than self exclusions */
337 bNonSelfExcl = TRUE;
339 if (bNonSelfExcl) {
340 npairs[g][i]=k;
341 srenew(pairs[g][i],npairs[g][i]);
342 } else {
343 /* Save a LOT of memory and some cpu cycles */
344 npairs[g][i]=-1;
345 sfree(pairs[g][i]);
347 } else {
348 npairs[g][i]=-1;
352 sfree(bExcl);
354 snew(x_i1,max_i);
355 nframes = 0;
356 invvol_sum = 0;
357 do {
358 /* Must init pbc every step because of pressure coupling */
359 copy_mat(box,box_pbc);
360 if (bPBC) {
361 if (top != NULL)
362 rm_pbc(&top->idef,ePBC,natoms,box,x,x);
363 if (bXY) {
364 check_box_c(box);
365 clear_rvec(box_pbc[ZZ]);
367 set_pbc(&pbc,ePBCrdf,box_pbc);
369 if (bXY)
370 /* Set z-size to 1 so we get the surface iso the volume */
371 box_pbc[ZZ][ZZ] = 1;
373 invvol = 1/det(box_pbc);
374 invvol_sum += invvol;
376 if (bCM) {
377 /* Calculate center of mass of the whole group */
378 calc_comg(is[0],coi[0],index[0],TRUE ,atom,x,x0);
379 } else if (!bClose && rdft[0][0] != 'a') {
380 calc_comg(is[0],coi[0],index[0],rdft[0][6]=='m',atom,x,x0);
383 for(g=0; g<ng; g++) {
384 if (rdft[0][0] == 'a') {
385 /* Copy the indexed coordinates to a continuous array */
386 for(i=0; i<isize[g+1]; i++)
387 copy_rvec(x[index[g+1][i]],x_i1[i]);
388 } else {
389 /* Calculate the COMs/COGs and store in x_i1 */
390 calc_comg(is[g+1],coi[g+1],index[g+1],rdft[0][6]=='m',atom,x,x_i1);
393 for(i=0; i<isize0; i++) {
394 if (bClose) {
395 /* Special loop, since we need to determine the minimum distance
396 * over all selected atoms in the reference molecule/residue.
398 if (rdft[0][0] == 'a')
399 isize_g = isize[g+1];
400 else
401 isize_g = is[g+1];
402 for(j=0; j<isize_g; j++) {
403 r2 = 1e30;
404 /* Loop over the selected atoms in the reference molecule */
405 for(ii=coi[0][i]; ii<coi[0][i+1]; ii++) {
406 if (bPBC)
407 pbc_dx(&pbc,x[index[0][ii]],x_i1[j],dx);
408 else
409 rvec_sub(x[index[0][ii]],x_i1[j],dx);
410 if (bXY)
411 r2ii = dx[XX]*dx[XX] + dx[YY]*dx[YY];
412 else
413 r2ii = iprod(dx,dx);
414 if (r2ii < r2)
415 r2 = r2ii;
417 if (r2>cut2 && r2<=rmax2)
418 count[g][(int)(sqrt(r2)*invhbinw)]++;
420 } else {
421 /* Real rdf between points in space */
422 if (bCM || rdft[0][0] != 'a') {
423 copy_rvec(x0[i],xi);
424 } else {
425 copy_rvec(x[index[0][i]],xi);
427 if (rdft[0][0] == 'a' && npairs[g][i] >= 0) {
428 /* Expensive loop, because of indexing */
429 for(j=0; j<npairs[g][i]; j++) {
430 jx=pairs[g][i][j];
431 if (bPBC)
432 pbc_dx(&pbc,xi,x[jx],dx);
433 else
434 rvec_sub(xi,x[jx],dx);
436 if (bXY)
437 r2 = dx[XX]*dx[XX] + dx[YY]*dx[YY];
438 else
439 r2=iprod(dx,dx);
440 if (r2>cut2 && r2<=rmax2)
441 count[g][(int)(sqrt(r2)*invhbinw)]++;
443 } else {
444 /* Cheaper loop, no exclusions */
445 if (rdft[0][0] == 'a')
446 isize_g = isize[g+1];
447 else
448 isize_g = is[g+1];
449 for(j=0; j<isize_g; j++) {
450 if (bPBC)
451 pbc_dx(&pbc,xi,x_i1[j],dx);
452 else
453 rvec_sub(xi,x_i1[j],dx);
454 if (bXY)
455 r2 = dx[XX]*dx[XX] + dx[YY]*dx[YY];
456 else
457 r2=iprod(dx,dx);
458 if (r2>cut2 && r2<=rmax2)
459 count[g][(int)(sqrt(r2)*invhbinw)]++;
465 nframes++;
466 } while (read_next_x(oenv,status,&t,natoms,x,box));
467 fprintf(stderr,"\n");
469 close_trj(status);
471 sfree(x);
473 /* Average volume */
474 invvol = invvol_sum/nframes;
476 /* Calculate volume of sphere segments or length of circle segments */
477 snew(inv_segvol,(nbin+1)/2);
478 prev_spherevol=0;
479 for(i=0; (i<(nbin+1)/2); i++) {
480 r = (i + 0.5)*binwidth;
481 if (bXY) {
482 spherevol=M_PI*r*r;
483 } else {
484 spherevol=(4.0/3.0)*M_PI*r*r*r;
486 segvol=spherevol-prev_spherevol;
487 inv_segvol[i]=1.0/segvol;
488 prev_spherevol=spherevol;
491 snew(rdf,ng);
492 for(g=0; g<ng; g++) {
493 /* We have to normalize by dividing by the number of frames */
494 if (rdft[0][0] == 'a')
495 normfac = 1.0/(nframes*invvol*isize0*isize[g+1]);
496 else
497 normfac = 1.0/(nframes*invvol*isize0*is[g+1]);
499 /* Do the normalization */
500 nrdf = max((nbin+1)/2,1+2*fade/binwidth);
501 snew(rdf[g],nrdf);
502 for(i=0; i<(nbin+1)/2; i++) {
503 r = i*binwidth;
504 if (i == 0)
505 j = count[g][0];
506 else
507 j = count[g][i*2-1] + count[g][i*2];
508 if ((fade > 0) && (r >= fade))
509 rdf[g][i] = 1 + (j*inv_segvol[i]*normfac-1)*exp(-16*sqr(r/fade-1));
510 else {
511 if (bNormalize)
512 rdf[g][i] = j*inv_segvol[i]*normfac;
513 else
514 rdf[g][i] = j/(binwidth*isize0*nframes);
517 for( ; (i<nrdf); i++)
518 rdf[g][i] = 1.0;
521 if (rdft[0][0] == 'a') {
522 sprintf(gtitle,"Radial distribution");
523 } else {
524 sprintf(gtitle,"Radial distribution of %s %s",
525 rdft[0][0]=='m' ? "molecule" : "residue",
526 rdft[0][6]=='m' ? "COM" : "COG");
528 fp=xvgropen(fnRDF,gtitle,"r","",oenv);
529 if (bCM) {
530 sprintf(refgt," %s","COM");
531 } else if (bClose) {
532 sprintf(refgt," closest atom in %s.",close);
533 } else {
534 sprintf(refgt,"%s","");
536 if (ng==1) {
537 if (output_env_get_print_xvgr_codes(oenv))
538 fprintf(fp,"@ subtitle \"%s%s - %s\"\n",grpname[0],refgt,grpname[1]);
540 else {
541 if (output_env_get_print_xvgr_codes(oenv))
542 fprintf(fp,"@ subtitle \"reference %s%s\"\n",grpname[0],refgt);
543 xvgr_legend(fp,ng,grpname+1,oenv);
545 for(i=0; (i<nrdf); i++) {
546 fprintf(fp,"%10g",i*binwidth);
547 for(g=0; g<ng; g++)
548 fprintf(fp," %10g",rdf[g][i]);
549 fprintf(fp,"\n");
551 ffclose(fp);
553 do_view(oenv,fnRDF,NULL);
555 /* h(Q) function: fourier transform of rdf */
556 if (fnHQ) {
557 int nhq = 401;
558 real *hq,*integrand,Q;
560 /* Get a better number density later! */
561 rho = isize[1]*invvol;
562 snew(hq,nhq);
563 snew(integrand,nrdf);
564 for(i=0; (i<nhq); i++) {
565 Q = i*0.5;
566 integrand[0] = 0;
567 for(j=1; (j<nrdf); j++) {
568 r = j*binwidth;
569 integrand[j] = (Q == 0) ? 1.0 : sin(Q*r)/(Q*r);
570 integrand[j] *= 4.0*M_PI*rho*r*r*(rdf[0][j]-1.0);
572 hq[i] = print_and_integrate(debug,nrdf,binwidth,integrand,NULL,0);
574 fp=xvgropen(fnHQ,"h(Q)","Q(/nm)","h(Q)",oenv);
575 for(i=0; (i<nhq); i++)
576 fprintf(fp,"%10g %10g\n",i*0.5,hq[i]);
577 ffclose(fp);
578 do_view(oenv,fnHQ,NULL);
579 sfree(hq);
580 sfree(integrand);
583 if (fnCNRDF) {
584 normfac = 1.0/(isize0*nframes);
585 fp=xvgropen(fnCNRDF,"Cumulative Number RDF","r","number",oenv);
586 if (ng==1) {
587 if (output_env_get_print_xvgr_codes(oenv))
588 fprintf(fp,"@ subtitle \"%s-%s\"\n",grpname[0],grpname[1]);
590 else {
591 if (output_env_get_print_xvgr_codes(oenv))
592 fprintf(fp,"@ subtitle \"reference %s\"\n",grpname[0]);
593 xvgr_legend(fp,ng,grpname+1,oenv);
595 snew(sum,ng);
596 for(i=0; (i<=nbin/2); i++) {
597 fprintf(fp,"%10g",i*binwidth);
598 for(g=0; g<ng; g++) {
599 fprintf(fp," %10g",(real)((double)sum[g]*normfac));
600 if (i*2+1 < nbin)
601 sum[g] += count[g][i*2] + count[g][i*2+1];
603 fprintf(fp,"\n");
605 ffclose(fp);
606 sfree(sum);
608 do_view(oenv,fnCNRDF,NULL);
611 for(g=0; g<ng; g++)
612 sfree(rdf[g]);
613 sfree(rdf);
617 int gmx_rdf(int argc,char *argv[])
619 const char *desc[] = {
620 "The structure of liquids can be studied by either neutron or X-ray",
621 "scattering. The most common way to describe liquid structure is by a",
622 "radial distribution function. However, this is not easy to obtain from",
623 "a scattering experiment.[PAR]",
624 "g_rdf calculates radial distribution functions in different ways.",
625 "The normal method is around a (set of) particle(s), the other methods",
626 "are around the center of mass of a set of particles ([TT]-com[tt])",
627 "or to the closest particle in a set ([TT]-surf[tt]).",
628 "With all methods rdf's can also be calculated around axes parallel",
629 "to the z-axis with option [TT]-xy[tt].",
630 "With option [TT]-surf[tt] normalization can not be used.[PAR]"
631 "The option [TT]-rdf[tt] sets the type of rdf to be computed.",
632 "Default is for atoms or particles, but one can also select center",
633 "of mass or geometry of molecules or residues. In all cases only",
634 "the atoms in the index groups are taken into account.",
635 "For molecules and/or the center of mass option a run input file",
636 "is required.",
637 "Other weighting than COM or COG can currently only be achieved",
638 "by providing a run input file with different masses.",
639 "Options [TT]-com[tt] and [TT]-surf[tt] also work in conjunction",
640 "with [TT]-rdf[tt].[PAR]",
641 "If a run input file is supplied ([TT]-s[tt]) and [TT]-rdf[tt] is set",
642 "to [TT]atom[tt], exclusions defined",
643 "in that file are taken into account when calculating the rdf.",
644 "The option [TT]-cut[tt] is meant as an alternative way to avoid",
645 "intramolecular peaks in the rdf plot.",
646 "It is however better to supply a run input file with a higher number of",
647 "exclusions. For eg. benzene a topology with nrexcl set to 5",
648 "would eliminate all intramolecular contributions to the rdf.",
649 "Note that all atoms in the selected groups are used, also the ones",
650 "that don't have Lennard-Jones interactions.[PAR]",
651 "Option [TT]-cn[tt] produces the cumulative number rdf,",
652 "i.e. the average number of particles within a distance r.[PAR]",
653 "To bridge the gap between theory and experiment structure factors can",
654 "be computed (option [TT]-sq[tt]). The algorithm uses FFT, the grid"
655 "spacing of which is determined by option [TT]-grid[tt]."
657 static bool bCM=FALSE,bXY=FALSE,bPBC=TRUE,bNormalize=TRUE;
658 static real cutoff=0,binwidth=0.002,grid=0.05,fade=0.0,lambda=0.1,distance=10;
659 static int npixel=256,nlevel=20,ngroups=1;
660 static real start_q=0.0, end_q=60.0, energy=12.0;
662 static const char *closet[]= { NULL, "no", "mol", "res", NULL };
663 static const char *rdft[]={ NULL, "atom", "mol_com", "mol_cog", "res_com", "res_cog", NULL };
665 t_pargs pa[] = {
666 { "-bin", FALSE, etREAL, {&binwidth},
667 "Binwidth (nm)" },
668 { "-com", FALSE, etBOOL, {&bCM},
669 "RDF with respect to the center of mass of first group" },
670 { "-surf", FALSE, etENUM, {closet},
671 "RDF with respect to the surface of the first group" },
672 { "-rdf", FALSE, etENUM, {rdft},
673 "RDF type" },
674 { "-pbc", FALSE, etBOOL, {&bPBC},
675 "Use periodic boundary conditions for computing distances. Without PBC the maximum range will be three times the largest box edge." },
676 { "-norm", FALSE, etBOOL, {&bNormalize},
677 "Normalize for volume and density" },
678 { "-xy", FALSE, etBOOL, {&bXY},
679 "Use only the x and y components of the distance" },
680 { "-cut", FALSE, etREAL, {&cutoff},
681 "Shortest distance (nm) to be considered"},
682 { "-ng", FALSE, etINT, {&ngroups},
683 "Number of secondary groups to compute RDFs around a central group" },
684 { "-fade", FALSE, etREAL, {&fade},
685 "From this distance onwards the RDF is tranformed by g'(r) = 1 + [g(r)-1] exp(-(r/fade-1)^2 to make it go to 1 smoothly. If fade is 0.0 nothing is done." },
686 { "-grid", FALSE, etREAL, {&grid},
687 "[HIDDEN]Grid spacing (in nm) for FFTs when computing structure factors" },
688 { "-npixel", FALSE, etINT, {&npixel},
689 "[HIDDEN]# pixels per edge of the square detector plate" },
690 { "-nlevel", FALSE, etINT, {&nlevel},
691 "Number of different colors in the diffraction image" },
692 { "-distance", FALSE, etREAL, {&distance},
693 "[HIDDEN]Distance (in cm) from the sample to the detector" },
694 { "-wave", FALSE, etREAL, {&lambda},
695 "[HIDDEN]Wavelength for X-rays/Neutrons for scattering. 0.1 nm corresponds to roughly 12 keV" },
697 {"-startq", FALSE, etREAL, {&start_q},
698 "Starting q (1/nm) "},
699 {"-endq", FALSE, etREAL, {&end_q},
700 "Ending q (1/nm)"},
701 {"-energy", FALSE, etREAL, {&energy},
702 "Energy of the incoming X-ray (keV) "}
704 #define NPA asize(pa)
705 const char *fnTPS,*fnNDX,*fnDAT;
706 bool bSQ,bRDF;
707 output_env_t oenv;
709 t_filenm fnm[] = {
710 { efTRX, "-f", NULL, ffREAD },
711 { efTPS, NULL, NULL, ffOPTRD },
712 { efNDX, NULL, NULL, ffOPTRD },
713 { efDAT, "-d", "sfactor", ffOPTRD },
714 { efXVG, "-o", "rdf", ffOPTWR },
715 { efXVG, "-sq", "sq", ffOPTWR },
716 { efXVG, "-cn", "rdf_cn", ffOPTWR },
717 { efXVG, "-hq", "hq", ffOPTWR },
718 /* { efXPM, "-image", "sq", ffOPTWR }*/
720 #define NFILE asize(fnm)
722 CopyRight(stderr,argv[0]);
723 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
724 NFILE,fnm,NPA,pa,asize(desc),desc,0,NULL,&oenv);
726 bSQ = opt2bSet("-sq",NFILE,fnm);
727 bRDF = opt2bSet("-o",NFILE,fnm) || !bSQ;
728 if (bSQ || bCM || closet[0][0]!='n' || rdft[0][0]=='m' || rdft[0][6]=='m') {
729 fnTPS = ftp2fn(efTPS,NFILE,fnm);
730 fnDAT = ftp2fn(efDAT,NFILE,fnm);
731 } else {
732 fnTPS = ftp2fn_null(efTPS,NFILE,fnm);
734 fnNDX = ftp2fn_null(efNDX,NFILE,fnm);
736 if (closet[0][0] != 'n') {
737 if (bCM) {
738 gmx_fatal(FARGS,"Can not have both -com and -surf");
740 if (bNormalize) {
741 fprintf(stderr,"Turning of normalization because of option -surf\n");
742 bNormalize = FALSE;
746 if (!bSQ && (!fnTPS && !fnNDX))
747 gmx_fatal(FARGS,"Neither index file nor topology file specified\n"
748 "Nothing to do!");
750 if (bSQ)
751 do_scattering_intensity(fnTPS,fnNDX,opt2fn("-sq",NFILE,fnm),
752 ftp2fn(efTRX,NFILE,fnm),fnDAT,
753 start_q, end_q, energy, ngroups,oenv);
755 if (bRDF)
756 do_rdf(fnNDX,fnTPS,ftp2fn(efTRX,NFILE,fnm),
757 opt2fn("-o",NFILE,fnm),opt2fn_null("-cn",NFILE,fnm),
758 opt2fn_null("-hq",NFILE,fnm),
759 bCM,closet[0],rdft,bXY,bPBC,bNormalize,cutoff,binwidth,fade,ngroups,
760 oenv);
762 thanx(stderr);
764 return 0;