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
49 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/fileio/futil.h"
56 #include "gromacs/fileio/matio.h"
57 #include "gromacs/fileio/tpxio.h"
58 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/utility/gmxomp.h"
61 int gmx_sans(int argc, char *argv[])
63 const char *desc[] = {
64 "[THISMODULE] computes SANS spectra using Debye formula.",
65 "It currently uses topology file (since it need to assigne element for each atom).",
68 "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]",
69 "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]",
70 "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]",
71 "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]",
72 "[TT]-startq[tt] Starting q value in nm[PAR]",
73 "[TT]-endq[tt] Ending q value in nm[PAR]",
74 "[TT]-qstep[tt] Stepping in q space[PAR]",
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 "WARNING: If sq or pr specified this tool can produce large number of files! Up to two times larger than number of frames!"
80 static gmx_bool bPBC = TRUE;
81 static gmx_bool bNORM = FALSE;
82 static real binwidth = 0.2, grid = 0.05; /* bins shouldnt be smaller then smallest bond (~0.1nm) length */
83 static real start_q = 0.0, end_q = 2.0, q_step = 0.01;
84 static real mcover = -1;
85 static unsigned int seed = 0;
86 static int nthreads = -1;
88 static const char *emode[] = { NULL, "direct", "mc", NULL };
89 static const char *emethod[] = { NULL, "debye", "fft", NULL };
91 gmx_neutron_atomic_structurefactors_t *gnsf;
97 { "-bin", FALSE, etREAL, {&binwidth},
98 "[HIDDEN]Binwidth (nm)" },
99 { "-mode", FALSE, etENUM, {emode},
100 "Mode for sans spectra calculation" },
101 { "-mcover", FALSE, etREAL, {&mcover},
102 "Monte-Carlo coverage should be -1(default) or (0,1]"},
103 { "-method", FALSE, etENUM, {emethod},
104 "[HIDDEN]Method for sans spectra calculation" },
105 { "-pbc", FALSE, etBOOL, {&bPBC},
106 "Use periodic boundary conditions for computing distances" },
107 { "-grid", FALSE, etREAL, {&grid},
108 "[HIDDEN]Grid spacing (in nm) for FFTs" },
109 {"-startq", FALSE, etREAL, {&start_q},
110 "Starting q (1/nm) "},
111 {"-endq", FALSE, etREAL, {&end_q},
113 { "-qstep", FALSE, etREAL, {&q_step},
114 "Stepping in q (1/nm)"},
115 { "-seed", FALSE, etINT, {&seed},
116 "Random seed for Monte-Carlo"},
118 { "-nt", FALSE, etINT, {&nthreads},
119 "Number of threads to start"},
123 const char *fnTPX, *fnNDX, *fnTRX, *fnDAT = NULL;
125 t_topology *top = NULL;
127 gmx_rmpbc_t gpbc = NULL;
129 gmx_bool bFFT = FALSE, bDEBYE = FALSE;
130 gmx_bool bMC = FALSE;
137 char **grpname = NULL;
138 atom_id *index = NULL;
143 t_filenm *fnmdup = NULL;
144 gmx_radial_distribution_histogram_t *prframecurrent = NULL, *pr = NULL;
145 gmx_static_structurefactor_t *sqframecurrent = NULL, *sq = NULL;
148 #define NFILE asize(fnm)
151 { efTPX, "-s", NULL, ffREAD },
152 { efTRX, "-f", NULL, ffREAD },
153 { efNDX, NULL, NULL, ffOPTRD },
154 { efDAT, "-d", "nsfactor", ffOPTRD },
155 { efXVG, "-pr", "pr", ffWRITE },
156 { efXVG, "-sq", "sq", ffWRITE },
157 { efXVG, "-prframe", "prframe", ffOPTWR },
158 { efXVG, "-sqframe", "sqframe", ffOPTWR }
161 nthreads = gmx_omp_get_max_threads();
163 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
164 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
169 /* check that binwidth not smaller than smallers distance */
170 check_binwidth(binwidth);
171 check_mcover(mcover);
173 /* setting number of omp threads globaly */
174 gmx_omp_set_num_threads(nthreads);
176 /* Now try to parse opts for modes */
177 switch (emethod[0][0])
204 fprintf(stderr, "Using Monte Carlo Debye method to calculate spectrum\n");
208 fprintf(stderr, "Using direct Debye method to calculate spectrum\n");
213 gmx_fatal(FARGS, "FFT method not implemented!");
217 gmx_fatal(FARGS, "Unknown combination for mode and method!");
220 /* Try to read files */
221 fnDAT = ftp2fn(efDAT, NFILE, fnm);
222 fnTPX = ftp2fn(efTPX, NFILE, fnm);
223 fnTRX = ftp2fn(efTRX, NFILE, fnm);
225 gnsf = gmx_neutronstructurefactors_init(fnDAT);
226 fprintf(stderr, "Read %d atom names from %s with neutron scattering parameters\n\n", gnsf->nratoms, fnDAT);
232 bTPX = read_tps_conf(fnTPX, title, top, &ePBC, &x, NULL, box, TRUE);
234 printf("\nPlease select group for SANS spectra calculation:\n");
235 get_index(&(top->atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, grpname);
237 gsans = gmx_sans_init(top, gnsf);
239 /* Prepare reference frame */
242 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
243 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
246 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
247 if (natoms != top->atoms.nr)
249 fprintf(stderr, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n", natoms, top->atoms.nr);
256 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
258 /* allocate memory for pr */
261 /* in case its first frame to read */
264 /* realy calc p(r) */
265 prframecurrent = calc_radial_distribution_histogram(gsans, x, box, index, isize, binwidth, bMC, bNORM, mcover, seed);
266 /* copy prframecurrent -> pr and summ up pr->gr[i] */
267 /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
270 /* check if we use pr->gr first time */
271 snew(pr->gr, prframecurrent->grn);
272 snew(pr->r, prframecurrent->grn);
276 /* resize pr->gr and pr->r if needed to preven overruns */
277 if (prframecurrent->grn > pr->grn)
279 srenew(pr->gr, prframecurrent->grn);
280 srenew(pr->r, prframecurrent->grn);
283 pr->grn = prframecurrent->grn;
284 pr->binwidth = prframecurrent->binwidth;
285 /* summ up gr and fill r */
286 for (i = 0; i < prframecurrent->grn; i++)
288 pr->gr[i] += prframecurrent->gr[i];
289 pr->r[i] = prframecurrent->r[i];
291 /* normalize histo */
292 normalize_probability(prframecurrent->grn, prframecurrent->gr);
293 /* convert p(r) to sq */
294 sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent, start_q, end_q, q_step);
295 /* print frame data if needed */
296 if (opt2fn_null("-prframe", NFILE, fnm))
299 snew(suffix, GMX_PATH_MAX);
301 sprintf(hdr, "g(r), t = %f", t);
302 /* prepare output filename */
303 fnmdup = dup_tfn(NFILE, fnm);
304 sprintf(suffix, "-t%.2f", t);
305 add_suffix_to_output_names(fnmdup, NFILE, suffix);
306 fp = xvgropen(opt2fn_null("-prframe", NFILE, fnmdup), hdr, "Distance (nm)", "Probability", oenv);
307 for (i = 0; i < prframecurrent->grn; i++)
309 fprintf(fp, "%10.6f%10.6f\n", prframecurrent->r[i], prframecurrent->gr[i]);
311 done_filenms(NFILE, fnmdup);
317 if (opt2fn_null("-sqframe", NFILE, fnm))
320 snew(suffix, GMX_PATH_MAX);
322 sprintf(hdr, "I(q), t = %f", t);
323 /* prepare output filename */
324 fnmdup = dup_tfn(NFILE, fnm);
325 sprintf(suffix, "-t%.2f", t);
326 add_suffix_to_output_names(fnmdup, NFILE, suffix);
327 fp = xvgropen(opt2fn_null("-sqframe", NFILE, fnmdup), hdr, "q (nm^-1)", "s(q)/s(0)", oenv);
328 for (i = 0; i < sqframecurrent->qn; i++)
330 fprintf(fp, "%10.6f%10.6f\n", sqframecurrent->q[i], sqframecurrent->s[i]);
332 done_filenms(NFILE, fnmdup);
338 /* free pr structure */
339 sfree(prframecurrent->gr);
340 sfree(prframecurrent->r);
341 sfree(prframecurrent);
342 /* free sq structure */
343 sfree(sqframecurrent->q);
344 sfree(sqframecurrent->s);
345 sfree(sqframecurrent);
347 while (read_next_x(oenv, status, &t, x, box));
350 /* normalize histo */
351 normalize_probability(pr->grn, pr->gr);
352 sq = convert_histogram_to_intensity_curve(pr, start_q, end_q, q_step);
354 fp = xvgropen(opt2fn_null("-pr", NFILE, fnm), "G(r)", "Distance (nm)", "Probability", oenv);
355 for (i = 0; i < pr->grn; i++)
357 fprintf(fp, "%10.6f%10.6f\n", pr->r[i], pr->gr[i]);
362 fp = xvgropen(opt2fn_null("-sq", NFILE, fnm), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv);
363 for (i = 0; i < sq->qn; i++)
365 fprintf(fp, "%10.6f%10.6f\n", sq->q[i], sq->s[i]);
378 please_cite(stdout, "Garmay2012");