EMAN2
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
EMAN::nnSSNR_ctfReconstructor Class Reference

#include <reconstructor.h>

Inheritance diagram for EMAN::nnSSNR_ctfReconstructor:
Inheritance graph
[legend]
Collaboration diagram for EMAN::nnSSNR_ctfReconstructor:
Collaboration graph
[legend]

Public Member Functions

 nnSSNR_ctfReconstructor ()
 
 nnSSNR_ctfReconstructor (const string &symmetry, int size, int npad, float snr, int sign)
 
 ~nnSSNR_ctfReconstructor ()
 
virtual void setup ()
 Initialize the reconstructor. More...
 
virtual int insert_slice (const EMData *const slice, const Transform &euler, const float weight)
 Insert a slice into a 3D volume, in a given orientation. More...
 
virtual EMDatafinish (bool doift=true)
 Finish reconstruction and return the complete model. More...
 
virtual string get_name () const
 Get the unique name of this class (especially for factory based instantiation access) More...
 
virtual string get_desc () const
 Get a clear, concise description of this class. More...
 
TypeDict get_param_types () const
 
void setup (const string &symmetry, int size, int npad, float snr, int sign)
 
int insert_padfft_slice (EMData *padded, const Transform &trans, float mult=1)
 
- Public Member Functions inherited from EMAN::Reconstructor
 Reconstructor ()
 
virtual ~Reconstructor ()
 
virtual void setup_seed (EMData *seed, float seed_weight)
 Initialize the reconstructor with a seed volume. More...
 
virtual void setup_seedandweights (EMData *seed, EMData *weight)
 Initialize the reconstructor with a seed volume, as above. More...
 
virtual EMDatapreprocess_slice (const EMData *const slice, const Transform &t=Transform())
 While you can just insert unprocessed slices, if you call preprocess_slice yourself, and insert the returned slice instead, repeatedly, it can save a fair bit of computation. More...
 
int insert_slice (const EMData *const slice, const Transform &euler)
 
virtual int determine_slice_agreement (EMData *slice, const Transform &euler, const float weight=1.0, bool sub=true)
 Compares a slice to the current reconstruction volume and computes a normalization factor and quality. More...
 
virtual EMDataprojection (const Transform &euler, int ret_fourier=1)
 This will create a projection from the current reconstruction. More...
 
virtual void clear ()
 set the volume and tmp_volume data to zero, for use in Monte Carlo reconstructors More...
 
void print_params () const
 Print the current parameters to std::out. More...
 
EMObjectoperator[] (const string &key)
 
- Public Member Functions inherited from EMAN::FactoryBase
 FactoryBase ()
 
virtual ~FactoryBase ()
 
Dict get_params () const
 get a copy of the parameters of this class More...
 
void set_params (const Dict &new_params)
 Set new parameters. More...
 
void set_param (const string key, const EMObject val)
 
void insert_params (const Dict &new_params)
 Insert parameters. More...
 
Dict copy_relevant_params (const FactoryBase *const that) const
 

Static Public Member Functions

static ReconstructorNEW ()
 

Static Public Attributes

static const string NAME = "nnSSNR_ctf"
 

Private Member Functions

void buildFFTVolume ()
 
void buildNormVolume ()
 
void buildNorm2Volume ()
 
void buildNorm3Volume ()
 

Private Attributes

EMDatam_volume
 
EMDatam_wptr
 
EMDatam_wptr2
 
EMDatam_wptr3
 
string m_symmetry
 
int m_weighting
 
int m_vnx
 
int m_vny
 
int m_vnz
 
int m_npad
 
int m_nsym
 
int m_vnzp
 
int m_vnyp
 
int m_vnxp
 
int m_vnzc
 
int m_vnyc
 
int m_vnxc
 
float m_wghta
 
float m_wghtb
 
int m_sign
 
float m_snr
 
int wiener
 

Additional Inherited Members

- Protected Attributes inherited from EMAN::FactoryBase
Dict params
 This is the dictionary the stores the parameters of the object. More...
 

Detailed Description

Definition at line 1642 of file reconstructor.h.

Constructor & Destructor Documentation

◆ nnSSNR_ctfReconstructor() [1/2]

