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) {
51 real smallest_bin=0.1;
52 if (binwidth<smallest_bin)
53 gmx_fatal(FARGS,"Binwidth shouldnt be smaller then smallest bond length (H-H bond ~0.1nm) in a box");
56 void check_mcover(real mcover) {
58 gmx_fatal(FARGS,"mcover should be -1 or (0,1]");
59 } else if ((mcover<0)&(mcover != -1)) {
60 gmx_fatal(FARGS,"mcover should be -1 or (0,1]");
66 void normalize_probability(int n,double *a){
69 for (i=0;i<n;i++) norm +=a[i];
70 for (i=0;i<n;i++) a[i]/=norm;
73 gmx_nentron_atomic_structurefactors_t *gmx_neutronstructurefactors_init(const char *datfn) {
74 /* read nsfactor.dat */
82 gmx_nentron_atomic_structurefactors_t *gnsf;
86 /* allocate memory for structure */
88 snew(gnsf->atomnm,nralloc);
89 snew(gnsf->p,nralloc);
90 snew(gnsf->n,nralloc);
91 snew(gnsf->slength,nralloc);
93 gnsf->nratoms=line_no;
95 while(get_a_line(fp,line,STRLEN)) {
97 if (sscanf(line,"%s %d %d %lf",atomnm,&p,&n,&slength) == 4) {
98 gnsf->atomnm[i]=strdup(atomnm);
101 gnsf->slength[i]=slength;
103 gnsf->nratoms=line_no;
104 if (line_no==nralloc){
106 srenew(gnsf->atomnm,nralloc);
107 srenew(gnsf->p,nralloc);
108 srenew(gnsf->n,nralloc);
109 srenew(gnsf->slength,nralloc);
112 fprintf(stderr,"WARNING: Error in file %s at line %d ignored\n",
115 srenew(gnsf->atomnm,gnsf->nratoms);
116 srenew(gnsf->p,gnsf->nratoms);
117 srenew(gnsf->n,gnsf->nratoms);
118 srenew(gnsf->slength,gnsf->nratoms);
122 return (gmx_nentron_atomic_structurefactors_t *) gnsf;
125 gmx_sans_t *gmx_sans_init (t_topology *top, gmx_nentron_atomic_structurefactors_t *gnsf) {
126 gmx_sans_t *gsans=NULL;
128 /* Try to assing scattering length from nsfactor.dat */
130 snew(gsans->slength,top->atoms.nr);
131 /* copy topology data */
133 for(i=0;i<top->atoms.nr;i++) {
134 for(j=0;j<gnsf->nratoms;j++) {
135 if(top->atoms.atom[i].atomnumber == gnsf->p[j]) {
136 /* we need special case for H and D */
137 if(top->atoms.atom[i].atomnumber == 1) {
138 if(top->atoms.atom[i].m == 1.008000) {
139 gsans->slength[i] = gnsf->slength[0];
141 gsans->slength[i] = gnsf->slength[1];
143 gsans->slength[i] = gnsf->slength[j];
148 return (gmx_sans_t *) gsans;
151 gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
161 gmx_radial_distribution_histogram_t *pr=NULL;
169 gmx_rng_t *trng=NULL;
171 gmx_large_int_t mc=0,max;
174 /* allocate memory for pr */
176 /* set some fields */
177 pr->binwidth=binwidth;
180 * create max dist rvec
181 * dist = box[xx] + box[yy] + box[zz]
183 rvec_add(box[XX],box[YY],dist);
184 rvec_add(box[ZZ],dist,dist);
188 pr->grn=(int)floor(rmax/pr->binwidth)+1;
189 rmax=pr->grn*pr->binwidth;
191 snew(pr->gr,pr->grn);
194 /* Special case for setting automaticaly number of mc iterations to 1% of total number of direct iterations */
196 max=(gmx_large_int_t)floor(0.5*0.01*isize*(isize-1));
198 max=(gmx_large_int_t)floor(0.5*mcover*isize*(isize-1));
200 rng=gmx_rng_init(seed);
202 nthreads = gmx_omp_get_max_threads();
205 for(i=0;i<nthreads;i++){
206 snew(tgr[i],pr->grn);
207 trng[i]=gmx_rng_init(gmx_rng_uniform_uint32(rng));
209 #pragma omp parallel shared(tgr,trng,mc) private(tid,i,j)
211 tid = gmx_omp_get_thread_num();
212 /* now starting parallel threads */
214 for(mc=0;mc<max;mc++) {
215 i=(int)floor(gmx_rng_uniform_real(trng[tid])*isize);
216 j=(int)floor(gmx_rng_uniform_real(trng[tid])*isize);
218 tgr[tid][(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
222 /* collecting data from threads */
223 for(i=0;i<pr->grn;i++) {
224 for(j=0;j<nthreads;j++) {
225 pr->gr[i] += tgr[j][i];
228 /* freeing memory for tgr and destroying trng */
229 for(i=0;i<nthreads;i++) {
231 gmx_rng_destroy(trng[i]);
236 for(mc=0;mc<max;mc++) {
237 i=(int)floor(gmx_rng_uniform_real(rng)*isize);
238 j=(int)floor(gmx_rng_uniform_real(rng)*isize);
240 pr->gr[(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
243 gmx_rng_destroy(rng);
246 nthreads = gmx_omp_get_max_threads();
247 /* Allocating memory for tgr arrays */
249 for(i=0;i<nthreads;i++) {
250 snew(tgr[i],pr->grn);
252 #pragma omp parallel shared(tgr) private(tid,i,j)
254 tid = gmx_omp_get_thread_num();
255 /* starting parallel threads */
257 for(i=0;i<isize;i++) {
259 tgr[tid][(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
263 /* collecating data for pr->gr */
264 for(i=0;i<pr->grn;i++) {
265 for(j=0;j<nthreads;j++) {
266 pr->gr[i] += tgr[j][i];
269 /* freeing memory for tgr */
270 for(i=0;i<nthreads;i++) {
275 for(i=0;i<isize;i++) {
277 pr->gr[(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
284 normalize_probability(pr->grn,pr->gr);
286 for(i=0;i<pr->grn;i++)
287 pr->r[i]=(pr->binwidth*i+pr->binwidth*0.5);
289 return (gmx_radial_distribution_histogram_t *) pr;
292 gmx_static_structurefator_t *convert_histogram_to_intensity_curve (gmx_radial_distribution_histogram_t *pr, double start_q, double end_q, double q_step) {
293 gmx_static_structurefator_t *sq=NULL;
297 sq->qn=(int)floor((end_q-start_q)/q_step);
300 for(i=0;i<sq->qn;i++)
301 sq->q[i]=start_q+i*q_step;
305 for(i=1;i<sq->qn;i++) {
306 for(j=0;j<pr->grn;j++)
307 sq->s[i]+=(pr->gr[j]/pr->r[j])*sin(sq->q[i]*pr->r[j]);
308 sq->s[i] /= sq->q[i];
311 for(i=0;i<sq->qn;i++) {
312 for(j=0;j<pr->grn;j++)
313 sq->s[i]+=(pr->gr[j]/pr->r[j])*sin(sq->q[i]*pr->r[j]);
314 sq->s[i] /= sq->q[i];
318 return (gmx_static_structurefator_t *) sq;