Functions and types of ThermovisorImages
ThermovisorImages.ThermovisorImages
— ModuleThermovisorData is a package designed to process thermal images stored as matrices.
Each element of thermal image represents a temperature value. The package enables users to load images from files, calculate temperature distributions, compute statistical analyses for temperatures along specified lines and recalculate the temperature by taking into account the emissivity and spectral range of thermal imaging camera.
Thermal image can be loaded using read_temperature_file
. User defined matrix can be wrapped in RescaledImage
which simple maps the values to [0,1] interval. After that, distinct areas (relative to their surroundings,like the most heated regions within the scene) can be fitted to regions of interest (ROIs) using fit_all_patterns
. Temperature distribution along the specified lined within each ROI can be obtained, and average radial and There is also possible to evaluate the averaged statistics over the fitted ROI's using CentredObjCollectionStat
and recalculate the temperature for a new emissivity value of the whole image or it's region using recalculate_with_new_emissivity!
Core.Type
— Method(::Type{T})() where T<:CentredObj
Empty object constructor
ThermovisorImages.CentredObj
— Type`CentredObj` is a sort of region of interest (ROI) marker object with coordinates and dimensions.
CentredObj
has centre coordinates and dimensions, object's center can be anywhere with respect to the image indices. ROI also has one or more size parameters (in pixels) coordinates of centre are equal to CartesianIndices, first element is the row index, the second element is the column index!! (y and x coordinate) This is opposite to the ImageDraw, where first Point coordinate corresponds to the column index and the second one to the row index. CentredObj
can also be used as for indexing image[c] - returns all elements of image within c, image[c]=x
sets all elements of image to the values of x, x firstindex should be 1. CentredObj
can also be used to set all image points within the ROI to a single value. e.g. image[c] = 30
By default we have three CentredObj
types in package:
CircleObj
, SquareObj
and RectangleObj
To impement CentredObj
abstraction one needs to implement:
side
- returns what is assumed to be the single number characterisation of the dimensions (by default the maximum value of dimensions(c) is taken)
area
- returns the CentredObj area, used in statistics calculation and in default two objs compirason see isless(c1::CentredObj,c2::CentredObj)
perimeter
- Object's perimeter
Base.size
- should always return the minimum size of frame which wrappes the centred obj, e.g. for circle it is [diameter,diemater]
is_within
- function to check if inds are within the CentredObj
line_within_mask
- function to check if all line points are within the CentredObj
fill_x0!
- function to fill the optimization starting vector during CentredObj
fitting the image
convert_to_drawable
fucntion to convert the CentredObj
to a drawable obj for ImageDraw
parnumber
function which returns the number of parameters
ThermovisorImages.CentredObjCollectionStat
— TypeCentredObjCollectionStat
Calculates the following statistics on CentredObj
s collection:
minimal, maximal, mean and standard deviation of sides,areas and perimeters of all CentredObj
objects in the collection v
.
ThermovisorImages.CircleObj
— TypeCircle object with defined diameter
ThermovisorImages.DistributionStatistics
— TypeDistributionStatistics
Basic descriptive statistics (mean and standard deviation) on temeprature distribution of value for a given coordinate
, D
is the matrix of distribution, where each column corresponds to sample, thus rows number of D
should be the same as the length of coordinate
, columns number is the number of samples.
ThermovisorImages.FilteredImage
— Type Type to store image with filtered temperature region
Fields:
full
- filtered rescaled image of the same size as the input with all pixels which are not the part of the pattern with label value
region_indices
- cartesian indices of the pattern in the input image
reduced
- image of reduced size where all not-inpatter pixels removed (the scaling of this image is the same as of the input imag.initial
see RescaledImage
type )
reduced_flag
- bitmatrix or Matrix{Bool} version of (reduced image)
ThermovisorImages.MarkeredImage
— TypeSpecial type to work with segmentated image, stored the matrix of indices and a vector of vectors of
patterns coordinates. Pattern can be sorted by area. And accessed by direct indexing into the MarkeredImage
object.
ThermovisorImages.RectangleObj
— TypeRectangular object with defined two sides
ThermovisorImages.RescaledImage
— Type RescaledImage - structure stores the image data mapped to region [0,1]
Fields:
initial - initial image before rescaling
sz - size of the image
min - minimum value
max - maximum value
im - image with all values from 0 to 1
ThermovisorImages.SquareObj
— TypeSquare with defined center and side
ThermovisorImages.StatDataPlot
— TypeType for custom plot recipe
Base.getindex
— MethodBase.getindex(img::AbstractMatrix,c::CentredObj)
CentredObj
can be used for matrix indexing, image[centred_object]
- returns the vector of temperatures of all points of image lying within the centred_object
of CentredObj
Base.isless
— MethodBase.isless(c1::CentredObj,c2::CentredObj)
By default centred objects are compared by their areas, thus the vector of centred objects can be sorted
Base.length
— MethodBase.length(m::MarkeredImage)
Return the total number of patterns in the markerred image MarkeredImage
Base.length
— MethodBase.length(c::CentredObj)
Total number of values needed to create CentredObj
of specified type
Base.setindex!
— MethodBase.setindex!(img::Matrix,x::Array,c::CentredObj)
img[c]=x assignes all x elements to the elements of img
with indices lying within the CentredObj
c
Base.setindex!
— MethodBase.setindex!(img::Matrix{T},x::Number,c::CentredObj) where T
Setting all elements within the CentredObj
to a single value
Base.size
— MethodBase.size(m::MarkeredImage)
The image size in pixels
ThermovisorImages._inbounds_flag
— Method_inbounds_flag(L,D,max_length,min_length)
Unsafe version of check! number of rows in D
should be the same as the number of elements in L
Returns Bool flag of all row not containing NaN's and lying within the minlength to maxlength range
ThermovisorImages._to_rgb
— Method_to_rgb(image::Matrix{Float64};
color_scheme::Symbol=:none)
Version without rescaling
ThermovisorImages.add_distrib_point!
— Methodadd_distrib_point!(points,distrib,point,value)
Internal fucntion to add the point to distribution
ThermovisorImages.along_line_distribution
— Methodalong_line_distribution(img::AbstractMatrix{T},x0,y0,x1,y1) where T
Function evaluates matrix values distribution along the line specified by two coordinates, img - input image returns the tuple of two vectors: coordinates and values see ImageDraw.bresenham
for details of finding the points of the line
returns points - vector of coordinates along the line distrib - distribution
ThermovisorImages.along_line_distribution_plot
— Methodalong_line_distribution_plot(along_line_length,along_line_distribution;
length_scaler::Float64=1.0,
is_centered::Bool=true,
kwargs...)
Plots temperature distribution along the line along_line_length
- coordinates, along_line_distribution
- values of temperature, length_scaler
- length scaler (can be used to convert pixels to the actual length) is_centered
- the line length is converted to the coordinates with zero value in the centre of the CentredObj
ThermovisorImages.along_line_distribution_xiaolin_wu
— Methodalong_line_distribution_xiaolin_wu(img::AbstractMatrix{T}, y0, x0, y1, x1) where T
Evaluates the value matrix content along the line with endpoint coordinates x0,y0,y1,x1, returns indices of all points. As far as Wu's algorithm returns two adjacent points the value is evaluated as an average of two point obtained with Wu's algorithm
see xiaolin_wu
function from ImageDraw
ThermovisorImages.along_mask_line_distribution
— Functionalong_mask_line_distribution(imag::AbstractMatrix,c::CentredObj,direction_angle=0.0,line_length=10.0;
length_per_pixel=1.0,
use_wu::Bool=false)
The same as within_mask_line_points_distribution
but returns the line length along the coordinates within the image.
line_length
- the length of line in the same units as length_per_pixel
. The calibration using mm_per_pixel
`, returns calibrated length along the line
ThermovisorImages.angular_distribution_statistics
— Methodangular_distribution_statistics(angles,along_length_coordinate,distrib;
max_length=-1.0,min_length=-1.0)
Fills DistributionStatistics
object for angular
Optional:
max_length
- maximal value of alonglengthcoordinate to be includet in to the statistics evaluation
min_length
- minimal value of alonglengthcoordinate to be includet in to the statistics evaluation
ThermovisorImages.angular_distribution_statistics_plot
— Methodangular_distribution_statistics_plot(ds::DistributionStatistics;
show_lower_bound::Bool=true,
show_upper_bound::Bool=true,
show_std::Bool=true,
probability::Float64=0.95,
is_use_student::Bool=true,
length_scaler::Float64=1.0,
label=nothing,
minorgrid=true,
gridlinewidth=2,
title="Average temperature angular distribution",framestyle = :box,
dpi=600,xlabel = "Angle , °", ylabel="Temperature, °C",
kwargs...)
The same as radial_distribution_statistics_plot
but plots averaged angular distribution
ThermovisorImages.area
— Methodarea(c::CentredObj)
Ealuates the surface area in pixels
ThermovisorImages.areas
— Methodareas(m::MarkeredImage)
Vector of patterns areas (number of pixels within the pattern)
ThermovisorImages.c_view
— Methodc_view(img::AbstractMatrix,c::CentredObj)
Returns the view of img
elements which are within the CentredObj
ThermovisorImages.cent_to_flag
— Methodcent_to_flag(c::CentredObj,sz::Tuple{Int,Int};external=false)
Converts CentredObj to bitmatrix of size sz
ThermovisorImages.cent_to_flag
— Methodcent_to_flag(::Type{T},c::CentredObj,sz::Tuple{Int,Int};external=false) where T<:FlagMatrix
Converts centred obj to BitMatrix or the Matrix of bool see `FlagMatrix
ThermovisorImages.center
— Methodcenter(c::CentredObj)
Vector of object center coordinates (MVector)
ThermovisorImages.change_default_colorscheme
— Methodchange_default_colorscheme(new_scheme::Symbol)
Change colorschee which is used by default to convert matrices to rgb images
ThermovisorImages.change_default_roi_color
— Methodchange_default_roi_color(color::RGB{Float64})
Changes default color to visualize the CentredObj
ThermovisorImages.convert_temperature_to_emissivity
— Methodconvert_temperature_to_emissivity(image::AbstractMatrix,
real_temperature::Float64,
image_emissivity::Float64;
λ_left::Union{Float64,Nothing}=14.0,
λ_right::Union{Float64,Nothing}=14.5,
is_in_kelvins::Bool=false)
Evaluates the emissivity, assuming the whole image
to be of the same temperature real_temperature
. If both λleft and λright are not equal and none of them is nothing
the emissivity ϵ
is calculated using:
ϵ = image_emissivity
⋅∫iᵦ(λ,Tᵣ)dλ/∫iᵦ(λ,Tₘ)dλ
Here iᵦ
is the Planck spectral radiance, Tₘ is the measured temperature (value of pixel temperature), Tᵣ is the real temperature of the surface. The integration is preformed over the [λ_left,λ_right]
spectral range with the help of band_power
function provied by the PlanckFunctions
package.
If one of (but not both) λleft and λright is set to nothing
or both of them are equal, than the following equaion is solved:
ϵ = image_emissivity
⋅iᵦ(λ,Tᵣ)/iᵦ(λ,Tₘ)
ThermovisorImages.convert_to_drawable
— Methodconvert_to_drawable(::CentredObj)
Converts CentredObj to a drawable structure appropriate to the ImageDraw
draw function, polygon,ellipse see [ImageDraw.draw
] function
ThermovisorImages.copyobj
— Methodcopyobj(c::T) where T<:CentreObj
Copies the CentredObj
creating new instance
ThermovisorImages.count_separate_patterns
— Methodcount_separate_patterns(markers::Matrix{Int})
This function takes matrix of markers see marker_image
and calculates the number of separate patterns
ThermovisorImages.diagonal_points
— Methoddiagonal_points(c::Union{SquareObj,CircleObj})
Returns diagonal points in row-column coordinates
ThermovisorImages.draw!
— Methoddraw!(rgbim::Matrix{RGB{Float64}},
c::CentredObj;
fill=false,
thickness::Int=55,
roi_color::RGB{Float64}=DefRoiColor[],
show_cross::Bool=true,
kwargs...)
Draws c
inside the image rgbim
(modified)
Keyword arguments:
fill - if true fills the object
thickness - thickness of c
boundary
(this two arguments are transfered to convert_to_drawable
)
roi_color - RGB color of ROI boundary and filling
show_cross - if true the cross is displayed
kwargs - Other keyword arguments are transfered directly to the ImageDraw.draw!
function
ThermovisorImages.draw
— Methoddraw(c::CentredObj;kwargs...)
Returns CentredObj
image of minimal possible size
ThermovisorImages.draw
— Methoddraw(image::Matrix{Float64},c::CentredObj;fill=false,thickness::Int=-1,
roi_color::RGB{Float64}=DefRoiColor[],
color_scheme::Symbol=:none,
show_cross=true,kwargs...)
Converts image
to RGB (if not all values of the image
are within [0,1] - rescales)
And draws CentreObj
c
inside the image
(does not affect the values of image
)
For kwargs
see draw!
ThermovisorImages.draw
— Methoddraw(image::RescaledImage,c::CentredObj)
Converts image to RGB and draws CentredObj
in it. For kwargs
see draw!
ThermovisorImages.draw
— Methoddraw(image::RescaledImage,c::Vector{T};kwargs...) where T<:CentredObj
Converts image
to RGB and draws multiple objects
ThermovisorImages.draw_line_within_mask
— Methoddraw_line_within_mask(image::Matrix{Float64},c::CentredObj,
ang,length;thickness::Int=55,
roi_color::RGB{Float64}=DefRoiColor[],
color_scheme::Symbol=:none,kwargs...)
Draws straight line on the image
converted to RGB, all points are located within the CentredObj
. This line goes through the centre of c
and oriented with the angle ang
in degrees with positive direction - counterclockwise.
ThermovisorImages.eval_bounds
— Methodeval_bounds(DS::DistributionStatistics;is_use_student::Bool=true)
Evaluates confidence bounds
ThermovisorImages.external_flag
— Methodexternal_flag(m::MarkeredImage,i::Int,::Type{T}=Matrix{Bool}) where T<:FlagMatrix
Inversed version of flag
(by default returns matrix of bool)
ThermovisorImages.fill_from_vect!
— Methodfill_from_vect!(c::CentredObj, v::AbstractVector)
Fills CentreObj parameters from the vector [centerindex1,centerindex2,dimension1,dimension2,...]
ThermovisorImages.fill_im!
— Methodfill_im!(img,c::CentredObj)
Fills bitmatrix or the matrix of Bool img
in a way that all pixels which are within the CentredObj
are true and false otherwise.
ThermovisorImages.fill_im_external!
— Methodfill_im_external!(img::FlagMatrix,c::CentredObj)
Fills image matrix img
in a way that all pixels which are not within the CentreObj set to true. See also is_within
ThermovisorImages.fill_vect!
— Methodfill_vect!(x::AbstractVector, c::CentredObj)
Converts CentredObj
to vector
ThermovisorImages.fill_x0!
— Methodfill_x0!(x0,im_bin::AbstractMatrix,c::CentredObj)
Fills the optimization starting vector by seraching the centre of the image im_bin
ThermovisorImages.fill_x0!
— Methodfill_x0!(x0,im_bin::FlagMatrix,::CircleObj)
Fills starting vector for the optimization of CentredObj
ThermovisorImages.filter_image!
— Methodfilter_image!(imag::AbstractMatrix,flag::BitMatrix)
Returns FilteredImage
taking all elements of imag which are not externalregionflag
ThermovisorImages.filter_image!
— Methodfilter_image!(imag::RescaledImage{Float64},external_region_flag::FlagMatrix)::FilteredImage
In-place filtering of RescaledImage
, filtered object is wrapped around the input RescaledImage
ThermovisorImages.filter_image
— Methodfilter_image(imag::AbstractMatrix,c::CentredObj;external=false)
Filters image according to centered object creating new image if external is true than as a filtering flag the inverse of centered object image is taken
ThermovisorImages.filter_image
— Methodfilter_image(imag::RescaledImage,c::CentredObj;external=false)
Filters image according
ThermovisorImages.filter_image
— Methodfilter_image(imag::RescaledImage,markers;label=0)
Funtion zeroes all pixels of the image, except those belonging to the specified pattern. image
- rescaled image (see RescaledImage
type) markers
- the matrix of the same size as the input image, each element of this matrix has unique value-label associated with some pattern. Function label_components
returns the markers matrix. (optional) - the value of the label to be selected as a pattern marker
Function returns FilteredImage
object
ThermovisorImages.find_temperature_files
— Function `find_temperature_files(folder::AbstractString)`
Searchs the folder for thermal images files using is_temperature_file
Returns dictionary Dict{String,Pair{Float64,String}}
with keys parts of files matched using is_temperature_file
, values - are temperature pairs of Float64
=> full-file-name
When file name contains "BB" it supposed to be the blackbody themperature distribution
ThermovisorImages.fit_all_patterns!
— Methodfit_all_patterns(img::RescaledImage,::Type{T}=CircleObj;
level_threshold::Float64=-1.0,
distance_threshold::Float64=-15.0,
max_centred_objs::Int=200,
sort_by_area::Bool = false,
is_descend::Bool = true,
optimizer::Optim.ZerothOrderOptimizer = NelderMead(),
options::Optim.Options=DEFAULT_FITTING_OPTIONS[]) where T<:CentredObj
Function fits all patterns of the image img
to the vector of CentredObj
ROI objects. The type of ROI should be provided as a second arguments (by default it is a CircleObj
)
img - input image of RescaledImage
type
For other input arguments see marker_image
and fit_centred_obj!
ThermovisorImages.fit_all_patterns!
— Methodfit_all_patterns!(c_vect::Vector{T},
markers::MarkeredImage,
optimizer::Optim.ZerothOrderOptimizer = NelderMead(),
options::Optim.Options=DEFAULT_FITTING_OPTIONS[]) where T<:CentredObj
Fits markered image pattern and fills precreated avector of Centerdobj
s See fit_all_patterns!
ThermovisorImages.fit_all_patterns
— Methodfit_all_patterns(img::RescaledImage,::Type{T}=CircleObj;
level_threshold::Float64=-1.0,
distance_threshold::Float64=-15.0,
max_centred_objs::Int=200,
sort_by_area::Bool = false,
is_descend::Bool = true,
optimizer::Optim.ZerothOrderOptimizer = NelderMead(),
options::Optim.Options=DEFAULT_FITTING_OPTIONS[]) where T<:CentredObj
Markers RescaledImage
and fits all patterns See fit_all_patterns!
ThermovisorImages.fit_centred_obj!
— Functionfit_centred_obj!(c::CentredObj,image::FilteredImage;kwargs...)
Fits CentredObj
(modified) to filtered image (not modified) fit_reduced
flag (default=true) indicates what version of the image should be fitted if true - reduced otherwise - full image. For other input arguments see fit_centred_obj!
ThermovisorImages.fit_centred_obj!
— Methodfit_centred_obj!(c::CentredObj,im_bin::FlagMatrix;
starting_point::Union{Nothing,Vector{Float64}}=nothing,
optimizer::Optim.ZerothOrderOptimizer = Optim.NelderMead(),
options::Optim.Options=DEFAULT_FITTING_OPTIONS[],refit::Bool = true)
Fits CentredObj
to binary image pattern (region of units) by adjusting centre coordinates and dimensions using zeroth-order optimizers from Optim.jl
package.
Input variables:
c - CentredObj
(modified)
im_bin - binarized image (BitMatrix or Matrix{Bool})
(optional)
startingpoint - staring vector (uses [`fillx0!`](@ref) function to fill starting point by default)
optimizer - zeroth-order optimizer from Optim.jl
package
options - optimization options from Optim.jl
package
refit - if true the starting point of the optimization recalculated otherwise it is taken from the current state vector of ROI
ThermovisorImages.flag!
— Methodflag!(fl::FlagMatrix,m::MarkeredImage,i::Int;negate::Bool=false)
Fills the fl matrix (Bitmatrix or Matrix{Bool}) with the same size as the entire image with all elements set to zero, except the pixels of i
'th pattern
ThermovisorImages.flag
— Methodflag(m::MarkeredImage,i::Int,::Type{T}=Matrix{Bool}) where T<:FlagMatrix
Creates matrix of Bool or Bitmatrix with the same size as the whole markers matrx with trues on the location of i'th pattern
ThermovisorImages.full_image_flag
— Methodfull_image_flag(filtered_im::FilteredImage)
Returns the BitMatrix flag of filtered pattern in the whole image.
Can be used as index matrix in the full image:
filtered_image.full.initial[full_image_flag(filtered_image)]
will return all elements which belong to the pattern
ThermovisorImages.generate_random_objs!
— Methodgenerate_random_objs!(im::Matrix{Float64},::Type{T},obj_number::Int,dimension_range::R) where {T<:CentredObj,R<:StepRange{Int,Int}}
Generates random objects on the image im
see generate_random_objs
ThermovisorImages.generate_random_objs
— Methodgenerate_random_objs(::Type{T},centers_range::NTuple{2,R},obj_number::Int,dimension_range::R) where {T<:CentredObj,R<:StepRange{Int,Int}}
Generates obj_number
of objects of CentredObj
subtype specified by the first argument, Second argument centers_range
is a tuple of StepRange
s for x and y coordinates of center to be taken from randomly, and dimension_range
is the range for dimensions
ThermovisorImages.image_discr
— Methodimage_discr(im1,im2)
Calculates the scalar distance between two matrices by checking the equality of their elements
ThermovisorImages.image_fill_discr
— Methodimage_fill_discr(image::FlagMatrix,c::CentredObj)
Function returns the function to evaluate the discrepancy between CentredObj
and the matrix, this function is used during the fitting procedure
ThermovisorImages.is_temperature_file
— Method is_temperature_file(file_name::AbstractString)
Checks if the file with file_name
has an appropriate name for thermovisor temperature distribution file
ThermovisorImages.is_within
— Methodis_within(c::CentredObj,_)
Function to check if indices are within CentredObj
ThermovisorImages.is_within
— Methodis_within(c::CentredObj,i::CartesianIndex)
CartesianIndex
support
ThermovisorImages.is_within_iterator
— Methodis_within_iterator(img::AbstractMatrix,c::CentredObj)
Iterator over all CartesianIndices
within the img
which are within the CentredObj c
ThermovisorImages.line_points_to_along_length
— Methodline_points_to_along_length(along_line_points::Vector{T},line_points) where T
Converts Cartesian indices of along_line_points
to the length along line
ThermovisorImages.line_within_mask
— Methodline_within_mask(c::CentredObj,ang::Float64,line_length::Int)
Function returns endpoint of the line lying fully within the mask - tuple of four point which can be directly splatted to the alonglinedistribution
ang - angle in degrees
line_length - the length of line
ThermovisorImages.line_within_mask
— Methodline_within_mask(c::CircleObj,ang,line_length)
Returns two endpoints of the line lying totally inside the CentredObj
ThermovisorImages.m_view
— Methodm_view(im::AbstractMatrix,m::MarkeredImage,i::Int)
Returns the view of matrix im
at points corresponding to the i
th pattern of MarkeredImage
image m
ThermovisorImages.m_view
— Methodm_view(m::MarkeredImage,i::Int)
Returns the view of specified marker in m.markers
ThermovisorImages.marker_image
— Methodmarker_image(rescaled::RescaledImage,level_threshold::Float64,distance_threshold::Float64=1e-3)
Markers image patterns, input umage is of RescaledImage
image type, levelthreshold - should be between 0.0 and 1.0 distancethreshold - criterium of image binarization after distance transform, should be less than 1
returns MarkeredImage
object
ThermovisorImages.mean_within_mask
— Methodmean_within_mask(img::AbstractMatrix,c::CentredObj)
Evaluates the average temperature of all points within the CentredObj
marker
ThermovisorImages.obj_from_vect
— Methodobj_from_vect(::Type{CentredObj},v::AbstractVector)
Creates object from parameters vector, first two arguments are center point other are dimensions [center[1],center[2],dimensions[1],...]
ThermovisorImages.parnumber
— Methodparnumber(::Type{T}) where T<:CentredObj
Returns total number of parameters needed to create new object
ThermovisorImages.perimeter
— Methodperimeter(c::CentredObj)
Returns perimeter of the object
ThermovisorImages.points_within_line!
— Methodpoints_within_line!(imag::AbstractMatrix,line_points::AbstractVector)
Forces all line points to lie within the possible region according to the image size
ThermovisorImages.radial_distribution
— Methodradial_distribution(imag::AbstractMatrix,c::CentredObj,angles_range::AbstractRange,line_length;mm_per_pixel=1.0)
Calls along_mask_line_distribution
on lines oriented with some angles range and puts the resulting distribution into one matrix j'th column of this matrix corresponds to the distribution along the line oriented with ang[j] angle. The length of line is line_length
, if it is less than 0.0 or greater than the smallest dimension of c
, than is is set to the smallest dimentsion of c
ThermovisorImages.radial_distribution_statistics
— Methodradial_distribution_statistics(along_length_coordinate::AbstractVector,
distrib::AbstractVecOrMat;
max_length=-1.0,min_length=-1.0)
Fills DistributionStatistics
object for radial distribution
Optional:
max_length
- maximal value of alonglengthcoordinate to be includet in to the statistics evaluation
min_length
- minimal value of alonglengthcoordinate to be includet in to the statistics evaluation
ThermovisorImages.radial_distribution_statistics_plot
— Methodradial_distribution_statistics_plot(ds::DistributionStatistics;
show_lower_bound::Bool=false,
show_upper_bound::Bool=false,
show_std::Bool=true,
is_use_student::Bool=true,
probability::Float64=0.95,
length_scaler::Float64=1.0,
is_centered::Bool=true,
bound_color::Symbol=:red,
# plot kwargs
label::Union{AbstractString,Nothing} =nothing,
minorgrid::Bool=true,
gridlinewidth=2,
title::Union{AbstractString,Nothing} ="Average temperature radial distribution",
framestyle::Symbol = :box,
dpi::Int=600,
xlabel::Union{AbstractString,Nothing} = "Distance across the sample ,mm",
ylabel::Union{AbstractString,Nothing} = "Temperature °C",
kwargs...)
Recipe to plot the distribution statistics DistributionStatistics
object. Returns structure of StatDataPlot
type, which has attached recipe to plot radial ditribution averaged value, lower and upper confidence bounds as <d> ± std (if show_std
is true) and confidence bounds multiplied by the Student's coefficient (if showlowerbound and showupperbound are true) calculated for the pobability
value, if length_scaler
is provied all coordinates are multiplied by this value (can be used to convert pixels to actual units) If is_centered
is true coordinate goes from -L/2 to +L/2 where L is the maximum of coordinates. Other key-word arguments are the same as for the plot functions, additional keyword arguments are transfered directly to the plot function Usage example:
using Plots
plot(radial_distribution_statistics_plot(distribution_statistics))
This code will plot the averaged distribution statistics together with confidence bounds
ThermovisorImages.read_temperature_file
— Methodread_temperature_file(f_name::AbstractString;inverse_intensity::Bool=false)
Reads temeprature file f_name
is a full file name, inverse_intensity
is true if the intensities in the loaded file should be inverted
ThermovisorImages.rearranged_diagonal
— Methodrearranged_diagonal(c::Union{SquareObj,CircleObj})
Returns diagonal points in x-y coordinates
ThermovisorImages.recalculate_with_new_emissivity!
— Methodrecalculate_with_new_emissivity!(image::AbstractArray,c::CentredObj,new_emissivity::Float64;
kwargs...)
Recalculates temperarture of each pixel within the CentredObj
with new_emissivity
see recalculate_with_new_emissivity!
for keyword arguments
ThermovisorImages.recalculate_with_new_emissivity!
— Methodrecalculate_with_new_emissivity!(image::AbstractArray,new_emissivity::Float64,
image_emissivity::Float64;
λ_left::Union{Float64,Nothing}=14.0,
λ_right::Union{Float64,Nothing}=14.5,
is_in_kelvins::Bool=false,
rel_tol::Float64=1e-3)
Function recalculates all temperatures in image
assuming the temperatures were measured for the surface with emissivity image_emissivity
; new_emissivity
is a new value of emissivity, λ_left
and λ_right
are left and right wavelength boundaries of infrared camera spectral range in μm. If is_in_kelvins
is false (default) all temperatures supposed to be in Celsius and Kelvins otherwise.
To find a new temperature value T
the function solves a non-linear equation. If both λleft and λright are not equal and none of them is nothing
, the following equation is solved:
new_emissivity
⋅∫iᵦ(λ,T)dλ = image_emissivity
⋅∫iᵦ(λ,Tₘ)dλ
Here iᵦ
is the Planck spectral radiance, Tₘ is the measured temperature (value of pixel intensity). The integration is preformed over the [λ_left,λ_right]
spectral range with the help of band_power
function provied by the PlanckFunctions
package.
If one of (but not both) λleft and λright is set to nothing
or both of them are equal, than the following equaion is solved:
new_emissivity
⋅iᵦ(λ,T) = image_emissivity
⋅iᵦ(λ,Tₘ)
In both cases the equation is solved by minimizing the univariate function of least-square difference of left and right parts of the equation. Second version of equation is much faster, but less precise.
ThermovisorImages.recalculate_with_new_emissivity!
— Methodrecalculate_with_new_emissivity!(image::AbstractArray,marker::MarkeredImage,
label::Int,new_emissivity::Float64,image_emissivity::Float64;
kwargs...)
Recalculates the temperarture of each pixel within image pattern label
of MarkeredImage
image marker
with new_emissivity
assuming image_emissivity
be the emissvity settled during measurements. Both marker
and image
should be of the same size. See recalculate_with_new_emissivity!
ThermovisorImages.revcentre
— Methodrevcentre(c::CentredObj)
Object coordinates in reverse order
ThermovisorImages.shift!
— Methodshift!(c::CentredObj,ind::CartesianIndex{2})
Relative shift of centred object center
ThermovisorImages.shrinked_flag
— Methodreduced_flag(m::MarkeredImage,i::Int)
returns flag matrix shrinked to the size of the pattern
ThermovisorImages.side
— Methodside(c::CentredObj)
Returns the side of CentredObj
ThermovisorImages.sort_markers!
— Functionsort_markers!(m::MarkeredImage,descending::Bool=true)
Sorts markers, see sort_reduce!
ThermovisorImages.sort_reduce!
— Methodsort_reduce!(m::MarkeredImage;total_number::Int = -1,descending::Bool=true)
Sorts markers by total area and reduces the total number of patterns if the maximum label is less then total_number
value
Input arguments: m - MarkeredImage
total_number - number of patterns retained in the image
ThermovisorImages.std_within_mask
— Methodstd_within_mask(img::AbstractMatrix, c::CentredObj)
Evaluates standard deviation of temperature for all points within the CentredObj
marker
ThermovisorImages.student_coefficient
— Methodstudent_coefficient(degrees_of_freedom::Int, probability; digits::Int = 3, side::Int = 2)
Evaluates Student's distribution coefficient
ThermovisorImages.to_rgb
— Methodto_rgb(image::Matrix{Float64};color_scheme::String="")
Converts matrix to rgb - martix by applyting the color scheme from ColorSchemes
package.
ThermovisorImages.within_mask_line_points_distribution
— Functionwithin_mask_line_points_distribution(imag::AbstractMatrix,c::CentredObj,direction_angle=0.0,line_length=10.0;use_wu::Bool=false)
Evaluates the distribution of values in imag
matrix along the line with length line_length
in pixels oriented with the angle direction_angle
in degrees with respect to the posistive direction of oX (column index increase), this line lies within the roi (CentreObj
) and goes through its center.
Function returns: a Tuple (points,distrib,linepoints)
points - vector of CartesianIndex
of image's points lying on the line
distrib - distribution of values
linepoints - endpoints of line the Tupple of (leftx,leftY,rightx,righty)