nnSSNR_ctfReconstructor::nnSSNR_ctfReconstructor ( )

Definition at line 5981 of file reconstructor.cpp.

5982{
5983 m_volume = NULL;
5984 m_wptr = NULL;
5985 m_wptr2 = NULL;
5986 m_wptr3 = NULL;
5987}

References m_volume, m_wptr, m_wptr2, and m_wptr3.

Referenced by NEW().

◆ nnSSNR_ctfReconstructor() [2/2]

nnSSNR_ctfReconstructor::nnSSNR_ctfReconstructor ( const string &  symmetry,
int  size,
int  npad,
float  snr,
int  sign 
)

Definition at line 5989 of file reconstructor.cpp.

5990{
5991 m_volume = NULL;
5992 m_wptr = NULL;
5993 m_wptr2 = NULL;
5994 m_wptr3 = NULL;
5995
5996 setup( symmetry, size, npad, snr, sign );
5997}
virtual void setup()
Initialize the reconstructor.

References m_volume, m_wptr, m_wptr2, m_wptr3, and setup().

◆ ~nnSSNR_ctfReconstructor()

nnSSNR_ctfReconstructor::~nnSSNR_ctfReconstructor ( )

Definition at line 5999 of file reconstructor.cpp.

6000{
6001
6002 //if( m_delete_volume ) checked_delete(m_volume);
6003 //if( m_delete_weight ) checked_delete( m_wptr );
6004 //if( m_delete_weight2 ) checked_delete( m_wptr2 );
6005 //if( m_delete_weight3 ) checked_delete( m_wptr3 );
6006 //checked_delete( m_result );
6007}

Member Function Documentation

◆ buildFFTVolume()

void nnSSNR_ctfReconstructor::buildFFTVolume ( )
private

Definition at line 6054 of file reconstructor.cpp.

6054 {
6055
6056 int offset = 2 - m_vnxp%2;
6057 m_volume = params["fftvol"];
6058
6059 m_volume->set_size(m_vnxp+offset,m_vnyp,m_vnzp);
6060 m_volume->to_zero();
6061
6062 m_volume->set_fftodd(m_vnxp % 2);
6063
6064 m_volume->set_nxc(m_vnxc);
6065 m_volume->set_complex(true);
6066 m_volume->set_ri(true); //(real, imaginary) instead of polar coordinate
6067 m_volume->set_fftpad(true);
6068 m_volume->set_attr("npad", m_npad);
6069 m_volume->set_array_offsets(0,1,1);
6070}
Dict params
This is the dictionary the stores the parameters of the object.
Definition: emobject.h:953

References m_npad, m_vnxc, m_vnxp, m_vnyp, m_vnzp, m_volume, and EMAN::FactoryBase::params.

Referenced by setup().

◆ buildNorm2Volume()

void nnSSNR_ctfReconstructor::buildNorm2Volume ( )
private

Definition at line 6081 of file reconstructor.cpp.

6081 {
6082
6083 m_wptr2 = params["weight2"];
6084 m_wptr2->set_size(m_vnxc+1,m_vnyp,m_vnzp);
6085 m_wptr2->to_zero();
6086 m_wptr2->set_array_offsets(0,1,1);
6087}

References m_vnxc, m_vnyp, m_vnzp, m_wptr2, and EMAN::FactoryBase::params.

Referenced by setup().

◆ buildNorm3Volume()

void nnSSNR_ctfReconstructor::buildNorm3Volume ( )
private

Definition at line 6089 of file reconstructor.cpp.

6089 {
6090
6091 m_wptr3 = params["weight3"];
6092 m_wptr3->set_size(m_vnxc+1,m_vnyp,m_vnzp);
6093 m_wptr3->to_zero();
6094 m_wptr3->set_array_offsets(0,1,1);
6095}

References m_vnxc, m_vnyp, m_vnzp, m_wptr3, and EMAN::FactoryBase::params.

Referenced by setup().

◆ buildNormVolume()

void nnSSNR_ctfReconstructor::buildNormVolume ( )
private

Definition at line 6073 of file reconstructor.cpp.

