Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / gmxana / gmx_sans.cpp
blobb9b608dbb92d3112991a68d978bd4eff71b4746c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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.
36 #include "gmxpre.h"
38 #include "config.h"
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/fileio/confio.h"
42 #include "gromacs/fileio/trxio.h"
43 #include "gromacs/fileio/xvgr.h"
44 #include "gromacs/gmxana/gmx_ana.h"
45 #include "gromacs/gmxana/gstat.h"
46 #include "gromacs/gmxana/nsfactor.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/pbcutil/rmpbc.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/arraysize.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/gmxomp.h"
57 #include "gromacs/utility/pleasecite.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).",
65 "[PAR]",
66 "Parameters:[PAR]"
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.",
76 "[PAR]",
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 shouldn't 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[] = { nullptr, "direct", "mc", nullptr };
88 static const char *emethod[] = { nullptr, "debye", "fft", nullptr };
90 gmx_neutron_atomic_structurefactors_t *gnsf;
91 gmx_sans_t *gsans;
93 #define NPA asize(pa)
95 t_pargs pa[] = {
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},
111 "Ending q (1/nm)"},
112 { "-qstep", FALSE, etREAL, {&q_step},
113 "Stepping in q (1/nm)"},
114 { "-seed", FALSE, etINT, {&seed},
115 "Random seed for Monte-Carlo"},
116 #if GMX_OPENMP
117 { "-nt", FALSE, etINT, {&nthreads},
118 "Number of threads to start"},
119 #endif
121 FILE *fp;
122 const char *fnTPX, *fnTRX, *fnDAT = nullptr;
123 t_trxstatus *status;
124 t_topology *top = nullptr;
125 gmx_rmpbc_t gpbc = nullptr;
126 gmx_bool bFFT = FALSE, bDEBYE = FALSE;
127 gmx_bool bMC = FALSE;
128 int ePBC = -1;
129 matrix box;
130 rvec *x;
131 int natoms;
132 real t;
133 char **grpname = nullptr;
134 int *index = nullptr;
135 int isize;
136 int i;
137 char *hdr = nullptr;
138 char *suffix = nullptr;
139 t_filenm *fnmdup = 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 #define NFILE asize(fnm)
146 t_filenm fnm[] = {
147 { efTPR, "-s", nullptr, ffREAD },
148 { efTRX, "-f", nullptr, ffREAD },
149 { efNDX, nullptr, nullptr, ffOPTRD },
150 { efDAT, "-d", "nsfactor", ffOPTRD },
151 { efXVG, "-pr", "pr", ffWRITE },
152 { efXVG, "-sq", "sq", ffWRITE },
153 { efXVG, "-prframe", "prframe", ffOPTWR },
154 { efXVG, "-sqframe", "sqframe", ffOPTWR }
157 nthreads = gmx_omp_get_max_threads();
159 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT,
160 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
162 return 0;
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])
176 case 'd':
177 bDEBYE = TRUE;
178 switch (emode[0][0])
180 case 'd':
181 bMC = FALSE;
182 break;
183 case 'm':
184 bMC = TRUE;
185 break;
186 default:
187 break;
189 break;
190 case 'f':
191 bFFT = TRUE;
192 break;
193 default:
194 break;
197 if (bDEBYE)
199 if (bMC)
201 fprintf(stderr, "Using Monte Carlo Debye method to calculate spectrum\n");
203 else
205 fprintf(stderr, "Using direct Debye method to calculate spectrum\n");
208 else if (bFFT)
210 gmx_fatal(FARGS, "FFT method not implemented!");
212 else
214 gmx_fatal(FARGS, "Unknown combination for mode and method!");
217 /* Try to read files */
218 fnDAT = ftp2fn(efDAT, NFILE, fnm);
219 fnTPX = ftp2fn(efTPR, NFILE, fnm);
220 fnTRX = ftp2fn(efTRX, NFILE, fnm);
222 gnsf = gmx_neutronstructurefactors_init(fnDAT);
223 fprintf(stderr, "Read %d atom names from %s with neutron scattering parameters\n\n", gnsf->nratoms, fnDAT);
225 snew(top, 1);
226 snew(grpname, 1);
227 snew(index, 1);
229 read_tps_conf(fnTPX, top, &ePBC, &x, nullptr, box, TRUE);
231 printf("\nPlease select group for SANS spectra calculation:\n");
232 get_index(&(top->atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, grpname);
234 gsans = gmx_sans_init(top, gnsf);
236 /* Prepare reference frame */
237 if (bPBC)
239 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
240 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
243 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
244 if (natoms != top->atoms.nr)
246 fprintf(stderr, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n", natoms, top->atoms.nr);
251 if (bPBC)
253 gmx_rmpbc(gpbc, top->atoms.nr, box, x);
255 /* allocate memory for pr */
256 if (pr == nullptr)
258 /* in case its first frame to read */
259 snew(pr, 1);
261 /* realy calc p(r) */
262 prframecurrent = calc_radial_distribution_histogram(gsans, x, box, index, isize, binwidth, bMC, bNORM, mcover, seed);
263 /* copy prframecurrent -> pr and summ up pr->gr[i] */
264 /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
265 if (pr->gr == nullptr)
267 /* check if we use pr->gr first time */
268 snew(pr->gr, prframecurrent->grn);
269 snew(pr->r, prframecurrent->grn);
271 else
273 /* resize pr->gr and pr->r if needed to preven overruns */
274 if (prframecurrent->grn > pr->grn)
276 srenew(pr->gr, prframecurrent->grn);
277 srenew(pr->r, prframecurrent->grn);
280 pr->grn = prframecurrent->grn;
281 pr->binwidth = prframecurrent->binwidth;
282 /* summ up gr and fill r */
283 for (i = 0; i < prframecurrent->grn; i++)
285 pr->gr[i] += prframecurrent->gr[i];
286 pr->r[i] = prframecurrent->r[i];
288 /* normalize histo */
289 normalize_probability(prframecurrent->grn, prframecurrent->gr);
290 /* convert p(r) to sq */
291 sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent, start_q, end_q, q_step);
292 /* print frame data if needed */
293 if (opt2fn_null("-prframe", NFILE, fnm))
295 snew(hdr, 25);
296 snew(suffix, GMX_PATH_MAX);
297 /* prepare header */
298 sprintf(hdr, "g(r), t = %f", t);
299 /* prepare output filename */
300 fnmdup = dup_tfn(NFILE, fnm);
301 sprintf(suffix, "-t%.2f", t);
302 add_suffix_to_output_names(fnmdup, NFILE, suffix);
303 fp = xvgropen(opt2fn_null("-prframe", NFILE, fnmdup), hdr, "Distance (nm)", "Probability", oenv);
304 for (i = 0; i < prframecurrent->grn; i++)
306 fprintf(fp, "%10.6f%10.6f\n", prframecurrent->r[i], prframecurrent->gr[i]);
308 done_filenms(NFILE, fnmdup);
309 xvgrclose(fp);
310 sfree(hdr);
311 sfree(suffix);
312 sfree(fnmdup);
314 if (opt2fn_null("-sqframe", NFILE, fnm))
316 snew(hdr, 25);
317 snew(suffix, GMX_PATH_MAX);
318 /* prepare header */
319 sprintf(hdr, "I(q), t = %f", t);
320 /* prepare output filename */
321 fnmdup = dup_tfn(NFILE, fnm);
322 sprintf(suffix, "-t%.2f", t);
323 add_suffix_to_output_names(fnmdup, NFILE, suffix);
324 fp = xvgropen(opt2fn_null("-sqframe", NFILE, fnmdup), hdr, "q (nm^-1)", "s(q)/s(0)", oenv);
325 for (i = 0; i < sqframecurrent->qn; i++)
327 fprintf(fp, "%10.6f%10.6f\n", sqframecurrent->q[i], sqframecurrent->s[i]);
329 done_filenms(NFILE, fnmdup);
330 xvgrclose(fp);
331 sfree(hdr);
332 sfree(suffix);
333 sfree(fnmdup);
335 /* free pr structure */
336 sfree(prframecurrent->gr);
337 sfree(prframecurrent->r);
338 sfree(prframecurrent);
339 /* free sq structure */
340 sfree(sqframecurrent->q);
341 sfree(sqframecurrent->s);
342 sfree(sqframecurrent);
344 while (read_next_x(oenv, status, &t, x, box));
345 close_trx(status);
347 /* normalize histo */
348 normalize_probability(pr->grn, pr->gr);
349 sq = convert_histogram_to_intensity_curve(pr, start_q, end_q, q_step);
350 /* prepare pr.xvg */
351 fp = xvgropen(opt2fn_null("-pr", NFILE, fnm), "G(r)", "Distance (nm)", "Probability", oenv);
352 for (i = 0; i < pr->grn; i++)
354 fprintf(fp, "%10.6f%10.6f\n", pr->r[i], pr->gr[i]);
356 xvgrclose(fp);
358 /* prepare sq.xvg */
359 fp = xvgropen(opt2fn_null("-sq", NFILE, fnm), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv);
360 for (i = 0; i < sq->qn; i++)
362 fprintf(fp, "%10.6f%10.6f\n", sq->q[i], sq->s[i]);
364 xvgrclose(fp);
366 * Clean up memory
368 sfree(pr->gr);
369 sfree(pr->r);
370 sfree(pr);
371 sfree(sq->q);
372 sfree(sq->s);
373 sfree(sq);
375 please_cite(stdout, "Garmay2012");
377 return 0;