#include <analyzer_sparx.h>
Inheritance diagram for EMAN::PCAlarge:
Public Member Functions | |
PCAlarge () | |
virtual int | insert_image (EMData *image) |
insert a image to the list of input images | |
virtual vector< EMData * > | analyze () |
main function for Analyzer, analyze input images and create output images | |
string | get_name () const |
Get the Analyzer's name. | |
string | get_desc () const |
Get the Analyzer's description. | |
void | set_params (const Dict &new_params) |
Set the Analyzer parameters using a key/value dictionary. | |
int | Lanczos (const string &maskedimages, int *kstep, float *diag, float *subdiag, float *V, float *beta) |
TypeDict | get_param_types () const |
Get Analyzer parameter information in a dictionary. | |
Static Public Member Functions | |
Analyzer * | NEW () |
Static Public Attributes | |
const string | NAME = "pca_large" |
Protected Attributes | |
EMData * | mask |
int | nvec |
Private Attributes | |
int | ncov |
int | nimages |
float * | eigval |
|
Definition at line 109 of file analyzer_sparx.h.
|
|
main function for Analyzer, analyze input images and create output images
Implements EMAN::Analyzer. Definition at line 390 of file analyzer.cpp. References diag, eigval, eigvec, EMDeletePtr(), EMAN::EMData::get_data(), Lanczos(), mask, ncov, nimages, nvec, nx, qmat, rdata, reconstitute_image_mask(), EMAN::EMData::set_attr(), EMAN::EMData::set_size(), sgemv_(), sstevd_(), status, and subdiag. 00391 { 00392 int status = 0; 00393 int ione = 1; 00394 float one = 1.0, zero = 0.0; 00395 char trans; 00396 float *eigvec; 00397 string scratchfile = string("maskedimages.scratch"); 00398 char command[100]; 00399 00400 printf("start analyzing..., ncov = %d\n", ncov); 00401 00402 float resnrm = 0.0; 00403 00404 if ( nvec > nimages || nvec ==0 ) nvec = nimages; 00405 int nx = ncov; 00406 00407 // the definition of kstep is purely a heuristic for right now 00408 int kstep = nvec + 20; 00409 if (kstep > nimages) kstep = nimages; 00410 00411 float *diag = new float[kstep]; 00412 float *subdiag = new float[kstep-1]; 00413 float *vmat = new float[nx*kstep]; 00414 00415 // run kstep-step Lanczos factorization 00416 status = Lanczos(scratchfile, &kstep, diag, subdiag, 00417 vmat, &resnrm); 00418 00419 // remove scratch file 00420 sprintf(command,"rm -f %s\n", scratchfile.c_str()); 00421 status = system(command); 00422 if (status != 0) { 00423 fprintf(stderr,"PCAlarge: cannot remove scratchfile\n"); 00424 } 00425 00426 char jobz[2] = "V"; 00427 float *qmat = new float[kstep*kstep]; 00428 // workspace size will be optimized later 00429 int lwork = 100 + 4*kstep + kstep*kstep; 00430 int liwork = 3+5*kstep; 00431 00432 float *work = new float[lwork]; 00433 int *iwork = new int[liwork]; 00434 int info = 0; 00435 00436 // call LAPACK tridiagonal eigensolver 00437 sstevd_(jobz, &kstep, diag, subdiag, qmat, &kstep, work, &lwork, 00438 iwork, &liwork, &info); 00439 00440 // store eigenvalues 00441 eigval = (float*)calloc(ncov,sizeof(float)); 00442 eigvec = (float*)calloc(ncov*nvec,sizeof(float)); 00443 00444 for (int j = 0; j < nvec; j++) { 00445 eigval[j] = diag(kstep-j); 00446 } 00447 00448 // for (int i=0; i<nvec; i++) printf("eigval = %11.4e\n", 00449 // eigval[i]); 00450 00451 // compute eigenvectors 00452 for (int j=1; j<=nvec; j++) { 00453 trans = 'N'; 00454 sgemv_(&trans, &nx, &kstep, &one, vmat, &nx, &qmat(1,kstep-j+1), 00455 &ione, &zero, &eigvec(1,j), &ione); 00456 } 00457 00458 // pack eigenvectors into the return imagelist 00459 EMData *eigenimage = new EMData(); 00460 eigenimage->set_size(ncov,1,1); 00461 float *rdata = eigenimage->get_data(); 00462 for (int j = 1; j<= nvec; j++) { 00463 for (int i = 1; i <= ncov; i++) 00464 rdata(i) = eigvec(i,j); 00465 00466 EMData* recons_eigvec = Util::reconstitute_image_mask(eigenimage,mask); 00467 00468 recons_eigvec->set_attr( "eigval", eigval[j-1] ); 00469 00470 images.push_back( recons_eigvec ); 00471 } 00472 00473 free(eigvec); 00474 EMDeletePtr(eigenimage); 00475 00476 return images; 00477 }
|
|
Get the Analyzer's description.
Implements EMAN::Analyzer. Definition at line 120 of file analyzer_sparx.h. 00121 { 00122 return "Principal component analysis"; 00123 }
|
|
Get the Analyzer's name. Each Analyzer is identified by a unique name.
Implements EMAN::Analyzer. Definition at line 115 of file analyzer_sparx.h. 00116 {
00117 return NAME;
00118 }
|
|
Get Analyzer parameter information in a dictionary. Each parameter has one record in the dictionary. Each record contains its name, data-type, and description.
Implements EMAN::Analyzer. Definition at line 137 of file analyzer_sparx.h. References EMAN::TypeDict::put(). 00138 { 00139 TypeDict d; 00140 d.put("mask", EMObject::EMDATA, "mask image"); 00141 d.put("nvec", EMObject::INT, "number of desired principal components"); 00142 return d; 00143 }
|
|
insert a image to the list of input images
Implements EMAN::Analyzer. Definition at line 334 of file analyzer.cpp. References compress_image_mask(), EMDeletePtr(), EMAN::EMData::get_data(), EMAN::EMData::get_xsize(), imgdata, mask, nimages, NullPointerException, and nx. 00335 { 00336 if(mask==0) 00337 throw NullPointerException("Null mask image pointer, set_params() first"); 00338 00339 EMData *maskedimage = Util::compress_image_mask(image,mask); 00340 00341 FILE *fp; 00342 string scratchfile = string("maskedimages.scratch"); 00343 00344 fp = fopen(scratchfile.c_str(),"a"); 00345 00346 int nx = maskedimage->get_xsize(); 00347 float *imgdata = maskedimage->get_data(); 00348 fwrite(imgdata, sizeof(float), nx, fp); 00349 nimages++; 00350 00351 fclose(fp); 00352 00353 EMDeletePtr(maskedimage); 00354 00355 return 0; 00356 }
|
|
Definition at line 492 of file analyzer.cpp. References Av, diag, EMDeleteArray(), hvec, imgdata, nimages, saxpy_(), sdot_(), sgemv_(), snrm2_(), status, subdiag, V, and v0. Referenced by analyze(). 00495 { 00496 /* 00497 Purpose: Compute a kstep-step Lanczos factorization 00498 on the covariant matrix X*trans(X), where 00499 X (imgstack) contains a set of images; 00500 00501 Input: 00502 imgstack (vector <EMData*>) a set of images on which PCA is 00503 to be performed; 00504 00505 kstep (int*) The maximum number of Lanczos iterations allowed. 00506 If Lanczos terminates before kstep steps 00507 is reached (an invariant subspace is found), 00508 kstep returns the number of steps taken; 00509 00510 Output: 00511 diag (float *) The projection of the covariant matrix into a 00512 Krylov subspace of dimension at most kstep. 00513 The projection is a tridiagonal matrix. The 00514 diagonal elements of this matrix is stored in 00515 the diag array. 00516 00517 subdiag (float*) The subdiagonal elements of the projection 00518 is stored here. 00519 00520 V (float *) an imgsize by kstep array that contains a 00521 set of orthonormal Lanczos basis vectors; 00522 00523 beta (float *) the residual norm of the factorization; 00524 */ 00525 int i, iter; 00526 00527 float alpha; 00528 int ione = 1; 00529 float zero = 0.0, one = 1.0, mone = -1.0; 00530 int status = 0; 00531 00532 char trans; 00533 int imgsize = 0; 00534 float *v0, *Av, *hvec, *htmp, *imgdata; 00535 FILE *fp=NULL; 00536 00537 if (nimages <= 0) { 00538 status = 2; // no image in the stack 00539 goto EXIT; 00540 } 00541 00542 imgsize = ncov; 00543 if (nimages <= 0) { 00544 status = 3; // no image in the stack 00545 goto EXIT; 00546 } 00547 00548 v0 = new float[imgsize]; 00549 Av = new float[imgsize]; 00550 hvec = new float[*kstep]; 00551 htmp = new float[*kstep]; 00552 imgdata = new float[imgsize]; 00553 00554 if (v0 == NULL || Av == NULL || hvec == NULL || 00555 htmp == NULL || imgdata == NULL) { 00556 fprintf(stderr, "Lanczos: failed to allocate v0,Av,hvec,htmp\n"); 00557 status = -1; 00558 goto EXIT; 00559 } 00560 00561 // may choose a random starting guess here 00562 for ( i = 1; i <= imgsize; i++) 00563 { 00564 v0(i) = 1.0; 00565 Av(i) = 0.0; 00566 } 00567 00568 // normalize the starting vector 00569 *beta = snrm2_(&imgsize, v0, &ione); 00570 for (i = 1; i<=imgsize; i++) 00571 V(i,1) = v0(i) / (*beta); 00572 00573 // do Av <-- A*v0, where A is a cov matrix 00574 fp = fopen(maskedimages.c_str(),"r"); 00575 if (fp==NULL) { 00576 fprintf(stderr,"Lanczos: cannot open %s\n", maskedimages.c_str()); 00577 } 00578 for (i = 0; i < nimages; i++) { 00579 fread(imgdata, sizeof(float), imgsize, fp); 00580 alpha = sdot_(&imgsize, imgdata, &ione, V, &ione); 00581 saxpy_(&imgsize, &alpha, imgdata, &ione, Av, &ione); 00582 } 00583 fclose(fp); 00584 00585 00586 // Av <--- Av - V(:,1)*V(:,1)'*Av 00587 diag(1) = sdot_(&imgsize, V, &ione, Av, &ione); 00588 alpha = -diag(1); 00589 saxpy_(&imgsize, &alpha, V, &ione, Av, &ione); 00590 00591 // main loop 00592 for ( iter = 2 ; iter <= *kstep ; iter++ ) { 00593 *beta = snrm2_(&imgsize, Av, &ione); 00594 00595 if (*beta < TOL) { 00596 // found an invariant subspace, exit 00597 *kstep = iter; 00598 break; 00599 } 00600 00601 subdiag(iter-1) = *beta; 00602 for ( i = 1 ; i <= imgsize ; i++ ) { 00603 V(i,iter) = Av(i) / (*beta); 00604 } 00605 00606 // do Av <-- A*V(:,iter), where A is a cov matrix 00607 for (i = 0; i < imgsize; i++) Av[i] = 0; 00608 fp = fopen(maskedimages.c_str(),"r"); 00609 for (i = 0; i < nimages; i++) { 00610 fread(imgdata, sizeof(float), imgsize, fp); 00611 alpha = sdot_(&imgsize, imgdata, &ione, &V(1,iter), &ione); 00612 saxpy_(&imgsize, &alpha, imgdata, &ione, Av, &ione); 00613 } 00614 fclose(fp); 00615 00616 // f <--- Av - V(:,1:iter)*V(:,1:iter)'*Av 00617 trans = 'T'; 00618 status = sgemv_(&trans, &imgsize, &iter, &one, V, &imgsize, Av, &ione, 00619 &zero , hvec , &ione); 00620 trans = 'N'; 00621 status = sgemv_(&trans, &imgsize, &iter, &mone, V, &imgsize, hvec, 00622 &ione , &one , Av, &ione); 00623 00624 // one step of reorthogonalization 00625 trans = 'T'; 00626 status = sgemv_(&trans, &imgsize, &iter, &one, V, &imgsize, Av, &ione, 00627 &zero , htmp , &ione); 00628 saxpy_(&iter, &one, htmp, &ione, hvec, &ione); 00629 trans = 'N'; 00630 status = sgemv_(&trans, &imgsize, &iter, &mone, V, &imgsize, htmp, 00631 &ione , &one , Av, &ione); 00632 diag(iter) = hvec(iter); 00633 } 00634 00635 EMDeleteArray(v0); 00636 EMDeleteArray(Av); 00637 EMDeleteArray(hvec); 00638 EMDeleteArray(htmp); 00639 EMDeleteArray(imgdata); 00640 00641 EXIT: 00642 return status; 00643 00644 }
|
|
Definition at line 125 of file analyzer_sparx.h. 00126 { 00127 return new PCAlarge(); 00128 }
|
|
Set the Analyzer parameters using a key/value dictionary.
Reimplemented from EMAN::Analyzer. Definition at line 358 of file analyzer.cpp. References compress_image_mask(), EMDeletePtr(), EMAN::EMData::get_xsize(), EMAN::EMData::get_ysize(), EMAN::EMData::get_zsize(), mask, ncov, nimages, nvec, nx, ny, and EMAN::EMData::set_size(). 00359 { 00360 params = new_params; 00361 mask = params["mask"]; 00362 nvec = params["nvec"]; 00363 00364 // count the number of pixels under the mask 00365 // (this is really ugly!!!) 00366 EMData *dummy = new EMData(); 00367 00368 int nx = mask->get_xsize(); 00369 int ny = mask->get_ysize(); 00370 int nz = mask->get_zsize(); 00371 00372 dummy->set_size(nx,ny,nz); 00373 00374 EMData *dummy1d = Util::compress_image_mask(dummy,mask); 00375 00376 ncov = dummy1d->get_xsize(); 00377 00378 EMDeletePtr(dummy); 00379 EMDeletePtr(dummy1d); 00380 // no need to allocate the covariance matrix 00381 nimages = 0; 00382 }
|
|
Definition at line 154 of file analyzer_sparx.h. Referenced by analyze(). |
|
Definition at line 148 of file analyzer_sparx.h. Referenced by analyze(), insert_image(), and set_params(). |
|
Definition at line 52 of file analyzer.cpp. |
|
Definition at line 152 of file analyzer_sparx.h. Referenced by analyze(), and set_params(). |
|
Definition at line 153 of file analyzer_sparx.h. Referenced by analyze(), insert_image(), Lanczos(), and set_params(). |
|
Definition at line 149 of file analyzer_sparx.h. Referenced by analyze(), and set_params(). |