6074{
6075 m_wptr = params["weight"];
6076 m_wptr->set_size(m_vnxc+1,m_vnyp,m_vnzp);
6077 m_wptr->to_zero();
6078 m_wptr->set_array_offsets(0,1,1);
6079}

References m_vnxc, m_vnyp, m_vnzp, m_wptr, and EMAN::FactoryBase::params.

Referenced by setup().

◆ finish()

EMData * nnSSNR_ctfReconstructor::finish ( bool  doift = true)
virtual

Finish reconstruction and return the complete model.

Parameters
doiftA flag indicating whether the returned object should be guaranteed to be in real-space (true) or should be left in whatever space the reconstructor generated
Returns
The result 3D model.

Reimplemented from EMAN::Reconstructor.

Definition at line 6127 of file reconstructor.cpp.

6128{
6129/*
6130 I changed the code on 05/15/2008 so it only returns variance.
6131 Lines commented out are marked by //#
6132 The version prior to the currect changes is r1.190
6133 PAP
6134*/
6135 /***
6136 m_volume ctf*(P^2D->3D(F^3D))
6137 m_wptr ctf^2
6138 m_wptr2 |P^2D->3D(F^3D)|^2
6139 m_wptr3 Kn
6140 nominator = sum_rot [ wght*signal ]
6141 denominator = sum_rot[ wght*variance ]
6142 ***/
6143 int kz, ky;
6144 int box = 7;
6145 int kc = (box-1)/2;
6146 float alpha = 0.0;
6147 float argx, argy, argz;
6148 vector< float > pow_a( 3*kc+1, 1.0 );
6149 float w = params["w"];
6150 float dx2 = 1.0f/float(m_vnxc)/float(m_vnxc);
6151 float dy2 = 1.0f/float(m_vnyc)/float(m_vnyc);
6152 float dz2 = 1.0f/Util::get_max(float(m_vnzc),1.0f)/Util::get_max(float(m_vnzc),1.0f);
6153 int inc = Util::round(float(Util::get_max(m_vnxc,m_vnyc,m_vnzc))/w);
6154
6155 EMData* SSNR = params["SSNR"];
6156 SSNR->set_size(inc+1,4,1);
6157 //#EMData* vol_ssnr = new EMData();
6158 //#vol_ssnr->set_size(m_vnxp, m_vnyp, m_vnzp);
6159 //#vol_ssnr->to_zero();
6160 //# new linea follow
6161 EMData* vol_ssnr = params["vol_ssnr"];
6162 vol_ssnr->set_size(m_vnxp+ 2 - m_vnxp%2, m_vnyp ,m_vnzp);
6163 vol_ssnr->to_zero();
6164 if ( m_vnxp % 2 == 0 ) vol_ssnr->set_fftodd(0);
6165 else vol_ssnr->set_fftodd(1);
6166 vol_ssnr->set_nxc(m_vnxc);
6167 vol_ssnr->set_complex(true);
6168 vol_ssnr->set_ri(true);
6169 vol_ssnr->set_fftpad(false);
6170 //#
6171 float *nom = new float[inc+1];
6172 float *denom = new float[inc+1];
6173 int *ka = new int[inc+1];
6174 int *nn = new int[inc+1];
6175 float wght = 1.f;
6176 for (int i = 0; i <= inc; i++) {
6177 nom[i] = 0.0f;
6178 denom[i] = 0.0f;
6179 nn[i] = 0;
6180 ka[i] = 0;
6181 }
6182 m_volume->symplane2(m_wptr, m_wptr2, m_wptr3);
6183 if ( m_weighting == ESTIMATE ) {
6184 int vol = box*box*box;
6185 for( unsigned int i=1; i < pow_a.size(); ++i ) pow_a[i] = pow_a[i-1] * exp(m_wghta);
6186 pow_a[3*kc] = 0.0;
6187 float max = max3d( kc, pow_a );
6188 alpha = ( 1.0f - 1.0f/(float)vol ) / max;
6189 }
6190 for (int iz = 1; iz <= m_vnzp; iz++) {
6191 if ( iz-1 > m_vnzc ) kz = iz-1-m_vnzp; else kz = iz-1;
6192 argz = float(kz*kz)*dz2;
6193 for (int iy = 1; iy <= m_vnyp; iy++) {
6194 if ( iy-1 > m_vnyc ) ky = iy-1-m_vnyp; else ky = iy-1;
6195 argy = argz + float(ky*ky)*dy2;
6196 for (int ix = 0; ix <= m_vnxc; ix++) {
6197 float Kn = (*m_wptr3)(ix,iy,iz);
6198 argx = std::sqrt(argy + float(ix*ix)*dx2);
6199 int r = Util::round(float(inc)*argx);
6200 if ( r >= 0 && Kn > 4.5f ) {
6201 if ( m_weighting == ESTIMATE ) {
6202 int cx = ix;
6203 int cy = (iy<=m_vnyc) ? iy - 1 : iy - 1 - m_vnyp;
6204 int cz = (iz<=m_vnzc) ? iz - 1 : iz - 1 - m_vnzp;
6205 float sum = 0.0;
6206 for( int ii = -kc; ii <= kc; ++ii ) {
6207 int nbrcx = cx + ii;
6208 if( nbrcx >= m_vnxc ) continue;
6209 for ( int jj= -kc; jj <= kc; ++jj ) {
6210 int nbrcy = cy + jj;
6211 if( nbrcy <= -m_vnyc || nbrcy >= m_vnyc ) continue;
6212 for( int kk = -kc; kk <= kc; ++kk ) {
6213 int nbrcz = cz + jj;
6214 if ( nbrcz <= -m_vnyc || nbrcz >= m_vnyc ) continue;
6215 if( nbrcx < 0 ) {
6216 nbrcx = -nbrcx;
6217 nbrcy = -nbrcy;
6218 nbrcz = -nbrcz;
6219 }
6220 int nbrix = nbrcx;
6221 int nbriy = nbrcy >= 0 ? nbrcy + 1 : nbrcy + 1 + m_vnyp;
6222 int nbriz = nbrcz >= 0 ? nbrcz + 1 : nbrcz + 1 + m_vnzp;
6223 if( (*m_wptr)( nbrix, nbriy, nbriz ) == 0 ) {
6224 int c = 3*kc+1 - std::abs(ii) - std::abs(jj) - std::abs(kk);
6225 sum = sum + pow_a[c];
6226 }
6227 }
6228 }
6229 }
6230// int r = std::abs(cx) + std::abs(cy) + std::abs(cz);
6231 wght = 1.0f / ( 1.0f - alpha * sum );
6232 } // end of ( m_weighting == ESTIMATE )
6233 float nominator = std::norm(m_volume->cmplx(ix,iy,iz))/(*m_wptr)(ix,iy,iz);
6234 float denominator = ((*m_wptr2)(ix,iy,iz)-std::norm(m_volume->cmplx(ix,iy,iz))/(*m_wptr)(ix,iy,iz))/(Kn-1.0f);
6235 // Skip Friedel related values
6236 if( (ix>0 || (kz>=0 && (ky>=0 || kz!=0)))) {
6237 if ( r <= inc ) {
6238 nom[r] += nominator*wght;
6239 denom[r] += denominator/(*m_wptr)(ix,iy,iz)*wght;
6240 nn[r] += 2;
6241 ka[r] += int(Kn);
6242 }
6243/*
6244 float tmp = Util::get_max(nominator/denominator/(*m_wptr)(ix,iy,iz)-1.0f,0.0f);
6245 // Create SSNR as a 3D array (-n/2:n/2+n%2-1)
6246 int iix = m_vnxc + ix; int iiy = m_vnyc + ky; int iiz = m_vnzc + kz;
6247 if( iix >= 0 && iix < m_vnxp && iiy >= 0 && iiy < m_vnyp && iiz >= 0 && iiz < m_vnzp )
6248 (*vol_ssnr)(iix, iiy, iiz) = tmp;
6249 // Friedel part
6250 iix = m_vnxc - ix; iiy = m_vnyc - ky; iiz = m_vnzc - kz;
6251 if( iix >= 0 && iix < m_vnxp && iiy >= 0 && iiy < m_vnyp && iiz >= 0 && iiz < m_vnzp )
6252 (*vol_ssnr)(iix, iiy, iiz) = tmp;
6253*/
6254 }
6255 (*vol_ssnr)(2*ix, iy-1, iz-1) = denominator*wght;
6256 } // end of Kn>4.5 or whatever
6257 }
6258 }
6259 }
6260 for (int i = 0; i <= inc; i++) {
6261 (*SSNR)(i,0,0) = nom[i];
6262 (*SSNR)(i,1,0) = denom[i];
6263 (*SSNR)(i,2,0) = static_cast<float>(nn[i]);
6264 (*SSNR)(i,3,0) = static_cast<float>(ka[i]);
6265 }
6266 vol_ssnr->update();
6267
6268 delete[] nom;
6269 delete[] denom;
6270 delete[] nn;
6271 delete[] ka;
6272
6273 return 0;
6274}
EMData stores an image's data and defines core image processing routines.
Definition: emdata.h:82
static int round(float x)
Get ceiling round of a float number x.
Definition: util.h:465
static float get_max(float f1, float f2)
Get the maximum of 2 numbers.
Definition: util.h:998
EMData * sqrt() const
return square root of current image
float max3d(int kc, const vector< float > &pow_a)
@ ESTIMATE

