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