#include <processor.h>
Inheritance diagram for EMAN::GradientPlaneRemoverProcessor:


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 | |
| Processor * | NEW () |
Static Public Attributes | |
| const string | NAME = "filter.gradientPlaneRemover" |
| 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 |
Definition at line 3579 of file processor.h.
|
|
Get the descrition of this specific processor. This function must be overwritten by a subclass.
Implements EMAN::Processor. Definition at line 3593 of file processor.h. 03594 {
03595 return "Remove gradient by least square plane fit";
03596 }
|
|
|
Get the processor's name. Each processor is identified by a unique name.
Implements EMAN::Processor. Definition at line 3584 of file processor.h. Referenced by process_inplace(). 03585 {
03586 return NAME;
03587 }
|
|
|
Get processor parameter information in a dictionary. Each parameter has one record in the dictionary. Each record contains its name, data-type, and description.
Reimplemented from EMAN::Processor. Definition at line 3598 of file processor.h. References EMAN::TypeDict::put(). 03599 {
03600 TypeDict d;
03601 d.put("mask", EMObject::EMDATA, "mask object: nonzero pixel positions will be used to fit plane. default = 0");
03602 d.put("changeZero", EMObject::INT, "if zero pixels are modified when removing gradient. default = 0");
03603 d.put("planeParam", EMObject::FLOATARRAY, "fitted plane parameters output");
03604 return d;
03605 }
|
|
|
Definition at line 3588 of file processor.h. 03589 {
03590 return new GradientPlaneRemoverProcessor();
03591 }
|
|
|
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.
Implements EMAN::Processor. Definition at line 2740 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, nx, ny, EMAN::EMData::update(), and V. 02741 {
02742 if (!image) {
02743 LOGWARN("NULL Image");
02744 return;
02745 }
02746
02747 int nz = image->get_zsize();
02748 if (nz > 1) {
02749 LOGERR("%s Processor doesn't support 3D model", get_name().c_str());
02750 throw ImageDimensionException("3D map not supported");
02751 }
02752
02753 int nx = image->get_xsize();
02754 int ny = image->get_ysize();
02755 float *d = image->get_data();
02756 EMData *mask = 0;
02757 float *dm = 0;
02758 if (params.has_key("mask")) {
02759 mask = params["mask"];
02760 if (nx!=mask->get_xsize() || ny!=mask->get_ysize()) {
02761 LOGERR("%s Processor requires same size mask image", get_name().c_str());
02762 throw ImageDimensionException("wrong size mask image");
02763 }
02764 dm = mask->get_data();
02765 }
02766 int count = 0;
02767 if (dm) {
02768 for(int i=0; i<nx*ny; i++) {
02769 if(dm[i]) count++;
02770 }
02771 }
02772 else {
02773 count = nx * ny;
02774 }
02775 if(count<3) {
02776 LOGERR("%s Processor requires at least 3 pixels to fit a plane", get_name().c_str());
02777 throw ImageDimensionException("too few usable pixels to fit a plane");
02778 }
02779 // Allocate the working space
02780 gsl_vector *S=gsl_vector_calloc(3);
02781 gsl_matrix *A=gsl_matrix_calloc(count,3);
02782 gsl_matrix *V=gsl_matrix_calloc(3,3);
02783
02784 double m[3] = {0, 0, 0};
02785 int index=0;
02786 if (dm) {
02787 for(int j=0; j<ny; j++){
02788 for(int i=0; i<nx; i++){
02789 int ij=j*nx+i;
02790 if(dm[ij]) {
02791 m[0]+=i; // x
02792 m[1]+=j; // y
02793 m[2]+=d[ij]; // z
02794 /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
02795 index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
02796 index++;
02797 }
02798 }
02799 }
02800 }
02801 else {
02802 for(int j=0; j<ny; j++){
02803 for(int i=0; i<nx; i++){
02804 int ij=j*nx+i;
02805 m[0]+=i; // x
02806 m[1]+=j; // y
02807 m[2]+=d[ij]; // z
02808 /*printf("index=%d/%d\ti,j=%d,%d\tval=%g\txm,ym,zm=%g,%g,%g\n", \
02809 index,count,i,j,d[ij],m[0]/(index+1),m[1]/(index+1),m[2]/(index+1));*/
02810 index++;
02811 }
02812 }
02813 }
02814
02815 for(int i=0; i<3; i++) m[i]/=count; // compute center of the plane
02816
02817 index=0;
02818 if (dm) {
02819 for(int j=0; j<ny; j++){
02820 for(int i=0; i<nx; i++){
02821 int ij=j*nx+i;
02822 if(dm[ij]) {
02823 //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
02824 gsl_matrix_set(A,index,0,i-m[0]);
02825 gsl_matrix_set(A,index,1,j-m[1]);
02826 gsl_matrix_set(A,index,2,d[ij]-m[2]);
02827 index++;
02828 }
02829 }
02830 }
02831 mask->update();
02832 }
02833 else {
02834 for(int j=0; j<ny; j++){
02835 for(int i=0; i<nx; i++){
02836 int ij=j*nx+i;
02837 //printf("index=%d/%d\ti,j=%d,%d\tval=%g\n",index,count,i,j,d[index]);
02838 gsl_matrix_set(A,index,0,i-m[0]);
02839 gsl_matrix_set(A,index,1,j-m[1]);
02840 gsl_matrix_set(A,index,2,d[ij]-m[2]);
02841 index++;
02842 }
02843 }
02844 }
02845
02846 // SVD decomposition and use the V vector associated with smallest singular value as the plan normal
02847 gsl_linalg_SV_decomp_jacobi(A, V, S);
02848
02849 double n[3];
02850 for(int i=0; i<3; i++) n[i] = gsl_matrix_get(V, i, 2);
02851
02852 #ifdef DEBUG
02853 printf("S=%g,%g,%g\n",gsl_vector_get(S,0), gsl_vector_get(S,1), gsl_vector_get(S,2));
02854 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));
02855 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));
02856 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));
02857 printf("Fitted plane: p0=%g,%g,%g\tn=%g,%g,%g\n",m[0],m[1],m[2],n[0],n[1],n[2]);
02858 #endif
02859
02860 int changeZero = 0;
02861 if (params.has_key("changeZero")) changeZero = params["changeZero"];
02862 if (changeZero) {
02863 for(int j=0; j<nx; j++){
02864 for(int i=0; i<ny; i++){
02865 int ij = j*nx+i;
02866 d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
02867 }
02868 }
02869 }
02870 else {
02871 for(int j=0; j<nx; j++){
02872 for(int i=0; i<ny; i++){
02873 int ij = j*nx+i;
02874 if(d[ij]) d[ij]-=static_cast<float>(-((i-m[0])*n[0]+(j-m[1])*n[1])/n[2]+m[2]);
02875 }
02876 }
02877 }
02878 image->update();
02879 // set return plane parameters
02880 vector< float > planeParam;
02881 planeParam.resize(6);
02882 for(int i=0; i<3; i++) planeParam[i] = static_cast<float>(n[i]);
02883 for(int i=0; i<3; i++) planeParam[i+3] = static_cast<float>(m[i]);
02884 params["planeParam"]=EMObject(planeParam);
02885 }
|
|
|
Definition at line 129 of file processor.cpp. |
1.3.9.1