References ESTIMATE, EMAN::Util::get_max(), m_vnxc, m_vnxp, m_vnyc, m_vnyp, m_vnzc, m_vnzp, m_volume, m_weighting, m_wghta, m_wptr, m_wptr2, m_wptr3, max3d(), EMAN::FactoryBase::params, EMAN::Util::round(), and sqrt().

◆ get_desc()

virtual string EMAN::nnSSNR_ctfReconstructor::get_desc ( ) const
inlinevirtual

Get a clear, concise description of this class.

Returns
a clear, concise description of this class

Implements EMAN::FactoryBase.

Definition at line 1673 of file reconstructor.h.

1674 {
1675 return "Reconstruction by nearest neighbor with 3D SSNR with CTF";
1676 }

◆ get_name()

virtual string EMAN::nnSSNR_ctfReconstructor::get_name ( ) const
inlinevirtual

Get the unique name of this class (especially for factory based instantiation access)

Returns
the unique name of this class

Implements EMAN::FactoryBase.

Definition at line 1668 of file reconstructor.h.

1669 {
1670 return NAME;
1671 }

References NAME.

◆ get_param_types()

TypeDict EMAN::nnSSNR_ctfReconstructor::get_param_types ( ) const
inlinevirtual
Returns
a TypeDict defining and describing the feasible parameters of this class

