Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / nsfactor.cpp
blobe56b568d756502b9ce707536d2b2b89d6a50ff6a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016, 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.
36 #include "gmxpre.h"
38 #include "nsfactor.h"
40 #include "config.h"
42 #include <cmath>
43 #include <cstring>
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/basedefinitions.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/exceptions.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/gmxomp.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/strdb.h"
58 void check_binwidth(real binwidth)
60 real smallest_bin = 0.1;
61 if (binwidth < smallest_bin)
63 gmx_fatal(FARGS,
64 "Binwidth shouldn't be smaller then smallest bond length (H-H bond ~0.1nm) in a "
65 "box");
69 void check_mcover(real mcover)
71 if (mcover > 1.0)
73 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
75 else if ((mcover < 0) && (mcover != -1))
77 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
79 else
81 return;
85 void normalize_probability(int n, double* a)
87 int i;
88 double norm = 0.0;
89 for (i = 0; i < n; i++)
91 norm += a[i];
93 for (i = 0; i < n; i++)
95 a[i] /= norm;
99 gmx_neutron_atomic_structurefactors_t* gmx_neutronstructurefactors_init(const char* datfn)
101 /* read nsfactor.dat */
102 char line[STRLEN];
103 int nralloc = 10;
104 int n, p;
105 int i, line_no;
106 char atomnm[8];
107 double slength;
108 gmx_neutron_atomic_structurefactors_t* gnsf;
110 gmx::FilePtr fp = gmx::openLibraryFile(datfn);
111 line_no = 0;
112 /* allocate memory for structure */
113 snew(gnsf, nralloc);
114 snew(gnsf->atomnm, nralloc);
115 snew(gnsf->p, nralloc);
116 snew(gnsf->n, nralloc);
117 snew(gnsf->slength, nralloc);
119 gnsf->nratoms = line_no;
121 while (get_a_line(fp.get(), line, STRLEN))
123 i = line_no;
124 if (sscanf(line, "%s %d %d %lf", atomnm, &p, &n, &slength) == 4)
126 gnsf->atomnm[i] = gmx_strdup(atomnm);
127 gnsf->n[i] = n;
128 gnsf->p[i] = p;
129 gnsf->slength[i] = slength;
130 line_no++;
131 gnsf->nratoms = line_no;
132 if (line_no == nralloc)
134 nralloc++;
135 srenew(gnsf->atomnm, nralloc);
136 srenew(gnsf->p, nralloc);
137 srenew(gnsf->n, nralloc);
138 srenew(gnsf->slength, nralloc);
141 else
143 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n", datfn, line_no);
146 srenew(gnsf->atomnm, gnsf->nratoms);
147 srenew(gnsf->p, gnsf->nratoms);
148 srenew(gnsf->n, gnsf->nratoms);
149 srenew(gnsf->slength, gnsf->nratoms);
151 return gnsf;
154 gmx_sans_t* gmx_sans_init(const t_topology* top, gmx_neutron_atomic_structurefactors_t* gnsf)
156 gmx_sans_t* gsans = nullptr;
157 int i, j;
158 /* Try to assing scattering length from nsfactor.dat */
159 snew(gsans, 1);
160 snew(gsans->slength, top->atoms.nr);
161 /* copy topology data */
162 gsans->top = top;
163 for (i = 0; i < top->atoms.nr; i++)
165 for (j = 0; j < gnsf->nratoms; j++)
167 if (top->atoms.atom[i].atomnumber == gnsf->p[j])
169 /* we need special case for H and D */
170 if (top->atoms.atom[i].atomnumber == 1)
172 if (top->atoms.atom[i].m == 1.008000)
174 gsans->slength[i] = gnsf->slength[0];
176 else
178 gsans->slength[i] = gnsf->slength[1];
181 else
183 gsans->slength[i] = gnsf->slength[j];
189 return gsans;
192 gmx_radial_distribution_histogram_t* calc_radial_distribution_histogram(gmx_sans_t* gsans,
193 rvec* x,
194 matrix box,
195 const int* index,
196 int isize,
197 double binwidth,
198 gmx_bool bMC,
199 gmx_bool bNORM,
200 real mcover,
201 unsigned int seed)
203 gmx_radial_distribution_histogram_t* pr = nullptr;
204 rvec dist;
205 double rmax;
206 int i, j;
207 #if GMX_OPENMP
208 double** tgr;
209 int tid;
210 int nthreads;
211 gmx::DefaultRandomEngine* trng = nullptr;
212 #endif
213 int64_t mc_max;
214 gmx::DefaultRandomEngine rng(seed);
216 /* allocate memory for pr */
217 snew(pr, 1);
218 /* set some fields */
219 pr->binwidth = binwidth;
222 * create max dist rvec
223 * dist = box[xx] + box[yy] + box[zz]
225 rvec_add(box[XX], box[YY], dist);
226 rvec_add(box[ZZ], dist, dist);
228 rmax = norm(dist);
230 pr->grn = static_cast<int>(std::floor(rmax / pr->binwidth) + 1);
232 snew(pr->gr, pr->grn);
234 if (bMC)
236 /* Special case for setting automaticaly number of mc iterations to 1% of total number of direct iterations */
237 if (mcover == -1)
239 mc_max = static_cast<int64_t>(std::floor(0.5 * 0.01 * isize * (isize - 1)));
241 else
243 mc_max = static_cast<int64_t>(std::floor(0.5 * mcover * isize * (isize - 1)));
245 #if GMX_OPENMP
246 nthreads = gmx_omp_get_max_threads();
247 snew(tgr, nthreads);
248 trng = new gmx::DefaultRandomEngine[nthreads];
249 for (i = 0; i < nthreads; i++)
251 snew(tgr[i], pr->grn);
252 trng[i].seed(rng());
254 # pragma omp parallel shared(tgr, trng) private(tid, i, j)
256 gmx::UniformIntDistribution<int> tdist(0, isize - 1);
257 tid = gmx_omp_get_thread_num();
258 /* now starting parallel threads */
259 INTEL_DIAGNOSTIC_IGNORE(593)
260 # pragma omp for
261 for (int64_t mc = 0; mc < mc_max; mc++)
265 i = tdist(trng[tid]); // [0,isize-1]
266 j = tdist(trng[tid]); // [0,isize-1]
267 if (i != j)
269 tgr[tid][static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
270 gsans->slength[index[i]] * gsans->slength[index[j]];
273 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
276 INTEL_DIAGNOSTIC_RESET
277 /* collecting data from threads */
278 for (i = 0; i < pr->grn; i++)
280 for (j = 0; j < nthreads; j++)
282 pr->gr[i] += tgr[j][i];
285 /* freeing memory for tgr and destroying trng */
286 for (i = 0; i < nthreads; i++)
288 sfree(tgr[i]);
290 sfree(tgr);
291 delete[] trng;
292 #else
293 gmx::UniformIntDistribution<int> dist(0, isize - 1);
294 for (int64_t mc = 0; mc < mc_max; mc++)
296 i = dist(rng); // [0,isize-1]
297 j = dist(rng); // [0,isize-1]
298 if (i != j)
300 pr->gr[static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
301 gsans->slength[index[i]] * gsans->slength[index[j]];
304 #endif
306 else
308 #if GMX_OPENMP
309 nthreads = gmx_omp_get_max_threads();
310 /* Allocating memory for tgr arrays */
311 snew(tgr, nthreads);
312 for (i = 0; i < nthreads; i++)
314 snew(tgr[i], pr->grn);
316 # pragma omp parallel shared(tgr) private(tid, i, j)
318 tid = gmx_omp_get_thread_num();
319 /* starting parallel threads */
320 # pragma omp for
321 for (i = 0; i < isize; i++)
325 for (j = 0; j < i; j++)
327 tgr[tid][static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
328 gsans->slength[index[i]] * gsans->slength[index[j]];
331 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
334 /* collecating data for pr->gr */
335 for (i = 0; i < pr->grn; i++)
337 for (j = 0; j < nthreads; j++)
339 pr->gr[i] += tgr[j][i];
342 /* freeing memory for tgr */
343 for (i = 0; i < nthreads; i++)
345 sfree(tgr[i]);
347 sfree(tgr);
348 #else
349 for (i = 0; i < isize; i++)
351 for (j = 0; j < i; j++)
353 pr->gr[static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]])) / binwidth))] +=
354 gsans->slength[index[i]] * gsans->slength[index[j]];
357 #endif
360 /* normalize if needed */
361 if (bNORM)
363 normalize_probability(pr->grn, pr->gr);
366 snew(pr->r, pr->grn);
367 for (i = 0; i < pr->grn; i++)
369 pr->r[i] = (pr->binwidth * i + pr->binwidth * 0.5);
372 return pr;
375 gmx_static_structurefactor_t* convert_histogram_to_intensity_curve(gmx_radial_distribution_histogram_t* pr,
376 double start_q,
377 double end_q,
378 double q_step)
380 gmx_static_structurefactor_t* sq = nullptr;
381 int i, j;
382 /* init data */
383 snew(sq, 1);
384 sq->qn = static_cast<int>(std::floor((end_q - start_q) / q_step));
385 snew(sq->q, sq->qn);
386 snew(sq->s, sq->qn);
387 for (i = 0; i < sq->qn; i++)
389 sq->q[i] = start_q + i * q_step;
392 if (start_q == 0.0)
394 sq->s[0] = 1.0;
395 for (i = 1; i < sq->qn; i++)
397 for (j = 0; j < pr->grn; j++)
399 sq->s[i] += (pr->gr[j] / pr->r[j]) * std::sin(sq->q[i] * pr->r[j]);
401 sq->s[i] /= sq->q[i];
404 else
406 for (i = 0; i < sq->qn; i++)
408 for (j = 0; j < pr->grn; j++)
410 sq->s[i] += (pr->gr[j] / pr->r[j]) * std::sin(sq->q[i] * pr->r[j]);
412 sq->s[i] /= sq->q[i];
416 return sq;