2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gmx_random.h"
53 void check_binwidth(real binwidth) {
54 real smallest_bin=0.1;
55 if (binwidth<smallest_bin)
56 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) {
61 gmx_fatal(FARGS,"mcover should be -1 or (0,1]");
62 } else if ((mcover<0)&(mcover != -1)) {
63 gmx_fatal(FARGS,"mcover should be -1 or (0,1]");
69 void normalize_probability(int n,double *a) {
72 for (i=0;i<n;i++) norm +=a[i];
73 for (i=0;i<n;i++) a[i]/=norm;
76 gmx_neutron_atomic_structurefactors_t *gmx_neutronstructurefactors_init(const char *datfn) {
77 /* read nsfactor.dat */
85 gmx_neutron_atomic_structurefactors_t *gnsf;
89 /* allocate memory for structure */
91 snew(gnsf->atomnm,nralloc);
92 snew(gnsf->p,nralloc);
93 snew(gnsf->n,nralloc);
94 snew(gnsf->slength,nralloc);
96 gnsf->nratoms=line_no;
98 while(get_a_line(fp,line,STRLEN)) {
100 if (sscanf(line,"%s %d %d %lf",atomnm,&p,&n,&slength) == 4) {
101 gnsf->atomnm[i]=strdup(atomnm);
104 gnsf->slength[i]=slength;
106 gnsf->nratoms=line_no;
107 if (line_no==nralloc){
109 srenew(gnsf->atomnm,nralloc);
110 srenew(gnsf->p,nralloc);
111 srenew(gnsf->n,nralloc);
112 srenew(gnsf->slength,nralloc);
115 fprintf(stderr,"WARNING: Error in file %s at line %d ignored\n",
118 srenew(gnsf->atomnm,gnsf->nratoms);
119 srenew(gnsf->p,gnsf->nratoms);
120 srenew(gnsf->n,gnsf->nratoms);
121 srenew(gnsf->slength,gnsf->nratoms);
125 return (gmx_neutron_atomic_structurefactors_t *) gnsf;
128 gmx_sans_t *gmx_sans_init (t_topology *top, gmx_neutron_atomic_structurefactors_t *gnsf) {
129 gmx_sans_t *gsans=NULL;
131 /* Try to assing scattering length from nsfactor.dat */
133 snew(gsans->slength,top->atoms.nr);
134 /* copy topology data */
136 for(i=0;i<top->atoms.nr;i++) {
137 for(j=0;j<gnsf->nratoms;j++) {
138 if(top->atoms.atom[i].atomnumber == gnsf->p[j]) {
139 /* we need special case for H and D */
140 if(top->atoms.atom[i].atomnumber == 1) {
141 if(top->atoms.atom[i].m == 1.008000) {
142 gsans->slength[i] = gnsf->slength[0];
144 gsans->slength[i] = gnsf->slength[1];
146 gsans->slength[i] = gnsf->slength[j];
151 return (gmx_sans_t *) gsans;
154 gmx_radial_distribution_histogram_t *calc_radial_distribution_histogram (
165 gmx_radial_distribution_histogram_t *pr=NULL;
173 gmx_rng_t *trng=NULL;
175 gmx_large_int_t mc=0,max;
178 /* allocate memory for pr */
180 /* set some fields */
181 pr->binwidth=binwidth;
184 * create max dist rvec
185 * dist = box[xx] + box[yy] + box[zz]
187 rvec_add(box[XX],box[YY],dist);
188 rvec_add(box[ZZ],dist,dist);
192 pr->grn=(int)floor(rmax/pr->binwidth)+1;
193 rmax=pr->grn*pr->binwidth;
195 snew(pr->gr,pr->grn);
198 /* Special case for setting automaticaly number of mc iterations to 1% of total number of direct iterations */
200 max=(gmx_large_int_t)floor(0.5*0.01*isize*(isize-1));
202 max=(gmx_large_int_t)floor(0.5*mcover*isize*(isize-1));
204 rng=gmx_rng_init(seed);
206 nthreads = gmx_omp_get_max_threads();
209 for(i=0;i<nthreads;i++){
210 snew(tgr[i],pr->grn);
211 trng[i]=gmx_rng_init(gmx_rng_uniform_uint32(rng));
213 #pragma omp parallel shared(tgr,trng,mc) private(tid,i,j)
215 tid = gmx_omp_get_thread_num();
216 /* now starting parallel threads */
218 for(mc=0;mc<max;mc++) {
219 i=(int)floor(gmx_rng_uniform_real(trng[tid])*isize);
220 j=(int)floor(gmx_rng_uniform_real(trng[tid])*isize);
222 tgr[tid][(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
226 /* collecting data from threads */
227 for(i=0;i<pr->grn;i++) {
228 for(j=0;j<nthreads;j++) {
229 pr->gr[i] += tgr[j][i];
232 /* freeing memory for tgr and destroying trng */
233 for(i=0;i<nthreads;i++) {
235 gmx_rng_destroy(trng[i]);
240 for(mc=0;mc<max;mc++) {
241 i=(int)floor(gmx_rng_uniform_real(rng)*isize);
242 j=(int)floor(gmx_rng_uniform_real(rng)*isize);
244 pr->gr[(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
247 gmx_rng_destroy(rng);
250 nthreads = gmx_omp_get_max_threads();
251 /* Allocating memory for tgr arrays */
253 for(i=0;i<nthreads;i++) {
254 snew(tgr[i],pr->grn);
256 #pragma omp parallel shared(tgr) private(tid,i,j)
258 tid = gmx_omp_get_thread_num();
259 /* starting parallel threads */
261 for(i=0;i<isize;i++) {
263 tgr[tid][(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
267 /* collecating data for pr->gr */
268 for(i=0;i<pr->grn;i++) {
269 for(j=0;j<nthreads;j++) {
270 pr->gr[i] += tgr[j][i];
273 /* freeing memory for tgr */
274 for(i=0;i<nthreads;i++) {
279 for(i=0;i<isize;i++) {
281 pr->gr[(int)floor(sqrt(distance2(x[index[i]],x[index[j]]))/binwidth)]+=gsans->slength[index[i]]*gsans->slength[index[j]];
287 /* normalize if needed */
289 normalize_probability(pr->grn,pr->gr);
293 for(i=0;i<pr->grn;i++)
294 pr->r[i]=(pr->binwidth*i+pr->binwidth*0.5);
296 return (gmx_radial_distribution_histogram_t *) pr;
299 gmx_static_structurefactor_t *convert_histogram_to_intensity_curve (gmx_radial_distribution_histogram_t *pr, double start_q, double end_q, double q_step) {
300 gmx_static_structurefactor_t *sq=NULL;
304 sq->qn=(int)floor((end_q-start_q)/q_step);
307 for(i=0;i<sq->qn;i++)
308 sq->q[i]=start_q+i*q_step;
312 for(i=1;i<sq->qn;i++) {
313 for(j=0;j<pr->grn;j++)
314 sq->s[i]+=(pr->gr[j]/pr->r[j])*sin(sq->q[i]*pr->r[j]);
315 sq->s[i] /= sq->q[i];
318 for(i=0;i<sq->qn;i++) {
319 for(j=0;j<pr->grn;j++)
320 sq->s[i]+=(pr->gr[j]/pr->r[j])*sin(sq->q[i]*pr->r[j]);
321 sq->s[i] /= sq->q[i];
325 return (gmx_static_structurefactor_t *) sq;