EMAN2
analyzer.h
Go to the documentation of this file.
1/*
2 * Author: Steven Ludtke, 04/10/2003 (sludtke@bcm.edu)
3 * Copyright (c) 2000-2006 Baylor College of Medicine
4 *
5 * This software is issued under a joint BSD/GNU license. You may use the
6 * source code in this file under either license. However, note that the
7 * complete EMAN2 and SPARX software packages have some GPL dependencies,
8 * so you are responsible for compliance with the licenses of these packages
9 * if you opt to use BSD licensing. The warranty disclaimer below holds
10 * in either instance.
11 *
12 * This complete copyright notice must be included in any revised version of the
13 * source code. Additional authorship citations may be added, but existing
14 * author citations must be preserved.
15 *
16 * This program is free software; you can redistribute it and/or modify
17 * it under the terms of the GNU General Public License as published by
18 * the Free Software Foundation; either version 2 of the License, or
19 * (at your option) any later version.
20 *
21 * This program is distributed in the hope that it will be useful,
22 * but WITHOUT ANY WARRANTY; without even the implied warranty of
23 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 * GNU General Public License for more details.
25 *
26 * You should have received a copy of the GNU General Public License
27 * along with this program; if not, write to the Free Software
28 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
29 *
30 * */
31
32#ifndef eman_analyzer_h__
33#define eman_analyzer_h__
34
35#include "emobject.h"
36#include <gsl/gsl_linalg.h>
37using std::vector;
38
39namespace EMAN
40{
41 class EMData;
42
57 {
58 public:
60
61 virtual ~Analyzer()
62 {}
63
68 virtual int insert_image(EMData * image) = 0;
69
74 virtual int insert_images_list(vector<EMData *> image_list);
75
79 virtual vector<EMData*> analyze() = 0;
80
84 virtual string get_name() const = 0;
85
89 virtual string get_desc() const = 0;
90
94 virtual void set_params(const Dict & new_params)
95 {
96 params = new_params;
97 }
98
102 virtual Dict get_params() const
103 {
104 return params;
105 }
106
113 virtual TypeDict get_param_types() const = 0;
114
115 protected:
116 mutable Dict params;
117 vector<EMData *> images;
118 };
119
128 {
129 public:
131
132 virtual int insert_image(EMData *image) {
133 images.push_back(image);
134 if (images.size()>1) { printf("InertiaMatrixAnalyzer takes only a single image\n"); return 1; }
135 return 0;
136 }
137
138 virtual vector<EMData*> analyze();
139
140 string get_name() const
141 {
142 return NAME;
143 }
144
145 string get_desc() const
146 {
147 return "Compute Inertia matrix for a volume";
148 }
149
150 static Analyzer * NEW()
151 {
152 return new InertiaMatrixAnalyzer();
153 }
154
156 {
157 TypeDict d;
158 d.put("verbose", EMObject::INT, "Display progress if set, more detail with larger numbers");
159 return d;
160 }
161
162 static const string NAME;
163
164 protected:
166 vector<EMData *> ret; // This will contain only a single image
167 };
168
178 {
179 public:
181
182 virtual int insert_image(EMData *image) {
183 images.push_back(image);
184 if (images.size()>1) { printf("ShapeAnalyzer takes only a single image\n"); return 1; }
185 return 0;
186 }
187
188 virtual vector<EMData*> analyze();
189
190 string get_name() const
191 {
192 return NAME;
193 }
194
195 string get_desc() const
196 {
197 return "Experimental. Computes a set of values characterizing a 3-D volume. Returns a 3x2x1 image containing X, Y and Z axial distributions using axis squared and axis linear weighting.";
198 }
199
200 static Analyzer * NEW()
201 {
202 return new ShapeAnalyzer();
203 }
204
206 {
207 TypeDict d;
208 d.put("verbose", EMObject::INT, "Display progress if set, more detail with larger numbers");
209 return d;
210 }
211
212 static const string NAME;
213
214 protected:
216 vector<EMData *> ret; // This will contain only a single image
217 };
218
219
220
237 {
238 public:
240
241 virtual int insert_image(EMData *image) {
242 images.push_back(image);
243 return 0;
244 }
245
246 virtual vector<EMData*> analyze();
247
248 string get_name() const
249 {
250 return NAME;
251 }
252
253 string get_desc() const
254 {
255 return "k-means classification";
256 }
257
258 static Analyzer * NEW()
259 {
260 return new KMeansAnalyzer();
261 }
262
263 void set_params(const Dict & new_params);
264
266 {
267 TypeDict d;
268 d.put("verbose", EMObject::INT, "Display progress if set, more detail with larger numbers (9 max)");
269 d.put("seedmode",EMObject::INT, "How to generate initial seeds. 0 - random element (default), 1 - max sum, min sum, linear");
270 d.put("ncls", EMObject::INT, "number of desired classes");
271 d.put("maxiter", EMObject::INT, "maximum number of iterations (default=100)");
272 d.put("minchange", EMObject::INT, "Terminate if fewer than minchange members move in an iteration");
273 d.put("mininclass", EMObject::INT, "Minumum number of particles to keep a class as good (not enforced at termination");
274 d.put("slowseed",EMObject::INT, "Instead of seeding all classes at once, it will gradually increase the number of classes by adding new seeds in groups with large standard deviations");
275 d.put("outlierclass",EMObject::INT, "The last class will be reserved for outliers. Any class containing fewer than n particles will be permanently moved to the outlier group. default = disabled");
276 d.put("calcsigmamean",EMObject::INT, "Computes standard deviation of the mean image for each class-average (center), and returns them at the end of the list of centers");
277 return d;
278 }
279
280 static const string NAME;
281
282 protected:
283 void update_centers(int sigmas=0);
284 void reclassify();
285 void reseed();
286 void resort();
287
288 vector<EMData *> centers;
289 int ncls; //number of current classes
290 int nclstot; //number of desired classes
299
300 };
301
307 class SVDAnalyzer : public Analyzer
308 {
309 public:
310 SVDAnalyzer() : mask(0), nvec(0), nimg(0), A(NULL) {}
311
312 virtual int insert_image(EMData * image);
313
314 virtual int insert_images_list(vector<EMData *> image_list) {
315 vector<EMData*>::const_iterator iter;
316 for(iter=image_list.begin(); iter!=image_list.end(); ++iter) {
317 images.push_back(*iter);
318 }
319 return 0;
320 }
321
322 virtual vector<EMData*> analyze();
323
324 string get_name() const
325 {
326 return NAME;
327 }
328
329 string get_desc() const
330 {
331 return "Singular Value Decomposition from GSL. Comparable to pca";
332 }
333
334 static Analyzer * NEW()
335 {
336 return new SVDAnalyzer();
337 }
338
339 void set_params(const Dict & new_params);
340
342 {
343 TypeDict d;
344 d.put("mask", EMObject::EMDATA, "mask image");
345 d.put("nvec", EMObject::INT, "number of desired basis vectors");
346 d.put("nimg", EMObject::INT, "total number of input images, required even with insert_image()");
347 return d;
348 }
349
350 static const string NAME;
351
352 protected:
354 int nvec; //number of desired principal components
355 int pixels; // pixels under the mask
356 int nimg; // number of input images
357
358 private:
360 gsl_matrix *A;
361 };
362
363
369 {
370 public:
372
373 virtual int insert_image(EMData *image) {
374 images.push_back(image);
375 return 0;
376 }
377
378 virtual vector<EMData*> analyze();
379
380 string get_name() const
381 {
382 return NAME;
383 }
384
385 string get_desc() const
386 {
387 return "Calculate the circular average around the center in real space";
388 }
389
390 static Analyzer * NEW()
391 {
392 return new CircularAverageAnalyzer();
393 }
394
396 {
397 TypeDict d;
398 d.put("verbose", EMObject::INT, "Display progress if set, more detail with larger numbers");
399 d.put("maxr", EMObject::INT, "Maximum radius.");
400 d.put("step", EMObject::INT, "Thickness of the ring.");
401 return d;
402 }
403
404 static const string NAME;
405
406 protected:
408 vector<EMData *> ret; // This will contain only a single image
409 };
410
412
413 void dump_analyzers();
414 map<string, vector<string> > dump_analyzers_list();
415}
416
417#endif //eman_analyzer_h__
Analyzer class defines a way to take a List of images as input, and returns a new List of images.
Definition: analyzer.h:57
virtual int insert_image(EMData *image)=0
insert a image to the list of input images
virtual ~Analyzer()
Definition: analyzer.h:61
virtual TypeDict get_param_types() const =0
Get Analyzer parameter information in a dictionary.
virtual vector< EMData * > analyze()=0
main function for Analyzer, analyze input images and create output images
virtual Dict get_params() const
Get the Reconstructor's parameters in a key/value dictionary.
Definition: analyzer.h:102
virtual void set_params(const Dict &new_params)
Set the Analyzer parameters using a key/value dictionary.
Definition: analyzer.h:94
virtual string get_desc() const =0
Get the Analyzer's description.
virtual string get_name() const =0
Get the Analyzer's name.
vector< EMData * > images
Definition: analyzer.h:117
virtual int insert_images_list(vector< EMData * > image_list)
insert a list of images to the list of input images
Definition: analyzer.cpp:69
Calculate the circular average around the center in real space.
Definition: analyzer.h:369
TypeDict get_param_types() const
Get Analyzer parameter information in a dictionary.
Definition: analyzer.h:395
virtual vector< EMData * > analyze()
main function for Analyzer, analyze input images and create output images
Definition: analyzer.cpp:1034
virtual int insert_image(EMData *image)
insert a image to the list of input images
Definition: analyzer.h:373
string get_desc() const
Get the Analyzer's description.
Definition: analyzer.h:385
vector< EMData * > ret
Definition: analyzer.h:408
static Analyzer * NEW()
Definition: analyzer.h:390
static const string NAME
Definition: analyzer.h:404
string get_name() const
Get the Analyzer's name.
Definition: analyzer.h:380
Dict is a dictionary to store <string, EMObject> pair.
Definition: emobject.h:385
EMData stores an image's data and defines core image processing routines.
Definition: emdata.h:82
Factory is used to store objects to create new instances.
Definition: emobject.h:725
Inertia Matrix computer Computes the Inertia Matrix for a 3-D volume.
Definition: analyzer.h:128
static const string NAME
Definition: analyzer.h:162
virtual int insert_image(EMData *image)
insert a image to the list of input images
Definition: analyzer.h:132
string get_desc() const
Get the Analyzer's description.
Definition: analyzer.h:145
string get_name() const
Get the Analyzer's name.
Definition: analyzer.h:140
static Analyzer * NEW()
Definition: analyzer.h:150
TypeDict get_param_types() const
Get Analyzer parameter information in a dictionary.
Definition: analyzer.h:155
vector< EMData * > ret
Definition: analyzer.h:166
virtual vector< EMData * > analyze()
main function for Analyzer, analyze input images and create output images
Definition: analyzer.cpp:78
KMeansAnalyzer Performs k-means classification on a set of input images (shape/size arbitrary) return...
Definition: analyzer.h:237
void set_params(const Dict &new_params)
Set the Analyzer parameters using a key/value dictionary.
Definition: analyzer.cpp:161
string get_desc() const
Get the Analyzer's description.
Definition: analyzer.h:253
vector< EMData * > centers
Definition: analyzer.h:288
void update_centers(int sigmas=0)
Definition: analyzer.cpp:255
static Analyzer * NEW()
Definition: analyzer.h:258
virtual int insert_image(EMData *image)
insert a image to the list of input images
Definition: analyzer.h:241
string get_name() const
Get the Analyzer's name.
Definition: analyzer.h:248
virtual vector< EMData * > analyze()
main function for Analyzer, analyze input images and create output images
Definition: analyzer.cpp:174
TypeDict get_param_types() const
Get Analyzer parameter information in a dictionary.
Definition: analyzer.h:265
static const string NAME
Definition: analyzer.h:280
Singular Value Decomposition from GSL.
Definition: analyzer.h:308
static Analyzer * NEW()
Definition: analyzer.h:334
virtual int insert_image(EMData *image)
insert a image to the list of input images
Definition: analyzer.cpp:937
virtual vector< EMData * > analyze()
main function for Analyzer, analyze input images and create output images
Definition: analyzer.cpp:959
string get_desc() const
Get the Analyzer's description.
Definition: analyzer.h:329
virtual int insert_images_list(vector< EMData * > image_list)
insert a list of images to the list of input images
Definition: analyzer.h:314
void set_params(const Dict &new_params)
Set the Analyzer parameters using a key/value dictionary.
Definition: analyzer.cpp:1003
gsl_matrix * A
Definition: analyzer.h:360
TypeDict get_param_types() const
Get Analyzer parameter information in a dictionary.
Definition: analyzer.h:341
static const string NAME
Definition: analyzer.h:350
string get_name() const
Get the Analyzer's name.
Definition: analyzer.h:324
EMData * mask
Definition: analyzer.h:353
Shape characterization Computes a set of values characteristic of the shape of a volume.
Definition: analyzer.h:178
static Analyzer * NEW()
Definition: analyzer.h:200
string get_desc() const
Get the Analyzer's description.
Definition: analyzer.h:195
virtual int insert_image(EMData *image)
insert a image to the list of input images
Definition: analyzer.h:182
vector< EMData * > ret
Definition: analyzer.h:216
static const string NAME
Definition: analyzer.h:212
string get_name() const
Get the Analyzer's name.
Definition: analyzer.h:190
TypeDict get_param_types() const
Get Analyzer parameter information in a dictionary.
Definition: analyzer.h:205
virtual vector< EMData * > analyze()
main function for Analyzer, analyze input images and create output images
Definition: analyzer.cpp:122
TypeDict is a dictionary to store <string, EMObject::ObjectType> pair.
Definition: emobject.h:305
void put(const string &key, EMObject::ObjectType o, const string &desc="")
Definition: emobject.h:330
E2Exception class.
Definition: aligner.h:40
map< string, vector< string > > dump_analyzers_list()
Definition: analyzer.cpp:1027
void dump_analyzers()
Definition: analyzer.cpp:1022