Implements EMAN::FactoryBase.

Definition at line 1683 of file reconstructor.h.

1684 {
1685 TypeDict d;
1686 d.put("size", EMObject::INT);
1687 d.put("npad", EMObject::INT);
1688 d.put("symmetry", EMObject::STRING);
1689 d.put("fftvol", EMObject::EMDATA);
1690 d.put("fftwvol", EMObject::EMDATA);
1691 d.put("weight", EMObject::EMDATA);
1692 d.put("weight2", EMObject::EMDATA);
1693 d.put("weight3", EMObject::EMDATA);
1694 d.put("SSNR", EMObject::EMDATA);
1695 d.put("vol_ssnr", EMObject::EMDATA);
1696 d.put("w", EMObject::FLOAT);
1697 d.put("sign", EMObject::INT);
1698 d.put("snr", EMObject::FLOAT);
1699 return d;
1700 }

References EMAN::EMObject::EMDATA, EMAN::EMObject::FLOAT, EMAN::EMObject::INT, EMAN::TypeDict::put(), and EMAN::EMObject::STRING.

◆ insert_padfft_slice()

int nnSSNR_ctfReconstructor::insert_padfft_slice ( EMData padded,
const Transform trans,
float  mult = 1 
)

Definition at line 6118 of file reconstructor.cpp.

6119{
6120
6121 // insert slice for all symmetry related positions
6122 vector<Transform> tsym = t.get_sym_proj(m_symmetry);
6123 for (int isym=0; isym < m_nsym; isym++) m_volume->nn_SSNR_ctf(m_wptr, m_wptr2, m_wptr3, padfft, tsym[isym], weight);
6124 return 0;
6125}

References EMAN::Transform::get_sym_proj(), m_nsym, m_symmetry, m_volume, m_wptr, m_wptr2, and m_wptr3.

