Convert support files in gmxana to C++
[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
43 #include <algorithm>
44
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/strdb.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/names.h"
51 #include "gromacs/legacyheaders/oenv.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
60
61
62 typedef struct gmx_structurefactors {
63     int    nratoms;
64     int   *p;      /* proton number */
65     int   *n;      /* neutron number */
66     /* Parameters for the Cromer Mann fit */
67     real **a;      /* parameter a */
68     real **b;      /* parameter b */
69     real  *c;      /* parameter c */
70     char **atomnm; /* atomname */
71
72 } gmx_structurefactors;
73
74 typedef struct reduced_atom{
75     rvec x;
76     int  t;
77 } reduced_atom;
78
79
80 typedef struct structure_factor
81 {
82     int       n_angles;
83     int       n_groups;
84     double    lambda;
85     double    energy;
86     double    momentum;
87     double    ref_k;
88     double  **F;
89     int       nSteps;
90     int       total_n_atoms;
91 } structure_factor;
92
93
94 extern int * create_indexed_atom_type (reduced_atom_t * atm, int size)
95 {
96 /*
97  * create an index of the atom types found in a  group
98  * i.e.: for water index_atp[0]=type_number_of_O and
99  *                 index_atp[1]=type_number_of_H
100  *
101  * the last element is set to 0
102  */
103     int          *index_atp, i, i_tmp, j;
104
105     reduced_atom *att = (reduced_atom *)atm;
106
107     snew (index_atp, 1);
108     i_tmp        = 1;
109     index_atp[0] = att[0].t;
110     for (i = 1; i < size; i++)
111     {
112         for (j = 0; j < i_tmp; j++)
113         {
114             if (att[i].t == index_atp[j])
115             {
116                 break;
117             }
118         }
119         if (j == i_tmp) /* i.e. no indexed atom type is  == to atm[i].t */
120         {
121             i_tmp++;
122             srenew (index_atp, i_tmp * sizeof (int));
123             index_atp[i_tmp - 1] = att[i].t;
124         }
125     }
126     i_tmp++;
127     srenew (index_atp, i_tmp * sizeof (int));
128     index_atp[i_tmp - 1] = 0;
129     return index_atp;
130 }
131
132
133
134 extern t_complex *** rc_tensor_allocation(int x, int y, int z)
135 {
136     t_complex ***t;
137     int          i, j;
138
139     snew(t, x);
140     snew(t[0], x*y);
141     snew(t[0][0], x*y*z);
142
143     for (j = 1; j < y; j++)
144     {
145         t[0][j] = t[0][j-1] + z;
146     }
147     for (i = 1; i < x; i++)
148     {
149         t[i]    = t[i-1] + y;
150         t[i][0] = t[i-1][0] + y*z;
151         for (j = 1; j < y; j++)
152         {
153             t[i][j] = t[i][j-1] + z;
154         }
155     }
156     return t;
157 }
158
159
160 extern void compute_structure_factor (structure_factor_t * sft, matrix box,
161                                       reduced_atom_t * red, int isize, real start_q,
162                                       real end_q, int group, real **sf_table)
163 {
164     structure_factor *sf   = (structure_factor *)sft;
165     reduced_atom     *redt = (reduced_atom *)red;
166
167     t_complex      ***tmpSF;
168     rvec              k_factor;
169     real              kdotx, asf, kx, ky, kz, krr;
170     int               kr, maxkx, maxky, maxkz, i, j, k, p, *counter;
171
172
173     k_factor[XX] = 2 * M_PI / box[XX][XX];
174     k_factor[YY] = 2 * M_PI / box[YY][YY];
175     k_factor[ZZ] = 2 * M_PI / box[ZZ][ZZ];
176
177     maxkx = static_cast<int>(end_q / k_factor[XX] + 0.5);
178     maxky = static_cast<int>(end_q / k_factor[YY] + 0.5);
179     maxkz = static_cast<int>(end_q / k_factor[ZZ] + 0.5);
180
181     snew (counter, sf->n_angles);
182
183     tmpSF = rc_tensor_allocation(maxkx, maxky, maxkz);
184 /*
185  * The big loop...
186  * compute real and imaginary part of the structure factor for every
187  * (kx,ky,kz))
188  */
189     fprintf(stderr, "\n");
190     for (i = 0; i < maxkx; i++)
191     {
192         fprintf (stderr, "\rdone %3.1f%%     ", (100.0*(i+1))/maxkx);
193         kx = i * k_factor[XX];
194         for (j = 0; j < maxky; j++)
195         {
196             ky = j * k_factor[YY];
197             for (k = 0; k < maxkz; k++)
198             {
199                 if (i != 0 || j != 0 || k != 0)
200                 {
201                     kz  = k * k_factor[ZZ];
202                     krr = std::sqrt (sqr (kx) + sqr (ky) + sqr (kz));
203                     if (krr >= start_q && krr <= end_q)
204                     {
205                         kr = static_cast<int>(krr/sf->ref_k + 0.5);
206                         if (kr < sf->n_angles)
207                         {
208                             counter[kr]++; /* will be used for the copmutation
209                                               of the average*/
210                             for (p = 0; p < isize; p++)
211                             {
212                                 asf = sf_table[redt[p].t][kr];
213
214                                 kdotx = kx * redt[p].x[XX] +
215                                     ky * redt[p].x[YY] + kz * redt[p].x[ZZ];
216
217                                 tmpSF[i][j][k].re += std::cos(kdotx) * asf;
218                                 tmpSF[i][j][k].im += std::sin(kdotx) * asf;
219                             }
220                         }
221                     }
222                 }
223             }
224         }
225     }               /* end loop on i */
226 /*
227  *  compute the square modulus of the structure factor, averaging on the surface
228  *  kx*kx + ky*ky + kz*kz = krr*krr
229  *  note that this is correct only for a (on the macroscopic scale)
230  *  isotropic system.
231  */
232     for (i = 0; i < maxkx; i++)
233     {
234         kx = i * k_factor[XX]; for (j = 0; j < maxky; j++)
235         {
236             ky = j * k_factor[YY]; for (k = 0; k < maxkz; k++)
237             {
238                 kz = k * k_factor[ZZ]; krr = std::sqrt (sqr (kx) + sqr (ky)
239                                                         + sqr (kz)); if (krr >= start_q && krr <= end_q)
240                 {
241                     kr = static_cast<int>(krr / sf->ref_k + 0.5);
242                     if (kr < sf->n_angles && counter[kr] != 0)
243                     {
244                         sf->F[group][kr] +=
245                             (sqr (tmpSF[i][j][k].re) +
246                              sqr (tmpSF[i][j][k].im))/ counter[kr];
247                     }
248                 }
249             }
250         }
251     }
252     sfree(counter);
253     sfree(tmpSF[0][0]);
254     sfree(tmpSF[0]);
255     sfree(tmpSF);
256 }
257
258
259 extern gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn)
260 {
261
262     /* Read the database for the structure factor of the different atoms */
263
264     FILE                 *fp;
265     char                  line[STRLEN];
266     gmx_structurefactors *gsf;
267     double                a1, a2, a3, a4, b1, b2, b3, b4, c;
268     int                   p;
269     int                   i;
270     int                   nralloc = 10;
271     int                   line_no;
272     char                  atomn[32];
273     fp      = libopen(datfn);
274     line_no = 0;
275     snew(gsf, 1);
276
277     snew(gsf->atomnm, nralloc);
278     snew(gsf->a, nralloc);
279     snew(gsf->b, nralloc);
280     snew(gsf->c, nralloc);
281     snew(gsf->p, nralloc);
282     gsf->n       = NULL;
283     gsf->nratoms = line_no;
284     while (get_a_line(fp, line, STRLEN))
285     {
286         i = line_no;
287         if (sscanf(line, "%s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf",
288                    atomn, &p, &a1, &a2, &a3, &a4, &b1, &b2, &b3, &b4, &c) == 11)
289         {
290             gsf->atomnm[i] = gmx_strdup(atomn);
291             gsf->p[i]      = p;
292             snew(gsf->a[i], 4);
293             snew(gsf->b[i], 4);
294             gsf->a[i][0] = a1;
295             gsf->a[i][1] = a2;
296             gsf->a[i][2] = a3;
297             gsf->a[i][3] = a4;
298             gsf->b[i][0] = b1;
299             gsf->b[i][1] = b2;
300             gsf->b[i][2] = b3;
301             gsf->b[i][3] = b4;
302             gsf->c[i]    = c;
303             line_no++;
304             gsf->nratoms = line_no;
305             if (line_no == nralloc)
306             {
307                 nralloc += 10;
308                 srenew(gsf->atomnm, nralloc);
309                 srenew(gsf->a, nralloc);
310                 srenew(gsf->b, nralloc);
311                 srenew(gsf->c, nralloc);
312                 srenew(gsf->p, nralloc);
313             }
314         }
315         else
316         {
317             fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n",
318                     datfn, line_no);
319         }
320     }
321
322     srenew(gsf->atomnm, gsf->nratoms);
323     srenew(gsf->a, gsf->nratoms);
324     srenew(gsf->b, gsf->nratoms);
325     srenew(gsf->c, gsf->nratoms);
326     srenew(gsf->p, gsf->nratoms);
327
328     fclose(fp);
329
330     return (gmx_structurefactors_t *) gsf;
331
332 }
333
334
335 extern void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, atom_id * index,
336                              int isize, t_topology * top, gmx_bool flag, gmx_structurefactors_t *gsf)
337 /* given the group's index, return the (continuous) array of atoms */
338 {
339     int           i;
340
341     reduced_atom *pos = (reduced_atom *)positions;
342
343     if (flag)
344     {
345         for (i = 0; i < isize; i++)
346         {
347             pos[i].t =
348                 return_atom_type (*(top->atoms.atomname[index[i]]), gsf);
349         }
350     }
351     for (i = 0; i < isize; i++)
352     {
353         copy_rvec (fr->x[index[i]], pos[i].x);
354     }
355 }
356
357
358 extern int return_atom_type (const char *name, gmx_structurefactors_t *gsf)
359 {
360     typedef struct {
361         const char *name;
362         int         nh;
363     } t_united_h;
364     t_united_h            uh[] = {
365         { "CH1", 1 }, { "CH2", 2 }, { "CH3", 3 },
366         { "CS1", 1 }, { "CS2", 2 }, { "CS3", 3 },
367         { "CP1", 1 }, { "CP2", 2 }, { "CP3", 3 }
368     };
369     int                   i, cnt = 0;
370     int                  *tndx;
371     int                   nrc;
372     int                   fndx = 0;
373     int                   NCMT;
374
375     gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
376
377     NCMT = gsft->nratoms;
378
379     snew(tndx, NCMT);
380
381     for (i = 0; (i < asize(uh)); i++)
382     {
383         if (strcmp(name, uh[i].name) == 0)
384         {
385             return NCMT-1+uh[i].nh;
386         }
387     }
388
389     for (i = 0; (i < NCMT); i++)
390     {
391         if (strncmp (name, gsft->atomnm[i], strlen(gsft->atomnm[i])) == 0)
392         {
393             tndx[cnt] = i;
394             cnt++;
395         }
396     }
397
398     if (cnt == 0)
399     {
400         gmx_fatal(FARGS, "\nError: atom (%s) not in list (%d types checked)!\n",
401                   name, i);
402     }
403     else
404     {
405         nrc = 0;
406         for (i = 0; i < cnt; i++)
407         {
408             if (strlen(gsft->atomnm[tndx[i]]) > (size_t)nrc)
409             {
410                 nrc  = strlen(gsft->atomnm[tndx[i]]);
411                 fndx = tndx[i];
412             }
413         }
414
415         return fndx;
416     }
417
418     return 0;
419 }
420
421 extern int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c)
422 {
423
424     int                   success;
425     int                   i;
426     gmx_structurefactors *gsft = (gmx_structurefactors *)gsf;
427     success = 0;
428
429     for (i = 0; i < 4; i++)
430     {
431         a[i] = gsft->a[elem][i];
432         b[i] = gsft->b[elem][i];
433         *c   = gsft->c[elem];
434     }
435
436     success += 1;
437     return success;
438 }
439
440 extern int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
441                                     const char* fnXVG, const char *fnTRX,
442                                     const char* fnDAT,
443                                     real start_q, real end_q,
444                                     real energy, int ng, const output_env_t oenv)
445 {
446     int                     i, *isize, flags = TRX_READ_X, **index_atp;
447     t_trxstatus            *status;
448     char                  **grpname, title[STRLEN];
449     atom_id               **index;
450     t_topology              top;
451     int                     ePBC;
452     t_trxframe              fr;
453     reduced_atom_t        **red;
454     structure_factor       *sf;
455     rvec                   *xtop;
456     real                  **sf_table;
457     matrix                  box;
458     real                    r_tmp;
459
460     gmx_structurefactors_t *gmx_sf;
461     real                   *a, *b, c;
462
463     snew(a, 4);
464     snew(b, 4);
465
466
467     gmx_sf = gmx_structurefactors_init(fnDAT);
468
469     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 = std::max(box[XX][XX], box[YY][YY]);
503     r_tmp = std::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 = static_cast<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                   = static_cast<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     xvgrclose (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;
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         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 = static_cast<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 = (static_cast<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 }