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
50 #include "gromacs/fileio/futil.h"
52 #include "gromacs/fileio/tpxio.h"
60 #include "gromacs/fileio/matio.h"
65 int gmx_saxs(int argc, char *argv[])
67 const char *desc[] = {
68 "[THISMODULE] calculates SAXS structure factors for given index",
69 "groups based on Cromer's method.",
70 "Both topology and trajectory files are required."
73 static real start_q = 0.0, end_q = 60.0, energy = 12.0;
74 static int ngroups = 1;
77 { "-ng", FALSE, etINT, {&ngroups},
78 "Number of groups to compute SAXS" },
79 {"-startq", FALSE, etREAL, {&start_q},
80 "Starting q (1/nm) "},
81 {"-endq", FALSE, etREAL, {&end_q},
83 {"-energy", FALSE, etREAL, {&energy},
84 "Energy of the incoming X-ray (keV) "}
87 const char *fnTPS, *fnTRX, *fnNDX, *fnDAT = NULL;
91 { efTRX, "-f", NULL, ffREAD },
92 { efTPS, NULL, NULL, ffREAD },
93 { efNDX, NULL, NULL, ffOPTRD },
94 { efDAT, "-d", "sfactor", ffOPTRD },
95 { efXVG, "-sq", "sq", ffWRITE },
97 #define NFILE asize(fnm)
98 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE,
99 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
104 fnTPS = ftp2fn(efTPS, NFILE, fnm);
105 fnTRX = ftp2fn(efTRX, NFILE, fnm);
106 fnDAT = ftp2fn(efDAT, NFILE, fnm);
107 fnNDX = ftp2fn_null(efNDX, NFILE, fnm);
109 do_scattering_intensity(fnTPS, fnNDX, opt2fn("-sq", NFILE, fnm),
111 start_q, end_q, energy, ngroups, oenv);
113 please_cite(stdout, "Cromer1968a");