EMAN2
Public Member Functions | Static Public Member Functions | Static Public Attributes
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]

List of all members.

Public Member Functions

void process_inplace (EMData *image)
 To process an image in-place.
string get_name () const
 Get the processor's name.
string get_desc () const
 Get the descrition of this specific processor.
TypeDict get_param_types () const
 Get processor parameter information in a dictionary.

Static Public Member Functions

static ProcessorNEW ()

Static Public Attributes

static const string NAME = "filter.gradientPlaneRemover"

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 3847 of file processor.h.


Member Function Documentation

string EMAN::GradientPlaneRemoverProcessor::get_desc ( ) const [inline, virtual]

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 3861 of file processor.h.

                {
                        return "Remove gradient by least square plane fit";
                }
string EMAN::GradientPlaneRemoverProcessor::get_name ( ) const [inline, virtual]

Get the processor's name.

Each processor is identified by a unique name.

Returns:
The processor's name.

Implements EMAN::Processor.

Definition at line 3852 of file processor.h.

References NAME.

                {
                        return NAME;
                }
TypeDict EMAN::GradientPlaneRemoverProcessor::get_param_types ( ) const [inline, virtual]

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 3866 of file processor.h.

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

                {
                        TypeDict d;
                        d.put("mask", EMObject::EMDATA, "mask object: nonzero pixel positions will be used to fit plane. default = 0");
                        d.put("changeZero", EMObject::INT, "if zero pixels are modified when removing gradient. default = 0");
                        d.put("planeParam", EMObject::FLOATARRAY, "fitted plane parameters output");
                        return d;
                }
static Processor* EMAN::GradientPlaneRemoverProcessor::NEW ( ) [inline, static]

Definition at line 3856 of file processor.h.

                {
                        return new GradientPlaneRemoverProcessor();
                }
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 2879 of file processor.cpp.

References dm, EMAN::EMData::get_data(), EMAN::CutoffBlockProcessor::get_name(), EMAN::EMData::get_xsize(), EMAN::EMData::get_ysize(), EMAN::EMData::get_zsize(), EMAN::Dict::has_key(), ImageDimensionException, LOGERR, LOGWARN, nx, ny, EMAN::Processor::params, EMAN::EMData::update(), and V.

{
        if (!image) {
                LOGWARN("NULL Image");
                return;
        }

        int nz = image->get_zsize();
        if (nz > 1) {
                LOGERR("%s Processor doesn't support 3D model", get_name().c_str());
                throw ImageDimensionException("3D map not supported");
        }

        int nx = image->get_xsize();
        int ny = image->get_ysize();
        float *d = image->get_data();
        EMData *mask = 0;
        float *dm = 0;
        if (params.has_key("mask")) {
                mask = params["mask"];
                if (nx!=mask->get_xsize() || ny!=mask->get_ysize()) {
                        LOGERR("%s Processor requires same size mask image", get_name().c_str());
                        throw ImageDimensionException("wrong size mask image");
                }
                dm = mask->get_data();
        }
        int count = 0;
        if (dm) {
                for(int i=0; i<nx*ny; i++) {
                        if(dm[i]) count++;
                }
        }
        else {
                count = nx * ny;
        }
        if(count<3) {
                LOGERR("%s Processor requires at least 3 pixels to fit a plane", get_name().c_str());
                throw ImageDimensionException("too few usable pixels to fit a plane");
        }
        // Allocate the working space
        gsl_vector *S=gsl_vector_calloc(3);
        gsl_matrix *A=gsl_matrix_calloc(count,3);
        gsl_matrix *V=gsl_matrix_calloc(3,3);

        double m[3] = {0, 0, 0};
        int index=0;
        if (dm) {
                for(int j=0; j<ny; j++){
                        for(int i=0; i<nx; i++){
                                int ij=j*nx+i;
                                if(dm[ij]) {
                                        m[0]+=i;        // x
                                        m[1]+=j;        // y
                                        m[2]+=d[ij];    // z
                                        /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
                                                index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
                                        index++;
                                }
                        }
                }
        }
        else {
                for(int j=0; j<ny; j++){
                        for(int i=0; i<nx; i++){
                                int ij=j*nx+i;
                                        m[0]+=i;        // x
                                        m[1]+=j;        // y
                                        m[2]+=d[ij];    // z
                                        /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
                                                index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
                                        index++;
                        }
                }
        }

        for(int i=0; i<3; i++) m[i]/=count;     // compute center of the plane

        index=0;
        if (dm) {
                for(int j=0; j<ny; j++){
                        for(int i=0; i<nx; i++){
                                int ij=j*nx+i;
                                if(dm[ij]) {
                                        //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
                                        gsl_matrix_set(A,index,0,i-m[0]);
                                        gsl_matrix_set(A,index,1,j-m[1]);
                                        gsl_matrix_set(A,index,2,d[ij]-m[2]);
                                        index++;
                                }
                        }
                }
                mask->update();
        }
        else {
                for(int j=0; j<ny; j++){
                        for(int i=0; i<nx; i++){
                                int ij=j*nx+i;
                                        //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
                                        gsl_matrix_set(A,index,0,i-m[0]);
                                        gsl_matrix_set(A,index,1,j-m[1]);
                                        gsl_matrix_set(A,index,2,d[ij]-m[2]);
                                        index++;
                        }
                }
        }

        // SVD decomposition and use the V vector associated with smallest singular value as the plan normal
        gsl_linalg_SV_decomp_jacobi(A, V, S);

        double n[3];
        for(int i=0; i<3; i++) n[i] = gsl_matrix_get(V, i, 2);

        #ifdef DEBUG
        printf("S=%g,%g,%g\n",gsl_vector_get(S,0), gsl_vector_get(S,1), gsl_vector_get(S,2));
        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));
        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));
        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));
        printf("Fitted plane: p0=%g,%g,%g\tn=%g,%g,%g\n",m[0],m[1],m[2],n[0],n[1],n[2]);
        #endif

        int changeZero = 0;
        if (params.has_key("changeZero")) changeZero = params["changeZero"];
        if (changeZero) {
                for(int j=0; j<nx; j++){
                        for(int i=0; i<ny; i++){
                                int ij = j*nx+i;
                                d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
                        }
                }
        }
        else {
                for(int j=0; j<nx; j++){
                        for(int i=0; i<ny; i++){
                                int ij = j*nx+i;
                                if(d[ij]) d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
                        }
                }
        }
        image->update();
        // set return plane parameters
        vector< float > planeParam;
        planeParam.resize(6);
        for(int i=0; i<3; i++) planeParam[i] = static_cast<float>(n[i]);
        for(int i=0; i<3; i++) planeParam[i+3] = static_cast<float>(m[i]);
        params["planeParam"]=EMObject(planeParam);
}

Member Data Documentation

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

Definition at line 3875 of file processor.h.

Referenced by get_name().


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