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)",
65 "[TT]-pr[tt] Computes normalized g(r) function",
67 "[TT]-sq[tt] Computes SANS intensity curve for needed q diapason",
69 "[TT]-startq[tt] Starting q value in nm",
71 "[TT]-endq[tt] Ending q value in nm",
73 "[TT]-qstep[tt] Stepping in q space",
75 "Note: When using Debye direct method computational cost increases as",
76 "1/2 * N * (N - 1) where N is atom number in group of interest"
78 static gmx_bool bPBC=TRUE;
79 static real binwidth=0.2,grid=0.05; /* bins shouldnt be smaller then bond (~0.1nm) length */
80 static real start_q=0.0, end_q=2.0, q_step=0.01;
81 static real mcover=-1;
82 static unsigned int seed=0;
83 static int nthreads=-1;
85 static const char *emode[]= { NULL, "direct", "mc", NULL };
86 static const char *emethod[]={ NULL, "debye", "fft", NULL };
88 gmx_nentron_atomic_structurefactors_t *gnsf;
94 { "-bin", FALSE, etREAL, {&binwidth},
95 "[HIDDEN]Binwidth (nm)" },
96 { "-mode", FALSE, etENUM, {emode},
97 "Mode for sans spectra calculation" },
98 { "-mcover", FALSE, etREAL, {&mcover},
99 "Monte-Carlo coverage should be -1(default) or (0,1]"},
100 { "-method", FALSE, etENUM, {emethod},
101 "[HIDDEN]Method for sans spectra calculation" },
102 { "-pbc", FALSE, etBOOL, {&bPBC},
103 "Use periodic boundary conditions for computing distances" },
104 { "-grid", FALSE, etREAL, {&grid},
105 "[HIDDEN]Grid spacing (in nm) for FFTs" },
106 {"-startq", FALSE, etREAL, {&start_q},
107 "Starting q (1/nm) "},
108 {"-endq", FALSE, etREAL, {&end_q},
110 { "-qstep", FALSE, etREAL, {&q_step},
111 "Stepping in q (1/nm)"},
112 { "-seed", FALSE, etINT, {&seed},
113 "Random seed for Monte-Carlo"},
115 { "-nt", FALSE, etINT, {&nthreads},
116 "Number of threads to start"},
120 const char *fnTPX,*fnNDX,*fnDAT=NULL;
122 t_topology *top=NULL;
124 gmx_rmpbc_t gpbc=NULL;
126 gmx_bool bFFT=FALSE, bDEBYE=FALSE;
138 gmx_radial_distribution_histogram_t *pr=NULL;
139 gmx_static_structurefator_t *sq=NULL;
142 #define NFILE asize(fnm)
145 { efTPX, "-s", NULL, ffREAD },
146 { efNDX, NULL, NULL, ffOPTRD },
147 { efDAT, "-d", "nsfactor", ffOPTRD },
148 { efXVG, "-sq", "sq", ffWRITE },
149 { efXVG, "-pr", "pr", ffWRITE }
152 nthreads = gmx_omp_get_max_threads();
154 CopyRight(stderr,argv[0]);
155 parse_common_args(&argc,argv,PCA_BE_NICE,
156 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
158 /* check that binwidth not smaller than smallers distance */
159 check_binwidth(binwidth);
160 check_mcover(mcover);
162 /* setting number of omp threads globaly */
163 gmx_omp_set_num_threads(nthreads);
165 /* Now try to parse opts for modes */
166 switch(emethod[0][0]) {
169 switch(emode[0][0]) {
187 if (!bDEBYE && !bFFT)
188 gmx_fatal(FARGS,"Unknown method. Set pr or fft!\n");
189 /* Try to read files */
190 fnDAT = ftp2fn(efDAT,NFILE,fnm);
191 fnTPX = ftp2fn(efTPX,NFILE,fnm);
193 gnsf = gmx_neutronstructurefactors_init(fnDAT);
194 fprintf(stderr,"Read %d atom names from %s with neutron scattering parameters\n\n",gnsf->nratoms,fnDAT);
200 bTPX=read_tps_conf(fnTPX,title,top,&ePBC,&x,NULL,box,TRUE);
202 printf("\nPlease select group for SANS spectra calculation:\n");
203 get_index(&(top->atoms),ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,grpname);
205 gsans = gmx_sans_init(top,gnsf);
207 /* Prepare reference frame */
209 gpbc = gmx_rmpbc_init(&top->idef,ePBC,top->atoms.nr,box);
210 gmx_rmpbc(gpbc,top->atoms.nr,box,x);
213 natoms=top->atoms.nr;
217 fprintf(stderr,"Using Monte Carlo Debye method to calculate spectrum\n");
219 fprintf(stderr,"Using direct Debye method to calculate spectrum\n");
222 gmx_fatal(FARGS,"Not implented!");
224 gmx_fatal(FARGS,"Whats this!");
227 /* realy calc p(r) */
228 pr = calc_radial_distribution_histogram(gsans,x,box,index,isize,binwidth,bMC,mcover,seed);
231 fp = xvgropen(opt2fn_null("-pr",NFILE,fnm),"G(r)","Distance (nm)","Probability",oenv);
232 for(i=0;i<pr->grn;i++)
233 fprintf(fp,"%10.6lf%10.6lf\n",pr->r[i],pr->gr[i]);
237 sq = convert_histogram_to_intensity_curve(pr,start_q,end_q,q_step);
238 fp = xvgropen(opt2fn_null("-sq",NFILE,fnm),"I(q)","q (nm^-1)","s(q)/s(0)",oenv);
239 for(i=0;i<sq->qn;i++) {
240 fprintf(fp,"%10.6lf%10.6lf\n",sq->q[i],sq->s[i]);
246 please_cite(stdout,"Garmay2012");