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.
45 #include "gromacs/math/vec.h"
46 #include "gromacs/random/threefry.h"
47 #include "gromacs/random/uniformintdistribution.h"
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/futil.h"
53 #include "gromacs/utility/gmxomp.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/strdb.h"
57 void check_binwidth(real binwidth)
59 real smallest_bin = 0.1;
60 if (binwidth < smallest_bin)
63 "Binwidth shouldn't be smaller then smallest bond length (H-H bond ~0.1nm) in a "
68 void check_mcover(real mcover)
72 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
74 else if ((mcover < 0) && (mcover != -1))
76 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
84 void normalize_probability(int n, double* a)
88 for (i = 0; i < n; i++)
92 for (i = 0; i < n; i++)
98 gmx_neutron_atomic_structurefactors_t* gmx_neutronstructurefactors_init(const char* datfn)
100 /* read nsfactor.dat */
107 gmx_neutron_atomic_structurefactors_t* gnsf;
109 gmx::FilePtr fp = gmx::openLibraryFile(datfn);
111 /* allocate memory for structure */
113 snew(gnsf->atomnm, nralloc);
114 snew(gnsf->p, nralloc);
115 snew(gnsf->n, nralloc);
116 snew(gnsf->slength, nralloc);
118 gnsf->nratoms = line_no;
120 while (get_a_line(fp.get(), line, STRLEN))
123 if (sscanf(line, "%s %d %d %lf", atomnm, &p, &n, &slength) == 4)
125 gnsf->atomnm[i] = gmx_strdup(atomnm);
128 gnsf->slength[i] = slength;
130 gnsf->nratoms = line_no;
131 if (line_no == nralloc)
134 srenew(gnsf->atomnm, nralloc);
135 srenew(gnsf->p, nralloc);
136 srenew(gnsf->n, nralloc);
137 srenew(gnsf->slength, nralloc);
142 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n", datfn, line_no);
145 srenew(gnsf->atomnm, gnsf->nratoms);
146 srenew(gnsf->p, gnsf->nratoms);
147 srenew(gnsf->n, gnsf->nratoms);
148 srenew(gnsf->slength, gnsf->nratoms);
153 gmx_sans_t* gmx_sans_init(const t_topology* top, gmx_neutron_atomic_structurefactors_t* gnsf)
155 gmx_sans_t* gsans = nullptr;
157 /* Try to assing scattering length from nsfactor.dat */
159 snew(gsans->slength, top->atoms.nr);
160 /* copy topology data */
162 for (i = 0; i < top->atoms.nr; i++)
164 for (j = 0; j < gnsf->nratoms; j++)
166 if (top->atoms.atom[i].atomnumber == gnsf->p[j])
168 /* we need special case for H and D */
169 if (top->atoms.atom[i].atomnumber == 1)
171 if (top->atoms.atom[i].m == 1.008000)
173 gsans->slength[i] = gnsf->slength[0];
177 gsans->slength[i] = gnsf->slength[1];
182 gsans->slength[i] = gnsf->slength[j];
191 gmx_radial_distribution_histogram_t* calc_radial_distribution_histogram(gmx_sans_t* gsans,
202 gmx_radial_distribution_histogram_t* pr = nullptr;
210 gmx::DefaultRandomEngine* trng = nullptr;
212 int64_t mc = 0, mc_max;
213 gmx::DefaultRandomEngine rng(seed);
215 /* allocate memory for pr */
217 /* set some fields */
218 pr->binwidth = binwidth;
221 * create max dist rvec
222 * dist = box[xx] + box[yy] + box[zz]
224 rvec_add(box[XX], box[YY], dist);
225 rvec_add(box[ZZ], dist, dist);
229 pr->grn = static_cast<int>(std::floor(rmax / pr->binwidth) + 1);
231 snew(pr->gr, pr->grn);
235 /* Special case for setting automaticaly number of mc iterations to 1% of total number of direct iterations */
238 mc_max = static_cast<int64_t>(std::floor(0.5 * 0.01 * isize * (isize - 1)));
242 mc_max = static_cast<int64_t>(std::floor(0.5 * mcover * isize * (isize - 1)));
245 nthreads = gmx_omp_get_max_threads();
247 trng = new gmx::DefaultRandomEngine[nthreads];
248 for (i = 0; i < nthreads; i++)
250 snew(tgr[i], pr->grn);
253 # pragma omp parallel shared(tgr, trng, mc) private(tid, i, j)
255 gmx::UniformIntDistribution<int> tdist(0, isize - 1);
256 tid = gmx_omp_get_thread_num();
257 /* now starting parallel threads */
259 for (mc = 0; mc < mc_max; mc++)
263 i = tdist(trng[tid]); // [0,isize-1]
264 j = tdist(trng[tid]); // [0,isize-1]
267 tgr[tid][static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
268 gsans->slength[index[i]] * gsans->slength[index[j]];
271 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
274 /* collecting data from threads */
275 for (i = 0; i < pr->grn; i++)
277 for (j = 0; j < nthreads; j++)
279 pr->gr[i] += tgr[j][i];
282 /* freeing memory for tgr and destroying trng */
283 for (i = 0; i < nthreads; i++)
290 gmx::UniformIntDistribution<int> dist(0, isize - 1);
291 for (mc = 0; mc < mc_max; mc++)
293 i = dist(rng); // [0,isize-1]
294 j = dist(rng); // [0,isize-1]
297 pr->gr[static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
298 gsans->slength[index[i]] * gsans->slength[index[j]];
306 nthreads = gmx_omp_get_max_threads();
307 /* Allocating memory for tgr arrays */
309 for (i = 0; i < nthreads; i++)
311 snew(tgr[i], pr->grn);
313 # pragma omp parallel shared(tgr) private(tid, i, j)
315 tid = gmx_omp_get_thread_num();
316 /* starting parallel threads */
318 for (i = 0; i < isize; i++)
322 for (j = 0; j < i; j++)
324 tgr[tid][static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
325 gsans->slength[index[i]] * gsans->slength[index[j]];
328 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
331 /* collecating data for pr->gr */
332 for (i = 0; i < pr->grn; i++)
334 for (j = 0; j < nthreads; j++)
336 pr->gr[i] += tgr[j][i];
339 /* freeing memory for tgr */
340 for (i = 0; i < nthreads; i++)
346 for (i = 0; i < isize; i++)
348 for (j = 0; j < i; j++)
350 pr->gr[static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
351 gsans->slength[index[i]] * gsans->slength[index[j]];
357 /* normalize if needed */
360 normalize_probability(pr->grn, pr->gr);
363 snew(pr->r, pr->grn);
364 for (i = 0; i < pr->grn; i++)
366 pr->r[i] = (pr->binwidth * i + pr->binwidth * 0.5);
372 gmx_static_structurefactor_t* convert_histogram_to_intensity_curve(gmx_radial_distribution_histogram_t* pr,
377 gmx_static_structurefactor_t* sq = nullptr;
381 sq->qn = static_cast<int>(std::floor((end_q - start_q) / q_step));
384 for (i = 0; i < sq->qn; i++)
386 sq->q[i] = start_q + i * q_step;
392 for (i = 1; i < sq->qn; i++)
394 for (j = 0; j < pr->grn; j++)
396 sq->s[i] += (pr->gr[j] / pr->r[j]) * std::sin(sq->q[i] * pr->r[j]);
398 sq->s[i] /= sq->q[i];
403 for (i = 0; i < sq->qn; i++)
405 for (j = 0; j < pr->grn; j++)
407 sq->s[i] += (pr->gr[j] / pr->r[j]) * std::sin(sq->q[i] * pr->r[j]);
409 sq->s[i] /= sq->q[i];