3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
42 #include "gmx_random.h"
50 void check_binwidth(real binwidth)
52 real smallest_bin = 0.1;
53 if (binwidth < smallest_bin)
55 gmx_fatal(FARGS, "Binwidth shouldnt be smaller then smallest bond length (H-H bond ~0.1nm) in a box");
59 void check_mcover(real mcover)
63 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
65 else if ((mcover < 0)&(mcover != -1))
67 gmx_fatal(FARGS, "mcover should be -1 or (0,1]");
75 void normalize_probability(int n, double *a)
79 for (i = 0; i < n; i++)
83 for (i = 0; i < n; i++)
89 gmx_neutron_atomic_structurefactors_t *gmx_neutronstructurefactors_init(const char *datfn)
91 /* read nsfactor.dat */
99 gmx_neutron_atomic_structurefactors_t *gnsf;
103 /* allocate memory for structure */
105 snew(gnsf->atomnm, nralloc);
106 snew(gnsf->p, nralloc);
107 snew(gnsf->n, nralloc);
108 snew(gnsf->slength, nralloc);
110 gnsf->nratoms = line_no;
112 while (get_a_line(fp, line, STRLEN))
115 if (sscanf(line, "%s %d %d %lf", atomnm, &p, &n, &slength) == 4)
117 gnsf->atomnm[i] = strdup(atomnm);
120 gnsf->slength[i] = slength;
122 gnsf->nratoms = line_no;
123 if (line_no == nralloc)
126 srenew(gnsf->atomnm, nralloc);
127 srenew(gnsf->p, nralloc);
128 srenew(gnsf->n, nralloc);
129 srenew(gnsf->slength, nralloc);
134 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
138 srenew(gnsf->atomnm, gnsf->nratoms);
139 srenew(gnsf->p, gnsf->nratoms);
140 srenew(gnsf->n, gnsf->nratoms);
141 srenew(gnsf->slength, gnsf->nratoms);
145 return (gmx_neutron_atomic_structurefactors_t *) gnsf;
148 gmx_sans_t *gmx_sans_init (t_topology *top, gmx_neutron_atomic_structurefactors_t *gnsf)
150 gmx_sans_t *gsans = NULL;
152 /* Try to assing scattering length from nsfactor.dat */
154 snew(gsans->slength, top->atoms.nr);
155 /* copy topology data */
157 for (i = 0; i < top->atoms.nr; i++)
159 for (j = 0; j < gnsf->nratoms; j++)
161 if (top->atoms.atom[i].atomnumber == gnsf->p[j])
163 /* we need special case for H and D */
164 if (top->atoms.atom[i].atomnumber == 1)
166 if (top->atoms.atom[i].m == 1.008000)
168 gsans->slength[i] = gnsf->slength[0];
172 gsans->slength[i] = gnsf->slength[1];
177 gsans->slength[i] = gnsf->slength[j];
183 return (gmx_sans_t *) gsans;
186 gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
198 gmx_radial_distribution_histogram_t *pr = NULL;
206 gmx_rng_t *trng = NULL;
208 gmx_large_int_t mc = 0, max;
209 gmx_rng_t rng = NULL;
211 /* allocate memory for pr */
213 /* set some fields */
214 pr->binwidth = binwidth;
217 * create max dist rvec
218 * dist = box[xx] + box[yy] + box[zz]
220 rvec_add(box[XX], box[YY], dist);
221 rvec_add(box[ZZ], dist, dist);
225 pr->grn = (int)floor(rmax/pr->binwidth)+1;
226 rmax = pr->grn*pr->binwidth;
228 snew(pr->gr, pr->grn);
232 /* Special case for setting automaticaly number of mc iterations to 1% of total number of direct iterations */
235 max = (gmx_large_int_t)floor(0.5*0.01*isize*(isize-1));
239 max = (gmx_large_int_t)floor(0.5*mcover*isize*(isize-1));
241 rng = gmx_rng_init(seed);
243 nthreads = gmx_omp_get_max_threads();
245 snew(trng, nthreads);
246 for (i = 0; i < nthreads; i++)
248 snew(tgr[i], pr->grn);
249 trng[i] = gmx_rng_init(gmx_rng_uniform_uint32(rng));
251 #pragma omp parallel shared(tgr,trng,mc) private(tid,i,j)
253 tid = gmx_omp_get_thread_num();
254 /* now starting parallel threads */
256 for (mc = 0; mc < max; mc++)
258 i = (int)floor(gmx_rng_uniform_real(trng[tid])*isize);
259 j = (int)floor(gmx_rng_uniform_real(trng[tid])*isize);
262 tgr[tid][(int)floor(sqrt(distance2(x[index[i]], x[index[j]]))/binwidth)] += gsans->slength[index[i]]*gsans->slength[index[j]];
266 /* collecting data from threads */
267 for (i = 0; i < pr->grn; i++)
269 for (j = 0; j < nthreads; j++)
271 pr->gr[i] += tgr[j][i];
274 /* freeing memory for tgr and destroying trng */
275 for (i = 0; i < nthreads; i++)
278 gmx_rng_destroy(trng[i]);
283 for (mc = 0; mc < max; mc++)
285 i = (int)floor(gmx_rng_uniform_real(rng)*isize);
286 j = (int)floor(gmx_rng_uniform_real(rng)*isize);
289 pr->gr[(int)floor(sqrt(distance2(x[index[i]], x[index[j]]))/binwidth)] += gsans->slength[index[i]]*gsans->slength[index[j]];
293 gmx_rng_destroy(rng);
298 nthreads = gmx_omp_get_max_threads();
299 /* Allocating memory for tgr arrays */
301 for (i = 0; i < nthreads; i++)
303 snew(tgr[i], pr->grn);
305 #pragma omp parallel shared(tgr) private(tid,i,j)
307 tid = gmx_omp_get_thread_num();
308 /* starting parallel threads */
310 for (i = 0; i < isize; i++)
312 for (j = 0; j < i; j++)
314 tgr[tid][(int)floor(sqrt(distance2(x[index[i]], x[index[j]]))/binwidth)] += gsans->slength[index[i]]*gsans->slength[index[j]];
318 /* collecating data for pr->gr */
319 for (i = 0; i < pr->grn; i++)
321 for (j = 0; j < nthreads; j++)
323 pr->gr[i] += tgr[j][i];
326 /* freeing memory for tgr */
327 for (i = 0; i < nthreads; i++)
333 for (i = 0; i < isize; i++)
335 for (j = 0; j < i; j++)
337 pr->gr[(int)floor(sqrt(distance2(x[index[i]], x[index[j]]))/binwidth)] += gsans->slength[index[i]]*gsans->slength[index[j]];
343 /* normalize if needed */
346 normalize_probability(pr->grn, pr->gr);
349 snew(pr->r, pr->grn);
350 for (i = 0; i < pr->grn; i++)
352 pr->r[i] = (pr->binwidth*i+pr->binwidth*0.5);
355 return (gmx_radial_distribution_histogram_t *) pr;
358 gmx_static_structurefactor_t *convert_histogram_to_intensity_curve (gmx_radial_distribution_histogram_t *pr, double start_q, double end_q, double q_step)
360 gmx_static_structurefactor_t *sq = NULL;
364 sq->qn = (int)floor((end_q-start_q)/q_step);
367 for (i = 0; i < sq->qn; i++)
369 sq->q[i] = start_q+i*q_step;
375 for (i = 1; i < sq->qn; i++)
377 for (j = 0; j < pr->grn; j++)
379 sq->s[i] += (pr->gr[j]/pr->r[j])*sin(sq->q[i]*pr->r[j]);
381 sq->s[i] /= sq->q[i];
386 for (i = 0; i < sq->qn; i++)
388 for (j = 0; j < pr->grn; j++)
390 sq->s[i] += (pr->gr[j]/pr->r[j])*sin(sq->q[i]*pr->r[j]);
392 sq->s[i] /= sq->q[i];
396 return (gmx_static_structurefactor_t *) sq;