Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members

EMAN::Util::KaiserBessel Class Reference

1-D Kaiser-Bessel window function class. More...

#include <util.h>

Inheritance diagram for EMAN::Util::KaiserBessel:

Inheritance graph
[legend]
List of all members.

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.

Detailed Description

1-D Kaiser-Bessel window function class.

(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.

See also:
P. A. Penczek, R. Renka, and H. Schomberg, J. Opt. Soc. Am. _21_, 449 (2004)

Definition at line 311 of file util.h.


Constructor & Destructor Documentation

Util::KaiserBessel::KaiserBessel float  alpha_,
int  K,
float  r_,
float  v_,
int  N_,
float  vtable_ = 0.f,
int  ntable_ = 5999
 

Definition at line 1988 of file util_sparx.cpp.

References alpha, alphar, build_I0table(), fac, facadj, K, v, v, vadjust, and vtable.

01990                 : alpha(alpha_), v(v_), r(r_), N(N_), K(K_), vtable(vtable_),
01991                   ntable(ntable_) {
01992         // Default values are alpha=1.25, K=6, r=0.5, v = K/2
01993         if (0.f == v) v = float(K)/2;
01994         if (0.f == vtable) vtable = v;
01995         alphar = alpha*r;
01996         fac = static_cast<float>(twopi)*alphar*v;
01997         vadjust = 1.0f*v;
01998         facadj = static_cast<float>(twopi)*alphar*vadjust;
01999         build_I0table();
02000 }

virtual EMAN::Util::KaiserBessel::~KaiserBessel  )  [inline, virtual]
 

Definition at line 331 of file util.h.

00451 {


Member Function Documentation

void Util::KaiserBessel::build_I0table  )  [protected, virtual]
 

2*pi*alpha*r*vadjust

Reimplemented in EMAN::Util::FakeKaiserBessel.

Definition at line 2011 of file util_sparx.cpp.

References facadj, fltb, i0table, K, ntable, EMAN::Util::round(), sqrt(), and vadjust.

Referenced by KaiserBessel().

02011                                      {
02012         i0table.resize(ntable+1); // i0table[0:ntable]
02013         int ltab = int(round(float(ntable)/1.25f));
02014         fltb = float(ltab)/(K/2);
02015         float val0 = static_cast<float>(gsl_sf_bessel_I0(facadj));
02016         for (int i=ltab+1; i <= ntable; i++) i0table[i] = 0.f;
02017         for (int i=0; i <= ltab; i++) {
02018                 float s = float(i)/fltb/N;
02019                 if (s < vadjust) {
02020                         float rt = sqrt(1.f - pow(s/vadjust, 2));
02021                         i0table[i] = static_cast<float>(gsl_sf_bessel_I0(facadj*rt))/val0;
02022                 } else {
02023                         i0table[i] = 0.f;
02024                 }
02025 //              cout << "  "<<s*N<<"  "<<i0table[i] <<endl;
02026         }
02027 }

vector<float> EMAN::Util::KaiserBessel::dump_table  )  [inline]
 

Definition at line 334 of file util.h.

00451                 {

kbi0_win EMAN::Util::KaiserBessel::get_kbi0_win  )  [inline]
 

I0 window function object factory.

Definition at line 381 of file util.h.

00451                 {

kbsinh_win EMAN::Util::KaiserBessel::get_kbsinh_win  )  [inline]
 

Sinh window function object factory.

Definition at line 367 of file util.h.

00451                 {

int EMAN::Util::KaiserBessel::get_window_size  )  const [inline]
 

Return the size of the I0 window.

Definition at line 355 of file util.h.

Referenced by EMAN::Util::get_pixel_conv_new(), and EMAN::Util::get_pixel_conv_new_background().

00451 {

float Util::KaiserBessel::I0table_maxerror  ) 
 

Compute the maximum error in the table.

Definition at line 2029 of file util_sparx.cpp.

References i0table.

02029                                          {
02030         float maxdiff = 0.f;
02031         for (int i = 1; i <= ntable; i++) {
02032                 float diff = fabs(i0table[i] - i0table[i-1]);
02033                 if (diff > maxdiff) maxdiff = diff;
02034         }
02035         return maxdiff;
02036 }

float Util::KaiserBessel::i0win float  x  )  const [virtual]
 

Kaiser-Bessel I0 window function.

Reimplemented in EMAN::Util::FakeKaiserBessel.

