00001
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037 #include <ctime>
00038 #include <memory>
00039 #include "emdata.h"
00040 #include "analyzer.h"
00041 #include "sparx/analyzer_sparx.h"
00042 #include "util.h"
00043 #include "cmp.h"
00044 #include "sparx/lapackblas.h"
00045 #include "sparx/varimax.h"
00046
00047 using namespace EMAN;
00048
00049 namespace EMAN {
00050
00051 const string PCAsmall::NAME = "pca";
00052 const string PCAlarge::NAME = "pca_large";
00053 const string varimax::NAME = "varimax";
00054 const string KMeansAnalyzer::NAME = "kmeans";
00055 const string SVDAnalyzer::NAME = "svd_gsl";
00056
00057 template <> Factory < Analyzer >::Factory()
00058 {
00059 force_add<PCAsmall>();
00060 force_add<PCAlarge>();
00061 force_add<varimax>();
00062 force_add<KMeansAnalyzer>();
00063 force_add<SVDAnalyzer>();
00064 }
00065
00066 }
00067
00068 int Analyzer::insert_images_list(vector<EMData *> image_list)
00069 {
00070 vector<EMData *>::const_iterator iter;
00071 for(iter=image_list.begin(); iter!=image_list.end(); ++iter) {
00072 insert_image(*iter);
00073 }
00074 return 0;
00075 }
00076
00077 void KMeansAnalyzer::set_params(const Dict & new_params)
00078 {
00079 params = new_params;
00080 if (params.has_key("ncls")) ncls = params["ncls"];
00081 if (params.has_key("maxiter"))maxiter = params["maxiter"];
00082 if (params.has_key("minchange"))minchange = params["minchange"];
00083 if (params.has_key("mininclass"))mininclass = params["mininclass"];
00084 if (params.has_key("slowseed"))slowseed = params["slowseed"];
00085 if (params.has_key("verbose"))verbose = params["verbose"];
00086 if (params.has_key("calcsigmamean")) calcsigmamean=params["calcsigmamean"];
00087
00088 }
00089
00090 vector<EMData *> KMeansAnalyzer::analyze()
00091 {
00092 if (ncls<=1) return vector<EMData *>();
00093
00094
00095
00096 int nptcl=images.size();
00097 int nclstot=ncls;
00098 if (calcsigmamean) centers.resize(nclstot*2);
00099 else centers.resize(nclstot);
00100 if (mininclass<1) mininclass=1;
00101
00102 for (int i=0; i<nptcl; i++) images[i]->set_attr("is_ok_center",(int)5);
00103
00104 if (slowseed) {
00105 if (ncls>25) slowseed=ncls/25+1;
00106
00107
00108 }
00109
00110 for (int i=0; i<ncls; i++) {
00111
00112 centers[i]=images[Util::get_irand(0,nptcl-1)]->copy();
00113
00114 }
00115
00116 if (calcsigmamean) {
00117 for (int i=nclstot; i<nclstot*2; i++) centers[i]=new EMData(images[0]->get_xsize(),images[0]->get_ysize(),images[0]->get_zsize());
00118 }
00119
00120
00121 for (int i=0; i<maxiter; i++) {
00122 nchanged=0;
00123 reclassify();
00124 if (verbose) printf("iter %d> %d (%d)\n",i,nchanged,ncls);
00125 if (nchanged<minchange && ncls==nclstot) break;
00126 update_centers();
00127
00128 if (slowseed && i%3==2 && ncls<nclstot) {
00129 for (int j=0; j<slowseed && ncls<nclstot; j++) {
00130 centers[ncls]=0;
00131 ncls++;
00132 }
00133 reseed();
00134 }
00135 }
00136 update_centers(calcsigmamean);
00137
00138 return centers;
00139 }
00140
00141 void KMeansAnalyzer::update_centers(int sigmas) {
00142 int nptcl=images.size();
00143
00144 int * repr = new int[ncls];
00145
00146 for (int i=0; i<ncls; i++) {
00147 centers[i]->to_zero();
00148 if (sigmas) centers[i+ncls]->to_zero();
00149 repr[i]=0;
00150 }
00151
00152
00153 for (int i=0; i<nptcl; i++) {
00154 int cid=images[i]->get_attr("class_id");
00155 if ((int)images[i]->get_attr("is_ok_center")>0) {
00156 centers[cid]->add(*images[i]);
00157 if (sigmas) centers[cid+ncls]->addsquare(*images[i]);
00158 repr[cid]++;
00159 }
00160 }
00161
00162 for (int i=0; i<ncls; i++) {
00163
00164 if (repr[i]<mininclass) {
00165
00166
00167 for (int j=0; j<nptcl; j++) {
00168 if ((int)images[j]->get_attr("class_id")==i) images[i]->set_attr("is_ok_center",(int)images[i]->get_attr("is_ok_center")-1);
00169 }
00170
00171 delete centers[i];
00172 centers[i]=0;
00173 repr[i]=0;
00174 }
00175
00176 else {
00177 centers[i]->mult((float)1.0/(float)(repr[i]));
00178 centers[i]->set_attr("ptcl_repr",repr[i]);
00179 if (sigmas) {
00180 centers[i+ncls]->mult((float)1.0/(float)(repr[i]));
00181 centers[i+ncls]->subsquare(*centers[i]);
00182 centers[i+ncls]->process("math.sqrt");
00183 centers[i+ncls]->mult((float)1.0/(float)sqrt((float)repr[i]));
00184 }
00185
00186 }
00187 if (verbose>1) printf("%d(%d)\t",i,(int)repr[i]);
00188 }
00189
00190 if (verbose>1) printf("\n");
00191
00192 reseed();
00193
00194 delete [] repr;
00195 }
00196
00197
00198 void KMeansAnalyzer::reseed() {
00199 int nptcl=images.size();
00200 int i,j;
00201
00202
00203 for (i=0; i<ncls; i++) {
00204 if (!centers[i]) break;
00205 }
00206 if (i==ncls) return;
00207
00208
00209 vector<int> goodcen;
00210 for (int i=0; i<nptcl; i++) if ((int)images[i]->get_attr("is_ok_center")>0) goodcen.push_back(i);
00211
00212 if (goodcen.size()==0) throw UnexpectedBehaviorException("Kmeans ran out of valid center particles with the provided parameters");
00213
00214
00215 for (i=0; i<ncls; i++) {
00216 if (centers[i]) continue;
00217 j=Util::get_irand(0,goodcen.size()-1);
00218 centers[i]=images[j]->copy();
00219 centers[i]->set_attr("ptcl_repr",1);
00220 printf("reseed %d -> %d\n",i,j);
00221 }
00222
00223
00224 }
00225
00226
00227 void KMeansAnalyzer::reclassify() {
00228 int nptcl=images.size();
00229
00230 Cmp *c = Factory < Cmp >::get("sqeuclidean");
00231 for (int i=0; i<nptcl; i++) {
00232 float best=1.0e38f;
00233 int bestn=0;
00234 for (int j=0; j<ncls; j++) {
00235 float d=c->cmp(images[i],centers[j]);
00236
00237 if (d<best) { best=d; bestn=j; }
00238 }
00239 int oldn=images[i]->get_attr_default("class_id",0);
00240 if (oldn!=bestn) nchanged++;
00241 images[i]->set_attr("class_id",bestn);
00242 }
00243 delete c;
00244 }
00245
00246 #define covmat(i,j) covmat[ ((j)-1)*nx + (i)-1 ]
00247 #define imgdata(i) imgdata[ (i)-1 ]
00248 int PCAsmall::insert_image(EMData * image)
00249 {
00250 if(mask==0)
00251 throw NullPointerException("Null mask image pointer, set_params() first");
00252
00253 EMData *maskedimage = Util::compress_image_mask(image,mask);
00254
00255 int nx = maskedimage->get_xsize();
00256 float *imgdata = maskedimage->get_data();
00257 if (nx != ncov) {
00258 fprintf(stderr,"insert_image: something is wrong...\n");
00259 exit(1);
00260 }
00261
00262
00263 nimages++;
00264 for (int j = 1; j <= nx; j++)
00265 for (int i = 1; i<=nx; i++) {
00266 covmat(i,j) += imgdata(i)*imgdata(j);
00267 }
00268
00269 EMDeletePtr(maskedimage);
00270 return 0;
00271 }
00272 #undef covmat
00273
00274 #define eigvec(i,j) eigvec[(j)*ncov + (i)]
00275 vector<EMData*> PCAsmall::analyze()
00276 {
00277 float *eigvec;
00278 int status = 0;
00279
00280 eigval = (float*)calloc(ncov,sizeof(float));
00281 eigvec = (float*)calloc(ncov*ncov,sizeof(float));
00282 status = Util::coveig(ncov, covmat, eigval, eigvec);
00283
00284
00285
00286
00287 EMData *eigenimage = new EMData();
00288 eigenimage->set_size(ncov,1,1);
00289 float *rdata = eigenimage->get_data();
00290 for (int j = 1; j<= nvec; j++) {
00291 for (int i = 0; i < ncov; i++) rdata[i] = eigvec(i,ncov-j);
00292
00293 EMData* recons_eigvec = Util::reconstitute_image_mask(eigenimage,mask);
00294 recons_eigvec->set_attr( "eigval", eigval[j-1] );
00295 images.push_back(recons_eigvec);
00296 }
00297
00298 free(eigvec);
00299 EMDeletePtr(eigenimage);
00300
00301 return images;
00302 }
00303 #undef eigvec
00304
00305 void PCAsmall::set_params(const Dict & new_params)
00306 {
00307 params = new_params;
00308 mask = params["mask"];
00309 nvec = params["nvec"];
00310
00311
00312
00313 EMData *dummy = new EMData();
00314
00315 int nx = mask->get_xsize();
00316 int ny = mask->get_ysize();
00317 int nz = mask->get_zsize();
00318
00319 dummy->set_size(nx,ny,nz);
00320
00321 EMData *dummy1d = Util::compress_image_mask(dummy,mask);
00322 ncov = dummy1d->get_xsize();
00323 EMDeletePtr(dummy);
00324 EMDeletePtr(dummy1d);
00325
00326
00327 nimages = 0;
00328 covmat = (float*)calloc(ncov*ncov,sizeof(float));
00329 }
00330
00331
00332
00333
00334 int PCAlarge::insert_image(EMData * image)
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 }
00357
00358 void PCAlarge::set_params(const Dict & new_params)
00359 {
00360 params = new_params;
00361 mask = params["mask"];
00362 nvec = params["nvec"];
00363
00364
00365
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
00381 nimages = 0;
00382 }
00383
00384 #define qmat(i,j) qmat[((j)-1)*kstep + (i) -1]
00385 #define diag(i) diag[(i)-1]
00386 #define rdata(i) rdata[(i)-1]
00387 #define eigvec(i,j) eigvec[((j)-1)*ncov + (i)-1]
00388 #define eigval(i) eigval[(i)-1]
00389
00390 vector<EMData*> PCAlarge::analyze()
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
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
00416 status = Lanczos(scratchfile, &kstep, diag, subdiag,
00417 vmat, &resnrm);
00418
00419
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
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
00437 sstevd_(jobz, &kstep, diag, subdiag, qmat, &kstep, work, &lwork,
00438 iwork, &liwork, &info);
00439
00440
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
00449
00450
00451
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
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 }
00478 #undef qmat
00479 #undef diag
00480 #undef rdata
00481 #undef eigvec
00482 #undef eigval
00483
00484 #define TOL 1e-7
00485 #define V(i,j) V[((j)-1)*imgsize + (i) - 1]
00486 #define v0(i) v0[(i)-1]
00487 #define Av(i) Av[(i)-1]
00488 #define subdiag(i) subdiag[(i)-1]
00489 #define diag(i) diag[(i)-1]
00490 #define hvec(i) hvec[(i)-1]
00491
00492 int PCAlarge::Lanczos(const string &maskedimages, int *kstep,
00493 float *diag, float *subdiag, float *V,
00494 float *beta)
00495 {
00496
00497
00498
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517
00518
00519
00520
00521
00522
00523
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;
00539 goto EXIT;
00540 }
00541
00542 imgsize = ncov;
00543 if (nimages <= 0) {
00544 status = 3;
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
00562 for ( i = 1; i <= imgsize; i++)
00563 {
00564 v0(i) = 1.0;
00565 Av(i) = 0.0;
00566 }
00567
00568
00569 *beta = snrm2_(&imgsize, v0, &ione);
00570 for (i = 1; i<=imgsize; i++)
00571 V(i,1) = v0(i) / (*beta);
00572
00573
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
00587 diag(1) = sdot_(&imgsize, V, &ione, Av, &ione);
00588 alpha = -diag(1);
00589 saxpy_(&imgsize, &alpha, V, &ione, Av, &ione);
00590
00591
00592 for ( iter = 2 ; iter <= *kstep ; iter++ ) {
00593 *beta = snrm2_(&imgsize, Av, &ione);
00594
00595 if (*beta < TOL) {
00596
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
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
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
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 }
00645 #undef v0
00646 #undef Av
00647 #undef V
00648 #undef hvec
00649 #undef diag
00650 #undef subdiag
00651 #undef TOL
00652
00653 void varimax::set_params(const Dict & new_params)
00654 {
00655 params = new_params;
00656 m_mask = params["mask"];
00657
00658
00659
00660 EMData *dummy = new EMData();
00661
00662 int nx = m_mask->get_xsize();
00663 int ny = m_mask->get_ysize();
00664 int nz = m_mask->get_zsize();
00665
00666 dummy->set_size(nx,ny,nz);
00667
00668 EMData *dummy1d = Util::compress_image_mask(dummy,m_mask);
00669
00670 m_nlen = dummy1d->get_xsize();
00671 m_nfac = 0;
00672
00673 EMDeletePtr(dummy);
00674 EMDeletePtr(dummy1d);
00675 }
00676
00677 int varimax::insert_image(EMData* image)
00678 {
00679 if(m_mask==0)
00680 throw NullPointerException("Null mask image pointer, set_params() first");
00681
00682 EMData* img1d = Util::compress_image_mask(image,m_mask);
00683
00684 m_data.insert( m_data.end(), img1d->get_data(), img1d->get_data() + m_nlen );
00685
00686 m_nfac++;
00687
00688 Assert( (int)m_data.size() == m_nfac*m_nlen);
00689
00690 return 0;
00691 }
00692
00693 vector<EMData*> varimax::analyze()
00694 {
00695 int itmax = 10000;
00696 float eps = 1e-4f;
00697 int verbose = 1;
00698 float params[4];
00699 params[0] = 1.0;
00700 varmx( &m_data[0], m_nlen, m_nfac, IVARIMAX, params, NULL, itmax, eps, verbose);
00701
00702 vector<EMData*> images;
00703
00704 EMData* img1d = new EMData();
00705 img1d->set_size(m_nlen, 1, 1);
00706 for( int i=0; i < m_nfac; ++i )
00707 {
00708 float* imgdata = img1d->get_data();
00709
00710 int offset = i * m_nlen;
00711 for( int i=0; i < m_nlen; ++i )
00712 {
00713 imgdata[i] = m_data[offset+i];
00714 }
00715
00716 EMData* img = Util::reconstitute_image_mask(img1d,m_mask);
00717 images.push_back(img);
00718 }
00719
00720 EMDeletePtr(img1d);
00721
00722 return images;
00723 }
00724
00725 int SVDAnalyzer::insert_image(EMData * image)
00726 {
00727 if (mask==0)
00728 throw NullPointerException("Null mask image pointer, set_params() first");
00729
00730
00731 size_t totpix=mask->get_xsize()*mask->get_ysize()*mask->get_zsize();
00732 float *d=image->get_data();
00733 float *md=mask ->get_data();
00734 for (size_t i=0,j=0; i<totpix; ++i) {
00735 if (md[i]) {
00736 gsl_matrix_set(A,j,nsofar,d[i]);
00737 j++;
00738 }
00739 }
00740 nsofar++;
00741
00742 return 0;
00743 }
00744 #undef covmat
00745
00746 #define eigvec(i,j) eigvec[(j)*ncov + (i)]
00747 vector<EMData*> SVDAnalyzer::analyze()
00748 {
00749
00750 gsl_vector *work=gsl_vector_alloc(nimg);
00751 gsl_vector *S=gsl_vector_alloc(nimg);
00752 gsl_matrix *V=gsl_matrix_alloc(nimg,nimg);
00753 gsl_matrix *X=gsl_matrix_alloc(nimg,nimg);
00754
00755
00756
00757 gsl_linalg_SV_decomp_mod (A,X, V, S, work);
00758
00759
00760 vector<EMData*> ret;
00761
00762 float *md=mask->get_data();
00763 size_t totpix=mask->get_xsize()*mask->get_ysize()*mask->get_zsize();
00764 for (int k=0; k<nvec; k++) {
00765 EMData *img = new EMData;
00766 img->set_size(mask->get_xsize(),mask->get_ysize(),mask->get_zsize());
00767
00768 float *d=img->get_data();
00769 for (size_t i=0,j=0; i<totpix; ++i) {
00770 if (md[i]) {
00771 d[i]=(float)gsl_matrix_get(A,j,k);
00772 j++;
00773 }
00774 }
00775 img->set_attr( "eigval", gsl_vector_get(S,k));
00776 ret.push_back(img);
00777 }
00778
00779 gsl_vector_free(work);
00780 gsl_vector_free(S);
00781 gsl_matrix_free(V);
00782 gsl_matrix_free(X);
00783
00784 gsl_matrix_free(A);
00785 A=NULL;
00786 mask=NULL;
00787
00788 return ret;
00789 }
00790
00791 void SVDAnalyzer::set_params(const Dict & new_params)
00792 {
00793 params = new_params;
00794 mask = params["mask"];
00795 nvec = params["nvec"];
00796 nimg = params["nimg"];
00797
00798
00799 pixels=0;
00800 size_t totpix=mask->get_xsize()*mask->get_ysize()*mask->get_zsize();
00801 float *d=mask->get_data();
00802 for (size_t i=0; i<totpix; ++i) if (d[i]) ++pixels;
00803
00804 printf("%d,%d\n",pixels,nimg);
00805 A=gsl_matrix_alloc(pixels,nimg);
00806 nsofar=0;
00807 }
00808
00809
00810 void EMAN::dump_analyzers()
00811 {
00812 dump_factory < Analyzer > ();
00813 }
00814
00815 map<string, vector<string> > EMAN::dump_analyzers_list()
00816 {
00817 return dump_factory_list < Analyzer > ();
00818 }
00819
00820
00821
00822
00823
00824
00825