Remove gmx custom fixed int (e.g. gmx_int64_t) types
[gromacs.git] / src / gromacs / gmxana / nsfactor.cpp
blob92e216526ea61d87ce244c4b974eedfd00dfd3cf
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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.
35 #include "gmxpre.h"
37 #include "nsfactor.h"
39 #include "config.h"
41 #include <cmath>
42 #include <cstring>
44 #include "gromacs/math/vec.h"
45 #include "gromacs/random/threefry.h"
46 #include "gromacs/random/uniformintdistribution.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/gmxomp.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/utility/strdb.h"
56 void check_binwidth(real binwidth)
58 real smallest_bin = 0.1;
59 if (binwidth < smallest_bin)
61 gmx_fatal(FARGS, "Binwidth shouldn't be smaller then smallest bond length (H-H bond ~0.1nm) in a box");
65 void check_mcover(real mcover)
67 if (mcover > 1.0)
69 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
71 else if ((mcover < 0)&(mcover != -1))
73 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
75 else
77 return;
81 void normalize_probability(int n, double *a)
83 int i;
84 double norm = 0.0;
85 for (i = 0; i < n; i++)
87 norm += a[i];
89 for (i = 0; i < n; i++)
91 a[i] /= norm;
95 gmx_neutron_atomic_structurefactors_t *gmx_neutronstructurefactors_init(const char *datfn)
97 /* read nsfactor.dat */
98 FILE *fp;
99 char line[STRLEN];
100 int nralloc = 10;
101 int n, p;
102 int i, line_no;
103 char atomnm[8];
104 double slength;
105 gmx_neutron_atomic_structurefactors_t *gnsf;
107 fp = libopen(datfn);
108 line_no = 0;
109 /* allocate memory for structure */
110 snew(gnsf, nralloc);
111 snew(gnsf->atomnm, nralloc);
112 snew(gnsf->p, nralloc);
113 snew(gnsf->n, nralloc);
114 snew(gnsf->slength, nralloc);
116 gnsf->nratoms = line_no;
118 while (get_a_line(fp, line, STRLEN))
120 i = line_no;
121 if (sscanf(line, "%s %d %d %lf", atomnm, &p, &n, &slength) == 4)
123 gnsf->atomnm[i] = gmx_strdup(atomnm);
124 gnsf->n[i] = n;
125 gnsf->p[i] = p;
126 gnsf->slength[i] = slength;
127 line_no++;
128 gnsf->nratoms = line_no;
129 if (line_no == nralloc)
131 nralloc++;
132 srenew(gnsf->atomnm, nralloc);
133 srenew(gnsf->p, nralloc);
134 srenew(gnsf->n, nralloc);
135 srenew(gnsf->slength, nralloc);
138 else
140 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
141 datfn, line_no);
144 srenew(gnsf->atomnm, gnsf->nratoms);
145 srenew(gnsf->p, gnsf->nratoms);
146 srenew(gnsf->n, gnsf->nratoms);
147 srenew(gnsf->slength, gnsf->nratoms);
149 fclose(fp);
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 (
193 gmx_sans_t *gsans,
194 rvec *x,
195 matrix box,
196 const int *index,
197 int isize,
198 double binwidth,
199 gmx_bool bMC,
200 gmx_bool bNORM,
201 real mcover,
202 unsigned int seed)
204 gmx_radial_distribution_histogram_t *pr = nullptr;
205 rvec dist;
206 double rmax;
207 int i, j;
208 #if GMX_OPENMP
209 double **tgr;
210 int tid;
211 int nthreads;
212 gmx::DefaultRandomEngine *trng = nullptr;
213 #endif
214 int64_t mc = 0, mc_max;
215 gmx::DefaultRandomEngine rng(seed);
217 /* allocate memory for pr */
218 snew(pr, 1);
219 /* set some fields */
220 pr->binwidth = binwidth;
223 * create max dist rvec
224 * dist = box[xx] + box[yy] + box[zz]
226 rvec_add(box[XX], box[YY], dist);
227 rvec_add(box[ZZ], dist, dist);
229 rmax = norm(dist);
231 pr->grn = static_cast<int>(std::floor(rmax/pr->binwidth)+1);
233 snew(pr->gr, pr->grn);
235 if (bMC)
237 /* Special case for setting automaticaly number of mc iterations to 1% of total number of direct iterations */
238 if (mcover == -1)
240 mc_max = static_cast<int64_t>(std::floor(0.5*0.01*isize*(isize-1)));
242 else
244 mc_max = static_cast<int64_t>(std::floor(0.5*mcover*isize*(isize-1)));
246 #if GMX_OPENMP
247 nthreads = gmx_omp_get_max_threads();
248 snew(tgr, nthreads);
249 trng = new gmx::DefaultRandomEngine[nthreads];
250 for (i = 0; i < nthreads; i++)
252 snew(tgr[i], pr->grn);
253 trng[i].seed(rng());
255 #pragma omp parallel shared(tgr,trng,mc) private(tid,i,j)
257 gmx::UniformIntDistribution<int> tdist(0, isize-1);
258 tid = gmx_omp_get_thread_num();
259 /* now starting parallel threads */
260 #pragma omp for
261 for (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))] += gsans->slength[index[i]]*gsans->slength[index[j]];
272 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
275 /* collecting data from threads */
276 for (i = 0; i < pr->grn; i++)
278 for (j = 0; j < nthreads; j++)
280 pr->gr[i] += tgr[j][i];
283 /* freeing memory for tgr and destroying trng */
284 for (i = 0; i < nthreads; i++)
286 sfree(tgr[i]);
288 sfree(tgr);
289 delete[] trng;
290 #else
291 gmx::UniformIntDistribution<int> dist(0, isize-1);
292 for (mc = 0; mc < mc_max; mc++)
294 i = dist(rng); // [0,isize-1]
295 j = dist(rng); // [0,isize-1]
296 if (i != j)
298 pr->gr[static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]]))/binwidth))] += gsans->slength[index[i]]*gsans->slength[index[j]];
301 #endif
303 else
305 #if GMX_OPENMP
306 nthreads = gmx_omp_get_max_threads();
307 /* Allocating memory for tgr arrays */
308 snew(tgr, nthreads);
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 */
317 #pragma omp for
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))] += gsans->slength[index[i]]*gsans->slength[index[j]];
327 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
330 /* collecating data for pr->gr */
331 for (i = 0; i < pr->grn; i++)
333 for (j = 0; j < nthreads; j++)
335 pr->gr[i] += tgr[j][i];
338 /* freeing memory for tgr */
339 for (i = 0; i < nthreads; i++)
341 sfree(tgr[i]);
343 sfree(tgr);
344 #else
345 for (i = 0; i < isize; i++)
347 for (j = 0; j < i; j++)
349 pr->gr[static_cast<int>(std::floor(std::sqrt(distance2(x[index[i]], x[index[j]]))/binwidth))] += gsans->slength[index[i]]*gsans->slength[index[j]];
352 #endif
355 /* normalize if needed */
356 if (bNORM)
358 normalize_probability(pr->grn, pr->gr);
361 snew(pr->r, pr->grn);
362 for (i = 0; i < pr->grn; i++)
364 pr->r[i] = (pr->binwidth*i+pr->binwidth*0.5);
367 return pr;
370 gmx_static_structurefactor_t *convert_histogram_to_intensity_curve (gmx_radial_distribution_histogram_t *pr, double start_q, double end_q, double q_step)
372 gmx_static_structurefactor_t *sq = nullptr;
373 int i, j;
374 /* init data */
375 snew(sq, 1);
376 sq->qn = static_cast<int>(std::floor((end_q-start_q)/q_step));
377 snew(sq->q, sq->qn);
378 snew(sq->s, sq->qn);
379 for (i = 0; i < sq->qn; i++)
381 sq->q[i] = start_q+i*q_step;
384 if (start_q == 0.0)
386 sq->s[0] = 1.0;
387 for (i = 1; i < sq->qn; i++)
389 for (j = 0; j < pr->grn; j++)
391 sq->s[i] += (pr->gr[j]/pr->r[j])*std::sin(sq->q[i]*pr->r[j]);
393 sq->s[i] /= sq->q[i];
396 else
398 for (i = 0; i < sq->qn; i++)
400 for (j = 0; j < pr->grn; j++)
402 sq->s[i] += (pr->gr[j]/pr->r[j])*std::sin(sq->q[i]*pr->r[j]);
404 sq->s[i] /= sq->q[i];
408 return sq;