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