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.
39 #include "gromacs/utility/smalloc.h"
40 #include "gromacs/utility/futil.h"
41 #include "gromacs/math/utilities.h"
42 #include "gromacs/utility/fatalerror.h"
43 #include "gromacs/math/vec.h"
45 #include "gromacs/topology/index.h"
46 #include "gromacs/fileio/strdb.h"
47 #include "gromacs/fileio/tpxio.h"
48 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/fileio/xvgr.h"
57 typedef struct gmx_structurefactors {
59 int *p; /* proton number */
60 int *n; /* neutron number */
61 /* Parameters for the Cromer Mann fit */
62 real **a; /* parameter a */
63 real **b; /* parameter b */
64 real *c; /* parameter c */
65 char **atomnm; /* atomname */
67 } gmx_structurefactors;
69 typedef struct reduced_atom{
75 typedef struct structure_factor
89 extern int * create_indexed_atom_type (reduced_atom_t * atm, int size)
92 * create an index of the atom types found in a group
93 * i.e.: for water index_atp[0]=type_number_of_O and
94 * index_atp[1]=type_number_of_H
96 * the last element is set to 0
98 int *index_atp, i, i_tmp, j;
100 reduced_atom *att = (reduced_atom *)atm;
104 index_atp[0] = att[0].t;
105 for (i = 1; i < size; i++)
107 for (j = 0; j < i_tmp; j++)
109 if (att[i].t == index_atp[j])
114 if (j == i_tmp) /* i.e. no indexed atom type is == to atm[i].t */
117 srenew (index_atp, i_tmp * sizeof (int));
118 index_atp[i_tmp - 1] = att[i].t;
122 srenew (index_atp, i_tmp * sizeof (int));
123 index_atp[i_tmp - 1] = 0;
129 extern t_complex *** rc_tensor_allocation(int x, int y, int z)
136 snew(t[0][0], x*y*z);
138 for (j = 1; j < y; j++)
140 t[0][j] = t[0][j-1] + z;
142 for (i = 1; i < x; i++)
145 t[i][0] = t[i-1][0] + y*z;
146 for (j = 1; j < y; j++)
148 t[i][j] = t[i][j-1] + z;
155 extern void compute_structure_factor (structure_factor_t * sft, matrix box,
156 reduced_atom_t * red, int isize, real start_q,
157 real end_q, int group, real **sf_table)
159 structure_factor *sf = (structure_factor *)sft;
160 reduced_atom *redt = (reduced_atom *)red;
164 real kdotx, asf, kx, ky, kz, krr;
165 int kr, maxkx, maxky, maxkz, i, j, k, p, *counter;
168 k_factor[XX] = 2 * M_PI / box[XX][XX];
169 k_factor[YY] = 2 * M_PI / box[YY][YY];
170 k_factor[ZZ] = 2 * M_PI / box[ZZ][ZZ];
172 maxkx = (int) (end_q / k_factor[XX] + 0.5);
173 maxky = (int) (end_q / k_factor[YY] + 0.5);
174 maxkz = (int) (end_q / k_factor[ZZ] + 0.5);
176 snew (counter, sf->n_angles);
178 tmpSF = rc_tensor_allocation(maxkx, maxky, maxkz);
181 * compute real and imaginary part of the structure factor for every
184 fprintf(stderr, "\n");
185 for (i = 0; i < maxkx; i++)
187 fprintf (stderr, "\rdone %3.1f%% ", (double)(100.0*(i+1))/maxkx);
188 kx = i * k_factor[XX];
189 for (j = 0; j < maxky; j++)
191 ky = j * k_factor[YY];
192 for (k = 0; k < maxkz; k++)
194 if (i != 0 || j != 0 || k != 0)
196 kz = k * k_factor[ZZ];
197 krr = sqrt (sqr (kx) + sqr (ky) + sqr (kz));
198 if (krr >= start_q && krr <= end_q)
200 kr = (int) (krr/sf->ref_k + 0.5);
201 if (kr < sf->n_angles)
203 counter[kr]++; /* will be used for the copmutation
205 for (p = 0; p < isize; p++)
207 asf = sf_table[redt[p].t][kr];
209 kdotx = kx * redt[p].x[XX] +
210 ky * redt[p].x[YY] + kz * redt[p].x[ZZ];
212 tmpSF[i][j][k].re += cos (kdotx) * asf;
213 tmpSF[i][j][k].im += sin (kdotx) * asf;
220 } /* end loop on i */
222 * compute the square modulus of the structure factor, averaging on the surface
223 * kx*kx + ky*ky + kz*kz = krr*krr
224 * note that this is correct only for a (on the macroscopic scale)
227 for (i = 0; i < maxkx; i++)
229 kx = i * k_factor[XX]; for (j = 0; j < maxky; j++)
231 ky = j * k_factor[YY]; for (k = 0; k < maxkz; k++)
233 kz = k * k_factor[ZZ]; krr = sqrt (sqr (kx) + sqr (ky)
234 + sqr (kz)); if (krr >= start_q && krr <= end_q)
236 kr = (int) (krr / sf->ref_k + 0.5);
237 if (kr < sf->n_angles && counter[kr] != 0)
240 (sqr (tmpSF[i][j][k].re) +
241 sqr (tmpSF[i][j][k].im))/ counter[kr];
254 extern gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn)
257 /* Read the database for the structure factor of the different atoms */
261 gmx_structurefactors *gsf;
262 double a1, a2, a3, a4, b1, b2, b3, b4, c;
272 snew(gsf->atomnm, nralloc);
273 snew(gsf->a, nralloc);
274 snew(gsf->b, nralloc);
275 snew(gsf->c, nralloc);
276 snew(gsf->p, nralloc);
278 gsf->nratoms = line_no;
279 while (get_a_line(fp, line, STRLEN))
282 if (sscanf(line, "%s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf",
283 atomn, &p, &a1, &a2, &a3, &a4, &b1, &b2, &b3, &b4, &c) == 11)
285 gsf->atomnm[i] = gmx_strdup(atomn);
299 gsf->nratoms = line_no;
300 if (line_no == nralloc)
303 srenew(gsf->atomnm, nralloc);
304 srenew(gsf->a, nralloc);
305 srenew(gsf->b, nralloc);
306 srenew(gsf->c, nralloc);
307 srenew(gsf->p, nralloc);
312 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
317 srenew(gsf->atomnm, gsf->nratoms);
318 srenew(gsf->a, gsf->nratoms);
319 srenew(gsf->b, gsf->nratoms);
320 srenew(gsf->c, gsf->nratoms);
321 srenew(gsf->p, gsf->nratoms);
325 return (gmx_structurefactors_t *) gsf;
330 extern void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, atom_id * index,
331 int isize, t_topology * top, gmx_bool flag, gmx_structurefactors_t *gsf)
332 /* given the group's index, return the (continuous) array of atoms */
336 reduced_atom *pos = (reduced_atom *)positions;
340 for (i = 0; i < isize; i++)
343 return_atom_type (*(top->atoms.atomname[index[i]]), gsf);
346 for (i = 0; i < isize; i++)
348 copy_rvec (fr->x[index[i]], pos[i].x);
353 extern int return_atom_type (const char *name, gmx_structurefactors_t *gsf)
360 { "CH1", 1 }, { "CH2", 2 }, { "CH3", 3 },
361 { "CS1", 1 }, { "CS2", 2 }, { "CS3", 3 },
362 { "CP1", 1 }, { "CP2", 2 }, { "CP3", 3 }
370 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
372 NCMT = gsft->nratoms;
376 for (i = 0; (i < asize(uh)); i++)
378 if (strcmp(name, uh[i].name) == 0)
380 return NCMT-1+uh[i].nh;
384 for (i = 0; (i < NCMT); i++)
386 if (strncmp (name, gsft->atomnm[i], strlen(gsft->atomnm[i])) == 0)
395 gmx_fatal(FARGS, "\nError: atom (%s) not in list (%d types checked)!\n",
401 for (i = 0; i < cnt; i++)
403 if (strlen(gsft->atomnm[tndx[i]]) > (size_t)nrc)
405 nrc = strlen(gsft->atomnm[tndx[i]]);
416 extern int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c)
421 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
424 for (i = 0; i < 4; i++)
426 a[i] = gsft->a[elem][i];
427 b[i] = gsft->b[elem][i];
435 extern int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
436 const char* fnXVG, const char *fnTRX,
438 real start_q, real end_q,
439 real energy, int ng, const output_env_t oenv)
441 int i, *isize, flags = TRX_READ_X, **index_atp;
443 char **grpname, title[STRLEN];
448 reduced_atom_t **red;
449 structure_factor *sf;
456 gmx_structurefactors_t *gmx_sf;
464 gmx_sf = gmx_structurefactors_init(fnDAT);
466 success = gmx_structurefactors_get_sf(gmx_sf, 0, a, b, &c);
471 /* Read the topology informations */
472 read_tps_conf (fnTPS, title, &top, &ePBC, &xtop, NULL, box, TRUE);
475 /* groups stuff... */
480 fprintf (stderr, "\nSelect %d group%s\n", ng,
484 get_index (&top.atoms, fnNDX, ng, isize, index, grpname);
488 rd_index (fnNDX, ng, isize, index, grpname);
491 /* The first time we read data is a little special */
492 read_first_frame (oenv, &status, fnTRX, &fr, flags);
494 sf->total_n_atoms = fr.natoms;
497 snew (index_atp, ng);
499 r_tmp = max (box[XX][XX], box[YY][YY]);
500 r_tmp = (double) max (box[ZZ][ZZ], r_tmp);
502 sf->ref_k = (2.0 * M_PI) / (r_tmp);
503 /* ref_k will be the reference momentum unit */
504 sf->n_angles = (int) (end_q / sf->ref_k + 0.5);
507 for (i = 0; i < ng; i++)
509 snew (sf->F[i], sf->n_angles);
511 for (i = 0; i < ng; i++)
513 snew (red[i], isize[i]);
514 rearrange_atoms (red[i], &fr, index[i], isize[i], &top, TRUE, gmx_sf);
515 index_atp[i] = create_indexed_atom_type (red[i], isize[i]);
518 sf_table = compute_scattering_factor_table (gmx_sf, (structure_factor_t *)sf);
521 /* This is the main loop over frames */
526 for (i = 0; i < ng; i++)
528 rearrange_atoms (red[i], &fr, index[i], isize[i], &top, FALSE, gmx_sf);
530 compute_structure_factor ((structure_factor_t *)sf, box, red[i], isize[i],
531 start_q, end_q, i, sf_table);
535 while (read_next_frame (oenv, status, &fr));
537 save_data ((structure_factor_t *)sf, fnXVG, ng, start_q, end_q, oenv);
543 gmx_structurefactors_done(gmx_sf);
549 extern void save_data (structure_factor_t *sft, const char *file, int ngrps,
550 real start_q, real end_q, const output_env_t oenv)
555 double *tmp, polarization_factor, A;
557 structure_factor *sf = (structure_factor *)sft;
559 fp = xvgropen (file, "Scattering Intensity", "q (1/nm)",
560 "Intensity (a.u.)", oenv);
564 for (g = 0; g < ngrps; g++)
566 for (i = 0; i < sf->n_angles; i++)
570 * theta is half the angle between incoming and scattered vectors.
572 * polar. fact. = 0.5*(1+cos^2(2*theta)) = 1 - 0.5 * sin^2(2*theta)
574 * sin(theta) = q/(2k) := A -> sin^2(theta) = 4*A^2 (1-A^2) ->
575 * -> 0.5*(1+cos^2(2*theta)) = 1 - 2 A^2 (1-A^2)
577 A = (double) (i * sf->ref_k) / (2.0 * sf->momentum);
578 polarization_factor = 1 - 2.0 * sqr (A) * (1 - sqr (A));
579 sf->F[g][i] *= polarization_factor;
582 for (i = 0; i < sf->n_angles; i++)
584 if (i * sf->ref_k >= start_q && i * sf->ref_k <= end_q)
586 fprintf (fp, "%10.5f ", i * sf->ref_k);
587 for (g = 0; g < ngrps; g++)
589 fprintf (fp, " %10.5f ", (sf->F[g][i]) /( sf->total_n_atoms*
600 extern double CMSF (gmx_structurefactors_t *gsf, int type, int nh, double lambda, double sin_theta)
602 * return Cromer-Mann fit for the atomic scattering factor:
603 * sin_theta is the sine of half the angle between incoming and scattered
604 * vectors. See g_sq.h for a short description of CM fit.
608 double tmp = 0.0, k2;
617 * f0[k] = c + [SUM a_i*EXP(-b_i*(k^2)) ]
627 tmp = (CMSF (gsf, return_atom_type ("C", gsf), 0, lambda, sin_theta) +
628 nh*CMSF (gsf, return_atom_type ("H", gsf), 0, lambda, sin_theta));
633 k2 = (sqr (sin_theta) / sqr (10.0 * lambda));
634 success = gmx_structurefactors_get_sf(gsf, type, a, b, &c);
636 for (i = 0; (i < 4); i++)
638 tmp += a[i] * exp (-b[i] * k2);
646 extern real **gmx_structurefactors_table(gmx_structurefactors_t *gsf, real momentum, real ref_k, real lambda, int n_angles)
654 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
656 NCMT = gsft->nratoms;
659 snew (sf_table, nsftable);
660 for (i = 0; (i < nsftable); i++)
662 snew (sf_table[i], n_angles);
663 for (j = 0; j < n_angles; j++)
665 q = ((double) j * ref_k);
666 /* theta is half the angle between incoming
667 and scattered wavevectors. */
668 sin_theta = q / (2.0 * momentum);
671 sf_table[i][j] = CMSF (gsf, i, 0, lambda, sin_theta);
675 sf_table[i][j] = CMSF (gsf, i, i-NCMT+1, lambda, sin_theta);
682 extern void gmx_structurefactors_done(gmx_structurefactors_t *gsf)
686 gmx_structurefactors *sf;
687 sf = (gmx_structurefactors *) gsf;
689 for (i = 0; i < sf->nratoms; i++)
693 sfree(sf->atomnm[i]);
706 extern real **compute_scattering_factor_table (gmx_structurefactors_t *gsf, structure_factor_t *sft)
709 * this function build up a table of scattering factors for every atom
710 * type and for every scattering angle.
713 double hc = 1239.842;
716 structure_factor *sf = (structure_factor *)sft;
719 /* \hbar \omega \lambda = hc = 1239.842 eV * nm */
720 sf->momentum = ((double) (2. * 1000.0 * M_PI * sf->energy) / hc);
721 sf->lambda = hc / (1000.0 * sf->energy);
722 fprintf (stderr, "\nwavelenght = %f nm\n", sf->lambda);
724 sf_table = gmx_structurefactors_table(gsf, sf->momentum, sf->ref_k, sf->lambda, sf->n_angles);