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

Gradient removed by least square plane fit. More...

#include <processor.h>

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

Public Member Functions

void process_inplace (EMData *image)
 To process an image in-place. More...
 
string get_name () const
 Get the processor's name. More...
 
string get_desc () const
 Get the descrition of this specific processor. More...
 
TypeDict get_param_types () const
 Get processor parameter information in a dictionary. More...
 
- Public Member Functions inherited from EMAN::Processor
virtual ~Processor ()
 
virtual EMDataprocess (const EMData *const image)
 To proccess an image out-of-place. More...
 
virtual void process_list_inplace (vector< EMData * > &images)
 To process multiple images using the same algorithm. More...
 
virtual Dict get_params () const
 Get the processor parameters in a key/value dictionary. More...
 
virtual void set_params (const Dict &new_params)
 Set the processor parameters using a key/value dictionary. More...
 

Static Public Member Functions

static ProcessorNEW ()
 
- Static Public Member Functions inherited from EMAN::Processor
static string get_group_desc ()
 Get the description of this group of processors. More...
 
static void EMFourierFilterInPlace (EMData *fimage, Dict params)
 Compute a Fourier-filter processed image in place. More...
 
static EMDataEMFourierFilter (EMData *fimage, Dict params)
 Compute a Fourier-processor processed image without altering the original image. More...
 

Static Public Attributes

static const string NAME = "filter.gradientPlaneRemover"
 

Additional Inherited Members

- Public Types inherited from EMAN::Processor
enum  fourier_filter_types {
  TOP_HAT_LOW_PASS , TOP_HAT_HIGH_PASS , TOP_HAT_BAND_PASS , TOP_HOMOMORPHIC ,
  GAUSS_LOW_PASS , GAUSS_HIGH_PASS , GAUSS_BAND_PASS , GAUSS_INVERSE ,
  GAUSS_HOMOMORPHIC , BUTTERWORTH_LOW_PASS , BUTTERWORTH_HIGH_PASS , BUTTERWORTH_HOMOMORPHIC ,
  KAISER_I0 , KAISER_SINH , KAISER_I0_INVERSE , KAISER_SINH_INVERSE ,
  SHIFT , TANH_LOW_PASS , TANH_HIGH_PASS , TANH_HOMOMORPHIC ,
  TANH_BAND_PASS , RADIAL_TABLE , CTF_
}
 Fourier filter Processor type enum. More...
 
- Protected Attributes inherited from EMAN::Processor
Dict params
 

Detailed Description

Gradient removed by least square plane fit.

Parameters
mask[in]optional EMData object to mask the pixels used to fit the plane
changeZero[in]optional bool to specify if the zero value pixels are modified
planeParam[out]optional return parameters [nx, ny, nz, cx, cy, cz] for the fitted plane defined as (x-cx)*nx+(y-cy)*ny+(z-cz)*nz=0
Author
Wen Jiang
Date
2006/7/18

Definition at line 5203 of file processor.h.

Member Function Documentation

◆ get_desc()

string EMAN::GradientPlaneRemoverProcessor::get_desc ( ) const
inlinevirtual

Get the descrition of this specific processor.

This function must be overwritten by a subclass.

Returns
The description of this processor.

Implements EMAN::Processor.

Definition at line 5217 of file processor.h.

5218 {
5219 return "Remove gradient by least square plane fit";
5220 }

◆ get_name()

string EMAN::GradientPlaneRemoverProcessor::get_name ( ) const
inlinevirtual

Get the processor's name.

Each processor is identified by a unique name.

Returns
The processor's name.

Implements EMAN::Processor.

Definition at line 5208 of file processor.h.

5209 {
5210 return NAME;
5211 }

References NAME.

Referenced by process_inplace().

◆ get_param_types()

TypeDict EMAN::GradientPlaneRemoverProcessor::get_param_types ( ) const
inlinevirtual

Get processor parameter information in a dictionary.

Each parameter has one record in the dictionary. Each record contains its name, data-type, and description.

Returns
A dictionary containing the parameter info.

