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,2015,2016,2017, 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.
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/topology/index.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/trajectory/trajectoryframe.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/strdb.h"
64 typedef struct gmx_structurefactors {
66 int *p; /* proton number */
67 int *n; /* neutron number */
68 /* Parameters for the Cromer Mann fit */
69 real **a; /* parameter a */
70 real **b; /* parameter b */
71 real *c; /* parameter c */
72 char **atomnm; /* atomname */
74 } gmx_structurefactors;
76 typedef struct reduced_atom{
82 typedef struct structure_factor
96 extern int * create_indexed_atom_type (reduced_atom_t * atm, int size)
99 * create an index of the atom types found in a group
100 * i.e.: for water index_atp[0]=type_number_of_O and
101 * index_atp[1]=type_number_of_H
103 * the last element is set to 0
105 int *index_atp, i, i_tmp, j;
107 reduced_atom *att = (reduced_atom *)atm;
111 index_atp[0] = att[0].t;
112 for (i = 1; i < size; i++)
114 for (j = 0; j < i_tmp; j++)
116 if (att[i].t == index_atp[j])
121 if (j == i_tmp) /* i.e. no indexed atom type is == to atm[i].t */
124 srenew (index_atp, i_tmp * sizeof (int));
125 index_atp[i_tmp - 1] = att[i].t;
129 srenew (index_atp, i_tmp * sizeof (int));
130 index_atp[i_tmp - 1] = 0;
136 extern t_complex *** rc_tensor_allocation(int x, int y, int z)
143 snew(t[0][0], x*y*z);
145 for (j = 1; j < y; j++)
147 t[0][j] = t[0][j-1] + z;
149 for (i = 1; i < x; i++)
152 t[i][0] = t[i-1][0] + y*z;
153 for (j = 1; j < y; j++)
155 t[i][j] = t[i][j-1] + z;
162 extern void compute_structure_factor (structure_factor_t * sft, matrix box,
163 reduced_atom_t * red, int isize, real start_q,
164 real end_q, int group, real **sf_table)
166 structure_factor *sf = (structure_factor *)sft;
167 reduced_atom *redt = (reduced_atom *)red;
171 real kdotx, asf, kx, ky, kz, krr;
172 int kr, maxkx, maxky, maxkz, i, j, k, p, *counter;
175 k_factor[XX] = 2 * M_PI / box[XX][XX];
176 k_factor[YY] = 2 * M_PI / box[YY][YY];
177 k_factor[ZZ] = 2 * M_PI / box[ZZ][ZZ];
179 maxkx = static_cast<int>(end_q / k_factor[XX] + 0.5);
180 maxky = static_cast<int>(end_q / k_factor[YY] + 0.5);
181 maxkz = static_cast<int>(end_q / k_factor[ZZ] + 0.5);
183 snew (counter, sf->n_angles);
185 tmpSF = rc_tensor_allocation(maxkx, maxky, maxkz);
188 * compute real and imaginary part of the structure factor for every
191 fprintf(stderr, "\n");
192 for (i = 0; i < maxkx; i++)
194 fprintf (stderr, "\rdone %3.1f%% ", (100.0*(i+1))/maxkx);
196 kx = i * k_factor[XX];
197 for (j = 0; j < maxky; j++)
199 ky = j * k_factor[YY];
200 for (k = 0; k < maxkz; k++)
202 if (i != 0 || j != 0 || k != 0)
204 kz = k * k_factor[ZZ];
205 krr = std::sqrt (gmx::square(kx) + gmx::square(ky) + gmx::square(kz));
206 if (krr >= start_q && krr <= end_q)
208 kr = static_cast<int>(krr/sf->ref_k + 0.5);
209 if (kr < sf->n_angles)
211 counter[kr]++; /* will be used for the copmutation
213 for (p = 0; p < isize; p++)
215 asf = sf_table[redt[p].t][kr];
217 kdotx = kx * redt[p].x[XX] +
218 ky * redt[p].x[YY] + kz * redt[p].x[ZZ];
220 tmpSF[i][j][k].re += std::cos(kdotx) * asf;
221 tmpSF[i][j][k].im += std::sin(kdotx) * asf;
228 } /* end loop on i */
230 * compute the square modulus of the structure factor, averaging on the surface
231 * kx*kx + ky*ky + kz*kz = krr*krr
232 * note that this is correct only for a (on the macroscopic scale)
235 for (i = 0; i < maxkx; i++)
237 kx = i * k_factor[XX]; for (j = 0; j < maxky; j++)
239 ky = j * k_factor[YY]; for (k = 0; k < maxkz; k++)
241 kz = k * k_factor[ZZ]; krr = std::sqrt (gmx::square(kx) + gmx::square(ky)
242 + gmx::square(kz)); if (krr >= start_q && krr <= end_q)
244 kr = static_cast<int>(krr / sf->ref_k + 0.5);
245 if (kr < sf->n_angles && counter[kr] != 0)
248 (gmx::square(tmpSF[i][j][k].re) +
249 gmx::square(tmpSF[i][j][k].im))/ counter[kr];
262 extern gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn)
265 /* Read the database for the structure factor of the different atoms */
269 gmx_structurefactors *gsf;
270 double a1, a2, a3, a4, b1, b2, b3, b4, c;
280 snew(gsf->atomnm, nralloc);
281 snew(gsf->a, nralloc);
282 snew(gsf->b, nralloc);
283 snew(gsf->c, nralloc);
284 snew(gsf->p, nralloc);
286 gsf->nratoms = line_no;
287 while (get_a_line(fp, line, STRLEN))
290 if (sscanf(line, "%s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf",
291 atomn, &p, &a1, &a2, &a3, &a4, &b1, &b2, &b3, &b4, &c) == 11)
293 gsf->atomnm[i] = gmx_strdup(atomn);
307 gsf->nratoms = line_no;
308 if (line_no == nralloc)
311 srenew(gsf->atomnm, nralloc);
312 srenew(gsf->a, nralloc);
313 srenew(gsf->b, nralloc);
314 srenew(gsf->c, nralloc);
315 srenew(gsf->p, nralloc);
320 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
325 srenew(gsf->atomnm, gsf->nratoms);
326 srenew(gsf->a, gsf->nratoms);
327 srenew(gsf->b, gsf->nratoms);
328 srenew(gsf->c, gsf->nratoms);
329 srenew(gsf->p, gsf->nratoms);
333 return (gmx_structurefactors_t *) gsf;
338 extern void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, int * index,
339 int isize, const t_topology * top, gmx_bool flag, gmx_structurefactors_t *gsf)
340 /* given the group's index, return the (continuous) array of atoms */
344 reduced_atom *pos = (reduced_atom *)positions;
348 for (i = 0; i < isize; i++)
351 return_atom_type (*(top->atoms.atomname[index[i]]), gsf);
354 for (i = 0; i < isize; i++)
356 copy_rvec (fr->x[index[i]], pos[i].x);
361 extern int return_atom_type (const char *name, gmx_structurefactors_t *gsf)
368 { "CH1", 1 }, { "CH2", 2 }, { "CH3", 3 },
369 { "CS1", 1 }, { "CS2", 2 }, { "CS3", 3 },
370 { "CP1", 1 }, { "CP2", 2 }, { "CP3", 3 }
378 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
380 NCMT = gsft->nratoms;
384 for (i = 0; (i < asize(uh)); i++)
386 if (std::strcmp(name, uh[i].name) == 0)
388 return NCMT-1+uh[i].nh;
392 for (i = 0; (i < NCMT); i++)
394 if (std::strncmp (name, gsft->atomnm[i], std::strlen(gsft->atomnm[i])) == 0)
403 gmx_fatal(FARGS, "\nError: atom (%s) not in list (%d types checked)!\n",
409 for (i = 0; i < cnt; i++)
411 if (std::strlen(gsft->atomnm[tndx[i]]) > (size_t)nrc)
413 nrc = std::strlen(gsft->atomnm[tndx[i]]);
424 extern int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c)
429 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
432 for (i = 0; i < 4; i++)
434 a[i] = gsft->a[elem][i];
435 b[i] = gsft->b[elem][i];
443 extern int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
444 const char* fnXVG, const char *fnTRX,
446 real start_q, real end_q,
447 real energy, int ng, const gmx_output_env_t *oenv)
449 int i, *isize, flags = TRX_READ_X, **index_atp;
456 reduced_atom_t **red;
457 structure_factor *sf;
463 gmx_structurefactors_t *gmx_sf;
470 gmx_sf = gmx_structurefactors_init(fnDAT);
472 gmx_structurefactors_get_sf(gmx_sf, 0, a, b, &c);
477 /* Read the topology informations */
478 read_tps_conf (fnTPS, &top, &ePBC, &xtop, nullptr, box, TRUE);
481 /* groups stuff... */
486 fprintf (stderr, "\nSelect %d group%s\n", ng,
490 get_index (&top.atoms, fnNDX, ng, isize, index, grpname);
494 rd_index (fnNDX, ng, isize, index, grpname);
497 /* The first time we read data is a little special */
498 read_first_frame (oenv, &status, fnTRX, &fr, flags);
500 sf->total_n_atoms = fr.natoms;
503 snew (index_atp, ng);
505 r_tmp = std::max(box[XX][XX], box[YY][YY]);
506 r_tmp = std::max(box[ZZ][ZZ], r_tmp);
508 sf->ref_k = (2.0 * M_PI) / (r_tmp);
509 /* ref_k will be the reference momentum unit */
510 sf->n_angles = static_cast<int>(end_q / sf->ref_k + 0.5);
513 for (i = 0; i < ng; i++)
515 snew (sf->F[i], sf->n_angles);
517 for (i = 0; i < ng; i++)
519 snew (red[i], isize[i]);
520 rearrange_atoms (red[i], &fr, index[i], isize[i], &top, TRUE, gmx_sf);
521 index_atp[i] = create_indexed_atom_type (red[i], isize[i]);
524 sf_table = compute_scattering_factor_table (gmx_sf, (structure_factor_t *)sf);
527 /* This is the main loop over frames */
532 for (i = 0; i < ng; i++)
534 rearrange_atoms (red[i], &fr, index[i], isize[i], &top, FALSE, gmx_sf);
536 compute_structure_factor ((structure_factor_t *)sf, box, red[i], isize[i],
537 start_q, end_q, i, sf_table);
541 while (read_next_frame (oenv, status, &fr));
543 save_data ((structure_factor_t *)sf, fnXVG, ng, start_q, end_q, oenv);
549 gmx_structurefactors_done(gmx_sf);
555 extern void save_data (structure_factor_t *sft, const char *file, int ngrps,
556 real start_q, real end_q, const gmx_output_env_t *oenv)
561 double *tmp, polarization_factor, A;
563 structure_factor *sf = (structure_factor *)sft;
565 fp = xvgropen (file, "Scattering Intensity", "q (1/nm)",
566 "Intensity (a.u.)", oenv);
570 for (g = 0; g < ngrps; g++)
572 for (i = 0; i < sf->n_angles; i++)
576 * theta is half the angle between incoming and scattered vectors.
578 * polar. fact. = 0.5*(1+cos^2(2*theta)) = 1 - 0.5 * sin^2(2*theta)
580 * sin(theta) = q/(2k) := A -> sin^2(theta) = 4*A^2 (1-A^2) ->
581 * -> 0.5*(1+cos^2(2*theta)) = 1 - 2 A^2 (1-A^2)
583 A = static_cast<double>(i * sf->ref_k) / (2.0 * sf->momentum);
584 polarization_factor = 1 - 2.0 * gmx::square(A) * (1 - gmx::square(A));
585 sf->F[g][i] *= polarization_factor;
588 for (i = 0; i < sf->n_angles; i++)
590 if (i * sf->ref_k >= start_q && i * sf->ref_k <= end_q)
592 fprintf (fp, "%10.5f ", i * sf->ref_k);
593 for (g = 0; g < ngrps; g++)
595 fprintf (fp, " %10.5f ", (sf->F[g][i]) /( sf->total_n_atoms*
606 extern double CMSF (gmx_structurefactors_t *gsf, int type, int nh, double lambda, double sin_theta)
608 * return Cromer-Mann fit for the atomic scattering factor:
609 * sin_theta is the sine of half the angle between incoming and scattered
610 * vectors. See g_sq.h for a short description of CM fit.
614 double tmp = 0.0, k2;
623 * f0[k] = c + [SUM a_i*EXP(-b_i*(k^2)) ]
633 tmp = (CMSF (gsf, return_atom_type ("C", gsf), 0, lambda, sin_theta) +
634 nh*CMSF (gsf, return_atom_type ("H", gsf), 0, lambda, sin_theta));
639 k2 = (gmx::square(sin_theta) / gmx::square(10.0 * lambda));
640 gmx_structurefactors_get_sf(gsf, type, a, b, &c);
642 for (i = 0; (i < 4); i++)
644 tmp += a[i] * exp (-b[i] * k2);
652 extern real **gmx_structurefactors_table(gmx_structurefactors_t *gsf, real momentum, real ref_k, real lambda, int n_angles)
660 gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
662 NCMT = gsft->nratoms;
665 snew (sf_table, nsftable);
666 for (i = 0; (i < nsftable); i++)
668 snew (sf_table[i], n_angles);
669 for (j = 0; j < n_angles; j++)
671 q = static_cast<double>(j * ref_k);
672 /* theta is half the angle between incoming
673 and scattered wavevectors. */
674 sin_theta = q / (2.0 * momentum);
677 sf_table[i][j] = CMSF (gsf, i, 0, lambda, sin_theta);
681 sf_table[i][j] = CMSF (gsf, i, i-NCMT+1, lambda, sin_theta);
688 extern void gmx_structurefactors_done(gmx_structurefactors_t *gsf)
692 gmx_structurefactors *sf;
693 sf = (gmx_structurefactors *) gsf;
695 for (i = 0; i < sf->nratoms; i++)
699 sfree(sf->atomnm[i]);
712 extern real **compute_scattering_factor_table (gmx_structurefactors_t *gsf, structure_factor_t *sft)
715 * this function build up a table of scattering factors for every atom
716 * type and for every scattering angle.
719 double hc = 1239.842;
722 structure_factor *sf = (structure_factor *)sft;
725 /* \hbar \omega \lambda = hc = 1239.842 eV * nm */
726 sf->momentum = (static_cast<double>(2. * 1000.0 * M_PI * sf->energy) / hc);
727 sf->lambda = hc / (1000.0 * sf->energy);
728 fprintf (stderr, "\nwavelenght = %f nm\n", sf->lambda);
730 sf_table = gmx_structurefactors_table(gsf, sf->momentum, sf->ref_k, sf->lambda, sf->n_angles);