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
63 int gmx_sans(int argc,char *argv[])
65 const char *desc[] = {
66 "This is simple tool to compute SANS spectra using Debye formula",
67 "It currently uses topology file (since it need to assigne element for each atom)",
69 "[TT]-pr[tt] Computes normalized g(r) function",
71 "[TT]-sq[tt] Computes SANS intensity curve for needed q diapason",
73 "[TT]-startq[tt] Starting q value in nm",
75 "[TT]-endq[tt] Ending q value in nm",
77 "[TT]-qstep[tt] Stepping in q space",
79 "Note: When using Debye direct method computational cost increases as",
80 "1/2 * N * (N - 1) where N is atom number in group of interest"
82 static gmx_bool bPBC=TRUE;
83 static real binwidth=0.2,grid=0.05; /* bins shouldnt be smaller then bond (~0.1nm) length */
84 static real start_q=0.0, end_q=2.0, q_step=0.01;
85 static real mcover=-1;
86 static unsigned int seed=0;
87 static int nthreads=-1;
89 static const char *emode[]= { NULL, "direct", "mc", NULL };
90 static const char *emethod[]={ NULL, "debye", "fft", NULL };
92 gmx_nentron_atomic_structurefactors_t *gnsf;
98 { "-bin", FALSE, etREAL, {&binwidth},
99 "[HIDDEN]Binwidth (nm)" },
100 { "-mode", FALSE, etENUM, {emode},
101 "Mode for sans spectra calculation" },
102 { "-mcover", FALSE, etREAL, {&mcover},
103 "Monte-Carlo coverage should be -1(default) or (0,1]"},
104 { "-method", FALSE, etENUM, {emethod},
105 "[HIDDEN]Method for sans spectra calculation" },
106 { "-pbc", FALSE, etBOOL, {&bPBC},
107 "Use periodic boundary conditions for computing distances" },
108 { "-grid", FALSE, etREAL, {&grid},
109 "[HIDDEN]Grid spacing (in nm) for FFTs" },
110 {"-startq", FALSE, etREAL, {&start_q},
111 "Starting q (1/nm) "},
112 {"-endq", FALSE, etREAL, {&end_q},
114 { "-qstep", FALSE, etREAL, {&q_step},
115 "Stepping in q (1/nm)"},
116 { "-seed", FALSE, etINT, {&seed},
117 "Random seed for Monte-Carlo"},
119 { "-nt", FALSE, etINT, {&nthreads},
120 "Number of threads to start"},
124 const char *fnTPX,*fnNDX,*fnDAT=NULL;
126 t_topology *top=NULL;
128 gmx_rmpbc_t gpbc=NULL;
130 gmx_bool bFFT=FALSE, bDEBYE=FALSE;
142 gmx_radial_distribution_histogram_t *pr=NULL;
143 gmx_static_structurefator_t *sq=NULL;
146 #define NFILE asize(fnm)
149 { efTPX, "-s", NULL, ffREAD },
150 { efNDX, NULL, NULL, ffOPTRD },
151 { efDAT, "-d", "nsfactor", ffOPTRD },
152 { efXVG, "-sq", "sq", ffWRITE },
153 { efXVG, "-pr", "pr", ffWRITE }
157 nthreads = omp_get_max_threads();
160 CopyRight(stderr,argv[0]);
161 parse_common_args(&argc,argv,PCA_BE_NICE,
162 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
164 /* check that binwidth not smaller than smallers distance */
165 check_binwidth(binwidth);
166 check_mcover(mcover);
168 /* setting number of omp threads globaly */
170 omp_set_num_threads(nthreads);
172 /* Now try to parse opts for modes */
173 switch(emethod[0][0]) {
176 switch(emode[0][0]) {
194 if (!bDEBYE && !bFFT)
195 gmx_fatal(FARGS,"Unknown method. Set pr or fft!\n");
196 /* Try to read files */
197 fnDAT = ftp2fn(efDAT,NFILE,fnm);
198 fnTPX = ftp2fn(efTPX,NFILE,fnm);
200 gnsf = gmx_neutronstructurefactors_init(fnDAT);
201 fprintf(stderr,"Read %d atom names from %s with neutron scattering parameters\n\n",gnsf->nratoms,fnDAT);
207 bTPX=read_tps_conf(fnTPX,title,top,&ePBC,&x,NULL,box,TRUE);
209 printf("\nPlease select group for SANS spectra calculation:\n");
210 get_index(&(top->atoms),ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,grpname);
212 gsans = gmx_sans_init(top,gnsf);
214 /* Prepare reference frame */
216 gpbc = gmx_rmpbc_init(&top->idef,ePBC,top->atoms.nr,box);
217 gmx_rmpbc(gpbc,top->atoms.nr,box,x);
220 natoms=top->atoms.nr;
224 fprintf(stderr,"Using Monte Carlo Debye method to calculate spectrum\n");
226 fprintf(stderr,"Using direct Debye method to calculate spectrum\n");
229 gmx_fatal(FARGS,"Not implented!");
231 gmx_fatal(FARGS,"Whats this!");
234 /* realy calc p(r) */
235 pr = calc_radial_distribution_histogram(gsans,x,box,index,isize,binwidth,bMC,mcover,seed);
238 fp = xvgropen(opt2fn_null("-pr",NFILE,fnm),"G(r)","Distance (nm)","Probability",oenv);
239 for(i=0;i<pr->grn;i++)
240 fprintf(fp,"%10.6lf%10.6lf\n",pr->r[i],pr->gr[i]);
244 sq = convert_histogram_to_intensity_curve(pr,start_q,end_q,q_step);
245 fp = xvgropen(opt2fn_null("-sq",NFILE,fnm),"I(q)","q (nm^-1)","s(q)/s(0)",oenv);
246 for(i=0;i<sq->qn;i++) {
247 fprintf(fp,"%10.6lf%10.6lf\n",sq->q[i],sq->s[i]);
253 please_cite(stdout,"Garmay2012");