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

03628                 {
03629                         return "Remove gradient by least square plane fit";
03630                 }

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

References NAME.

Referenced by process_inplace().

03619                 {
03620                         return NAME;
03621                 }

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

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

03633                 {
03634                         TypeDict d;
03635                         d.put("mask", EMObject::EMDATA, "mask object: nonzero pixel positions will be used to fit plane. default = 0");
03636                         d.put("changeZero", EMObject::INT, "if zero pixels are modified when removing gradient. default = 0");
03637                         d.put("planeParam", EMObject::FLOATARRAY, "fitted plane parameters output");
03638                         return d;
03639                 }

static Processor* EMAN::GradientPlaneRemoverProcessor::NEW (  )  [inline, static]

Definition at line 3622 of file processor.h.

03623                 {
03624                         return new GradientPlaneRemoverProcessor();
03625                 }

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:
image The image to be processed.

Implements EMAN::Processor.

Definition at line 2787 of file processor.cpp.

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

02788 {
02789         if (!image) {
02790                 LOGWARN("NULL Image");
02791                 return;
02792         }
02793 
02794         int nz = image->get_zsize();
02795         if (nz > 1) {
02796                 LOGERR("%s Processor doesn't support 3D model", get_name().c_str());
02797                 throw ImageDimensionException("3D map not supported");
02798         }
02799 
02800         int nx = image->get_xsize();
02801         int ny = image->get_ysize();
02802         float *d = image->get_data();
02803         EMData *mask = 0;
02804         float *dm = 0;
02805         if (params.has_key("mask")) {
02806                 mask = params["mask"];
02807                 if (nx!=mask->get_xsize() || ny!=mask->get_ysize()) {
02808                         LOGERR("%s Processor requires same size mask image", get_name().c_str());
02809                         throw ImageDimensionException("wrong size mask image");
02810                 }
02811                 dm = mask->get_data();
02812         }
02813         int count = 0;
02814         if (dm) {
02815                 for(int i=0; i<nx*ny; i++) {
02816                         if(dm[i]) count++;
02817                 }
02818         }
02819         else {
02820                 count = nx * ny;
02821         }
02822         if(count<3) {
02823                 LOGERR("%s Processor requires at least 3 pixels to fit a plane", get_name().c_str());
02824                 throw ImageDimensionException("too few usable pixels to fit a plane");
02825         }
02826         // Allocate the working space
02827         gsl_vector *S=gsl_vector_calloc(3);
02828         gsl_matrix *A=gsl_matrix_calloc(count,3);
02829         gsl_matrix *V=gsl_matrix_calloc(3,3);
02830 
02831         double m[3] = {0, 0, 0};
02832         int index=0;
02833         if (dm) {
02834                 for(int j=0; j<ny; j++){
02835                         for(int i=0; i<nx; i++){
02836                                 int ij=j*nx+i;
02837                                 if(dm[ij]) {
02838                                         m[0]+=i;        // x
02839                                         m[1]+=j;        // y
02840                                         m[2]+=d[ij];    // z
02841                                         /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
02842                                                 index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
02843                                         index++;
02844                                 }
02845                         }
02846                 }
02847         }
02848         else {
02849                 for(int j=0; j<ny; j++){
02850                         for(int i=0; i<nx; i++){
02851                                 int ij=j*nx+i;
02852                                         m[0]+=i;        // x
02853                                         m[1]+=j;        // y
02854                                         m[2]+=d[ij];    // z
02855                                         /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
02856                                                 index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
02857                                         index++;
02858                         }
02859                 }
02860         }
02861 
02862         for(int i=0; i<3; i++) m[i]/=count;     // compute center of the plane
02863 
02864         index=0;
02865         if (dm) {
02866                 for(int j=0; j<ny; j++){
02867                         for(int i=0; i<nx; i++){
02868                                 int ij=j*nx+i;
02869                                 if(dm[ij]) {
02870                                         //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
02871                                         gsl_matrix_set(A,index,0,i-m[0]);
02872                                         gsl_matrix_set(A,index,1,j-m[1]);
02873                                         gsl_matrix_set(A,index,2,d[ij]-m[2]);
02874                                         index++;
02875                                 }
02876                         }
02877                 }
02878                 mask->update();
02879         }
02880         else {
02881                 for(int j=0; j<ny; j++){
02882                         for(int i=0; i<nx; i++){
02883                                 int ij=j*nx+i;
02884                                         //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
02885                                         gsl_matrix_set(A,index,0,i-m[0]);
02886                                         gsl_matrix_set(A,index,1,j-m[1]);
02887                                         gsl_matrix_set(A,index,2,d[ij]-m[2]);
02888                                         index++;
02889                         }
02890                 }
02891         }
02892 
02893         // SVD decomposition and use the V vector associated with smallest singular value as the plan normal
02894         gsl_linalg_SV_decomp_jacobi(A, V, S);
02895 
02896         double n[3];
02897         for(int i=0; i<3; i++) n[i] = gsl_matrix_get(V, i, 2);
02898 
02899         #ifdef DEBUG
02900         printf("S=%g,%g,%g\n",gsl_vector_get(S,0), gsl_vector_get(S,1), gsl_vector_get(S,2));
02901         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));
02902         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));
02903         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));
02904         printf("Fitted plane: p0=%g,%g,%g\tn=%g,%g,%g\n",m[0],m[1],m[2],n[0],n[1],n[2]);
02905         #endif
02906 
02907         int changeZero = 0;
02908         if (params.has_key("changeZero")) changeZero = params["changeZero"];
02909         if (changeZero) {
02910                 for(int j=0; j<nx; j++){
02911                         for(int i=0; i<ny; i++){
02912                                 int ij = j*nx+i;
02913                                 d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
02914                         }
02915                 }
02916         }
02917         else {
02918                 for(int j=0; j<nx; j++){
02919                         for(int i=0; i<ny; i++){
02920                                 int ij = j*nx+i;
02921                                 if(d[ij]) d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
02922                         }
02923                 }
02924         }
02925         image->update();
02926         // set return plane parameters
02927         vector< float > planeParam;
02928         planeParam.resize(6);
02929         for(int i=0; i<3; i++) planeParam[i] = static_cast<float>(n[i]);
02930         for(int i=0; i<3; i++) planeParam[i+3] = static_cast<float>(m[i]);
02931         params["planeParam"]=EMObject(planeParam);
02932 }


Member Data Documentation

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

Definition at line 3641 of file processor.h.

Referenced by get_name().


The documentation for this class was generated from the following files:
Generated on Thu Nov 17 12:46:24 2011 for EMAN2 by  doxygen 1.4.7