2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/gmxana/nsfactor.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/pbcutil/rmpbc.h"
52 #include "gromacs/topology/index.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/arraysize.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/gmxomp.h"
60 #include "gromacs/utility/pleasecite.h"
61 #include "gromacs/utility/smalloc.h"
63 int gmx_sans(int argc
, char* argv
[])
65 const char* desc
[] = {
66 "[THISMODULE] computes SANS spectra using Debye formula.",
67 "It currently uses topology file (since it need to assigne element for each atom).",
70 "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]",
71 "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]",
72 "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]",
73 "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]",
74 "[TT]-startq[tt] Starting q value in nm[PAR]",
75 "[TT]-endq[tt] Ending q value in nm[PAR]",
76 "[TT]-qstep[tt] Stepping in q space[PAR]",
77 "Note: When using Debye direct method computational cost increases as",
78 "1/2 * N * (N - 1) where N is atom number in group of interest.",
80 "WARNING: If sq or pr specified this tool can produce large number of files! Up to ",
81 "two times larger than number of frames!"
83 static gmx_bool bPBC
= TRUE
;
84 static gmx_bool bNORM
= FALSE
;
85 static real binwidth
= 0.2,
86 grid
= 0.05; /* bins shouldn't be smaller then smallest bond (~0.1nm) length */
87 static real start_q
= 0.0, end_q
= 2.0, q_step
= 0.01;
88 static real mcover
= -1;
89 static unsigned int seed
= 0;
90 static int nthreads
= -1;
92 static const char* emode
[] = { nullptr, "direct", "mc", nullptr };
93 static const char* emethod
[] = { nullptr, "debye", "fft", nullptr };
95 gmx_neutron_atomic_structurefactors_t
* gnsf
;
101 { "-bin", FALSE
, etREAL
, { &binwidth
}, "[HIDDEN]Binwidth (nm)" },
102 { "-mode", FALSE
, etENUM
, { emode
}, "Mode for sans spectra calculation" },
107 "Monte-Carlo coverage should be -1(default) or (0,1]" },
108 { "-method", FALSE
, etENUM
, { emethod
}, "[HIDDEN]Method for sans spectra calculation" },
113 "Use periodic boundary conditions for computing distances" },
114 { "-grid", FALSE
, etREAL
, { &grid
}, "[HIDDEN]Grid spacing (in nm) for FFTs" },
115 { "-startq", FALSE
, etREAL
, { &start_q
}, "Starting q (1/nm) " },
116 { "-endq", FALSE
, etREAL
, { &end_q
}, "Ending q (1/nm)" },
117 { "-qstep", FALSE
, etREAL
, { &q_step
}, "Stepping in q (1/nm)" },
118 { "-seed", FALSE
, etINT
, { &seed
}, "Random seed for Monte-Carlo" },
120 { "-nt", FALSE
, etINT
, { &nthreads
}, "Number of threads to start" },
124 const char * fnTPX
, *fnTRX
, *fnDAT
= nullptr;
126 t_topology
* top
= nullptr;
127 gmx_rmpbc_t gpbc
= nullptr;
128 gmx_bool bFFT
= FALSE
, bDEBYE
= FALSE
;
129 gmx_bool bMC
= FALSE
;
135 char** grpname
= nullptr;
136 int* index
= nullptr;
140 char* suffix
= nullptr;
141 gmx_radial_distribution_histogram_t
*prframecurrent
= nullptr, *pr
= nullptr;
142 gmx_static_structurefactor_t
* sqframecurrent
= nullptr, *sq
= nullptr;
143 gmx_output_env_t
* oenv
;
145 std::array
<t_filenm
, 8> filenames
= { { { efTPR
, "-s", nullptr, ffREAD
},
146 { efTRX
, "-f", nullptr, ffREAD
},
147 { efNDX
, nullptr, nullptr, ffOPTRD
},
148 { efDAT
, "-d", "nsfactor", ffOPTRD
},
149 { efXVG
, "-pr", "pr", ffWRITE
},
150 { efXVG
, "-sq", "sq", ffWRITE
},
151 { efXVG
, "-prframe", "prframe", ffOPTWR
},
152 { efXVG
, "-sqframe", "sqframe", ffOPTWR
} } };
153 t_filenm
* fnm
= filenames
.data();
155 const auto NFILE
= filenames
.size();
157 nthreads
= gmx_omp_get_max_threads();
159 if (!parse_common_args(&argc
, argv
, PCA_CAN_TIME
| PCA_TIME_UNIT
, NFILE
, fnm
, asize(pa
), pa
,
160 asize(desc
), desc
, 0, nullptr, &oenv
))
165 /* check that binwidth not smaller than smallers distance */
166 check_binwidth(binwidth
);
167 check_mcover(mcover
);
169 /* setting number of omp threads globaly */
170 gmx_omp_set_num_threads(nthreads
);
172 /* Now try to parse opts for modes */
173 GMX_RELEASE_ASSERT(emethod
[0] != nullptr, "Options inconsistency; emethod[0] is NULL");
174 switch (emethod
[0][0])
180 case 'd': bMC
= FALSE
; break;
181 case 'm': bMC
= TRUE
; break;
185 case 'f': bFFT
= TRUE
; break;
193 fprintf(stderr
, "Using Monte Carlo Debye method to calculate spectrum\n");
197 fprintf(stderr
, "Using direct Debye method to calculate spectrum\n");
202 gmx_fatal(FARGS
, "FFT method not implemented!");
206 gmx_fatal(FARGS
, "Unknown combination for mode and method!");
209 /* Try to read files */
210 fnDAT
= ftp2fn(efDAT
, NFILE
, fnm
);
211 fnTPX
= ftp2fn(efTPR
, NFILE
, fnm
);
212 fnTRX
= ftp2fn(efTRX
, NFILE
, fnm
);
214 gnsf
= gmx_neutronstructurefactors_init(fnDAT
);
215 fprintf(stderr
, "Read %d atom names from %s with neutron scattering parameters\n\n",
216 gnsf
->nratoms
, fnDAT
);
222 read_tps_conf(fnTPX
, top
, &ePBC
, &x
, nullptr, box
, TRUE
);
224 printf("\nPlease select group for SANS spectra calculation:\n");
225 get_index(&(top
->atoms
), ftp2fn_null(efNDX
, NFILE
, fnm
), 1, &isize
, &index
, grpname
);
227 gsans
= gmx_sans_init(top
, gnsf
);
229 /* Prepare reference frame */
232 gpbc
= gmx_rmpbc_init(&top
->idef
, ePBC
, top
->atoms
.nr
);
233 gmx_rmpbc(gpbc
, top
->atoms
.nr
, box
, x
);
236 natoms
= read_first_x(oenv
, &status
, fnTRX
, &t
, &x
, box
);
237 if (natoms
!= top
->atoms
.nr
)
239 fprintf(stderr
, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n",
240 natoms
, top
->atoms
.nr
);
247 gmx_rmpbc(gpbc
, top
->atoms
.nr
, box
, x
);
249 /* allocate memory for pr */
252 /* in case its first frame to read */
255 /* realy calc p(r) */
256 prframecurrent
= calc_radial_distribution_histogram(gsans
, x
, box
, index
, isize
, binwidth
,
257 bMC
, bNORM
, mcover
, seed
);
258 /* copy prframecurrent -> pr and summ up pr->gr[i] */
259 /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
260 if (pr
->gr
== nullptr)
262 /* check if we use pr->gr first time */
263 snew(pr
->gr
, prframecurrent
->grn
);
264 snew(pr
->r
, prframecurrent
->grn
);
268 /* resize pr->gr and pr->r if needed to preven overruns */
269 if (prframecurrent
->grn
> pr
->grn
)
271 srenew(pr
->gr
, prframecurrent
->grn
);
272 srenew(pr
->r
, prframecurrent
->grn
);
275 pr
->grn
= prframecurrent
->grn
;
276 pr
->binwidth
= prframecurrent
->binwidth
;
277 /* summ up gr and fill r */
278 for (i
= 0; i
< prframecurrent
->grn
; i
++)
280 pr
->gr
[i
] += prframecurrent
->gr
[i
];
281 pr
->r
[i
] = prframecurrent
->r
[i
];
283 /* normalize histo */
284 normalize_probability(prframecurrent
->grn
, prframecurrent
->gr
);
285 /* convert p(r) to sq */
286 sqframecurrent
= convert_histogram_to_intensity_curve(prframecurrent
, start_q
, end_q
, q_step
);
287 /* print frame data if needed */
288 if (opt2fn_null("-prframe", NFILE
, fnm
))
291 snew(suffix
, GMX_PATH_MAX
);
293 sprintf(hdr
, "g(r), t = %f", t
);
294 /* prepare output filename */
295 auto fnmdup
= filenames
;
296 sprintf(suffix
, "-t%.2f", t
);
297 add_suffix_to_output_names(fnmdup
.data(), NFILE
, suffix
);
298 fp
= xvgropen(opt2fn_null("-prframe", NFILE
, fnmdup
.data()), hdr
, "Distance (nm)",
299 "Probability", oenv
);
300 for (i
= 0; i
< prframecurrent
->grn
; i
++)
302 fprintf(fp
, "%10.6f%10.6f\n", prframecurrent
->r
[i
], prframecurrent
->gr
[i
]);
308 if (opt2fn_null("-sqframe", NFILE
, fnm
))
311 snew(suffix
, GMX_PATH_MAX
);
313 sprintf(hdr
, "I(q), t = %f", t
);
314 /* prepare output filename */
315 auto fnmdup
= filenames
;
316 sprintf(suffix
, "-t%.2f", t
);
317 add_suffix_to_output_names(fnmdup
.data(), NFILE
, suffix
);
318 fp
= xvgropen(opt2fn_null("-sqframe", NFILE
, fnmdup
.data()), hdr
, "q (nm^-1)",
320 for (i
= 0; i
< sqframecurrent
->qn
; i
++)
322 fprintf(fp
, "%10.6f%10.6f\n", sqframecurrent
->q
[i
], sqframecurrent
->s
[i
]);
328 /* free pr structure */
329 sfree(prframecurrent
->gr
);
330 sfree(prframecurrent
->r
);
331 sfree(prframecurrent
);
332 /* free sq structure */
333 sfree(sqframecurrent
->q
);
334 sfree(sqframecurrent
->s
);
335 sfree(sqframecurrent
);
336 } while (read_next_x(oenv
, status
, &t
, x
, box
));
339 /* normalize histo */
340 normalize_probability(pr
->grn
, pr
->gr
);
341 sq
= convert_histogram_to_intensity_curve(pr
, start_q
, end_q
, q_step
);
343 fp
= xvgropen(opt2fn_null("-pr", NFILE
, fnm
), "G(r)", "Distance (nm)", "Probability", oenv
);
344 for (i
= 0; i
< pr
->grn
; i
++)
346 fprintf(fp
, "%10.6f%10.6f\n", pr
->r
[i
], pr
->gr
[i
]);
351 fp
= xvgropen(opt2fn_null("-sq", NFILE
, fnm
), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv
);
352 for (i
= 0; i
< sq
->qn
; i
++)
354 fprintf(fp
, "%10.6f%10.6f\n", sq
->q
[i
], sq
->s
[i
]);
367 please_cite(stdout
, "Garmay2012");