3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Green Red Orange Magenta Azure Cyan Skyblue
59 int gmx_sans(int argc,char *argv[])
61 const char *desc[] = {
62 "This is simple tool to compute SANS spectra using Debye formula",
63 "It currently uses topology file (since it need to assigne element for each atom)",
66 "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]",
67 "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]",
68 "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]",
69 "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]",
70 "[TT]-startq[tt] Starting q value in nm[PAR]",
71 "[TT]-endq[tt] Ending q value in nm[PAR]",
72 "[TT]-qstep[tt] Stepping in q space[PAR]",
73 "Note: When using Debye direct method computational cost increases as",
74 "1/2 * N * (N - 1) where N is atom number in group of interest",
76 "WARNING: If sq or pr specified this tool can produce large number of files! Up to two times larger than number of frames!"
78 static gmx_bool bPBC=TRUE;
79 static gmx_bool bNORM=FALSE;
80 static real binwidth=0.2,grid=0.05; /* bins shouldnt be smaller then smallest bond (~0.1nm) length */
81 static real start_q=0.0, end_q=2.0, q_step=0.01;
82 static real mcover=-1;
83 static unsigned int seed=0;
84 static int nthreads=-1;
86 static const char *emode[]= { NULL, "direct", "mc", NULL };
87 static const char *emethod[]={ NULL, "debye", "fft", NULL };
89 gmx_neutron_atomic_structurefactors_t *gnsf;
95 { "-bin", FALSE, etREAL, {&binwidth},
96 "[HIDDEN]Binwidth (nm)" },
97 { "-mode", FALSE, etENUM, {emode},
98 "Mode for sans spectra calculation" },
99 { "-mcover", FALSE, etREAL, {&mcover},
100 "Monte-Carlo coverage should be -1(default) or (0,1]"},
101 { "-method", FALSE, etENUM, {emethod},
102 "[HIDDEN]Method for sans spectra calculation" },
103 { "-pbc", FALSE, etBOOL, {&bPBC},
104 "Use periodic boundary conditions for computing distances" },
105 { "-grid", FALSE, etREAL, {&grid},
106 "[HIDDEN]Grid spacing (in nm) for FFTs" },
107 {"-startq", FALSE, etREAL, {&start_q},
108 "Starting q (1/nm) "},
109 {"-endq", FALSE, etREAL, {&end_q},
111 { "-qstep", FALSE, etREAL, {&q_step},
112 "Stepping in q (1/nm)"},
113 { "-seed", FALSE, etINT, {&seed},
114 "Random seed for Monte-Carlo"},
116 { "-nt", FALSE, etINT, {&nthreads},
117 "Number of threads to start"},
121 const char *fnTPX,*fnNDX,*fnTRX,*fnDAT=NULL;
123 t_topology *top=NULL;
125 gmx_rmpbc_t gpbc=NULL;
127 gmx_bool bFFT=FALSE, bDEBYE=FALSE;
141 t_filenm *fnmdup=NULL;
142 gmx_radial_distribution_histogram_t *prframecurrent=NULL, *pr=NULL;
143 gmx_static_structurefactor_t *sqframecurrent=NULL, *sq=NULL;
146 #define NFILE asize(fnm)
149 { efTPX, "-s", NULL, ffREAD },
150 { efTRX, "-f", NULL, ffREAD },
151 { efNDX, NULL, NULL, ffOPTRD },
152 { efDAT, "-d", "nsfactor", ffOPTRD },
153 { efXVG, "-pr", "pr", ffWRITE },
154 { efXVG, "-sq", "sq", ffWRITE },
155 { efXVG, "-prframe", "prframe", ffOPTWR },
156 { efXVG, "-sqframe", "sqframe", ffOPTWR }
159 nthreads = gmx_omp_get_max_threads();
161 CopyRight(stderr,argv[0]);
162 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
163 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
165 /* check that binwidth not smaller than smallers distance */
166 check_binwidth(binwidth);
167 check_mcover(mcover);
169 /* setting number of omp threads globaly */
170 gmx_omp_set_num_threads(nthreads);
172 /* Now try to parse opts for modes */
173 switch(emethod[0][0]) {
176 switch(emode[0][0]) {
196 fprintf(stderr,"Using Monte Carlo Debye method to calculate spectrum\n");
198 fprintf(stderr,"Using direct Debye method to calculate spectrum\n");
201 gmx_fatal(FARGS,"FFT method not implemented!");
203 gmx_fatal(FARGS,"Unknown combination for mode and method!");
206 /* Try to read files */
207 fnDAT = ftp2fn(efDAT,NFILE,fnm);
208 fnTPX = ftp2fn(efTPX,NFILE,fnm);
209 fnTRX = ftp2fn(efTRX,NFILE,fnm);
211 gnsf = gmx_neutronstructurefactors_init(fnDAT);
212 fprintf(stderr,"Read %d atom names from %s with neutron scattering parameters\n\n",gnsf->nratoms,fnDAT);
218 bTPX=read_tps_conf(fnTPX,title,top,&ePBC,&x,NULL,box,TRUE);
220 printf("\nPlease select group for SANS spectra calculation:\n");
221 get_index(&(top->atoms),ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,grpname);
223 gsans = gmx_sans_init(top,gnsf);
225 /* Prepare reference frame */
227 gpbc = gmx_rmpbc_init(&top->idef,ePBC,top->atoms.nr,box);
228 gmx_rmpbc(gpbc,top->atoms.nr,box,x);
231 natoms=read_first_x(oenv,&status,fnTRX,&t,&x,box);
232 if (natoms != top->atoms.nr) {
233 fprintf(stderr,"\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n",natoms,top->atoms.nr);
238 gmx_rmpbc(gpbc,top->atoms.nr,box,x);
240 /* allocate memory for pr */
242 /* in case its first frame to read */
245 /* realy calc p(r) */
246 prframecurrent = calc_radial_distribution_histogram(gsans,x,box,index,isize,binwidth,bMC,bNORM,mcover,seed);
247 /* copy prframecurrent -> pr and summ up pr->gr[i] */
248 /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
249 if (pr->gr == NULL) {
250 /* check if we use pr->gr first time */
251 snew(pr->gr,prframecurrent->grn);
252 snew(pr->r,prframecurrent->grn);
254 /* resize pr->gr and pr->r if needed to preven overruns */
255 if(prframecurrent->grn > pr->grn) {
256 srenew(pr->gr,prframecurrent->grn);
257 srenew(pr->r,prframecurrent->grn);
260 pr->grn = prframecurrent->grn;
261 pr->binwidth = prframecurrent->binwidth;
262 /* summ up gr and fill r */
263 for(i=0;i<prframecurrent->grn;i++) {
264 pr->gr[i] += prframecurrent->gr[i];
265 pr->r[i] = prframecurrent->r[i];
267 /* normalize histo */
268 normalize_probability(prframecurrent->grn,prframecurrent->gr);
269 /* convert p(r) to sq */
270 sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent,start_q,end_q,q_step);
271 /* print frame data if needed */
272 if(opt2fn_null("-prframe",NFILE,fnm)) {
274 snew(suffix,GMX_PATH_MAX);
276 sprintf(hdr,"g(r), t = %f",t);
277 /* prepare output filename */
278 fnmdup = dup_tfn(NFILE,fnm);
279 sprintf(suffix,"-t%.2f",t);
280 add_suffix_to_output_names(fnmdup,NFILE,suffix);
281 fp = xvgropen(opt2fn_null("-prframe",NFILE,fnmdup),hdr,"Distance (nm)","Probability",oenv);
282 for(i=0;i<prframecurrent->grn;i++) {
283 fprintf(fp,"%10.6f%10.6f\n",prframecurrent->r[i],prframecurrent->gr[i]);
285 done_filenms(NFILE,fnmdup);
291 if(opt2fn_null("-sqframe",NFILE,fnm)) {
293 snew(suffix,GMX_PATH_MAX);
295 sprintf(hdr,"I(q), t = %f",t);
296 /* prepare output filename */
297 fnmdup = dup_tfn(NFILE,fnm);
298 sprintf(suffix,"-t%.2f",t);
299 add_suffix_to_output_names(fnmdup,NFILE,suffix);
300 fp = xvgropen(opt2fn_null("-sqframe",NFILE,fnmdup),hdr,"q (nm^-1)","s(q)/s(0)",oenv);
301 for(i=0;i<sqframecurrent->qn;i++) {
302 fprintf(fp,"%10.6f%10.6f\n",sqframecurrent->q[i],sqframecurrent->s[i]);
304 done_filenms(NFILE,fnmdup);
310 /* free pr structure */
311 sfree(prframecurrent->gr);
312 sfree(prframecurrent->r);
313 sfree(prframecurrent);
314 /* free sq structure */
315 sfree(sqframecurrent->q);
316 sfree(sqframecurrent->s);
317 sfree(sqframecurrent);
318 } while (read_next_x(oenv,status,&t,natoms,x,box));
321 /* normalize histo */
322 normalize_probability(pr->grn,pr->gr);
323 sq = convert_histogram_to_intensity_curve(pr,start_q,end_q,q_step);
325 fp = xvgropen(opt2fn_null("-pr",NFILE,fnm),"G(r)","Distance (nm)","Probability",oenv);
326 for(i=0;i<pr->grn;i++)
327 fprintf(fp,"%10.6f%10.6f\n",pr->r[i],pr->gr[i]);
331 fp = xvgropen(opt2fn_null("-sq",NFILE,fnm),"I(q)","q (nm^-1)","s(q)/s(0)",oenv);
332 for(i=0;i<sq->qn;i++) {
333 fprintf(fp,"%10.6f%10.6f\n",sq->q[i],sq->s[i]);
346 please_cite(stdout,"Garmay2012");