Reimplemented from EMAN::Processor.

Definition at line 5222 of file processor.h.

5223 {
5224 TypeDict d;
5225 d.put("mask", EMObject::EMDATA, "mask object: nonzero pixel positions will be used to fit plane. default = 0");
5226 d.put("changeZero", EMObject::INT, "if zero pixels are modified when removing gradient. default = 0");
5227 d.put("planeParam", EMObject::FLOATARRAY, "fitted plane parameters output");
5228 return d;
5229 }
TypeDict is a dictionary to store <string, EMObject::ObjectType> pair.
Definition: emobject.h:305
void put(const string &key, EMObject::ObjectType o, const string &desc="")
Definition: emobject.h:330

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

◆ NEW()

static Processor * EMAN::GradientPlaneRemoverProcessor::NEW ( )
inlinestatic

Definition at line 5212 of file processor.h.

5213 {
5214 return new GradientPlaneRemoverProcessor();
5215 }
Gradient removed by least square plane fit.
Definition: processor.h:5204

◆ process_inplace()

void GradientPlaneRemoverProcessor::process_inplace ( EMData image)
virtual

To process an image in-place.

For those processors which can only be processed out-of-place, override this function to just print out some error message to remind user call the out-of-place version.

Parameters
imageThe image to be processed.

Implements EMAN::Processor.

Definition at line 3809 of file processor.cpp.

