029f6e2c8c56b1febfe44590e1bb7d0a5e5007ef
[alexxy/gromacs.git] / src / gromacs / gmxana / sfactor.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "sfactor.h"
40
41 #include <cmath>
42 #include <cstring>
43
44 #include <algorithm>
45
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/strdb.h"
48 #include "gromacs/fileio/trxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/oenv.h"
53 #include "gromacs/legacyheaders/typedefs.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/topology/index.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
62
63 typedef struct gmx_structurefactors {
64     int    nratoms;
65     int   *p;      /* proton number */
66     int   *n;      /* neutron number */
67     /* Parameters for the Cromer Mann fit */
68     real **a;      /* parameter a */
69     real **b;      /* parameter b */
70     real  *c;      /* parameter c */
71     char **atomnm; /* atomname */
72
73 } gmx_structurefactors;
74
75 typedef struct reduced_atom{
76     rvec x;
77     int  t;
78 } reduced_atom;
79
80
81 typedef struct structure_factor
82 {
83     int       n_angles;
84     int       n_groups;
85     double    lambda;
86     double    energy;
87     double    momentum;
88     double    ref_k;
89     double  **F;
90     int       nSteps;
91     int       total_n_atoms;
92 } structure_factor;
93
94
95 extern int * create_indexed_atom_type (reduced_atom_t * atm, int size)
96 {
97 /*
98  * create an index of the atom types found in a  group
99  * i.e.: for water index_atp[0]=type_number_of_O and
100  *                 index_atp[1]=type_number_of_H
101  *
102  * the last element is set to 0
103  */
104     int          *index_atp, i, i_tmp, j;
105
106     reduced_atom *att = (reduced_atom *)atm;
107
108     snew (index_atp, 1);
109     i_tmp        = 1;
110     index_atp[0] = att[0].t;
111     for (i = 1; i < size; i++)
112     {
113         for (j = 0; j < i_tmp; j++)
114         {
115             if (att[i].t == index_atp[j])
116             {
117                 break;
118             }
119         }
120         if (j == i_tmp) /* i.e. no indexed atom type is  == to atm[i].t */
121         {
122             i_tmp++;
123             srenew (index_atp, i_tmp * sizeof (int));
124             index_atp[i_tmp - 1] = att[i].t;
125         }
126     }
127     i_tmp++;
128     srenew (index_atp, i_tmp * sizeof (int));
129     index_atp[i_tmp - 1] = 0;
130     return index_atp;
131 }
132
133
134
135 extern t_complex *** rc_tensor_allocation(int x, int y, int z)
136 {
137     t_complex ***t;
138     int          i, j;
139
140     snew(t, x);
141     snew(t[0], x*y);
142     snew(t[0][0], x*y*z);
143
144     for (j = 1; j < y; j++)
145     {
146         t[0][j] = t[0][j-1] + z;
147     }
148     for (i = 1; i < x; i++)
149     {
150         t[i]    = t[i-1] + y;
151         t[i][0] = t[i-1][0] + y*z;
152         for (j = 1; j < y; j++)
153         {
154             t[i][j] = t[i][j-1] + z;
155         }
156     }
157     return t;
158 }
159
160
161 extern void compute_structure_factor (structure_factor_t * sft, matrix box,
162                                       reduced_atom_t * red, int isize, real start_q,
163                                       real end_q, int group, real **sf_table)
164 {
165     structure_factor *sf   = (structure_factor *)sft;
166     reduced_atom     *redt = (reduced_atom *)red;
167
168     t_complex      ***tmpSF;
169     rvec              k_factor;
170     real              kdotx, asf, kx, ky, kz, krr;
171     int               kr, maxkx, maxky, maxkz, i, j, k, p, *counter;
172
173
174     k_factor[XX] = 2 * M_PI / box[XX][XX];
175     k_factor[YY] = 2 * M_PI / box[YY][YY];
176     k_factor[ZZ] = 2 * M_PI / box[ZZ][ZZ];
177
178     maxkx = static_cast<int>(end_q / k_factor[XX] + 0.5);
179     maxky = static_cast<int>(end_q / k_factor[YY] + 0.5);
180     maxkz = static_cast<int>(end_q / k_factor[ZZ] + 0.5);
181
182     snew (counter, sf->n_angles);
183
184     tmpSF = rc_tensor_allocation(maxkx, maxky, maxkz);
185 /*
186  * The big loop...
187  * compute real and imaginary part of the structure factor for every
188  * (kx,ky,kz))
189  */
190     fprintf(stderr, "\n");
191     for (i = 0; i < maxkx; i++)
192     {
193         fprintf (stderr, "\rdone %3.1f%%     ", (100.0*(i+1))/maxkx);
194         kx = i * k_factor[XX];
195         for (j = 0; j < maxky; j++)
196         {
197             ky = j * k_factor[YY];
198             for (k = 0; k < maxkz; k++)
199             {
200                 if (i != 0 || j != 0 || k != 0)
201                 {
202                     kz  = k * k_factor[ZZ];
203                     krr = std::sqrt (sqr (kx) + sqr (ky) + sqr (kz));
204                     if (krr >= start_q && krr <= end_q)
205                     {
206                         kr = static_cast<int>(krr/sf->ref_k + 0.5);
207                         if (kr < sf->n_angles)
208                         {
209                             counter[kr]++; /* will be used for the copmutation
210                                               of the average*/
211                             for (p = 0; p < isize; p++)
212                             {
213                                 asf = sf_table[redt[p].t][kr];
214
215                                 kdotx = kx * redt[p].x[XX] +
216                                     ky * redt[p].x[YY] + kz * redt[p].x[ZZ];
217
218                                 tmpSF[i][j][k].re += std::cos(kdotx) * asf;
219                                 tmpSF[i][j][k].im += std::sin(kdotx) * asf;
220                             }
221                         }
222                     }
223                 }
224             }
225         }
226     }               /* end loop on i */
227 /*
228  *  compute the square modulus of the structure factor, averaging on the surface
229  *  kx*kx + ky*ky + kz*kz = krr*krr
230  *  note that this is correct only for a (on the macroscopic scale)
231  *  isotropic system.
232  */
233     for (i = 0; i < maxkx; i++)
234     {
235         kx = i * k_factor[XX]; for (j = 0; j < maxky; j++)
236         {
237             ky = j * k_factor[YY]; for (k = 0; k < maxkz; k++)
238             {
239                 kz = k * k_factor[ZZ]; krr = std::sqrt (sqr (kx) + sqr (ky)
240                                                         + sqr (kz)); if (krr >= start_q && krr <= end_q)
241                 {
242                     kr = static_cast<int>(krr / sf->ref_k + 0.5);
243                     if (kr < sf->n_angles && counter[kr] != 0)
244                     {
245                         sf->F[group][kr] +=
246                             (sqr (tmpSF[i][j][k].re) +
247                              sqr (tmpSF[i][j][k].im))/ counter[kr];
248                     }
249                 }
250             }
251         }
252     }
253     sfree(counter);
254     sfree(tmpSF[0][0]);
255     sfree(tmpSF[0]);
256     sfree(tmpSF);
257 }
258
259
260 extern gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn)
261 {
262
263     /* Read the database for the structure factor of the different atoms */
264
265     FILE                 *fp;
266     char                  line[STRLEN];
267     gmx_structurefactors *gsf;
268     double                a1, a2, a3, a4, b1, b2, b3, b4, c;
269     int                   p;
270     int                   i;
271     int                   nralloc = 10;
272     int                   line_no;
273     char                  atomn[32];
274     fp      = libopen(datfn);
275     line_no = 0;
276     snew(gsf, 1);
277
278     snew(gsf->atomnm, nralloc);
279     snew(gsf->a, nralloc);
280     snew(gsf->b, nralloc);
281     snew(gsf->c, nralloc);
282     snew(gsf->p, nralloc);
283     gsf->n       = NULL;
284     gsf->nratoms = line_no;
285     while (get_a_line(fp, line, STRLEN))
286     {
287         i = line_no;
288         if (sscanf(line, "%s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf",
289                    atomn, &p, &a1, &a2, &a3, &a4, &b1, &b2, &b3, &b4, &c) == 11)
290         {
291             gsf->atomnm[i] = gmx_strdup(atomn);
292             gsf->p[i]      = p;
293             snew(gsf->a[i], 4);
294             snew(gsf->b[i], 4);
295             gsf->a[i][0] = a1;
296             gsf->a[i][1] = a2;
297             gsf->a[i][2] = a3;
298             gsf->a[i][3] = a4;
299             gsf->b[i][0] = b1;
300             gsf->b[i][1] = b2;
301             gsf->b[i][2] = b3;
302             gsf->b[i][3] = b4;
303             gsf->c[i]    = c;
304             line_no++;
305             gsf->nratoms = line_no;
306             if (line_no == nralloc)
307             {
308                 nralloc += 10;
309                 srenew(gsf->atomnm, nralloc);
310                 srenew(gsf->a, nralloc);
311                 srenew(gsf->b, nralloc);
312                 srenew(gsf->c, nralloc);
313                 srenew(gsf->p, nralloc);
314             }
315         }
316         else
317         {
318             fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
319                     datfn, line_no);
320         }
321     }
322
323     srenew(gsf->atomnm, gsf->nratoms);
324     srenew(gsf->a, gsf->nratoms);
325     srenew(gsf->b, gsf->nratoms);
326     srenew(gsf->c, gsf->nratoms);
327     srenew(gsf->p, gsf->nratoms);
328
329     fclose(fp);
330
331     return (gmx_structurefactors_t *) gsf;
332
333 }
334
335
336 extern void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, atom_id * index,
337                              int isize, t_topology * top, gmx_bool flag, gmx_structurefactors_t *gsf)
338 /* given the group's index, return the (continuous) array of atoms */
339 {
340     int           i;
341
342     reduced_atom *pos = (reduced_atom *)positions;
343
344     if (flag)
345     {
346         for (i = 0; i < isize; i++)
347         {
348             pos[i].t =
349                 return_atom_type (*(top->atoms.atomname[index[i]]), gsf);
350         }
351     }
352     for (i = 0; i < isize; i++)
353     {
354         copy_rvec (fr->x[index[i]], pos[i].x);
355     }
356 }
357
358
359 extern int return_atom_type (const char *name, gmx_structurefactors_t *gsf)
360 {
361     typedef struct {
362         const char *name;
363         int         nh;
364     } t_united_h;
365     t_united_h            uh[] = {
366         { "CH1", 1 }, { "CH2", 2 }, { "CH3", 3 },
367         { "CS1", 1 }, { "CS2", 2 }, { "CS3", 3 },
368         { "CP1", 1 }, { "CP2", 2 }, { "CP3", 3 }
369     };
370     int                   i, cnt = 0;
371     int                  *tndx;
372     int                   nrc;
373     int                   fndx = 0;
374     int                   NCMT;
375
376     gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
377
378     NCMT = gsft->nratoms;
379
380     snew(tndx, NCMT);
381
382     for (i = 0; (i < asize(uh)); i++)
383     {
384         if (std::strcmp(name, uh[i].name) == 0)
385         {
386             return NCMT-1+uh[i].nh;
387         }
388     }
389
390     for (i = 0; (i < NCMT); i++)
391     {
392         if (std::strncmp (name, gsft->atomnm[i], std::strlen(gsft->atomnm[i])) == 0)
393         {
394             tndx[cnt] = i;
395             cnt++;
396         }
397     }
398
399     if (cnt == 0)
400     {
401         gmx_fatal(FARGS, "\nError: atom (%s) not in list (%d types checked)!\n",
402                   name, i);
403     }
404     else
405     {
406         nrc = 0;
407         for (i = 0; i < cnt; i++)
408         {
409             if (std::strlen(gsft->atomnm[tndx[i]]) > (size_t)nrc)
410             {
411                 nrc  = std::strlen(gsft->atomnm[tndx[i]]);
412                 fndx = tndx[i];
413             }
414         }
415
416         return fndx;
417     }
418
419     return 0;
420 }
421
422 extern int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c)
423 {
424
425     int                   success;
426     int                   i;
427     gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
428     success = 0;
429
430     for (i = 0; i < 4; i++)
431     {
432         a[i] = gsft->a[elem][i];
433         b[i] = gsft->b[elem][i];
434         *c   = gsft->c[elem];
435     }
436
437     success += 1;
438     return success;
439 }
440
441 extern int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
442                                     const char* fnXVG, const char *fnTRX,
443                                     const char* fnDAT,
444                                     real start_q, real end_q,
445                                     real energy, int ng, const output_env_t oenv)
446 {
447     int                     i, *isize, flags = TRX_READ_X, **index_atp;
448     t_trxstatus            *status;
449     char                  **grpname, title[STRLEN];
450     atom_id               **index;
451     t_topology              top;
452     int                     ePBC;
453     t_trxframe              fr;
454     reduced_atom_t        **red;
455     structure_factor       *sf;
456     rvec                   *xtop;
457     real                  **sf_table;
458     matrix                  box;
459     real                    r_tmp;
460
461     gmx_structurefactors_t *gmx_sf;
462     real                   *a, *b, c;
463
464     snew(a, 4);
465     snew(b, 4);
466
467
468     gmx_sf = gmx_structurefactors_init(fnDAT);
469
470     gmx_structurefactors_get_sf(gmx_sf, 0, a, b, &c);
471
472     snew (sf, 1);
473     sf->energy = energy;
474
475     /* Read the topology informations */
476     read_tps_conf (fnTPS, title, &top, &ePBC, &xtop, NULL, box, TRUE);
477     sfree (xtop);
478
479     /* groups stuff... */
480     snew (isize, ng);
481     snew (index, ng);
482     snew (grpname, ng);
483
484     fprintf (stderr, "\nSelect %d group%s\n", ng,
485              ng == 1 ? "" : "s");
486     if (fnTPS)
487     {
488         get_index (&top.atoms, fnNDX, ng, isize, index, grpname);
489     }
490     else
491     {
492         rd_index (fnNDX, ng, isize, index, grpname);
493     }
494
495     /* The first time we read data is a little special */
496     read_first_frame (oenv, &status, fnTRX, &fr, flags);
497
498     sf->total_n_atoms = fr.natoms;
499
500     snew (red, ng);
501     snew (index_atp, ng);
502
503     r_tmp = std::max(box[XX][XX], box[YY][YY]);
504     r_tmp = std::max(box[ZZ][ZZ], r_tmp);
505
506     sf->ref_k = (2.0 * M_PI) / (r_tmp);
507     /* ref_k will be the reference momentum unit */
508     sf->n_angles = static_cast<int>(end_q / sf->ref_k + 0.5);
509
510     snew (sf->F, ng);
511     for (i = 0; i < ng; i++)
512     {
513         snew (sf->F[i], sf->n_angles);
514     }
515     for (i = 0; i < ng; i++)
516     {
517         snew (red[i], isize[i]);
518         rearrange_atoms (red[i], &fr, index[i], isize[i], &top, TRUE, gmx_sf);
519         index_atp[i] = create_indexed_atom_type (red[i], isize[i]);
520     }
521
522     sf_table = compute_scattering_factor_table (gmx_sf, (structure_factor_t *)sf);
523
524
525     /* This is the main loop over frames */
526
527     do
528     {
529         sf->nSteps++;
530         for (i = 0; i < ng; i++)
531         {
532             rearrange_atoms (red[i], &fr, index[i], isize[i], &top, FALSE, gmx_sf);
533
534             compute_structure_factor ((structure_factor_t *)sf, box, red[i], isize[i],
535                                       start_q, end_q, i, sf_table);
536         }
537     }
538
539     while (read_next_frame (oenv, status, &fr));
540
541     save_data ((structure_factor_t *)sf, fnXVG, ng, start_q, end_q, oenv);
542
543
544     sfree(a);
545     sfree(b);
546
547     gmx_structurefactors_done(gmx_sf);
548
549     return 0;
550 }
551
552
553 extern void save_data (structure_factor_t *sft, const char *file, int ngrps,
554                        real start_q, real end_q, const output_env_t oenv)
555 {
556
557     FILE             *fp;
558     int               i, g = 0;
559     double           *tmp, polarization_factor, A;
560
561     structure_factor *sf = (structure_factor *)sft;
562
563     fp = xvgropen (file, "Scattering Intensity", "q (1/nm)",
564                    "Intensity (a.u.)", oenv);
565
566     snew (tmp, ngrps);
567
568     for (g = 0; g < ngrps; g++)
569     {
570         for (i = 0; i < sf->n_angles; i++)
571         {
572
573 /*
574  *          theta is half the angle between incoming and scattered vectors.
575  *
576  *          polar. fact. = 0.5*(1+cos^2(2*theta)) = 1 - 0.5 * sin^2(2*theta)
577  *
578  *          sin(theta) = q/(2k) := A  ->  sin^2(theta) = 4*A^2 (1-A^2) ->
579  *          -> 0.5*(1+cos^2(2*theta)) = 1 - 2 A^2 (1-A^2)
580  */
581             A                   = static_cast<double>(i * sf->ref_k) / (2.0 * sf->momentum);
582             polarization_factor = 1 - 2.0 * sqr (A) * (1 - sqr (A));
583             sf->F[g][i]        *= polarization_factor;
584         }
585     }
586     for (i = 0; i < sf->n_angles; i++)
587     {
588         if (i * sf->ref_k >= start_q && i * sf->ref_k <= end_q)
589         {
590             fprintf (fp, "%10.5f  ", i * sf->ref_k);
591             for (g = 0; g < ngrps; g++)
592             {
593                 fprintf (fp, "  %10.5f ", (sf->F[g][i]) /( sf->total_n_atoms*
594                                                            sf->nSteps));
595             }
596             fprintf (fp, "\n");
597         }
598     }
599
600     xvgrclose (fp);
601 }
602
603
604 extern double CMSF (gmx_structurefactors_t *gsf, int type, int nh, double lambda, double sin_theta)
605 /*
606  * return Cromer-Mann fit for the atomic scattering factor:
607  * sin_theta is the sine of half the angle between incoming and scattered
608  * vectors. See g_sq.h for a short description of CM fit.
609  */
610 {
611     int    i;
612     double tmp = 0.0, k2;
613     real  *a, *b;
614     real   c;
615
616     snew(a, 4);
617     snew(b, 4);
618
619     /*
620      *
621      * f0[k] = c + [SUM a_i*EXP(-b_i*(k^2)) ]
622      *             i=1,4
623      */
624
625     /*
626      *  united atoms case
627      *  CH2 / CH3 groups
628      */
629     if (nh > 0)
630     {
631         tmp = (CMSF (gsf, return_atom_type ("C", gsf), 0, lambda, sin_theta) +
632                nh*CMSF (gsf, return_atom_type ("H", gsf), 0, lambda, sin_theta));
633     }
634     /* all atom case */
635     else
636     {
637         k2      = (sqr (sin_theta) / sqr (10.0 * lambda));
638         gmx_structurefactors_get_sf(gsf, type, a, b, &c);
639         tmp     = c;
640         for (i = 0; (i < 4); i++)
641         {
642             tmp += a[i] * exp (-b[i] * k2);
643         }
644     }
645     return tmp;
646 }
647
648
649
650 extern real **gmx_structurefactors_table(gmx_structurefactors_t *gsf, real momentum, real ref_k, real lambda, int n_angles)
651 {
652
653     int                   NCMT;
654     int                   nsftable;
655     int                   i, j;
656     double                q, sin_theta;
657     real                **sf_table;
658     gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
659
660     NCMT     = gsft->nratoms;
661     nsftable = NCMT+3;
662
663     snew (sf_table, nsftable);
664     for (i = 0; (i < nsftable); i++)
665     {
666         snew (sf_table[i], n_angles);
667         for (j = 0; j < n_angles; j++)
668         {
669             q = static_cast<double>(j * ref_k);
670             /* theta is half the angle between incoming
671                and scattered wavevectors. */
672             sin_theta = q / (2.0 * momentum);
673             if (i < NCMT)
674             {
675                 sf_table[i][j] = CMSF (gsf, i, 0, lambda, sin_theta);
676             }
677             else
678             {
679                 sf_table[i][j] = CMSF (gsf, i, i-NCMT+1, lambda, sin_theta);
680             }
681         }
682     }
683     return sf_table;
684 }
685
686 extern void gmx_structurefactors_done(gmx_structurefactors_t *gsf)
687 {
688
689     int                   i;
690     gmx_structurefactors *sf;
691     sf = (gmx_structurefactors *) gsf;
692
693     for (i = 0; i < sf->nratoms; i++)
694     {
695         sfree(sf->a[i]);
696         sfree(sf->b[i]);
697         sfree(sf->atomnm[i]);
698     }
699
700     sfree(sf->a);
701     sfree(sf->b);
702     sfree(sf->atomnm);
703     sfree(sf->p);
704     sfree(sf->c);
705
706     sfree(sf);
707
708 }
709
710 extern real **compute_scattering_factor_table (gmx_structurefactors_t *gsf, structure_factor_t *sft)
711 {
712 /*
713  *  this function build up a table of scattering factors for every atom
714  *  type and for every scattering angle.
715  */
716
717     double            hc = 1239.842;
718     real           ** sf_table;
719
720     structure_factor *sf = (structure_factor *)sft;
721
722
723     /* \hbar \omega \lambda = hc = 1239.842 eV * nm */
724     sf->momentum = (static_cast<double>(2. * 1000.0 * M_PI * sf->energy) / hc);
725     sf->lambda   = hc / (1000.0 * sf->energy);
726     fprintf (stderr, "\nwavelenght = %f nm\n", sf->lambda);
727
728     sf_table = gmx_structurefactors_table(gsf, sf->momentum, sf->ref_k, sf->lambda, sf->n_angles);
729
730     return sf_table;
731 }