Title: | Tools to analyze in vivo microscopy imaging data focused on tracking flowing blood cells |
---|---|
Description: | flowcatchR is a set of tools to analyze in vivo microscopy imaging data, focused on tracking flowing blood cells. It guides the steps from segmentation to calculation of features, filtering out particles not of interest, providing also a set of utilities to help checking the quality of the performed operations (e.g. how good the segmentation was). It allows investigating the issue of tracking flowing cells such as in blood vessels, to categorize the particles in flowing, rolling and adherent. This classification is applied in the study of phenomena such as hemostasis and study of thrombosis development. Moreover, flowcatchR presents an integrated workflow solution, based on the integration with a Shiny App and Jupyter notebooks, which is delivered alongside the package, and can enable fully reproducible bioimage analysis in the R environment. |
Authors: | Federico Marini [aut, cre] |
Maintainer: | Federico Marini <[email protected]> |
License: | BSD_3_clause + file LICENSE |
Version: | 1.41.0 |
Built: | 2024-11-18 03:24:22 UTC |
Source: | https://github.com/bioc/flowcatchR |
Frames
object
Creates a Frames
object containing raw information, combined with the segmented
images and the relative trajectory under analysisIf a TrajectorySet
is provided and mode is set to trajectories
, returns
a Frames
with all trajectories included in the IDs
vector painted accordingly.
If the mode is set to particles
, it will just plot the particles (all) on all frames.
If no TrajectorySet
is provided, it will be computed with default parameters.
If no binary.frames
is provided, it will be computed also with default parameters
add.contours( raw.frames, binary.frames = NULL, trajectoryset = NULL, trajIDs = NULL, mode = "particles", col = NULL, channel = NULL )
add.contours( raw.frames, binary.frames = NULL, trajectoryset = NULL, trajIDs = NULL, mode = "particles", col = NULL, channel = NULL )
raw.frames |
A |
binary.frames |
A |
trajectoryset |
A |
trajIDs |
Numeric vector, the ID(s) of the trajectory. |
mode |
A character string, can assume the values |
col |
A vector of color strings |
channel |
A character string, to select which channel to process |
A new Frames
object with contours of the objects added
Federico Marini, [email protected], 2014
data("MesenteriumSubset") ## Not run: paintedTrajectories <- add.contours(raw.frames = MesenteriumSubset, mode = "trajectories",channel="red") paintedParticles <- add.contours(raw.frames = MesenteriumSubset, mode = "particles",channel="red") inspect.Frames(paintedTrajectories) inspect.Frames(paintedParticles) ## End(Not run)
data("MesenteriumSubset") ## Not run: paintedTrajectories <- add.contours(raw.frames = MesenteriumSubset, mode = "trajectories",channel="red") paintedParticles <- add.contours(raw.frames = MesenteriumSubset, mode = "particles",channel="red") inspect.Frames(paintedTrajectories) inspect.Frames(paintedParticles) ## End(Not run)
Frames
object and the corresponding preprocessed oneAll objects are painted with a unique colour - for sake of speed
addParticles(raw.frames, binary.frames, col = NULL)
addParticles(raw.frames, binary.frames, col = NULL)
raw.frames |
A |
binary.frames |
A |
col |
A color character string, to select which color will be used for drawing the contours of the particles. If not specified, it will default according to the objects provided |
A Frames
object, whose images are the combination of the raw images with the segmented objects drawn on them
Federico Marini, [email protected], 2014
Auxiliary function to return the dimensions of the field of interest
axesInfo(frames)
axesInfo(frames)
frames |
A |
A list object, containing the extremes of the field of interest (x-y-z, where z is time)
Federico Marini, [email protected], 2014
ParticleSet
objectThe sample ParticleSet
object is constituted by the platelets identified from the MesenteriumSubset
data
Federico Marini, [email protected], 2014
objects
channel
channel.Frames(frames, mode)
channel.Frames(frames, mode)
frames |
A |
mode |
A character value specifying the target mode for conversion. |
A Frames
object with just the infotmation on the selected channel
data("MesenteriumSubset") channel.Frames(MesenteriumSubset,"red")
data("MesenteriumSubset") channel.Frames(MesenteriumSubset,"red")
Calculates the Mean Squared Displacement for a trajectory
computeMSD(sx, sy, until = 4)
computeMSD(sx, sy, until = 4)
sx |
x axis positions along the trajectory |
sy |
y axis positions along the trajectory |
until |
how many points should be included in the Mean Squared Displacement curve |
A numeric vector containing the values of the MSD
Federico Marini, [email protected], 2014
Frames
objectPerforms cropping on the Frames
object, selecting how many pixels should be cut on each side
crop.Frames( frames, cutLeft = 5, cutRight = 5, cutUp = 5, cutDown = 5, cutAll = 0, testing = FALSE, ... )
crop.Frames( frames, cutLeft = 5, cutRight = 5, cutUp = 5, cutDown = 5, cutAll = 0, testing = FALSE, ... )
frames |
An input |
cutLeft |
Amount of pixels to be cut at the side |
cutRight |
Amount of pixels to be cut at the side |
cutUp |
Amount of pixels to be cut at the side |
cutDown |
Amount of pixels to be cut at the side |
cutAll |
Amount of pixels to be cut at all sides. Overrides the single side values |
testing |
Logical, whether to just test the cropping or to actually perform it. Default set to |
... |
Arguments to be passed to |
Cropping can be performed with careful choice of all cutting sides, or cropping a single value from all sides
A Frames
object, with cropped frames in the image
slot
Federico Marini, [email protected], 2014
data("MesenteriumSubset") crop.Frames(MesenteriumSubset)
data("MesenteriumSubset") crop.Frames(MesenteriumSubset)
Frames
objectWrites the images contained in the image
slot of the Frames
object elements.
The images can be exported as single frames, or as a .gif image that is composed
by the single frames.
export.Frames( frames, dir = tempdir(), nameStub = "testExport", createGif = FALSE, removeAfterCreatingGif = TRUE )
export.Frames( frames, dir = tempdir(), nameStub = "testExport", createGif = FALSE, removeAfterCreatingGif = TRUE )
frames |
A |
dir |
The path of the folder where the image should be written |
nameStub |
The stub for the file name, that will be used as a prefix for the exported images |
createGif |
Logical, whether to create or not an animated .gif file |
removeAfterCreatingGif |
Logical, whether to remove the single exported .png images after creating the single .gif |
Image files are written in the desired location
Federico Marini, [email protected], 2014
data("MesenteriumSubset") ## Not run: export.Frames(MesenteriumSubset,nameStub="subset_export_", createGif=TRUE,removeAfterCreatingGif=FALSE) ## End(Not run)
data("MesenteriumSubset") ## Not run: export.Frames(MesenteriumSubset,nameStub="subset_export_", createGif=TRUE,removeAfterCreatingGif=FALSE) ## End(Not run)
ParticleSet
objectWrites the particles contained in the particles
data frame slot of the ParticleSet
object elements.
A track of the provenience of the particles is stored as a comment line above the header
export.particles( particleset, dir = tempdir(), nameStub = "testExport_particles" )
export.particles( particleset, dir = tempdir(), nameStub = "testExport_particles" )
particleset |
A |
dir |
The path of the folder where the particle sets should be written |
nameStub |
The stub for the file name, that will be used as a prefix for the exported particle sets |
Particle sets files are written in the desired location
Federico Marini, [email protected], 2014
data("candidate.platelets") ## Not run: export.particles(candidate.platelets)
data("candidate.platelets") ## Not run: export.particles(candidate.platelets)
The computed set of parameters include delta.x
, delta.t
and delta.v
(displacements and instantaneous velocity), totalTime
, totalDistance
,
distStartToEnd
, curvilinearVelocity
, straightLineVelocity
and
linearityForwardProgression
, Mean Squared Displacement, velocity
autocorrelation, and more
extractKinematics.traj( trajectoryset, trajectoryID, acquisitionFrequency = 30, scala = 50 )
extractKinematics.traj( trajectoryset, trajectoryID, acquisitionFrequency = 30, scala = 50 )
trajectoryset |
A |
trajectoryID |
The ID of a single trajectory |
acquisitionFrequency |
The frame rate of acquisition for the images, in milliseconds |
scala |
The value of micro(?)meters to which each single pixel corresponds |
A KinematicsFeatures
object
Federico Marini, [email protected], 2014
A toolset to analyze in vivo microscopy imaging data focused on tracking flowing blood cells.
flowcatchR is a set of tools to analyze in vivo microscopy imaging data, focused on tracking flowing blood cells. It guides the steps from segmentation to calculation of features, filtering out particles not of interest, providing also a set of utilities to help checking the quality of the performed operations (e.g. how good the segmentation was). The main novel contribution investigates the issue of tracking flowing cells such as in blood vessels, to categorize the particles in flowing, rolling and adherent. This classification is applied in the study of phenomena such as hemostasis and study of thrombosis development.
Federico Marini [email protected], Johanna Mazur [email protected], Harald Binder [email protected], 2015
Maintainer: Federico Marini [email protected]
Frames
objectConstructor for a Frames
object
Frames(x, channel)
Frames(x, channel)
x |
A multi-dimensional |
channel |
A character vector, can be 'red','green','blue' or 'all' (if in color mode) |
The created Frames
object.
Federico Marini, [email protected], 2014
data("MesenteriumSubset") inputImg <- Image(MesenteriumSubset) Frames(inputImg,"red")
data("MesenteriumSubset") inputImg <- Image(MesenteriumSubset) Frames(inputImg,"red")
S4 class for storing information on multiple images belonging to the same time-lapse experiment. It is designed as a subclass of the existing Image class from the EBImage package
channel
A character vector, can be 'red','green','blue' or 'all' (if in color mode)
ParticleSet
object for subsequent linking/trackingInitialize a ParticleSet
object for subsequent linking/tracking
initialize.LinkedParticleSet(particleset, linkrange = 1)
initialize.LinkedParticleSet(particleset, linkrange = 1)
particleset |
A |
linkrange |
The number of frames to look for candidate particles potentially belonging to the same track |
A ParticleSet
object with slots dedicated for the tracking pre-filled
Federico Marini, [email protected], 2014
Frames
The first frames of a Frames
are displayed in the browser, and are interactively navigable.
inspect.Frames( frames, nframes = NULL, display.method = "browser", verbose = FALSE )
inspect.Frames( frames, nframes = NULL, display.method = "browser", verbose = FALSE )
frames |
A |
nframes |
The number of frames to display (default value: |
display.method |
Method for displaying, can be either |
verbose |
Logical, whether to provide additional output on the command line alongside with the images themselves |
inspect.Frames
returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("MesenteriumSubset") ## Not run: inspect.Frames(MesenteriumSubset)
data("MesenteriumSubset") ## Not run: inspect.Frames(MesenteriumSubset)
TrajectorySet
object, or a single parameter, or from a single trajectory (all possible combinations)The computed set of parameters include delta.x
, delta.t
and delta.v
(displacements and instantaneous velocity), totalTime
, totalDistance
,
distStartToEnd
, curvilinearVelocity
, straightLineVelocity
and
linearityForwardProgression
, Mean Squared Displacement, velocity
autocorrelation, and more. If a single trajectory is specified, the computation is
performed for that trajectory alone. If a parameter is specified, only that
parameter is reported, either for one or all trajectories
kinematics( trajectoryset, trajectoryIDs = NULL, acquisitionFrequency = 30, scala = 50, feature = NULL )
kinematics( trajectoryset, trajectoryIDs = NULL, acquisitionFrequency = 30, scala = 50, feature = NULL )
trajectoryset |
A |
trajectoryIDs |
The ID of a single trajectory |
acquisitionFrequency |
The frame rate of acquisition for the images, in milliseconds |
scala |
The value of micro(?)meters to which each single pixel corresponds |
feature |
Character string, the name of the feature to be computed |
A KinematicsFeaturesSet
object, or a KinematicsFeatures
object,
or an atomic value, or a list(eventually coerced to a vector)
Federico Marini, [email protected], 2014
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) # for all trajectories, all features alltrajs.features <- kinematics(platelets.trajectories) # for one trajectory, all features traj11features <- kinematics(platelets.trajectories,trajectoryIDs = 11) # for all trajectories, one feature alltrajs.curvVel <- kinematics(platelets.trajectories,feature = "curvilinearVelocity")
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) # for all trajectories, all features alltrajs.features <- kinematics(platelets.trajectories) # for one trajectory, all features traj11features <- kinematics(platelets.trajectories,trajectoryIDs = 11) # for all trajectories, one feature alltrajs.curvVel <- kinematics(platelets.trajectories,feature = "curvilinearVelocity")
S4 class for storing information on all kinematics features identified for a single trajectory
.Data
A list storing the information for the kinematics features
S4 class for storing information on all kinematics features identified for all trajectories. Single
KinematicsFeatures
objects are the element of the main list
.Data
A list storing the information for the sets of kinematics features
Frames
objectCompute the length of render frames in a Frames
object
## S3 method for class 'Frames' length(x)
## S3 method for class 'Frames' length(x)
x |
A |
An integer number
Federico Marini, [email protected], 2014
data("MesenteriumSubset") length(MesenteriumSubset)
data("MesenteriumSubset") length(MesenteriumSubset)
ParticleSet
objectPerforms linking of the particles by tracking them through the frames
link.particles( particleset, L, R = 2, epsilon1 = 0.1, epsilon2 = 2, lambda1 = 1, lambda2 = 1, penaltyFunction = penaltyFunctionGenerator(), verboseOutput = FALSE, prog = FALSE, include.intensity = TRUE, include.area = FALSE )
link.particles( particleset, L, R = 2, epsilon1 = 0.1, epsilon2 = 2, lambda1 = 1, lambda2 = 1, penaltyFunction = penaltyFunctionGenerator(), verboseOutput = FALSE, prog = FALSE, include.intensity = TRUE, include.area = FALSE )
particleset |
A |
L |
Maximum number of pixels an object can move in two consecutive frames |
R |
Linkrange, i.e. the number of consecutive frames to search for potential candidate links |
epsilon1 |
A numeric value, to be used in the formula. Jitter for allowing angular displacements |
epsilon2 |
A numeric value, to be used in the formula. Jitter for allowing spatial displacements |
lambda1 |
A numeric value. Multiplicative factor for the penalty function |
lambda2 |
A numeric value. Multiplicative factor applied to the angular displacement |
penaltyFunction |
A function structured in such a way to be applied as penalty function in the linking |
verboseOutput |
Logical, whether the output should report additional intermediate steps. For debugging use mainly |
prog |
Logical, whether the a progress bar should be shown during the tracking phase |
include.intensity |
Logical, whether to include also intensity change of the particles in the cost function calculation |
include.area |
Logical, whether to include also area change of the particles in the cost function calculation |
A LinkedParticleSet
object
Federico Marini, [email protected], 2014
I F Sbalzarini and P Koumoutsakos."Feature point tracking and trajectory analysis for video imaging in cell biology." In: Journal of structural biology 151.2 (Aug. 2005), pp. 182-95. ISSN: 1047-8477. DOI: 10.1016/j.jsb.2005.06.002. URL: http://www.ncbi.nlm.nih.gov/pubmed/16043363
data("candidate.platelets") tracked.platelets <- link.particles(candidate.platelets, L= 40)
data("candidate.platelets") tracked.platelets <- link.particles(candidate.platelets, L= 40)
S4 class for storing information of particles after they have been tracked. It inherits the slots
from the ParticleSet
class.
tracking
A list storing all necessary information for the tracking algorithm to work, and for providing the information to the function to determine the trajectories
Match trajectories to the related particles in the TrajectorySet
and
ParticleSet
objects. This function returns a new ParticleSet
object that contains as additional column the trajectory ID that the particular
particle was assigned to. Used also by other routines, such as snap
matchTrajToParticles(particleset, trajectoryset)
matchTrajToParticles(particleset, trajectoryset)
particleset |
A |
trajectoryset |
A |
A ParticleSet
object with an additional column with the trajectory
IDs
Federico Marini, [email protected], 2015
data(candidate.platelets) trajs <- trajectories(candidate.platelets) matchTrajToParticles(candidate.platelets, trajs)
data(candidate.platelets) trajs <- trajectories(candidate.platelets) matchTrajToParticles(candidate.platelets, trajs)
Frames
objectThe sample Frames
object is constituted by a subset of a time-lapse intravital microscopy imaging dataset.
Green channel marks leukocytes, red channel focuses on blood platelets. 20 frames are provided in this subset.
Images are kindly provided by Sven Jaeckel ([email protected]).
Federico Marini, [email protected], 2014
Frames
objectApplies a transformation to the Frames
object in a way that the intensities
throughout the acquisition are normalized overall in term of pixel values sums.
It can be used to compensate for example a global change in the illumination values,
e.g. due to changed acquisition conditions in experiments that span long timescales.
normalizeFrames(frames, normFun = "median")
normalizeFrames(frames, normFun = "median")
frames |
A |
normFun |
The normalization function chosen. Can be one of |
A Frames
object with normalized pixel values.
Federico Marini, [email protected], 2014
data(MesenteriumSubset) normalizeFrames(MesenteriumSubset,normFun="median")
data(MesenteriumSubset) normalizeFrames(MesenteriumSubset,normFun="median")
Frames
object.Extracts particles from the images of a Frames
object.
particles( raw.frames, binary.frames = NULL, channel = NULL, BPPARAM = bpparam() )
particles( raw.frames, binary.frames = NULL, channel = NULL, BPPARAM = bpparam() )
raw.frames |
A |
binary.frames |
A |
channel |
Character string. The channel to perform the operations on. Can be |
BPPARAM |
a |
A ParticleSet
object, containing all detected particles for each frame
Federico Marini, [email protected], 2015
data("MesenteriumSubset")
data("MesenteriumSubset")
S4 class for storing information on particles detected in distinct frames.
.Data
A list storing the information for the particles
channel
A character vector, can be 'red','green', or 'blue'. It refers to which channel the particles were detected
A function to generate penalty functions to use while linking particles
penaltyFunctionGenerator( epsilon1 = 0.1, epsilon2 = 2, lambda1 = 1, lambda2 = 1 )
penaltyFunctionGenerator( epsilon1 = 0.1, epsilon2 = 2, lambda1 = 1, lambda2 = 1 )
epsilon1 |
A numeric value, to be used in the formula. Jitter for allowing angular displacements |
epsilon2 |
A numeric value, to be used in the formula. Jitter for allowing spatial displacements |
lambda1 |
A numeric value. Multiplicative factor for the penalty function |
lambda2 |
A numeric value. Multiplicative factor applied to the angular displacement |
A function object, to be used as penalty function
Federico Marini, [email protected], 2014
custom.function <- penaltyFunctionGenerator(epsilon1=0.1,epsilon2=6,lambda1=1.5,lambda2=0)
custom.function <- penaltyFunctionGenerator(epsilon1=0.1,epsilon2=6,lambda1=1.5,lambda2=0)
TrajectorySet
objectProvides a visual representation of a TrajectorySet
object
## S3 method for class 'TrajectorySet' plot(x, frames, verbose = FALSE, ...)
## S3 method for class 'TrajectorySet' plot(x, frames, verbose = FALSE, ...)
x |
A |
frames |
A |
verbose |
Logical, whether to provide additional output on the command line |
... |
Arguments to be passed to methods |
Based on the plotly
library, the function extracts the region of interests
from the dimensions of an image of the Frames
object,
and afterwards plots the x-y-time representation of the identified trajectories
plot.TrajectorySet
returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("MesenteriumSubset") data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) ## Not run: plot(platelets.trajectories,MesenteriumSubset) ## End(Not run)
data("MesenteriumSubset") data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) ## Not run: plot(platelets.trajectories,MesenteriumSubset) ## End(Not run)
TrajectorySet
objectProvides a bird's eye view of a TrajectorySet
object on a bidimensional space
plot2D.TrajectorySet( trajectoryset, frames, trajIDs = NULL, addGrid = FALSE, verbose = FALSE, ... )
plot2D.TrajectorySet( trajectoryset, frames, trajIDs = NULL, addGrid = FALSE, verbose = FALSE, ... )
trajectoryset |
A |
frames |
A |
trajIDs |
A vector containing the ids of the desired trajectories |
addGrid |
Logical, add an additional grid to the 2-dimensional plot (visual aid for backtracking trajectory point locations) |
verbose |
Logical, whether to provide additional output on the command line |
... |
Arguments to be passed to methods |
This function extracts the region of interests from the dimensions
of an image of the Frames
object,
and afterwards plots the x-y-time representation of the identified
trajectories on a 2d plane. It is possible to subset the TrajectorySet
object with the IDs of the desired trajectories
plot2D.TrajectorySet
returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("MesenteriumSubset") data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) plot2D.TrajectorySet(platelets.trajectories,MesenteriumSubset)
data("MesenteriumSubset") data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) plot2D.TrajectorySet(platelets.trajectories,MesenteriumSubset)
Frames
objectsFrames
objects are processed according to the chosen set of parameters. Many of them refer directly to
existing EBImage
functions, please see the corresponding help for additional information
preprocess.Frames( frames, brush.size = 3, brush.shape = "disc", at.offset = 0.15, at.wwidth = 10, at.wheight = 10, kern.size = 3, kern.shape = "disc", ws.tolerance = 1, ws.radius = 1, displayprocessing = FALSE, ... )
preprocess.Frames( frames, brush.size = 3, brush.shape = "disc", at.offset = 0.15, at.wwidth = 10, at.wheight = 10, kern.size = 3, kern.shape = "disc", ws.tolerance = 1, ws.radius = 1, displayprocessing = FALSE, ... )
frames |
A |
brush.size |
Size in pixels of the brush to be used for initial smoothing (low-pass filtering) |
brush.shape |
Shape of the brush to be used for initial smoothing (low-pass filtering) |
at.offset |
Offset to be used in the adaptive thresholding step - see also |
at.wwidth |
Width of the window for the adaptive thresholding step - see also |
at.wheight |
Height of the window for the adaptive thresholding step - see also |
kern.size |
Size in pixels of the kernel used for morphological operations - e.g., opening, which is an erosion followed by a dilation, and closing which is a dilation followed by an erosion - see also |
kern.shape |
Shape of the kernel used for morphological operations |
ws.tolerance |
Tolerance allowed in performing the watershed-based segmentation (see also |
ws.radius |
Radius for the watershed-based segmentation (see also |
displayprocessing |
Logical, whether to display intermediate steps while performing preprocessing. Dismissed currently, it could increase runtime a lot |
... |
Arguments to be passed to methods |
A Frames
object, whose frame images are the preprocessed versions of the input images
Federico Marini, [email protected], 2014
data("MesenteriumSubset") preprocess.Frames(channel.Frames(MesenteriumSubset,"red"))
data("MesenteriumSubset") preprocess.Frames(channel.Frames(MesenteriumSubset,"red"))
Frames
objectThis function is used to create a Frames
object from a vector of image files (or a folder specifying the directory
containing them).
The number of frames is also specified, as just a subset of the images can be used for this
read.Frames(image.files, nframes = NULL)
read.Frames(image.files, nframes = NULL)
image.files |
Vector of strings containing the locations where the (raw) images are to be found, or alternatively, the path to the folder |
nframes |
Number of frames that will constitute the |
An object of the Frames
class, which holds the info on a list of frames, specifying for each the following elements:
image |
The |
location |
The complete path to the location of the original image |
Federico Marini, [email protected], 2014
## see vignette ## Not run: fullData <- read.Frames(image.files = "/path/to/the/directory", nframes = 100)
## see vignette ## Not run: fullData <- read.Frames(image.files = "/path/to/the/directory", nframes = 100)
ParticleSet
objectThis function is used to create a ParticleSet
object from a vector/list of tab separated text files, each of one containing one line for each
particle in the related frame, alongside with its coordinates and if available, the computed features
The number of frames is also specified, as just a subset of the particle lists can be used for this
read.particles(particle.files, nframes = NULL)
read.particles(particle.files, nframes = NULL)
particle.files |
Vector of strings containing the locations where the particle coordinates are to be found, or alternatively, the path to the folder |
nframes |
Number of frames that will constitute the |
An object of the ParticleSet
class
Federico Marini, [email protected], 2014
## see vignette and export.particles
## see vignette and export.particles
A more flexible and stylish alternative to replicate the behaviour of the repmat function of MATLAB
repmat(a, n, m)
repmat(a, n, m)
a |
The matrix to copy |
n |
The n value for the tiling |
m |
The m value for the tiling |
Creates a large matrix consisting of an m-by-n tiling of copies of a.
Robin Hankin, 2001
http://cran.r-project.org/doc/contrib/R-and-octave.txt
Frames
objectRotation is performed exploiting the rotate function of the EBImage
package. Could be automated if support for coordinate/pixel interaction is included
rotate.Frames(frames, angle, testing = FALSE)
rotate.Frames(frames, angle, testing = FALSE)
frames |
A |
angle |
The rotation angle (clockwise) specified in degrees |
testing |
Logical, whether to just test the rotation or to actually perform it. Default set to |
A Frames
object containing the rotated frames
Federico Marini, [email protected], 2014
data("MesenteriumSubset") rotate.Frames(MesenteriumSubset,angle = 40)
data("MesenteriumSubset") rotate.Frames(MesenteriumSubset,angle = 40)
Frames
objectAn input Frames
object is subject to subsetting. This function is useful e.g. when the trajectory of interest
is presenting gaps (i.e. does not actually include a frame)
select.Frames(frames, framesToKeep = 1, ...)
select.Frames(frames, framesToKeep = 1, ...)
frames |
A |
framesToKeep |
A vector containing the indexes of the frames to keep in the selection |
... |
Arguments to be passed to methods |
A Frames
object, composed by the subset of frames of the input Frames
Federico Marini, [email protected], 2014
data("MesenteriumSubset") select.Frames(MesenteriumSubset, framesToKeep = c(1:10, 14:20))
data("MesenteriumSubset") select.Frames(MesenteriumSubset, framesToKeep = c(1:10, 14:20))
ParticleSet
objectAccording to parameters of interests, such as size, eccentricity/shape, filters out the particles that do not satisfy the indicated requirements
select.particles(particleset, min.area = 1, max.area = 1000)
select.particles(particleset, min.area = 1, max.area = 1000)
particleset |
A |
min.area |
Size in pixels of the minimum area needed to detect the object as a potential particle of interest |
max.area |
Size in pixels of the maximum area allowed to detect the object as a potential particle of interest |
A ParticleSet
object
Federico Marini, [email protected], 2014
data("candidate.platelets") selected.platelets <- select.particles(candidate.platelets, min.area = 5) selected.platelets
data("candidate.platelets") selected.platelets <- select.particles(candidate.platelets, min.area = 5) selected.platelets
flowcatchR
Launches a Shiny Web Application for interactive data exploration. Default data loaded are
the frames from the MesenteriumSubset
object, custom values can be inserted by typing
the location of the data stored in a local folder. The Application is structured in a variety
of tabs that mirror the steps in the usual workflow in time-lapse microscopy images. These can
allow the user to interactively explore the parameters and their effect in the reactive
framework provided by Shiny.
shinyFlow()
shinyFlow()
The Shiny Application is launched in the web browser
Federico Marini, [email protected], 2015
## Not run: shinyFlow()
## Not run: shinyFlow()
Frames
objectDisplay conveniently a Frames
object
## S4 method for signature 'Frames' show(object)
## S4 method for signature 'Frames' show(object)
object |
A |
This returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("MesenteriumSubset") print(MesenteriumSubset)
data("MesenteriumSubset") print(MesenteriumSubset)
KinematicsFeatures
objectDisplaying conveniently a KinematicsFeatures
object
## S4 method for signature 'KinematicsFeatures' show(object)
## S4 method for signature 'KinematicsFeatures' show(object)
object |
A |
This returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) traj11features <- kinematics(platelets.trajectories,trajectoryIDs = 11) print(traj11features)
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) traj11features <- kinematics(platelets.trajectories,trajectoryIDs = 11) print(traj11features)
KinematicsFeatureSet
objectDisplay conveniently a KinematicsFeatureSet
object
## S4 method for signature 'KinematicsFeaturesSet' show(object)
## S4 method for signature 'KinematicsFeaturesSet' show(object)
object |
A |
This returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) alltrajs.features <- kinematics(platelets.trajectories) print(alltrajs.features)
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) alltrajs.features <- kinematics(platelets.trajectories) print(alltrajs.features)
LinkedParticleSet
objectDisplay conveniently a LinkedParticleSet
object
## S4 method for signature 'LinkedParticleSet' show(object)
## S4 method for signature 'LinkedParticleSet' show(object)
object |
A |
This returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("candidate.platelets") linked.platelets <- link.particles(candidate.platelets,L=26,R=3,epsilon1=0, epsilon2=0,lambda1=1,lambda2=0,penaltyFunction=penaltyFunctionGenerator(), include.area=FALSE) print(linked.platelets)
data("candidate.platelets") linked.platelets <- link.particles(candidate.platelets,L=26,R=3,epsilon1=0, epsilon2=0,lambda1=1,lambda2=0,penaltyFunction=penaltyFunctionGenerator(), include.area=FALSE) print(linked.platelets)
ParticleSet
objectDisplay conveniently a ParticleSet
object
## S4 method for signature 'ParticleSet' show(object)
## S4 method for signature 'ParticleSet' show(object)
object |
A |
This returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("candidate.platelets") print(candidate.platelets)
data("candidate.platelets") print(candidate.platelets)
TrajectorySet
objectDisplay conveniently a TrajectorySet
object
## S4 method for signature 'TrajectorySet' show(object)
## S4 method for signature 'TrajectorySet' show(object)
object |
A |
This returns an invisible NULL
.
Federico Marini, [email protected], 2014
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) print(platelets.trajectories)
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets) print(platelets.trajectories)
This function combines all classes related to a single experiment in order to deliver a clickable feedback on one of the frames.
snap( raw.frames, binary.frames, particleset, trajectoryset, frameID = 1, infocol = "yellow", infocex = 1, showVelocity = FALSE )
snap( raw.frames, binary.frames, particleset, trajectoryset, frameID = 1, infocol = "yellow", infocex = 1, showVelocity = FALSE )
raw.frames |
A |
binary.frames |
A |
particleset |
A |
trajectoryset |
A |
frameID |
The ID of the frame to inspect |
infocol |
The color to use for plotting the contours and the information on the clicked particle |
infocex |
The numeric character expansion value as in |
showVelocity |
Logical, whether to display additional information on the instantaneous velocity of the particle |
An image of the selected frame, rendered in R native graphics, and additionally a list with the coordinates as well as the trajectory ID of the particle closest to the clicked location
Federico Marini, [email protected], 2015
## Not run: data(MesenteriumSubset) binary.frames <- preprocess.Frames(channel.Frames(MesenteriumSubset,"red")) particleset <- particles(MesenteriumSubset,binary.frames,"red") trajectoryset <- trajectories(particleset) snap(MesenteriumSubset,binary.frames,particleset,trajectoryset,frameID=1) ## End(Not run)
## Not run: data(MesenteriumSubset) binary.frames <- preprocess.Frames(channel.Frames(MesenteriumSubset,"red")) particleset <- particles(MesenteriumSubset,binary.frames,"red") trajectoryset <- trajectories(particleset) snap(MesenteriumSubset,binary.frames,particleset,trajectoryset,frameID=1) ## End(Not run)
Conversion from (radius,theta) to (x,y)
toCartesianCoords(Theta, Radius)
toCartesianCoords(Theta, Radius)
Theta |
The Theta angle |
Radius |
The radius value in polar coordinates |
A list containing Theta and Radius, as in polar coordinates
Federico Marini, [email protected], 2014
Conversion from (x,y) to (radius,theta)
toPolarCoords(x, y)
toPolarCoords(x, y)
x |
x coordinate |
y |
y coordinate |
A list containing Theta and Radius, as in polar coordinates
Federico Marini, [email protected], 2014
Generates a TrajectorySet
object from a (Linked
)ParticleSet
trajectories(particleset, verbose = FALSE, ...)
trajectories(particleset, verbose = FALSE, ...)
particleset |
A ( |
verbose |
Logical, currently not used - could be introduced for providing additional info on the trajectories |
... |
Arguments to be passed to methods |
A TrajectorySet
object
Federico Marini, [email protected], 2014
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets)
data("candidate.platelets") platelets.trajectories <- trajectories(candidate.platelets)
S4 class for storing information on the trajectories identified, including whether there were gaps, the number of points, and more
.Data
A list storing the information for the particles
channel
A character vector, can be 'red','green', or 'blue'. It refers to which channel the particles were detected