#include <util.h>


Public Member Functions | |
| KaiserBessel (float alpha_, int K, float r_, float v_, int N_, float vtable_=0.f, int ntable_=5999) | |
| virtual | ~KaiserBessel () |
| float | I0table_maxerror () |
| Compute the maximum error in the table. | |
| vector< float > | dump_table () |
| virtual float | sinhwin (float x) const |
| Kaiser-Bessel Sinh window function. | |
| virtual float | i0win (float x) const |
| Kaiser-Bessel I0 window function. | |
| float | i0win_tab (float x) const |
| Kaiser-Bessel I0 window function (uses table lookup). | |
| int | get_window_size () const |
| Return the size of the I0 window. | |
| kbsinh_win | get_kbsinh_win () |
| Sinh window function object factory. | |
| kbi0_win | get_kbi0_win () |
| I0 window function object factory. | |
Protected Member Functions | |
| virtual void | build_I0table () |
| 2*pi*alpha*r*vadjust | |
Protected Attributes | |
| float | alpha |
| float | v |
| float | r |
| int | N |
| Kaiser-Bessel parameters. | |
| int | K |
| size in Ix-space | |
| float | vtable |
| I0 window size. | |
| int | ntable |
| table I0 non-zero domain maximum | |
| vector< float > | i0table |
| float | dtable |
| float | alphar |
| table spacing | |
| float | fac |
| alpha*r | |
| float | vadjust |
| 2*pi*alpha*r*v | |
| float | facadj |
| float | fltb |
| Tabulate I0 window for speed. | |
Classes | |
| class | kbi0_win |
| I0 window function object. More... | |
| class | kbsinh_win |
| Sinh window function object. More... | |
(It's a class so that the windowing parameters may be instantiated and held in the instance object.)
The I0 version can be tabulated and interpolated upon demand, but the max error needs to be checked. The "vtable" parameter corresponds to the maximum value of x for which the I0 window is non-zero. Setting "vtable" different from "v" corresponds to a change in units of x. In practice, it is often handy to replace x in some sort of absolute units with x described in terms of grid intervals.
The get_kbsinh_win and get_kbi0_win functions return single-argument function objects, which is what a generic routine is likely to want.
Definition at line 271 of file util.h.
| Util::KaiserBessel::KaiserBessel | ( | float | alpha_, | |
| int | K, | |||
| float | r_, | |||
| float | v_, | |||
| int | N_, | |||
| float | vtable_ = 0.f, |
|||
| int | ntable_ = 5999 | |||
| ) |
Definition at line 2008 of file util_sparx.cpp.
References alpha, alphar, build_I0table(), fac, facadj, K, r, twopi, v, vadjust, and vtable.
02010 : alpha(alpha_), v(v_), r(r_), N(N_), K(K_), vtable(vtable_), 02011 ntable(ntable_) { 02012 // Default values are alpha=1.25, K=6, r=0.5, v = K/2 02013 if (0.f == v) v = float(K)/2; 02014 if (0.f == vtable) vtable = v; 02015 alphar = alpha*r; 02016 fac = static_cast<float>(twopi)*alphar*v; 02017 vadjust = 1.0f*v; 02018 facadj = static_cast<float>(twopi)*alphar*vadjust; 02019 build_I0table(); 02020 }
| virtual EMAN::Util::KaiserBessel::~KaiserBessel | ( | ) | [inline, virtual] |
| void Util::KaiserBessel::build_I0table | ( | ) | [protected, virtual] |
2*pi*alpha*r*vadjust
Reimplemented in EMAN::Util::FakeKaiserBessel.
Definition at line 2031 of file util_sparx.cpp.
References facadj, fltb, i0table, K, N, ntable, EMAN::Util::round(), sqrt(), and vadjust.
Referenced by KaiserBessel().
02031 { 02032 i0table.resize(ntable+1); // i0table[0:ntable] 02033 int ltab = int(round(float(ntable)/1.25f)); 02034 fltb = float(ltab)/(K/2); 02035 float val0 = static_cast<float>(gsl_sf_bessel_I0(facadj)); 02036 for (int i=ltab+1; i <= ntable; i++) i0table[i] = 0.f; 02037 for (int i=0; i <= ltab; i++) { 02038 float s = float(i)/fltb/N; 02039 if (s < vadjust) { 02040 float rt = sqrt(1.f - pow(s/vadjust, 2)); 02041 i0table[i] = static_cast<float>(gsl_sf_bessel_I0(facadj*rt))/val0; 02042 } else { 02043 i0table[i] = 0.f; 02044 } 02045 // cout << " "<<s*N<<" "<<i0table[i] <<endl; 02046 } 02047 }
| float Util::KaiserBessel::I0table_maxerror | ( | ) |
Compute the maximum error in the table.
Definition at line 2049 of file util_sparx.cpp.
References i0table, and ntable.
02049 { 02050 float maxdiff = 0.f; 02051 for (int i = 1; i <= ntable; i++) { 02052 float diff = fabs(i0table[i] - i0table[i-1]); 02053 if (diff > maxdiff) maxdiff = diff; 02054 } 02055 return maxdiff; 02056 }
| vector<float> EMAN::Util::KaiserBessel::dump_table | ( | ) | [inline] |
| float Util::KaiserBessel::sinhwin | ( | float | x | ) | const [virtual] |
Kaiser-Bessel Sinh window function.
Reimplemented in EMAN::Util::FakeKaiserBessel.
Definition at line 2058 of file util_sparx.cpp.
References alphar, fac, and sqrt().
Referenced by EMAN::EMData::divkbsinh(), and EMAN::Processor::EMFourierFilterFunc().
02058 { 02059 float val0 = sinh(fac)/fac; 02060 float absx = fabs(x); 02061 if (0.0 == x) { 02062 float res = 1.0f; 02063 return res; 02064 } else if (absx == alphar) { 02065 return 1.0f/val0; 02066 } else if (absx < alphar) { 02067 float rt = sqrt(1.0f - pow((x/alphar), 2)); 02068 float facrt = fac*rt; 02069 float res = (sinh(facrt)/facrt)/val0; 02070 return res; 02071 } else { 02072 float rt = sqrt(pow((x/alphar),2) - 1.f); 02073 float facrt = fac*rt; 02074 float res = (sin(facrt)/facrt)/val0; 02075 return res; 02076 } 02077 }
| float Util::KaiserBessel::i0win | ( | float | x | ) | const [virtual] |
Kaiser-Bessel I0 window function.
Reimplemented in EMAN::Util::FakeKaiserBessel.
Definition at line 2022 of file util_sparx.cpp.
References facadj, sqrt(), and vadjust.
Referenced by EMAN::Processor::EMFourierFilterFunc().
02022 { 02023 float val0 = float(gsl_sf_bessel_I0(facadj)); 02024 float absx = fabs(x); 02025 if (absx > vadjust) return 0.f; 02026 float rt = sqrt(1.f - pow(absx/vadjust, 2)); 02027 float res = static_cast<float>(gsl_sf_bessel_I0(facadj*rt))/val0; 02028 return res; 02029 }
| float EMAN::Util::KaiserBessel::i0win_tab | ( | float | x | ) | const [inline] |
Kaiser-Bessel I0 window function (uses table lookup).
Definition at line 302 of file util.h.
Referenced by EMAN::EMData::extract_plane(), EMAN::EMData::extractline(), EMAN::EMData::extractpoint(), EMAN::Util::extractpoint2(), EMAN::EMData::get_pixel_conv(), EMAN::Util::get_pixel_conv_new(), EMAN::Util::get_pixel_conv_new_background(), EMAN::EMData::rot_scale_conv(), and EMAN::EMData::rot_scale_conv7().
| int EMAN::Util::KaiserBessel::get_window_size | ( | ) | const [inline] |
Return the size of the I0 window.
Definition at line 315 of file util.h.
Referenced by EMAN::EMData::extract_plane(), EMAN::EMData::extractline(), EMAN::EMData::extractpoint(), EMAN::EMData::get_pixel_conv(), EMAN::Util::get_pixel_conv_new(), EMAN::Util::get_pixel_conv_new_background(), EMAN::EMData::rot_scale_conv(), EMAN::EMData::rot_scale_conv7(), EMAN::EMData::rot_scale_conv_new(), and EMAN::EMData::rot_scale_conv_new_background().
| kbsinh_win EMAN::Util::KaiserBessel::get_kbsinh_win | ( | ) | [inline] |
| kbi0_win EMAN::Util::KaiserBessel::get_kbi0_win | ( | ) | [inline] |
float EMAN::Util::KaiserBessel::alpha [protected] |
float EMAN::Util::KaiserBessel::v [protected] |
float EMAN::Util::KaiserBessel::r [protected] |
int EMAN::Util::KaiserBessel::N [protected] |
Kaiser-Bessel parameters.
Definition at line 275 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), and build_I0table().
int EMAN::Util::KaiserBessel::K [protected] |
size in Ix-space
Definition at line 276 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), build_I0table(), and KaiserBessel().
float EMAN::Util::KaiserBessel::vtable [protected] |
int EMAN::Util::KaiserBessel::ntable [protected] |
table I0 non-zero domain maximum
Definition at line 278 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), build_I0table(), and I0table_maxerror().
vector<float> EMAN::Util::KaiserBessel::i0table [protected] |
Definition at line 279 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), build_I0table(), and I0table_maxerror().
float EMAN::Util::KaiserBessel::dtable [protected] |
float EMAN::Util::KaiserBessel::alphar [protected] |
table spacing
Definition at line 281 of file util.h.
Referenced by KaiserBessel(), EMAN::Util::FakeKaiserBessel::sinhwin(), and sinhwin().
float EMAN::Util::KaiserBessel::fac [protected] |
alpha*r
Definition at line 282 of file util.h.
Referenced by KaiserBessel(), EMAN::Util::FakeKaiserBessel::sinhwin(), and sinhwin().
float EMAN::Util::KaiserBessel::vadjust [protected] |
2*pi*alpha*r*v
Definition at line 283 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), build_I0table(), EMAN::Util::FakeKaiserBessel::i0win(), i0win(), and KaiserBessel().
float EMAN::Util::KaiserBessel::facadj [protected] |
Definition at line 284 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), build_I0table(), EMAN::Util::FakeKaiserBessel::i0win(), i0win(), and KaiserBessel().
float EMAN::Util::KaiserBessel::fltb [protected] |
Tabulate I0 window for speed.
Definition at line 286 of file util.h.
Referenced by EMAN::Util::FakeKaiserBessel::build_I0table(), and build_I0table().
1.5.6