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