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 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gromacs/utility/smalloc.h"
43 #include "gromacs/fileio/futil.h"
44 #include "gromacs/math/utilities.h"
45 #include "gmx_fatal.h"
49 #include "gromacs/fileio/strdb.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/fileio/gmxfio.h"
56 #include "gromacs/fileio/matio.h"
61 typedef struct gmx_structurefactors {
63 int *p; /* proton number */
64 int *n; /* neutron number */
65 /* Parameters for the Cromer Mann fit */
66 real **a; /* parameter a */
67 real **b; /* parameter b */
68 real *c; /* parameter c */
69 char **atomnm; /* atomname */
71 } gmx_structurefactors;
73 typedef struct reduced_atom{
79 typedef struct structure_factor
93 extern int * create_indexed_atom_type (reduced_atom_t * atm, int size)
96 * create an index of the atom types found in a group
97 * i.e.: for water index_atp[0]=type_number_of_O and
98 * index_atp[1]=type_number_of_H
100 * the last element is set to 0
102 int *index_atp, i, i_tmp, j;
104 reduced_atom *att = (reduced_atom *)atm;
108 index_atp[0] = att[0].t;
109 for (i = 1; i < size; i++)
111 for (j = 0; j < i_tmp; j++)
113 if (att[i].t == index_atp[j])
118 if (j == i_tmp) /* i.e. no indexed atom type is == to atm[i].t */
121 srenew (index_atp, i_tmp * sizeof (int));
122 index_atp[i_tmp - 1] = att[i].t;
126 srenew (index_atp, i_tmp * sizeof (int));
127 index_atp[i_tmp - 1] = 0;
133 extern t_complex *** rc_tensor_allocation(int x, int y, int z)
138 t = (t_complex ***)calloc(x, sizeof(t_complex**));
141 exit(fprintf(stderr, "\nallocation error"));
143 t[0] = (t_complex **)calloc(x*y, sizeof(t_complex*));
146 exit(fprintf(stderr, "\nallocation error"));
148 t[0][0] = (t_complex *)calloc(x*y*z, sizeof(t_complex));
151 exit(fprintf(stderr, "\nallocation error"));
154 for (j = 1; j < y; j++)
156 t[0][j] = t[0][j-1] + z;
158 for (i = 1; i < x; i++)
161 t[i][0] = t[i-1][0] + y*z;
162 for (j = 1; j < y; j++)
164 t[i][j] = t[i][j-1] + z;
171 extern void compute_structure_factor (structure_factor_t * sft, matrix box,
172 reduced_atom_t * red, int isize, real start_q,
173 real end_q, int group, real **sf_table)
175 structure_factor *sf = (structure_factor *)sft;
176 reduced_atom *redt = (reduced_atom *)red;
180 real kdotx, asf, kx, ky, kz, krr;
181 int kr, maxkx, maxky, maxkz, i, j, k, p, *counter;
184 k_factor[XX] = 2 * M_PI / box[XX][XX];
185 k_factor[YY] = 2 * M_PI / box[YY][YY];
186 k_factor[ZZ] = 2 * M_PI / box[ZZ][ZZ];
188 maxkx = (int) (end_q / k_factor[XX] + 0.5);
189 maxky = (int) (end_q / k_factor[YY] + 0.5);
190 maxkz = (int) (end_q / k_factor[ZZ] + 0.5);
192 snew (counter, sf->n_angles);
194 tmpSF = rc_tensor_allocation(maxkx, maxky, maxkz);
197 * compute real and imaginary part of the structure factor for every
200 fprintf(stderr, "\n");
201 for (i = 0; i < maxkx; i++)
203 fprintf (stderr, "\rdone %3.1f%% ", (double)(100.0*(i+1))/maxkx);
204 kx = i * k_factor[XX];
205 for (j = 0; j < maxky; j++)
207 ky = j * k_factor[YY];
208 for (k = 0; k < maxkz; k++)
210 if (i != 0 || j != 0 || k != 0)
212 kz = k * k_factor[ZZ];
213 krr = sqrt (sqr (kx) + sqr (ky) + sqr (kz));
214 if (krr >= start_q && krr <= end_q)
216 kr = (int) (krr/sf->ref_k + 0.5);
217 if (kr < sf->n_angles)
219 counter[kr]++; /* will be used for the copmutation
221 for (p = 0; p < isize; p++)
223 asf = sf_table[redt[p].t][kr];
225 kdotx = kx * redt[p].x[XX] +
226 ky * redt[p].x[YY] + kz * redt[p].x[ZZ];
228 tmpSF[i][j][k].re += cos (kdotx) * asf;
229 tmpSF[i][j][k].im += sin (kdotx) * asf;
236 } /* end loop on i */
238 * compute the square modulus of the structure factor, averaging on the surface
239 * kx*kx + ky*ky + kz*kz = krr*krr
240 * note that this is correct only for a (on the macroscopic scale)
243 for (i = 0; i < maxkx; i++)
245 kx = i * k_factor[XX]; for (j = 0; j < maxky; j++)
247 ky = j * k_factor[YY]; for (k = 0; k < maxkz; k++)
249 kz = k * k_factor[ZZ]; krr = sqrt (sqr (kx) + sqr (ky)
250 + sqr (kz)); if (krr >= start_q && krr <= end_q)
252 kr = (int) (krr / sf->ref_k + 0.5);
253 if (kr < sf->n_angles && counter[kr] != 0)
256 (sqr (tmpSF[i][j][k].re) +
257 sqr (tmpSF[i][j][k].im))/ counter[kr];
263 sfree (counter); free(tmpSF[0][0]); free(tmpSF[0]); free(tmpSF);
267 extern gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn)
270 /* Read the database for the structure factor of the different atoms */
274 gmx_structurefactors *gsf;
275 double a1, a2, a3, a4, b1, b2, b3, b4, c;
285 snew(gsf->atomnm, nralloc);
286 snew(gsf->a, nralloc);
287 snew(gsf->b, nralloc);
288 snew(gsf->c, nralloc);
289 snew(gsf->p, nralloc);
291 gsf->nratoms = line_no;
292 while (get_a_line(fp, line, STRLEN))
295 if (sscanf(line, "%s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf",
296 atomn, &p, &a1, &a2, &a3, &a4, &b1, &b2, &b3, &b4, &c) == 11)
298 gsf->atomnm[i] = strdup(atomn);
312 gsf->nratoms = line_no;
313 if (line_no == nralloc)
316 srenew(gsf->atomnm, nralloc);
317 srenew(gsf->a, nralloc);
318 srenew(gsf->b, nralloc);
319 srenew(gsf->c, nralloc);
320 srenew(gsf->p, nralloc);
325 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
330 srenew(gsf->atomnm, gsf->nratoms);
331 srenew(gsf->a, gsf->nratoms);
332 srenew(gsf->b, gsf->nratoms);
333 srenew(gsf->c, gsf->nratoms);
334 srenew(gsf->p, gsf->nratoms);
338 return (gmx_structurefactors_t *) gsf;
343 extern void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, atom_id * index,
344 int isize, t_topology * top, gmx_bool flag, gmx_structurefactors_t *gsf)
345 /* given the group's index, return the (continuous) array of atoms */
349 reduced_atom *pos = (reduced_atom *)positions;
353 for (i = 0; i < isize; i++)
356 return_atom_type (*(top->atoms.atomname[index[i]]), gsf);
359 for (i = 0; i < isize; i++)
361 copy_rvec (fr->x[index[i]], pos[i].x);
364 positions = (reduced_atom_t *)pos;
368 extern int return_atom_type (const char *name, gmx_structurefactors_t *gsf)
375 { "CH1", 1 }, { "CH2", 2 }, { "CH3", 3 },
376 { "CS1", 1 }, { "CS2", 2 }, { "CS3", 3 },
377 { "CP1", 1 }, { "CP2", 2 }, { "CP3", 3 }
385 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
387 NCMT = gsft->nratoms;
391 for (i = 0; (i < asize(uh)); i++)
393 if (strcmp(name, uh[i].name) == 0)
395 return NCMT-1+uh[i].nh;
399 for (i = 0; (i < NCMT); i++)
401 if (strncmp (name, gsft->atomnm[i], strlen(gsft->atomnm[i])) == 0)
410 gmx_fatal(FARGS, "\nError: atom (%s) not in list (%d types checked)!\n",
416 for (i = 0; i < cnt; i++)
418 if (strlen(gsft->atomnm[tndx[i]]) > (size_t)nrc)
420 nrc = strlen(gsft->atomnm[tndx[i]]);
431 extern int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c)
436 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
439 for (i = 0; i < 4; i++)
441 a[i] = gsft->a[elem][i];
442 b[i] = gsft->b[elem][i];
450 extern int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
451 const char* fnXVG, const char *fnTRX,
453 real start_q, real end_q,
454 real energy, int ng, const output_env_t oenv)
456 int i, *isize, flags = TRX_READ_X, **index_atp;
458 char **grpname, title[STRLEN];
463 reduced_atom_t **red;
464 structure_factor *sf;
471 gmx_structurefactors_t *gmx_sf;
479 gmx_sf = gmx_structurefactors_init(fnDAT);
481 success = gmx_structurefactors_get_sf(gmx_sf, 0, a, b, &c);
486 /* Read the topology informations */
487 read_tps_conf (fnTPS, title, &top, &ePBC, &xtop, NULL, box, TRUE);
490 /* groups stuff... */
495 fprintf (stderr, "\nSelect %d group%s\n", ng,
499 get_index (&top.atoms, fnNDX, ng, isize, index, grpname);
503 rd_index (fnNDX, ng, isize, index, grpname);
506 /* The first time we read data is a little special */
507 read_first_frame (oenv, &status, fnTRX, &fr, flags);
509 sf->total_n_atoms = fr.natoms;
512 snew (index_atp, ng);
514 r_tmp = max (box[XX][XX], box[YY][YY]);
515 r_tmp = (double) max (box[ZZ][ZZ], r_tmp);
517 sf->ref_k = (2.0 * M_PI) / (r_tmp);
518 /* ref_k will be the reference momentum unit */
519 sf->n_angles = (int) (end_q / sf->ref_k + 0.5);
522 for (i = 0; i < ng; i++)
524 snew (sf->F[i], sf->n_angles);
526 for (i = 0; i < ng; i++)
528 snew (red[i], isize[i]);
529 rearrange_atoms (red[i], &fr, index[i], isize[i], &top, TRUE, gmx_sf);
530 index_atp[i] = create_indexed_atom_type (red[i], isize[i]);
533 sf_table = compute_scattering_factor_table (gmx_sf, (structure_factor_t *)sf);
536 /* This is the main loop over frames */
541 for (i = 0; i < ng; i++)
543 rearrange_atoms (red[i], &fr, index[i], isize[i], &top, FALSE, gmx_sf);
545 compute_structure_factor ((structure_factor_t *)sf, box, red[i], isize[i],
546 start_q, end_q, i, sf_table);
550 while (read_next_frame (oenv, status, &fr));
552 save_data ((structure_factor_t *)sf, fnXVG, ng, start_q, end_q, oenv);
558 gmx_structurefactors_done(gmx_sf);
564 extern void save_data (structure_factor_t *sft, const char *file, int ngrps,
565 real start_q, real end_q, const output_env_t oenv)
570 double *tmp, polarization_factor, A;
572 structure_factor *sf = (structure_factor *)sft;
574 fp = xvgropen (file, "Scattering Intensity", "q (1/nm)",
575 "Intensity (a.u.)", oenv);
579 for (g = 0; g < ngrps; g++)
581 for (i = 0; i < sf->n_angles; i++)
585 * theta is half the angle between incoming and scattered vectors.
587 * polar. fact. = 0.5*(1+cos^2(2*theta)) = 1 - 0.5 * sin^2(2*theta)
589 * sin(theta) = q/(2k) := A -> sin^2(theta) = 4*A^2 (1-A^2) ->
590 * -> 0.5*(1+cos^2(2*theta)) = 1 - 2 A^2 (1-A^2)
592 A = (double) (i * sf->ref_k) / (2.0 * sf->momentum);
593 polarization_factor = 1 - 2.0 * sqr (A) * (1 - sqr (A));
594 sf->F[g][i] *= polarization_factor;
597 for (i = 0; i < sf->n_angles; i++)
599 if (i * sf->ref_k >= start_q && i * sf->ref_k <= end_q)
601 fprintf (fp, "%10.5f ", i * sf->ref_k);
602 for (g = 0; g < ngrps; g++)
604 fprintf (fp, " %10.5f ", (sf->F[g][i]) /( sf->total_n_atoms*
615 extern double CMSF (gmx_structurefactors_t *gsf, int type, int nh, double lambda, double sin_theta)
617 * return Cromer-Mann fit for the atomic scattering factor:
618 * sin_theta is the sine of half the angle between incoming and scattered
619 * vectors. See g_sq.h for a short description of CM fit.
623 double tmp = 0.0, k2;
632 * f0[k] = c + [SUM a_i*EXP(-b_i*(k^2)) ]
642 tmp = (CMSF (gsf, return_atom_type ("C", gsf), 0, lambda, sin_theta) +
643 nh*CMSF (gsf, return_atom_type ("H", gsf), 0, lambda, sin_theta));
648 k2 = (sqr (sin_theta) / sqr (10.0 * lambda));
649 success = gmx_structurefactors_get_sf(gsf, type, a, b, &c);
651 for (i = 0; (i < 4); i++)
653 tmp += a[i] * exp (-b[i] * k2);
661 extern real **gmx_structurefactors_table(gmx_structurefactors_t *gsf, real momentum, real ref_k, real lambda, int n_angles)
669 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
671 NCMT = gsft->nratoms;
674 snew (sf_table, nsftable);
675 for (i = 0; (i < nsftable); i++)
677 snew (sf_table[i], n_angles);
678 for (j = 0; j < n_angles; j++)
680 q = ((double) j * ref_k);
681 /* theta is half the angle between incoming
682 and scattered wavevectors. */
683 sin_theta = q / (2.0 * momentum);
686 sf_table[i][j] = CMSF (gsf, i, 0, lambda, sin_theta);
690 sf_table[i][j] = CMSF (gsf, i, i-NCMT+1, lambda, sin_theta);
697 extern void gmx_structurefactors_done(gmx_structurefactors_t *gsf)
701 gmx_structurefactors *sf;
702 sf = (gmx_structurefactors *) gsf;
704 for (i = 0; i < sf->nratoms; i++)
708 sfree(sf->atomnm[i]);
721 extern real **compute_scattering_factor_table (gmx_structurefactors_t *gsf, structure_factor_t *sft)
724 * this function build up a table of scattering factors for every atom
725 * type and for every scattering angle.
728 double hc = 1239.842;
731 structure_factor *sf = (structure_factor *)sft;
734 /* \hbar \omega \lambda = hc = 1239.842 eV * nm */
735 sf->momentum = ((double) (2. * 1000.0 * M_PI * sf->energy) / hc);
736 sf->lambda = hc / (1000.0 * sf->energy);
737 fprintf (stderr, "\nwavelenght = %f nm\n", sf->lambda);
739 sf_table = gmx_structurefactors_table(gsf, sf->momentum, sf->ref_k, sf->lambda, sf->n_angles);