module #include "diplib/statistics.h"
Projection operators Operators that project the image data onto fewer spatial dimensions, computing image statistics.
Contents
- Reference
Functions
-
void dip::
Mean(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::String const& mode = "", dip::BooleanArray const& process = {}) - Calculates the (arithmetic) mean of the pixel values over all those dimensions which are specified by
process
. -
void dip::
Sum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the sum of the pixel values over all those dimensions which are specified by
process
. -
void dip::
GeometricMean(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the geometric mean of the pixel values over all those dimensions which are specified by
process
. -
void dip::
Product(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the product of the pixel values over all those dimensions which are specified by
process
. -
void dip::
MeanAbs(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the mean of the absolute pixel values over all those dimensions which are specified by
process
. -
void dip::
MeanModulus(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the mean of the modulus of the pixel values. Alias to
dip::MeanAbs
. -
void dip::
SumAbs(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the sum of the absolute pixel values over all those dimensions which are specified by
process
. -
void dip::
SumModulus(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the sum of the modulus of the pixel values. Alias to
dip::SumAbs
. -
void dip::
MeanSquare(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the mean of the square pixel values over all those dimensions which are specified by
process
. -
void dip::
SumSquare(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the sum of the square pixel values over all those dimensions which are specified by
process
. -
void dip::
MeanSquareModulus(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the mean of the square modulus of the pixel values over all those dimensions which are specified by
process
. -
void dip::
SumSquareModulus(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the sum of the square modulus of the pixel values over all those dimensions which are specified by
process
. -
void dip::
Variance(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::String mode = S::FAST, dip::BooleanArray const& process = {}) - Calculates the variance of the pixel values over all those dimensions which are specified by
process
. -
void dip::
StandardDeviation(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::String mode = S::FAST, dip::BooleanArray const& process = {}) - Calculates the standard deviation of the pixel values over all those dimensions which are specified by
process
. -
void dip::
Maximum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the maximum of the pixel values over all those dimensions which are specified by
process
. -
void dip::
Minimum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the minimum of the pixel values over all those dimensions which are specified by
process
. -
void dip::
MaximumAbs(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the maximum of the absolute pixel values over all those dimensions which are specified by
process
. -
void dip::
MinimumAbs(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the minimum of the absolute pixel values over all those dimensions which are specified by
process
. -
void dip::
Percentile(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::dfloat percentile = 50, dip::BooleanArray const& process = {}) - Calculates the percentile of the pixel values over all those dimensions which are specified by
process
. -
void dip::
Median(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Calculates the median of the pixel values over all those dimensions which are specified by
process
. -
void dip::
MedianAbsoluteDeviation(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Computes the median absolute deviation (MAD)
-
void dip::
All(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Determines if all pixels have non-zero values over all those dimensions which are specified by
process
. -
void dip::
Any(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::BooleanArray const& process = {}) - Determines if any pixel has a non-zero value over all those dimensions which are specified by
process
. -
void dip::
PositionMaximum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::uint dim = 0, dip::String const& mode = S::FIRST) - Calculates the position of the maximum of the pixel values in a single dimension specified by
dim
. -
void dip::
PositionMinimum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::uint dim = 0, dip::String const& mode = S::FIRST) - Calculates the position of the minimum of the pixel values in a single dimension specified by
dim
. -
void dip::
PositionPercentile(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::dfloat percentile = 50, dip::uint dim = 0, dip::String const& mode = S::FIRST) - Calculates the position of the percentile of the pixel values in a single dimension specified by
dim
. -
void dip::
PositionMedian(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::uint dim = 0, dip::String const& mode = S::FIRST) - Calculates the position of the median of the pixel values in a single dimension specified by
dim
. -
void dip::
RadialSum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::dfloat binSize = 1, dip::String const& maxRadius = S::OUTERRADIUS, dip::FloatArray const& center = {}) - Computes the radial projection of the sum of the pixel values of
in
. -
void dip::
RadialMean(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::dfloat binSize = 1, dip::String const& maxRadius = S::OUTERRADIUS, dip::FloatArray const& center = {}) - Computes the radial projection of the mean of the pixel values of
in
. -
void dip::
RadialMinimum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::dfloat binSize = 1, dip::String const& maxRadius = S::OUTERRADIUS, dip::FloatArray const& center = {}) - Computes the radial projection of the minimum of the pixel values of
in
. -
void dip::
RadialMaximum(dip::Image const& in, dip::Image const& mask, dip::Image& out, dip::dfloat binSize = 1, dip::String const& maxRadius = S::OUTERRADIUS, dip::FloatArray const& center = {}) - Computes the radial projection of the maximum of the pixel values of
in
.
Function documentation
void
dip:: Mean(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::String const& mode = "",
dip::BooleanArray const& process = {})
Calculates the (arithmetic) mean of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the mean pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the mean projection along the processing dimensions. To get the mean value of all pixels in the
image:
dip::Mean( img ).As< double >();
If mode
is "directional"
, the data in in
are assumed to be angles, and directional statistics are used.
If in
contains orientations, multiply it by 2 before applying this function, and divide the result by 2.
For directional statistics, the input must be floating point.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::Mean( img.At( mask ), mode )
is the same as dip::Mean( img, mask, mode )
.
void
dip:: Sum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the sum of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the sum of pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the sum projection along the processing dimensions. To get the sum of all pixel values in the
image:
dip::Sum( img ).As< double >();
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::Sum( img.At( mask ))
is the same as dip::Sum( img, mask )
.
void
dip:: GeometricMean(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the geometric mean of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the product of pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the product projection along the processing dimensions. To get the product of all pixel values in the
image:
dip::GeometricMean( img ).As< double >();
For tensor images, the geometric mean is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::GeometricMean( img.At( mask ))
is the same as dip::GeometricMean( img, mask )
.
void
dip:: Product(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the product of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the product of pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the product projection along the processing dimensions. To get the product of all pixel values in the
image:
dip::Product( img ).As< double >();
For tensor images, the product is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::Product( img.At( mask ))
is the same as dip::Product( img, mask )
.
void
dip:: MeanAbs(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the mean of the absolute pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the mean absolute pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the mean absolute projection along the processing dimensions. To get the mean absolute value of all pixels in the
image:
dip::MeanAbs( img ).As< double >();
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::MeanAbs( img.At( mask ))
is the same as dip::MeanAbs( img, mask )
.
void
dip:: SumAbs(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the sum of the absolute pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the sum absolute pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the sum absolute projection along the processing dimensions. To get the sum absolute value of all pixels in the
image:
dip::SumAbs( img ).As< double >();
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::SumAbs( img.At( mask ))
is the same as dip::SumAbs( img, mask )
.
void
dip:: MeanSquare(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the mean of the square pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the mean square pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the mean square projection along the processing dimensions. To get the mean square value of all pixels in the
image:
dip::MeanSquare( img ).As< double >();
For tensor images, the result is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::MeanSquare( img.At( mask ))
is the same as dip::MeanSquare( img, mask )
.
void
dip:: SumSquare(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the sum of the square pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the sum square pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the sum square projection along the processing dimensions. To get the sum square value of all pixels in the
image:
dip::SumSquare( img ).As< double >();
For tensor images, the result is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::SumSquare( img.At( mask ))
is the same as dip::SumSquare( img, mask )
.
void
dip:: MeanSquareModulus(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the mean of the square modulus of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the mean square pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the mean square projection along the processing dimensions. To get the mean square modulus value of all
pixels in the image:
dip::MeanSquareModulus( img ).As< double >();
For tensor images, the result is computed for each element independently. If img
is complex, out
is of the
corresponding floating-point type. For other input data types, this function is identical to dip::MeanSquare
.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::MeanSquareModulus( img.At( mask ))
is the same as dip::MeanSquareModulus( img, mask )
.
void
dip:: SumSquareModulus(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the sum of the square modulus of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the sum square pixel value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the sum square projection along the processing dimensions. To get the sum square modulus value of all
pixels in the image:
dip::SumSquareModulus( img ).As< double >();
For tensor images, the result is computed for each element independently. If img
is complex, out
is of the
corresponding floating-point type. For other input data types, this function is identical to dip::SumSquare
.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::SumSquareModulus( img.At( mask ))
is the same as dip::SumSquareModulus( img, mask )
.
void
dip:: Variance(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::String mode = S::FAST,
dip::BooleanArray const& process = {})
Calculates the variance of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the variance of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the variance projection along the processing dimensions. To get the variance of all pixels in the
image:
dip::Variance( img ).As< double >();
If mode
is "fast"
, a simplistic method to compute variance is used; this method can result in catastrophic
cancellation if the variance is very small with respect to the mean. If mode
is "stable"
, a stable algorithm
is used that avoids catastrophic cancellation, but is slower (see dip::VarianceAccumulator
and
dip::FastVarianceAccumulator
). For 8 and 16-bit integer images, the fast algorithm is always used.
If mode
is "directional"
, the data in in
are assumed to be angles, and directional statistics are used.
If in
contains orientations, multiply it by 2 before applying this function, and divide the result by 2.
For tensor images, the result is computed for each element independently. Input must be not complex, and for directional statistics it must be floating point.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::Variance( img.At( mask ), mode )
is the same as dip::Variance( img, mask, mode )
.
void
dip:: StandardDeviation(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::String mode = S::FAST,
dip::BooleanArray const& process = {})
Calculates the standard deviation of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the standard deviation of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the standard deviation projection along the processing dimensions. To get the standard deviation of all pixels in the
image:
dip::StandardDeviation( img ).As< double >();
If mode
is "fast"
, a simplistic method to compute standard deviation is used; this method can result in
catastrophic cancellation if the variance is very small with respect to the mean. If mode
is "stable"
, a
stable algorithm is used that avoids catastrophic cancellation, but is slower (see dip::VarianceAccumulator
and dip::FastVarianceAccumulator
). For 8 and 16-bit integer images, the fast algorithm is always used.
If mode
is "directional"
, the data in in
are assumed to be angles, and directional statistics are used.
If in
contains orientations, multiply it by 2 before applying this function, and divide the result by 2.
For tensor images, the result is computed for each element independently. Input must be not complex, and for directional statistics it must be floating point.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::StandardDeviation( img.At( mask ), mode )
is the same as dip::StandardDeviation( img, mask, mode )
.
void
dip:: Maximum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the maximum of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the maximum of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the maximum projection along the processing dimensions. To get the maximum of all pixels in the
image:
dip::Maximum( img ).As< double >();
For tensor images, the result is computed for each element independently. Input must be not complex.
If mask
is forged, only those pixels selected by the mask image are used.
To compute the sample-wise maximum over two or more images, use dip::Supremum
.
An alias is defined such that dip::Maximum( img.At( mask ))
is the same as dip::Maximum( img, mask )
.
void
dip:: Minimum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the minimum of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the minimum of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the minimum projection along the processing dimensions. To get the minimum of all pixels in the
image:
dip::Minimum( img ).As< double >();
For tensor images, the result is computed for each element independently. Input must be not complex.
If mask
is forged, only those pixels selected by the mask image are used.
To compute the sample-wise minimum over two or more images, use dip::Infimum
.
An alias is defined such that dip::Minimum( img.At( mask ))
is the same as dip::Minimum( img, mask )
.
void
dip:: MaximumAbs(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the maximum of the absolute pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the maximum of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the maximum projection along the processing dimensions. To get the maximum of all pixels in the
image:
dip::Maximum( img ).As< double >();
For tensor images, the result is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::MaximumAbs( img.At( mask ))
is the same as dip::MaximumAbs( img, mask )
.
void
dip:: MinimumAbs(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the minimum of the absolute pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the minimum of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the minimum projection along the processing dimensions. To get the minimum of all pixels in the
image:
dip::Minimum( img ).As< double >();
For tensor images, the result is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::MinimumAbs( img.At( mask ))
is the same as dip::MinimumAbs( img, mask )
.
void
dip:: Percentile(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::dfloat percentile = 50,
dip::BooleanArray const& process = {})
Calculates the percentile of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the percentile
percentile of the pixel values. Otherwise, the output has as many dimensions as elements
in process
that are false
, and equals the percentile projection along the processing dimensions.
To get the 30th percentile of all pixels in the image:
dip::Percentile( img, {}, 30.0 ).As< double >();
Note that the sample nearest the partition is picked, values are not interpolated if the partition falls in between samples.
For tensor images, the result is computed for each element independently. Input must be not complex.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::Percentile( img.At( mask ))
is the same as dip::Percentile( img, mask )
.
void
dip:: Median(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Calculates the median of the pixel values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the median (50th percentile) of the pixel values. Otherwise, the output has as many dimensions as elements
in process
that are false
, and equals the median projection along the processing dimensions. To get the
median of all pixels in the image:
dip::Median( img ).As< double >();
For tensor images, the result is computed for each element independently. Input must be not complex.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::Median( img.At( mask ))
is the same as dip::Median( img, mask )
.
void
dip:: MedianAbsoluteDeviation(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Computes the median absolute deviation (MAD)
The MAD is a measure of statistical dispersion. It can be used as a robust estimate of the standard deviation.
For normally distributed data, the standard deviation equals 1.4826 * MAD
. It is computed as if by
dip::Median( dip::Abs( img - dip::Median( img )));
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
the MAD of the pixel values. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the MAD projection along the processing dimensions. To get the MAD of all pixels in the
image:
dip::MedianAbsoluteDeviation( img ).As< double >();
For tensor images, the result is computed for each element independently. Input must be not complex.
If mask
is forged, only those pixels selected by the mask image are used.
An alias is defined such that dip::MedianAbsoluteDeviation( img.At( mask ))
is the same as dip::MedianAbsoluteDeviation( img, mask )
.
void
dip:: All(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Determines if all pixels have non-zero values over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
a boolean value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the “all” projection along the processing dimensions. To test if all the pixels in the image are
non-zero:
dip::All( img ).As< bool >();
For tensor images, the result is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
void
dip:: Any(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::BooleanArray const& process = {})
Determines if any pixel has a non-zero value over all those dimensions which are specified by process
.
If process
is an empty array, all dimensions are processed, and a 0D output image is generated containing
a boolean value. Otherwise, the output has as many dimensions as elements in process
that are false
,
and equals the “any” projection along the processing dimensions. To test if any pixel in the image is
non-zero:
dip::Any( img ).As< bool >();
For tensor images, the result is computed for each element independently.
If mask
is forged, only those pixels selected by the mask image are used.
void
dip:: PositionMaximum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::uint dim = 0,
dip::String const& mode = S::FIRST)
Calculates the position of the maximum of the pixel values in a single dimension specified by dim
.
The out
image has size 1 in the dim
dimension and is equally sized as in
in the other dimensions.
For each image line in the dim
dimension, the position of the maximum is computed
and its dim
-coordinate is stored in out
at the coordinates of that image line.
If mask
is forged, only those pixels selected by the mask image are used.
If mode
is "first"
, the first maximum is found, in linear index order.
If it is "last"
, the last one is found.
For tensor images, the result is computed for each element independently. Input must be not complex.
void
dip:: PositionMinimum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::uint dim = 0,
dip::String const& mode = S::FIRST)
Calculates the position of the minimum of the pixel values in a single dimension specified by dim
.
The out
image has size 1 in the dim
dimension and is equally sized as in
in the other dimensions.
For each image line in the dim
dimension, the position of the minimum is computed
and its dim
-coordinate is stored in out
at the coordinates of that image line.
If mask
is forged, only those pixels selected by the mask image are used.
If mode
is "first"
, the first minimum is found, in linear index order.
If it is "last"
, the last one is found.
For tensor images, the result is computed for each element independently. Input must be not complex.
void
dip:: PositionPercentile(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::dfloat percentile = 50,
dip::uint dim = 0,
dip::String const& mode = S::FIRST)
Calculates the position of the percentile of the pixel values in a single dimension specified by dim
.
The out
image has size 1 in the dim
dimension and is equally sized as in
in the other dimensions.
For each image line in the dim
dimension, the position of the percentile is computed
and its dim
-coordinate is stored in out
at the coordinates of that image line.
If mask
is forged, only those pixels selected by the mask image are used.
percentile
must be between 0.0 and 100.0.
If mode
is "first"
, the first percentile is found, in linear index order.
If it is "last"
, the last one is found.
For tensor images, the result is computed for each element independently. Input must be not complex.
A call to this function with percentile
set to 0.0 redirects to dip::PositionMinimum
and
a value of 100.0 redirects to dip::PositionMaximum
.
void
dip:: PositionMedian(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::uint dim = 0,
dip::String const& mode = S::FIRST)
Calculates the position of the median of the pixel values in a single dimension specified by dim
.
The out
image has size 1 in the dim
dimension and is equally sized as in
in the other dimensions.
For each image line in the dim
dimension, the position of the median is computed
and its dim
-coordinate is stored in out
at the coordinates of that image line.
If mask
is forged, only those pixels selected by the mask image are used.
If mode
is "first"
, the first percentile is found, in linear index order.
If it is "last"
, the last one is found.
For tensor images, the result is computed for each element independently. Input must be not complex.
This function redirects to dip::PositionPercentile
with percentile
set to 50.0.
void
dip:: RadialSum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::dfloat binSize = 1,
dip::String const& maxRadius = S::OUTERRADIUS,
dip::FloatArray const& center = {})
Computes the radial projection of the sum of the pixel values of in
.
If the radial distance of a pixel to center
is r
, then the sum of
the intensities of all pixels with n * binSize <= r < (n + 1) * binSize
is
stored at position n
in the radial dimension of out
.
binSize
sets the size of the bins (pixels) in the radial output dimension.
If maxRadius
is set to "inner radius"
, the maximum radius that is projected is equal to
the smallest distance from the given center
to any edge pixel of the image. Otherwise, when
maxRadius
is set to "outer radius"
, the maximum radius is set to largest distance from the given
center
to any edge pixel.
center
defaults to in.GetCenter()
.
The output data type is dip::DT_DFLOAT
for non-complex inputs and dip::DT_DCOMPLEX
for complex inputs.
void
dip:: RadialMean(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::dfloat binSize = 1,
dip::String const& maxRadius = S::OUTERRADIUS,
dip::FloatArray const& center = {})
Computes the radial projection of the mean of the pixel values of in
.
If the radial distance of a pixel to center
is r
, then the mean of
the intensities of all pixels with n * binSize <= r < (n + 1) * binSize
is
stored at position n
in the radial dimension of out
.
binSize
sets the size of the bins (pixels) in the radial output dimension.
If maxRadius
is set to "inner radius"
, the maximum radius that is projected is equal to
the smallest distance from the given center
to any edge pixel of the image. Otherwise, when
maxRadius
is set to "outer radius"
, the maximum radius is set to largest distance from the given
center
to any edge pixel.
center
defaults to in.GetCenter()
.
The output data type is dip::DT_DFLOAT
for non-complex inputs and dip::DT_DCOMPLEX
for complex inputs.
void
dip:: RadialMinimum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::dfloat binSize = 1,
dip::String const& maxRadius = S::OUTERRADIUS,
dip::FloatArray const& center = {})
Computes the radial projection of the minimum of the pixel values of in
.
If the radial distance of a pixel to center
is r
, then the minimum of
the intensities of all pixels with n * binSize <= r < (n + 1) * binSize
is
stored at position n
in the radial dimension of out
.
binSize
sets the size of the bins (pixels) in the radial output dimension.
If maxRadius
is set to "inner radius"
, the maximum radius that is projected is equal to
the smallest distance from the given center
to any edge pixel of the image. Otherwise, when
maxRadius
is set to "outer radius"
, the maximum radius is set to largest distance from the given
center
to any edge pixel.
center
defaults to in.GetCenter()
.
The output data type is equal to the input data type.
void
dip:: RadialMaximum(dip::Image const& in,
dip::Image const& mask,
dip::Image& out,
dip::dfloat binSize = 1,
dip::String const& maxRadius = S::OUTERRADIUS,
dip::FloatArray const& center = {})
Computes the radial projection of the maximum of the pixel values of in
.
If the radial distance of a pixel to center
is r
, then the maximum of
the intensities of all pixels with n * binSize <= r < (n + 1) * binSize
is
stored at position n
in the radial dimension of out
.
binSize
sets the size of the bins (pixels) in the radial output dimension.
If maxRadius
is set to "inner radius"
, the maximum radius that is projected is equal to
the smallest distance from the given center
to any edge pixel of the image. Otherwise, when
maxRadius
is set to "outer radius"
, the maximum radius is set to largest distance from the given
center
to any edge pixel.
center
defaults to in.GetCenter()
.
The output data type is equal to the input data type.