Referenced by insert_slice().

◆ insert_slice()

int nnSSNR_ctfReconstructor::insert_slice ( const EMData *const  slice,
const Transform euler,
const float  weight 
)
virtual

Insert a slice into a 3D volume, in a given orientation.

Returns
0 if successful, 1 otherwise
Parameters
slicethe image slice to be inserted into the 3D volume
eulerEuler angle of this image slice.
weightA weighting factor for this slice, generally the number of particles in a class-average. May be ignored by some reconstructors
Returns
0 if OK. 1 if error.
Exceptions
NullPointerExceptionif the input EMData pointer is null
ImageFormatExceptionif the image is complex as opposed to real

Reimplemented from EMAN::Reconstructor.

Definition at line 6097 of file reconstructor.cpp.

6097 {
6098 // sanity checks
6099 if (!slice) {
6100 LOGERR("try to insert NULL slice");
6101 return 1;
6102 }
6103 if( weight > 0.0f ) {
6104 int padffted= slice->get_attr_default("padffted", 0);
6105 if ( padffted==0 && (slice->get_xsize()!=slice->get_ysize() || slice->get_xsize()!=m_vnx) ) {
6106 // FIXME: Why doesn't this throw an exception?
6107 LOGERR("Tried to insert a slice that is the wrong size.");
6108 return 1;
6109 }
6110 EMData* padfft = padfft_slice( slice, t, m_npad );
6111
6112 insert_padfft_slice( padfft, t, weight );
6113
6114 checked_delete( padfft );
6115 }
6116 return 0;
6117}
int insert_padfft_slice(EMData *padded, const Transform &trans, float mult=1)
#define LOGERR
Definition: log.h:51
EMData * padfft_slice(const EMData *const slice, const Transform &t, int npad)
Direct Fourier inversion Reconstructor.
void checked_delete(T *&x)

References checked_delete(), insert_padfft_slice(), LOGERR, m_npad, m_vnx, and EMAN::padfft_slice().

◆ NEW()

static Reconstructor * EMAN::nnSSNR_ctfReconstructor::NEW ( )
inlinestatic

Definition at line 1678 of file reconstructor.h.

1679 {
1680 return new nnSSNR_ctfReconstructor();
1681 }

References nnSSNR_ctfReconstructor().

◆ setup() [1/2]

void nnSSNR_ctfReconstructor::setup ( )
virtual

Initialize the reconstructor.

Implements EMAN::Reconstructor.

Definition at line 6009 of file reconstructor.cpp.

6010{
6011 int size = params["size"];
6012 int npad = params["npad"];
6013 int sign = params["sign"];
6014 float snr = params["snr"];
6015 string symmetry;
6016 if( params.has_key("symmetry") ) symmetry = params["symmetry"].to_str();
6017 else symmetry = "c1";
6018
6019 setup( symmetry, size, npad, snr, sign );
6020}
bool has_key(const string &key) const
Ask the Dictionary if it as a particular key.
Definition: emobject.h:511

References EMAN::Dict::has_key(), EMAN::FactoryBase::params, and setup().

Referenced by nnSSNR_ctfReconstructor(), and setup().

◆ setup() [2/2]

void nnSSNR_ctfReconstructor::setup ( const string &  symmetry,
int  size,
int  npad,
float  snr,
int  sign 
)

Definition at line 6021 of file reconstructor.cpp.

6022{
6023
6025 m_wghta = 0.2f;
6026 m_wghtb = 0.004f;
6027 wiener = 1;
6028
6029 m_symmetry = symmetry;
6030 m_npad = npad;
6032
6033 m_sign = sign;
6034 m_snr = snr;
6035
6036 m_vnx = size;
6037 m_vny = size;
6038 m_vnz = size;
6039
6040 m_vnxp = size*npad;
6041 m_vnyp = size*npad;
6042 m_vnzp = size*npad;
6043
6044 m_vnxc = m_vnxp/2;
6045 m_vnyc = m_vnyp/2;
6046 m_vnzc = m_vnzp/2;
6047
6052}
static int get_nsym(const string &sym)
get the number of symmetries associated with the given symmetry name
Definition: transform.cpp:1570

