2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gromacs/math/vec.h"
50 #include "gromacs/commandline/pargs.h"
51 #include "gromacs/fileio/tpxio.h"
52 #include "gromacs/fileio/trxio.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/gmxomp.h"
58 #include "gromacs/utility/smalloc.h"
60 int gmx_sans(int argc, char *argv[])
62 const char *desc[] = {
63 "[THISMODULE] computes SANS spectra using Debye formula.",
64 "It currently uses topology file (since it need to assigne element for each atom).",
67 "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]",
68 "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]",
69 "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]",
70 "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]",
71 "[TT]-startq[tt] Starting q value in nm[PAR]",
72 "[TT]-endq[tt] Ending q value in nm[PAR]",
73 "[TT]-qstep[tt] Stepping in q space[PAR]",
74 "Note: When using Debye direct method computational cost increases as",
75 "1/2 * N * (N - 1) where N is atom number in group of interest.",
77 "WARNING: If sq or pr specified this tool can produce large number of files! Up to two times larger than number of frames!"
79 static gmx_bool bPBC = TRUE;
80 static gmx_bool bNORM = FALSE;
81 static real binwidth = 0.2, grid = 0.05; /* bins shouldnt be smaller then smallest bond (~0.1nm) length */
82 static real start_q = 0.0, end_q = 2.0, q_step = 0.01;
83 static real mcover = -1;
84 static unsigned int seed = 0;
85 static int nthreads = -1;
87 static const char *emode[] = { NULL, "direct", "mc", NULL };
88 static const char *emethod[] = { NULL, "debye", "fft", NULL };
90 gmx_neutron_atomic_structurefactors_t *gnsf;
96 { "-bin", FALSE, etREAL, {&binwidth},
97 "[HIDDEN]Binwidth (nm)" },
98 { "-mode", FALSE, etENUM, {emode},
99 "Mode for sans spectra calculation" },
100 { "-mcover", FALSE, etREAL, {&mcover},
101 "Monte-Carlo coverage should be -1(default) or (0,1]"},
102 { "-method", FALSE, etENUM, {emethod},
103 "[HIDDEN]Method for sans spectra calculation" },
104 { "-pbc", FALSE, etBOOL, {&bPBC},
105 "Use periodic boundary conditions for computing distances" },
106 { "-grid", FALSE, etREAL, {&grid},
107 "[HIDDEN]Grid spacing (in nm) for FFTs" },
108 {"-startq", FALSE, etREAL, {&start_q},
109 "Starting q (1/nm) "},
110 {"-endq", FALSE, etREAL, {&end_q},
112 { "-qstep", FALSE, etREAL, {&q_step},
113 "Stepping in q (1/nm)"},
114 { "-seed", FALSE, etINT, {&seed},
115 "Random seed for Monte-Carlo"},
117 { "-nt", FALSE, etINT, {&nthreads},
118 "Number of threads to start"},
122 const char *fnTPX, *fnNDX, *fnTRX, *fnDAT = NULL;
124 t_topology *top = NULL;
126 gmx_rmpbc_t gpbc = NULL;
128 gmx_bool bFFT = FALSE, bDEBYE = FALSE;
129 gmx_bool bMC = FALSE;
136 char **grpname = NULL;
137 atom_id *index = NULL;
142 t_filenm *fnmdup = NULL;
143 gmx_radial_distribution_histogram_t *prframecurrent = NULL, *pr = NULL;
144 gmx_static_structurefactor_t *sqframecurrent = NULL, *sq = NULL;
147 #define NFILE asize(fnm)
150 { efTPX, "-s", NULL, ffREAD },
151 { efTRX, "-f", NULL, ffREAD },
152 { efNDX, NULL, NULL, ffOPTRD },
153 { efDAT, "-d", "nsfactor", ffOPTRD },
154 { efXVG, "-pr", "pr", ffWRITE },
155 { efXVG, "-sq", "sq", ffWRITE },
156 { efXVG, "-prframe", "prframe", ffOPTWR },
157 { efXVG, "-sqframe", "sqframe", ffOPTWR }
160 nthreads = gmx_omp_get_max_threads();
162 if (!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))
168 /* check that binwidth not smaller than smallers distance */
169 check_binwidth(binwidth);
170 check_mcover(mcover);
172 /* setting number of omp threads globaly */
173 gmx_omp_set_num_threads(nthreads);
175 /* Now try to parse opts for modes */
176 switch (emethod[0][0])
203 fprintf(stderr, "Using Monte Carlo Debye method to calculate spectrum\n");
207 fprintf(stderr, "Using direct Debye method to calculate spectrum\n");
212 gmx_fatal(FARGS, "FFT method not implemented!");
216 gmx_fatal(FARGS, "Unknown combination for mode and method!");
219 /* Try to read files */
220 fnDAT = ftp2fn(efDAT, NFILE, fnm);
221 fnTPX = ftp2fn(efTPX, NFILE, fnm);
222 fnTRX = ftp2fn(efTRX, NFILE, fnm);
224 gnsf = gmx_neutronstructurefactors_init(fnDAT);
225 fprintf(stderr, "Read %d atom names from %s with neutron scattering parameters\n\n", gnsf->nratoms, fnDAT);
231 bTPX = read_tps_conf(fnTPX, title, top, &ePBC, &x, NULL, box, TRUE);
233 printf("\nPlease select group for SANS spectra calculation:\n");
234 get_index(&(top->atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, grpname);
236 gsans = gmx_sans_init(top, gnsf);
238 /* Prepare reference frame */
241 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
242 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
245 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
246 if (natoms != top->atoms.nr)
248 fprintf(stderr, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n", natoms, top->atoms.nr);
255 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
257 /* allocate memory for pr */
260 /* in case its first frame to read */
263 /* realy calc p(r) */
264 prframecurrent = calc_radial_distribution_histogram(gsans, x, box, index, isize, binwidth, bMC, bNORM, mcover, seed);
265 /* copy prframecurrent -> pr and summ up pr->gr[i] */
266 /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
269 /* check if we use pr->gr first time */
270 snew(pr->gr, prframecurrent->grn);
271 snew(pr->r, prframecurrent->grn);
275 /* resize pr->gr and pr->r if needed to preven overruns */
276 if (prframecurrent->grn > pr->grn)
278 srenew(pr->gr, prframecurrent->grn);
279 srenew(pr->r, prframecurrent->grn);
282 pr->grn = prframecurrent->grn;
283 pr->binwidth = prframecurrent->binwidth;
284 /* summ up gr and fill r */
285 for (i = 0; i < prframecurrent->grn; i++)
287 pr->gr[i] += prframecurrent->gr[i];
288 pr->r[i] = prframecurrent->r[i];
290 /* normalize histo */
291 normalize_probability(prframecurrent->grn, prframecurrent->gr);
292 /* convert p(r) to sq */
293 sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent, start_q, end_q, q_step);
294 /* print frame data if needed */
295 if (opt2fn_null("-prframe", NFILE, fnm))
298 snew(suffix, GMX_PATH_MAX);
300 sprintf(hdr, "g(r), t = %f", t);
301 /* prepare output filename */
302 fnmdup = dup_tfn(NFILE, fnm);
303 sprintf(suffix, "-t%.2f", t);
304 add_suffix_to_output_names(fnmdup, NFILE, suffix);
305 fp = xvgropen(opt2fn_null("-prframe", NFILE, fnmdup), hdr, "Distance (nm)", "Probability", oenv);
306 for (i = 0; i < prframecurrent->grn; i++)
308 fprintf(fp, "%10.6f%10.6f\n", prframecurrent->r[i], prframecurrent->gr[i]);
310 done_filenms(NFILE, fnmdup);
316 if (opt2fn_null("-sqframe", NFILE, fnm))
319 snew(suffix, GMX_PATH_MAX);
321 sprintf(hdr, "I(q), t = %f", t);
322 /* prepare output filename */
323 fnmdup = dup_tfn(NFILE, fnm);
324 sprintf(suffix, "-t%.2f", t);
325 add_suffix_to_output_names(fnmdup, NFILE, suffix);
326 fp = xvgropen(opt2fn_null("-sqframe", NFILE, fnmdup), hdr, "q (nm^-1)", "s(q)/s(0)", oenv);
327 for (i = 0; i < sqframecurrent->qn; i++)
329 fprintf(fp, "%10.6f%10.6f\n", sqframecurrent->q[i], sqframecurrent->s[i]);
331 done_filenms(NFILE, fnmdup);
337 /* free pr structure */
338 sfree(prframecurrent->gr);
339 sfree(prframecurrent->r);
340 sfree(prframecurrent);
341 /* free sq structure */
342 sfree(sqframecurrent->q);
343 sfree(sqframecurrent->s);
344 sfree(sqframecurrent);
346 while (read_next_x(oenv, status, &t, x, box));
349 /* normalize histo */
350 normalize_probability(pr->grn, pr->gr);
351 sq = convert_histogram_to_intensity_curve(pr, start_q, end_q, q_step);
353 fp = xvgropen(opt2fn_null("-pr", NFILE, fnm), "G(r)", "Distance (nm)", "Probability", oenv);
354 for (i = 0; i < pr->grn; i++)
356 fprintf(fp, "%10.6f%10.6f\n", pr->r[i], pr->gr[i]);
361 fp = xvgropen(opt2fn_null("-sq", NFILE, fnm), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv);
362 for (i = 0; i < sq->qn; i++)
364 fprintf(fp, "%10.6f%10.6f\n", sq->q[i], sq->s[i]);
377 please_cite(stdout, "Garmay2012");