3810{
3811 if (!image) {
3812 LOGWARN("NULL Image");
3813 return;
3814 }
3815
3816 int nz = image->get_zsize();
3817 if (nz > 1) {
3818 LOGERR("%s Processor doesn't support 3D model", get_name().c_str());
3819 throw ImageDimensionException("3D map not supported");
3820 }
3821
3822 int nx = image->get_xsize();
3823 int ny = image->get_ysize();
3824 float *d = image->get_data();
3825 EMData *mask = 0;
3826 float *dm = 0;
3827 if (params.has_key("mask")) {
3828 mask = params["mask"];
3829 if (nx!=mask->get_xsize() || ny!=mask->get_ysize()) {
3830 LOGERR("%s Processor requires same size mask image", get_name().c_str());
3831 throw ImageDimensionException("wrong size mask image");
3832 }
3833 dm = mask->get_data();
3834 }
3835 int count = 0;
3836 if (dm) {
3837 for(int i=0; i<nx*ny; i++) {
3838 if(dm[i]) count++;
3839 }
3840 }
3841 else {
3842 count = nx * ny;
3843 }
3844 if(count<3) {
3845 LOGERR("%s Processor requires at least 3 pixels to fit a plane", get_name().c_str());
3846 throw ImageDimensionException("too few usable pixels to fit a plane");
3847 }
3848 // Allocate the working space
3849 gsl_vector *S=gsl_vector_calloc(3);
3850 gsl_matrix *A=gsl_matrix_calloc(count,3);
3851 gsl_matrix *V=gsl_matrix_calloc(3,3);
3852
3853 double m[3] = {0, 0, 0};
3854 int index=0;
3855 if (dm) {
3856 for(int j=0; j<ny; j++){
3857 for(int i=0; i<nx; i++){
3858 int ij=j*nx+i;
3859 if(dm[ij]) {
3860 m[0]+=i; // x
3861 m[1]+=j; // y
3862 m[2]+=d[ij]; // z
3863 /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
3864 index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
3865 index++;
3866 }
3867 }
3868 }
3869 }
3870 else {
3871 for(int j=0; j<ny; j++){
3872 for(int i=0; i<nx; i++){
3873 int ij=j*nx+i;
3874 m[0]+=i; // x
3875 m[1]+=j; // y
3876 m[2]+=d[ij]; // z
3877 /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
3878 index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
3879 index++;
3880 }
3881 }
3882 }
3883
3884 for(int i=0; i<3; i++) m[i]/=count; // compute center of the plane
3885
3886 index=0;
3887 if (dm) {
3888 for(int j=0; j<ny; j++){
3889 for(int i=0; i<nx; i++){
3890 int ij=j*nx+i;
3891 if(dm[ij]) {
3892 //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
3893 gsl_matrix_set(A,index,0,i-m[0]);
3894 gsl_matrix_set(A,index,1,j-m[1]);
3895 gsl_matrix_set(A,index,2,d[ij]-m[2]);
3896 index++;
3897 }
3898 }
3899 }
3900 mask->update();
3901 }
3902 else {
3903 for(int j=0; j<ny; j++){
3904 for(int i=0; i<nx; i++){
3905 int ij=j*nx+i;
3906 //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
3907 gsl_matrix_set(A,index,0,i-m[0]);
3908 gsl_matrix_set(A,index,1,j-m[1]);
3909 gsl_matrix_set(A,index,2,d[ij]-m[2]);
3910 index++;
3911 }
3912 }
3913 }
3914
3915 // SVD decomposition and use the V vector associated with smallest singular value as the plan normal
3916 gsl_linalg_SV_decomp_jacobi(A, V, S);
3917
3918 double n[3];
3919 for(int i=0; i<3; i++) n[i] = gsl_matrix_get(V, i, 2);
3920
3921 #ifdef DEBUG
3922 printf("S=%g,%g,%g\n",gsl_vector_get(S,0), gsl_vector_get(S,1), gsl_vector_get(S,2));
3923 printf("V[0,:]=%g,%g,%g\n",gsl_matrix_get(V,0,0), gsl_matrix_get(V,0,1),gsl_matrix_get(V,0,2));
3924 printf("V[1,:]=%g,%g,%g\n",gsl_matrix_get(V,1,0), gsl_matrix_get(V,1,1),gsl_matrix_get(V,1,2));
3925 printf("V[2,:]=%g,%g,%g\n",gsl_matrix_get(V,2,0), gsl_matrix_get(V,2,1),gsl_matrix_get(V,2,2));
3926 printf("Fitted plane: p0=%g,%g,%g\tn=%g,%g,%g\n",m[0],m[1],m[2],n[0],n[1],n[2]);
3927 #endif
3928
3929 int changeZero = 0;
3930 if (params.has_key("changeZero")) changeZero = params["changeZero"];
3931 if (changeZero) {
3932 for(int j=0; j<nx; j++){
3933 for(int i=0; i<ny; i++){
3934 int ij = j*nx+i;
3935 d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
3936 }
3937 }
3938 }
3939 else {
3940 for(int j=0; j<nx; j++){
3941 for(int i=0; i<ny; i++){
3942 int ij = j*nx+i;
3943 if(d[ij]) d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
3944 }
3945 }
3946 }
3947 image->update();
3948 // set return plane parameters
3949 vector< float > planeParam;
3950 planeParam.resize(6);
3951 for(int i=0; i<3; i++) planeParam[i] = static_cast<float>(n[i]);
3952 for(int i=0; i<3; i++) planeParam[i+3] = static_cast<float>(m[i]);
3953 params["planeParam"]=EMObject(planeParam);
3954}
#define V(i, j)
Definition: analyzer.cpp:697
bool has_key(const string &key) const
Ask the Dictionary if it as a particular key.
Definition: emobject.h:511
EMData stores an image's data and defines core image processing routines.
Definition: emdata.h:82
EMObject is a wrapper class for types including int, float, double, etc as defined in ObjectType.
Definition: emobject.h:123
string get_name() const
Get the processor's name.
Definition: processor.h:5208
#define ImageDimensionException(desc)
Definition: exception.h:166
#define LOGWARN
Definition: log.h:53
#define LOGERR
Definition: log.h:51
#define dm(i)
Definition: projector.cpp:1606

References dm, get_name(), EMAN::Dict::has_key(), ImageDimensionException, LOGERR, LOGWARN, EMAN::Processor::params, and V.

Member Data Documentation

◆ NAME

const string GradientPlaneRemoverProcessor::NAME = "filter.gradientPlaneRemover"
static

Definition at line 5231 of file processor.h.

Referenced by get_name().


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