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