Definition at line 2002 of file util_sparx.cpp.

References facadj, sqrt(), vadjust, and x.

Referenced by EMAN::Processor::EMFourierFilterFunc().

02002                                            {
02003         float val0 = float(gsl_sf_bessel_I0(facadj));
02004         float absx = fabs(x);
02005         if (absx > vadjust) return 0.f;
02006         float rt = sqrt(1.f - pow(absx/vadjust, 2));
02007         float res = static_cast<float>(gsl_sf_bessel_I0(facadj*rt))/val0;
02008         return res;
02009 }

float EMAN::Util::KaiserBessel::i0win_tab float  x  )  const [inline]
 

Kaiser-Bessel I0 window function (uses table lookup).

Definition at line 342 of file util.h.

Referenced by EMAN::Util::extractpoint2(), EMAN::Util::get_pixel_conv_new(), and EMAN::Util::get_pixel_conv_new_background().

00451                 {

float Util::KaiserBessel::sinhwin float  x  )  const [virtual]
 

Kaiser-Bessel Sinh window function.

Reimplemented in EMAN::Util::FakeKaiserBessel.

Definition at line 2038 of file util_sparx.cpp.

References alphar, fac, sqrt(), and x.

Referenced by EMAN::Processor::EMFourierFilterFunc().

02038                                              {
02039         float val0 = sinh(fac)/fac;
02040         float absx = fabs(x);
02041         if (0.0 == x) {
02042                 float res = 1.0f;
02043                 return res;
02044         } else if (absx == alphar) {
02045                 return 1.0f/val0;
02046         } else if (absx < alphar) {
02047                 float rt = sqrt(1.0f - pow((x/alphar), 2));
02048                 float facrt = fac*rt;
02049                 float res = (sinh(facrt)/facrt)/val0;
02050                 return res;
02051         } else {
02052                 float rt = sqrt(pow((x/alphar),2) - 1.f);
02053                 float facrt = fac*rt;
02054                 float res = (sin(facrt)/facrt)/val0;
02055                 return res;
02056         }
02057 }


Member Data Documentation

float EMAN::Util::KaiserBessel::alpha [protected]
 

Definition at line 314 of file util.h.

Referenced by KaiserBessel().

float EMAN::Util::KaiserBessel::alphar [protected]
 

table spacing

Definition at line 321 of file util.h.

Referenced by KaiserBessel(), and sinhwin().

float EMAN::Util::KaiserBessel::dtable [protected]
 

Definition at line 320 of file util.h.

float EMAN::Util::KaiserBessel::fac [protected]
 

alpha*r

Definition at line 322 of file util.h.

Referenced by KaiserBessel(), and sinhwin().

float EMAN::Util::KaiserBessel::facadj [protected]
 

Definition at line 324 of file util.h.

Referenced by build_I0table(), i0win(), and KaiserBessel().

float EMAN::Util::KaiserBessel::fltb [protected]
 

Tabulate I0 window for speed.

Definition at line 326 of file util.h.

Referenced by build_I0table().

vector<float> EMAN::Util::KaiserBessel::i0table [protected]
 

Definition at line 319 of file util.h.

Referenced by build_I0table(), and I0table_maxerror().

int EMAN::Util::KaiserBessel::K [protected]
 

size in Ix-space

Definition at line 316 of file util.h.

Referenced by build_I0table(), and KaiserBessel().

int EMAN::Util::KaiserBessel::N [protected]
 

Kaiser-Bessel parameters.

Definition at line 315 of file util.h.

int EMAN::Util::KaiserBessel::ntable [protected]
 

table I0 non-zero domain maximum

Definition at line 318 of file util.h.

Referenced by build_I0table().

float EMAN::Util::KaiserBessel::r [protected]
 

Definition at line 314 of file util.h.

float EMAN::Util::KaiserBessel::v [protected]
 

Definition at line 314 of file util.h.

Referenced by KaiserBessel().

float EMAN::Util::KaiserBessel::vadjust [protected]
 

2*pi*alpha*r*v

Definition at line 323 of file util.h.

Referenced by build_I0table(), i0win(), and KaiserBessel().

float EMAN::Util::KaiserBessel::vtable [protected]
 

I0 window size.

Definition at line 317 of file util.h.

Referenced by KaiserBessel().


The documentation for this class was generated from the following files:
Generated on Tue Jun 11 13:49:52 2013 for EMAN2 by  doxygen 1.3.9.1