Xmipp
v3.23.11-Nereus
|
Classes | |
class | Mask |
class | ProgMask |
Functions | |
void | apply_geo_binary_2D_mask (MultidimArray< int > &mask, const Matrix2D< double > &A) |
void | apply_geo_cont_2D_mask (MultidimArray< double > &mask, const Matrix2D< double > &A) |
Variables | |
constexpr int | INNER_MASK = 1 |
constexpr int | OUTSIDE_MASK = 2 |
constexpr int | NO_ACTIVATE = 0 |
constexpr int | ACTIVATE = 1 |
Actual masks | |
void | RaisedCosineMask (MultidimArray< double > &mask, double r1, double r2, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | RaisedCrownMask (MultidimArray< double > &mask, double r1, double r2, double pix_width, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | KaiserMask (MultidimArray< double > &mask, double delta=0.01, double Deltaw=1.0/12.0) |
void | SincMask (MultidimArray< double > &mask, double omega, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | SincKaiserMask (MultidimArray< double > &mask, double omega, double delta=0.01, double Deltaw=1.0/12.0) |
void | BlackmanMask (MultidimArray< double > &mask, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | SincBlackmanMask (MultidimArray< double > &mask, double omega, double power_percentage, double x0=0, double y0=0, double z0=0) |
void | BinaryCircularMask (MultidimArray< int > &mask, double radius, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | BlobCircularMask (MultidimArray< double > &mask, double r1, blobtype blob, int mode, double x0=0, double y0=0, double z0=0) |
void | BinaryCrownMask (MultidimArray< int > &mask, double R1, double R2, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | BlobCrownMask (MultidimArray< double > &mask, double r1, double r2, blobtype blob, int mode, double x0=0, double y0=0, double z0=0) |
void | GaussianMask (MultidimArray< double > &mask, double sigma, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | BinaryDWTCircularMask2D (MultidimArray< int > &mask, double radius, int smin, int smax, const std::string &quadrant) |
void | SeparableSincKaiserMask2D (MultidimArray< double > &mask, double omega, double delta=0.01, double Deltaw=1.0/12.0) |
void | mask2D_4neig (MultidimArray< int > &mask, int value=1, int center=NO_ACTIVATE) |
void | mask2D_8neig (MultidimArray< int > &mask, int value1=1, int value2=1, int center=NO_ACTIVATE) |
void | BinaryDWTSphericalMask2D (MultidimArray< int > &mask, double radius, int smin, int smax, const std::string &quadrant) |
void | BinaryCylinderMask (MultidimArray< int > &mask, double R, double H, int mode=INNER_MASK, double x0=0, double y0=0, int z0=0) |
void | BinaryFrameMask (MultidimArray< int > &mask, int Xrect, int Yrect, int Zrect, int mode=INNER_MASK, double x0=0, double y0=0, double z0=0) |
void | BinaryConeMask (MultidimArray< int > &mask, double theta, int mode=INNER_MASK, bool centerOrigin=false) |
void | BinaryWedgeMask (MultidimArray< int > &mask, double theta0, double thetaF, const Matrix2D< double > &A, bool centerOrigin=false) |
void | mask3D_6neig (MultidimArray< int > &mask, int value=1, int center=NO_ACTIVATE) |
void | mask3D_18neig (MultidimArray< int > &mask, int value1=1, int value2=1, int center=NO_ACTIVATE) |
void | mask3D_26neig (MultidimArray< int > &mask, int value1=1, int value2=1, int value3=1, int center=NO_ACTIVATE) |
Mask Tools | |
All Mask tools work only in the overlapping area of the given image/volume and the mask in logical coordinates. Ie, if you have a mask defined from -2 to 2 and you apply it to an image defined from 0 to 63 then only those values of the mask between 0 and 2 will be applied. The rest of the image will remain untouched. This region where the mask is active within the overlapping area will be called in this documentation: active area. | |
constexpr int | COUNT_ABOVE = 1 |
constexpr int | COUNT_BELOW = 2 |
constexpr int | COUNT_BETWEEN = 3 |
template<typename T1 , typename T > | |
void | computeStats_within_binary_mask (const MultidimArray< T1 > &mask, const MultidimArray< T > &m, double &min_val, double &max_val, double &avg, double &stddev) |
void | computeStats_within_binary_mask (const MultidimArray< int > &mask, const MultidimArrayGeneric &m, double &min_val, double &max_val, double &avg, double &stddev) |
template<typename T > | |
void | apply_binary_mask (const MultidimArray< int > &mask, const MultidimArray< T > &m_in, MultidimArray< T > &m_out, T subs_val=(T) 0) |
template<typename T > | |
void | apply_cont_mask (const MultidimArray< double > &mask, const MultidimArray< T > &m_in, MultidimArray< T > &m_out) |
template<typename T > | |
void | compute_hist_within_binary_mask (const MultidimArray< int > &mask, MultidimArray< T > &v, Histogram1D &hist, int no_steps) |
template<typename T > | |
void | compute_hist_within_binary_mask (const MultidimArray< int > &mask, const MultidimArray< T > &v, Histogram1D &hist, T min, T max, int no_steps) |
template<typename T > | |
int | count_with_mask (const MultidimArray< int > &mask, const MultidimArray< T > &m, int mode, double th1, double th2) |
int | count_with_mask (const MultidimArray< int > &mask, const MultidimArray< std::complex< double > > &m, int mode, double th1, double th2) |
template<typename T > | |
void | invert_binary_mask (MultidimArray< T > &mask) |
void | rangeAdjust_within_mask (const MultidimArray< double > *mask, const MultidimArray< double > &m1, MultidimArray< double > &m2) |
#define | count_with_mask_above(mask, m, th) count_with_mask(mask, m, COUNT_ABOVE, th, 0); |
#define | count_with_mask_below(mask, m, th) count_with_mask(mask, m, COUNT_BELOW, th, 0); |
#define | count_with_mask_between(mask, m, th1, th2) count_with_mask(mask, m, COUNT_BETWEEN, th1, th2); |
#define count_with_mask_above | ( | mask, | |
m, | |||
th | |||
) | count_with_mask(mask, m, COUNT_ABOVE, th, 0); |
#define count_with_mask_below | ( | mask, | |
m, | |||
th | |||
) | count_with_mask(mask, m, COUNT_BELOW, th, 0); |
#define count_with_mask_between | ( | mask, | |
m, | |||
th1, | |||
th2 | |||
) | count_with_mask(mask, m, COUNT_BETWEEN, th1, th2); |
void apply_binary_mask | ( | const MultidimArray< int > & | mask, |
const MultidimArray< T > & | m_in, | ||
MultidimArray< T > & | m_out, | ||
T | subs_val = (T) 0 |
||
) |
Apply mask to a MultidimArray
The volume values for which the input mask is 0 are set to 0. The input and output volumes can be the same ones.
Definition at line 857 of file mask.h.
void apply_cont_mask | ( | const MultidimArray< double > & | mask, |
const MultidimArray< T > & | m_in, | ||
MultidimArray< T > & | m_out | ||
) |
Apply continuous mask to a MultidimArray
The image is multiplied by the mask. The input and output matrices can be the same ones. Only the overlapping values are affected by the mask.
Definition at line 881 of file mask.h.
void apply_geo_binary_2D_mask | ( | MultidimArray< int > & | mask, |
const Matrix2D< double > & | A | ||
) |
Apply geometric transformation to a binary (2D) mask
Definition at line 1682 of file mask.cpp.
void apply_geo_cont_2D_mask | ( | MultidimArray< double > & | mask, |
const Matrix2D< double > & | A | ||
) |
Apply geometric transformation to a continuous (2D) mask
Definition at line 1703 of file mask.cpp.
void BinaryCircularMask | ( | MultidimArray< int > & | mask, |
double | radius, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a circular mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0,z0), by default (0,0,0), is created with the radius indicated. The only two valid modes are INNER_MASK (by default) or OUTSIDE_MASK. When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 193 of file mask.cpp.
void BinaryConeMask | ( | MultidimArray< int > & | mask, |
double | theta, | ||
int | mode = INNER_MASK , |
||
bool | centerOrigin = false |
||
) |
Creates a 3D cone mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A cone with angle theta is placed parallel to the Z-axis and centered at (x0,y0,z0). The only two valid modes are INNER_MASK (by default, inside the crown is zero, outside is 1) or OUTSIDE_MASK (the negative of the crown).
When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 488 of file mask.cpp.
void BinaryCrownMask | ( | MultidimArray< int > & | mask, |
double | R1, | ||
double | R2, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a crown mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0,z0), by default (0,0,0), is created with the two radii indicated. The only two valid modes are INNER_MASK (by default, between the two radii) or OUTSIDE_MASK (the negative of the crown). It is supposed that R1 is smaller than R2.
When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 244 of file mask.cpp.
void BinaryCylinderMask | ( | MultidimArray< int > & | mask, |
double | R, | ||
double | H, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
int | z0 = 0 |
||
) |
Creates a 3D Cylinder mask for already sized masks
The mask is supposed to be already resized and with its logical origin defined. A cylinder placed logically at (z0,x0,y0), by default (0,0,0), is created with the radius and height indicated. The only two valid modes ar INNER_MASK (by default) or OUTSIDE_MASK. When entering the mask is initialiazed to 0 and then the mask is created.
void BinaryDWTCircularMask2D | ( | MultidimArray< int > & | mask, |
double | radius, | ||
int | smin, | ||
int | smax, | ||
const std::string & | quadrant | ||
) |
Binary Circular 2D mask in wavelet space
Definition at line 356 of file mask.cpp.
void BinaryDWTSphericalMask2D | ( | MultidimArray< int > & | mask, |
double | radius, | ||
int | smin, | ||
int | smax, | ||
const std::string & | quadrant | ||
) |
Creates a 3D DWT spherical for already sized masks
The mask size must be a power of 2. The radius must be expressed in pixel units corresponding to the size of the image. For instance, a 64x64x64 image might have a radius of 32 pixels to concentrate on the central part only.
If quadrant=xxx then 001,010,011,100,101,110 and 111 are generated together
void BinaryFrameMask | ( | MultidimArray< int > & | mask, |
int | Xrect, | ||
int | Yrect, | ||
int | Zrect, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a 3D frame mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A square placed logically at (x0,y0,z0), by default (0,0,0), is created with the two rectangular dimensions indicated. The only two valid modes are INNER_MASK (by default, between the two radii) or OUTSIDE_MASK (the negative of the crown). It is supposed that R1 is smaller than R2.
When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 309 of file mask.cpp.
void BinaryWedgeMask | ( | MultidimArray< int > & | mask, |
double | theta0, | ||
double | thetaF, | ||
const Matrix2D< double > & | A, | ||
bool | centerOrigin = false |
||
) |
Creates a 3D missing wedge mask
The mask is supposed to be resized and with its logical origin already set. theta0 and thetaF are the tilting angles (around the Y-axis) wherin between the data is supposed to be collected. In this region the mask will be one, outside (ie in the missing wedge) it will be zero. The mask is centered at (x0,y0,z0), and rotated with respect to euler angle matrix A.
Definition at line 524 of file mask.cpp.
void BlackmanMask | ( | MultidimArray< double > & | mask, |
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Blackman selfWindow
It receives no parameter.
Definition at line 151 of file mask.cpp.
void BlobCircularMask | ( | MultidimArray< double > & | mask, |
double | r1, | ||
blobtype | blob, | ||
int | mode, | ||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a circular mask with blob-shaped edges for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0,z0), by default (0,0,0), is created with the radius indicated. The only two valid modes are INNER_MASK (by default) or OUTSIDE_MASK. When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 219 of file mask.cpp.
void BlobCrownMask | ( | MultidimArray< double > & | mask, |
double | r1, | ||
double | r2, | ||
blobtype | blob, | ||
int | mode, | ||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a crown mask with blob-shaped edges for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0,z0), by default (0,0,0), is created with the two radii indicated. The only two valid modes are INNER_MASK (by default, between the two radii) or OUTSIDE_MASK (the negative of the crown). It is supposed that R1 is smaller than R2.
When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 278 of file mask.cpp.
void compute_hist_within_binary_mask | ( | const MultidimArray< int > & | mask, |
MultidimArray< T > & | v, | ||
Histogram1D & | hist, | ||
int | no_steps | ||
) |
Compute histogram inside mask within its minimum and maximum value (3D)
Given a volume as input, this function returns the histogram of values inside the mask within the minimum and maximum of the volume, in this way all the values in the volume are counted. The volume can be of any numerical type (short int, int, double, ...). The number of steps must always be given.
Definition at line 906 of file mask.h.
void compute_hist_within_binary_mask | ( | const MultidimArray< int > & | mask, |
const MultidimArray< T > & | v, | ||
Histogram1D & | hist, | ||
T | min, | ||
T | max, | ||
int | no_steps | ||
) |
Compute histogram inside mask within two values (3D)
Given a volume as input, this function returns the histogram of the values inside the mask within two values, the volume values outside this range are not counted. This can be used to avoid the effect of outliers which causes a "compression" in the histogram. The volume can be of any numerical type (short int, int, double, ...). The number of steps must always be given.
Definition at line 928 of file mask.h.
void computeStats_within_binary_mask | ( | const MultidimArray< T1 > & | mask, |
const MultidimArray< T > & | m, | ||
double & | min_val, | ||
double & | max_val, | ||
double & | avg, | ||
double & | stddev | ||
) |
Compute statistics in the active area
Only the statistics for values in the overlapping between the mask and the volume for those the mask is not 0 are computed.
Definition at line 799 of file mask.h.
|
inline |
int count_with_mask | ( | const MultidimArray< int > & | mask, |
const MultidimArray< T > & | m, | ||
int | mode, | ||
double | th1, | ||
double | th2 | ||
) |
Count voxels with mask and threshold.
This function returns the number of voxels in the ACTIVE area of an volume with a value:
COUNT_ABOVE: greater or equal than th1 COUNT_BELOW: smaller or equal than th1 COUNT_BETWEEN: smaller or equal than th1 and greater or equal than th2
For complex matrices the absolute value is compared.
Definition at line 983 of file mask.h.
int count_with_mask | ( | const MultidimArray< int > & | mask, |
const MultidimArray< std::complex< double > > & | m, | ||
int | mode, | ||
double | th1, | ||
double | th2 | ||
) |
Definition at line 1712 of file mask.cpp.
void GaussianMask | ( | MultidimArray< double > & | mask, |
double | sigma, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a gaussian mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0,z0), by default (0,0,0), is created with the radius indicated. The only two valid modes are INNER_MASK (by default) or OUTSIDE_MASK. Inner mask are normal sincs, and outside masks are 1 -gaussian.
When entering the mask is initialiazed to 0 and then the mask is created.
Definition at line 325 of file mask.cpp.
void invert_binary_mask | ( | MultidimArray< T > & | mask | ) |
Invert binary mask.
0's are converted in 1's and viceversa
Definition at line 1020 of file mask.h.
void KaiserMask | ( | MultidimArray< double > & | mask, |
double | delta = 0.01 , |
||
double | Deltaw = 1.0/12.0 |
||
) |
Kaiser selfWindow The mask is resized. delta=ripple (in natural units) in the pass band. Deltaw=transition bandwidth (normalized to 1.0).
Definition at line 83 of file mask.cpp.
void mask2D_4neig | ( | MultidimArray< int > & | mask, |
int | value = 1 , |
||
int | center = NO_ACTIVATE |
||
) |
Creates a 3x3 mask with value (1 by default) for those 4-neighbours of the central point (0 otherwise).
The parameter center controls whether the center pixel is set to 1 or not
Definition at line 409 of file mask.cpp.
void mask2D_8neig | ( | MultidimArray< int > & | mask, |
int | value1 = 1 , |
||
int | value2 = 1 , |
||
int | center = NO_ACTIVATE |
||
) |
Creates a 3x3 mask with value1 for those 4-neighbors of the central point and value2 for the 8 neighbours.
The parameter center controls whether the center pixel is set to 1 or not
Definition at line 417 of file mask.cpp.
void mask3D_18neig | ( | MultidimArray< int > & | mask, |
int | value1 = 1 , |
||
int | value2 = 1 , |
||
int | center = NO_ACTIVATE |
||
) |
Creates a 3x3x3 mask with value1 (1 by default) for those 6-neighbors and value2 for the 18 neighbors of the central point (0 otherwise).
The parameter center controls whether the center pixel is set to 1 or not
Definition at line 598 of file mask.cpp.
void mask3D_26neig | ( | MultidimArray< int > & | mask, |
int | value1 = 1 , |
||
int | value2 = 1 , |
||
int | value3 = 1 , |
||
int | center = NO_ACTIVATE |
||
) |
Creates a 3x3x3 mask with value1 (1 by default) for those 6-neighbors, value2 for the 18 neighbors and value3 for the 26 neighbors of the central point (0 otherwise).
The parameter center controls whether the center pixel is set to 1 or not
Definition at line 612 of file mask.cpp.
void mask3D_6neig | ( | MultidimArray< int > & | mask, |
int | value = 1 , |
||
int | center = NO_ACTIVATE |
||
) |
Creates a 3x3x3 mask with value (1 by default) for those 6-neighbors of the central point (0 otherwise).
The parameter center controls whether the center pixel is set to 1 or not
Definition at line 589 of file mask.cpp.
void RaisedCosineMask | ( | MultidimArray< double > & | mask, |
double | r1, | ||
double | r2, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a RaisedCosine mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0,z0), by default (0,0,0), is created with the radius indicated. The only two valid modes are INNER_MASK (by default) or OUTSIDE_MASK. Inner mask are normal RaisedCosines, and outside masks are 1 - RaisedCosine. When entering, the mask is initialiazed to 0 and then the mask is created.
Definition at line 36 of file mask.cpp.
void RaisedCrownMask | ( | MultidimArray< double > & | mask, |
double | r1, | ||
double | r2, | ||
double | pix_width, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a RaisedCrown mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0,y0), by default (0,0), is created within the two the radii indicated with an extra region of (pix_width) pixels. The only two valid modes are INNER_MASK (by default) or OUTSIDE_MASK. Inner mask are normal RaisedCrowns, and outside masks are 1 - RaisedCrowns. When entering, the mask is initialiazed to 0 and then the mask is created.
Definition at line 67 of file mask.cpp.
void rangeAdjust_within_mask | ( | const MultidimArray< double > * | mask, |
const MultidimArray< double > & | m1, | ||
MultidimArray< double > & | m2 | ||
) |
Range adjust within binary mask
Make the grey values of m2 fit, in L2 sense, with those in m1. Only the voxels within the mask are used to compute the linear transformation. If no mask is provided then all voxels are used.
Definition at line 1738 of file mask.cpp.
void SeparableSincKaiserMask2D | ( | MultidimArray< double > & | mask, |
double | omega, | ||
double | delta = 0.01 , |
||
double | Deltaw = 1.0/12.0 |
||
) |
Creates a 2D separable-sinc-kaiser mask, the mask is resized. This function returns a sinc mask windowed by a Kaiser selfWindow. delta=ripple (in natural units) in the pass band. Deltaw=transition bandwidth (normalized to 1). omega=low pass frequency (normalized to 1).
Definition at line 379 of file mask.cpp.
void SincBlackmanMask | ( | MultidimArray< double > & | mask, |
double | omega, | ||
double | power_percentage, | ||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a sinc-blackman mask, the mask is resized
This function returns a sinc mask windowed by a Blackman selfWindow. The selfWindow is designed to cover a certain power of the sinc
Definition at line 166 of file mask.cpp.
void SincKaiserMask | ( | MultidimArray< double > & | mask, |
double | omega, | ||
double | delta = 0.01 , |
||
double | Deltaw = 1.0/12.0 |
||
) |
Creates a radial-sinc-kaiser mask, the mask is resized. This function returns a sinc mask windowed by a Kaiser selfWindow. delta=ripple (in natural units) in the pass band. Deltaw=transition bandwidth (normalized to 1). omega=low pass frequency (normalized to 1).
Definition at line 139 of file mask.cpp.
void SincMask | ( | MultidimArray< double > & | mask, |
double | omega, | ||
int | mode = INNER_MASK , |
||
double | x0 = 0 , |
||
double | y0 = 0 , |
||
double | z0 = 0 |
||
) |
Creates a sinc mask for already sized masks
The mask is supposed to be resized and with its logical origin already set. A circle placed logically at (x0), by default (0), is created with the radius indicated. The only two valid modes are INNER_MASK (by default) or OUTSIDE_MASK. Inner mask are normal sincs, and outside masks are 1 - sinc. When entering the mask is initialiazed to 0 and then the mask is created.
Remind that sinc(w*n) is zero at n=1/w;
Definition at line 127 of file mask.cpp.