2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/gstat.h"
48 #include "gromacs/gmxana/nsfactor.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/rmpbc.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/arraysize.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/gmxomp.h"
59 #include "gromacs/utility/pleasecite.h"
60 #include "gromacs/utility/smalloc.h"
62 int gmx_sans(int argc, char* argv[])
64 const char* desc[] = {
65 "[THISMODULE] computes SANS spectra using Debye formula.",
66 "It currently uses topology file (since it need to assigne element for each atom).",
69 "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]",
70 "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]",
71 "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]",
72 "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]",
73 "[TT]-startq[tt] Starting q value in nm[PAR]",
74 "[TT]-endq[tt] Ending q value in nm[PAR]",
75 "[TT]-qstep[tt] Stepping in q space[PAR]",
76 "Note: When using Debye direct method computational cost increases as",
77 "1/2 * N * (N - 1) where N is atom number in group of interest.",
79 "WARNING: If sq or pr specified this tool can produce large number of files! Up to ",
80 "two times larger than number of frames!"
82 static gmx_bool bPBC = TRUE;
83 static gmx_bool bNORM = FALSE;
84 static real binwidth = 0.2,
85 grid = 0.05; /* bins shouldn't be smaller then smallest bond (~0.1nm) length */
86 static real start_q = 0.0, end_q = 2.0, q_step = 0.01;
87 static real mcover = -1;
88 static unsigned int seed = 0;
89 static int nthreads = -1;
91 static const char* emode[] = { nullptr, "direct", "mc", nullptr };
92 static const char* emethod[] = { nullptr, "debye", "fft", nullptr };
94 gmx_neutron_atomic_structurefactors_t* gnsf;
100 { "-bin", FALSE, etREAL, { &binwidth }, "[HIDDEN]Binwidth (nm)" },
101 { "-mode", FALSE, etENUM, { emode }, "Mode for sans spectra calculation" },
106 "Monte-Carlo coverage should be -1(default) or (0,1]" },
107 { "-method", FALSE, etENUM, { emethod }, "[HIDDEN]Method for sans spectra calculation" },
112 "Use periodic boundary conditions for computing distances" },
113 { "-grid", FALSE, etREAL, { &grid }, "[HIDDEN]Grid spacing (in nm) for FFTs" },
114 { "-startq", FALSE, etREAL, { &start_q }, "Starting q (1/nm) " },
115 { "-endq", FALSE, etREAL, { &end_q }, "Ending q (1/nm)" },
116 { "-qstep", FALSE, etREAL, { &q_step }, "Stepping in q (1/nm)" },
117 { "-seed", FALSE, etINT, { &seed }, "Random seed for Monte-Carlo" },
119 { "-nt", FALSE, etINT, { &nthreads }, "Number of threads to start" },
123 const char * fnTPX, *fnTRX, *fnDAT = nullptr;
125 t_topology* top = nullptr;
126 gmx_rmpbc_t gpbc = nullptr;
127 gmx_bool bFFT = FALSE, bDEBYE = FALSE;
128 gmx_bool bMC = FALSE;
134 char** grpname = nullptr;
135 int* index = nullptr;
139 char* suffix = nullptr;
140 gmx_radial_distribution_histogram_t *prframecurrent = nullptr, *pr = nullptr;
141 gmx_static_structurefactor_t * sqframecurrent = nullptr, *sq = nullptr;
142 gmx_output_env_t* oenv;
144 std::array<t_filenm, 8> filenames = { { { efTPR, "-s", nullptr, ffREAD },
145 { efTRX, "-f", nullptr, ffREAD },
146 { efNDX, nullptr, nullptr, ffOPTRD },
147 { efDAT, "-d", "nsfactor", ffOPTRD },
148 { efXVG, "-pr", "pr", ffWRITE },
149 { efXVG, "-sq", "sq", ffWRITE },
150 { efXVG, "-prframe", "prframe", ffOPTWR },
151 { efXVG, "-sqframe", "sqframe", ffOPTWR } } };
152 t_filenm* fnm = filenames.data();
154 const auto NFILE = filenames.size();
156 nthreads = gmx_omp_get_max_threads();
158 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT, NFILE, fnm, asize(pa), pa,
159 asize(desc), desc, 0, nullptr, &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 */
169 gmx_omp_set_num_threads(nthreads);
171 /* Now try to parse opts for modes */
172 GMX_RELEASE_ASSERT(emethod[0] != nullptr, "Options inconsistency; emethod[0] is NULL");
173 switch (emethod[0][0])
179 case 'd': bMC = FALSE; break;
180 case 'm': bMC = TRUE; break;
184 case 'f': bFFT = TRUE; break;
192 fprintf(stderr, "Using Monte Carlo Debye method to calculate spectrum\n");
196 fprintf(stderr, "Using direct Debye method to calculate spectrum\n");
201 gmx_fatal(FARGS, "FFT method not implemented!");
205 gmx_fatal(FARGS, "Unknown combination for mode and method!");
208 /* Try to read files */
209 fnDAT = ftp2fn(efDAT, NFILE, fnm);
210 fnTPX = ftp2fn(efTPR, NFILE, fnm);
211 fnTRX = ftp2fn(efTRX, NFILE, fnm);
213 gnsf = gmx_neutronstructurefactors_init(fnDAT);
214 fprintf(stderr, "Read %d atom names from %s with neutron scattering parameters\n\n",
215 gnsf->nratoms, fnDAT);
221 read_tps_conf(fnTPX, top, &ePBC, &x, nullptr, box, TRUE);
223 printf("\nPlease select group for SANS spectra calculation:\n");
224 get_index(&(top->atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, grpname);
226 gsans = gmx_sans_init(top, gnsf);
228 /* Prepare reference frame */
231 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
232 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
235 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
236 if (natoms != top->atoms.nr)
238 fprintf(stderr, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n",
239 natoms, top->atoms.nr);
246 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
248 /* allocate memory for pr */
251 /* in case its first frame to read */
254 /* realy calc p(r) */
255 prframecurrent = calc_radial_distribution_histogram(gsans, x, box, index, isize, binwidth,
256 bMC, bNORM, mcover, seed);
257 /* copy prframecurrent -> pr and summ up pr->gr[i] */
258 /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
259 if (pr->gr == nullptr)
261 /* check if we use pr->gr first time */
262 snew(pr->gr, prframecurrent->grn);
263 snew(pr->r, prframecurrent->grn);
267 /* resize pr->gr and pr->r if needed to preven overruns */
268 if (prframecurrent->grn > pr->grn)
270 srenew(pr->gr, prframecurrent->grn);
271 srenew(pr->r, prframecurrent->grn);
274 pr->grn = prframecurrent->grn;
275 pr->binwidth = prframecurrent->binwidth;
276 /* summ up gr and fill r */
277 for (i = 0; i < prframecurrent->grn; i++)
279 pr->gr[i] += prframecurrent->gr[i];
280 pr->r[i] = prframecurrent->r[i];
282 /* normalize histo */
283 normalize_probability(prframecurrent->grn, prframecurrent->gr);
284 /* convert p(r) to sq */
285 sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent, start_q, end_q, q_step);
286 /* print frame data if needed */
287 if (opt2fn_null("-prframe", NFILE, fnm))
290 snew(suffix, GMX_PATH_MAX);
292 sprintf(hdr, "g(r), t = %f", t);
293 /* prepare output filename */
294 auto fnmdup = filenames;
295 sprintf(suffix, "-t%.2f", t);
296 add_suffix_to_output_names(fnmdup.data(), NFILE, suffix);
297 fp = xvgropen(opt2fn_null("-prframe", NFILE, fnmdup.data()), hdr, "Distance (nm)",
298 "Probability", oenv);
299 for (i = 0; i < prframecurrent->grn; i++)
301 fprintf(fp, "%10.6f%10.6f\n", prframecurrent->r[i], prframecurrent->gr[i]);
307 if (opt2fn_null("-sqframe", NFILE, fnm))
310 snew(suffix, GMX_PATH_MAX);
312 sprintf(hdr, "I(q), t = %f", t);
313 /* prepare output filename */
314 auto fnmdup = filenames;
315 sprintf(suffix, "-t%.2f", t);
316 add_suffix_to_output_names(fnmdup.data(), NFILE, suffix);
317 fp = xvgropen(opt2fn_null("-sqframe", NFILE, fnmdup.data()), hdr, "q (nm^-1)",
319 for (i = 0; i < sqframecurrent->qn; i++)
321 fprintf(fp, "%10.6f%10.6f\n", sqframecurrent->q[i], sqframecurrent->s[i]);
327 /* free pr structure */
328 sfree(prframecurrent->gr);
329 sfree(prframecurrent->r);
330 sfree(prframecurrent);
331 /* free sq structure */
332 sfree(sqframecurrent->q);
333 sfree(sqframecurrent->s);
334 sfree(sqframecurrent);
335 } while (read_next_x(oenv, status, &t, x, box));
338 /* normalize histo */
339 normalize_probability(pr->grn, pr->gr);
340 sq = convert_histogram_to_intensity_curve(pr, start_q, end_q, q_step);
342 fp = xvgropen(opt2fn_null("-pr", NFILE, fnm), "G(r)", "Distance (nm)", "Probability", oenv);
343 for (i = 0; i < pr->grn; i++)
345 fprintf(fp, "%10.6f%10.6f\n", pr->r[i], pr->gr[i]);
350 fp = xvgropen(opt2fn_null("-sq", NFILE, fnm), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv);
351 for (i = 0; i < sq->qn; i++)
353 fprintf(fp, "%10.6f%10.6f\n", sq->q[i], sq->s[i]);
366 please_cite(stdout, "Garmay2012");