References buildFFTVolume(), buildNorm2Volume(), buildNorm3Volume(), buildNormVolume(), ESTIMATE, EMAN::Transform::get_nsym(), m_npad, m_nsym, m_sign, m_snr, m_symmetry, m_vnx, m_vnxc, m_vnxp, m_vny, m_vnyc, m_vnyp, m_vnz, m_vnzc, m_vnzp, m_weighting, m_wghta, m_wghtb, and wiener.

Member Data Documentation

◆ m_npad

int EMAN::nnSSNR_ctfReconstructor::m_npad
private

Definition at line 1715 of file reconstructor.h.

Referenced by buildFFTVolume(), insert_slice(), and setup().

◆ m_nsym

int EMAN::nnSSNR_ctfReconstructor::m_nsym
private

Definition at line 1716 of file reconstructor.h.

Referenced by insert_padfft_slice(), and setup().

◆ m_sign

int EMAN::nnSSNR_ctfReconstructor::m_sign
private

Definition at line 1725 of file reconstructor.h.

Referenced by setup().

◆ m_snr

float EMAN::nnSSNR_ctfReconstructor::m_snr
private

Definition at line 1726 of file reconstructor.h.

Referenced by setup().

◆ m_symmetry

string EMAN::nnSSNR_ctfReconstructor::m_symmetry
private

Definition at line 1712 of file reconstructor.h.

Referenced by insert_padfft_slice(), and setup().

◆ m_vnx

int EMAN::nnSSNR_ctfReconstructor::m_vnx
private

Definition at line 1714 of file reconstructor.h.

Referenced by insert_slice(), and setup().

◆ m_vnxc

int EMAN::nnSSNR_ctfReconstructor::m_vnxc
private

◆ m_vnxp

int EMAN::nnSSNR_ctfReconstructor::m_vnxp
private

Definition at line 1717 of file reconstructor.h.

Referenced by buildFFTVolume(), finish(), and setup().

◆ m_vny

int EMAN::nnSSNR_ctfReconstructor::m_vny
private

Definition at line 1714 of file reconstructor.h.

Referenced by setup().

◆ m_vnyc

int EMAN::nnSSNR_ctfReconstructor::m_vnyc
private

Definition at line 1718 of file reconstructor.h.

Referenced by finish(), and setup().

◆ m_vnyp

int EMAN::nnSSNR_ctfReconstructor::m_vnyp
private

◆ m_vnz

int EMAN::nnSSNR_ctfReconstructor::m_vnz
private

Definition at line 1714 of file reconstructor.h.

Referenced by setup().

◆ m_vnzc

int EMAN::nnSSNR_ctfReconstructor::m_vnzc
private

Definition at line 1718 of file reconstructor.h.

Referenced by finish(), and setup().

◆ m_vnzp

int EMAN::nnSSNR_ctfReconstructor::m_vnzp
private

◆ m_volume

EMData* EMAN::nnSSNR_ctfReconstructor::m_volume
private

◆ m_weighting

int EMAN::nnSSNR_ctfReconstructor::m_weighting
private

Definition at line 1713 of file reconstructor.h.

Referenced by finish(), and setup().

◆ m_wghta

float EMAN::nnSSNR_ctfReconstructor::m_wghta
private

Definition at line 1723 of file reconstructor.h.

Referenced by finish(), and setup().

◆ m_wghtb

float EMAN::nnSSNR_ctfReconstructor::m_wghtb
private

Definition at line 1724 of file reconstructor.h.

Referenced by setup().

◆ m_wptr

EMData* EMAN::nnSSNR_ctfReconstructor::m_wptr
private

◆ m_wptr2

EMData* EMAN::nnSSNR_ctfReconstructor::m_wptr2
private

◆ m_wptr3

EMData* EMAN::nnSSNR_ctfReconstructor::m_wptr3
private

◆ NAME

const string nnSSNR_ctfReconstructor::NAME = "nnSSNR_ctf"
static

Definition at line 1705 of file reconstructor.h.

Referenced by get_name().

◆ wiener

int EMAN::nnSSNR_ctfReconstructor::wiener
private

Definition at line 1727 of file reconstructor.h.

Referenced by setup().


The documentation for this class was generated from the following files: