Title: | A package for the analysis of anatomical tissue structures in spatial omics data |
---|---|
Description: | `sosta` (Spatial Omics STructure Analysis) is a package for analyzing spatial omics data to explore tissue organization at the anatomical structure level. It reconstructs morphologically relevant structures based on molecular features or cell types. It further calculates a range of structural and shape metrics to quantitatively describe tissue architecture. The package is designed to integrate with other packages for the analysis of spatial (omics) data. |
Authors: | Samuel Gunz [aut, cre] , Mark D. Robinson [aut, fnd] |
Maintainer: | Samuel Gunz <[email protected]> |
License: | GPL (>= 3) + file LICENSE |
Version: | 0.99.3 |
Built: | 2025-01-20 03:29:34 UTC |
Source: | https://github.com/bioc/sosta |
Function to estimate the intensity image of a point pattern
.intensityImage(ppp, mark_select = NULL, bndw = NULL, dim)
.intensityImage(ppp, mark_select = NULL, bndw = NULL, dim)
ppp |
point pattern object of class |
mark_select |
character; name of mark that is to be selected for the reconstruction |
bndw |
bandwidth of kernel density estimator |
dim |
numeric; x dimension of the final reconstruction. |
list; list with the intensity image and the bandwidth and dimension parameters
Function to estimate the intensity threshold for the reconstruction of spatial structures
.intensityThreshold(density_image, steps = 250)
.intensityThreshold(density_image, steps = 250)
density_image |
real-valued pixel image; output from the function |
steps |
numeric; value used to filter the density estimates, where only
densities greater than the maximum value divided by |
numeric; estimated threshold
Converts a binary matrix to an sf polygon
binaryImageToSF(binaryMatrix, xmin, xmax, ymin, ymax)
binaryImageToSF(binaryMatrix, xmin, xmax, ymin, ymax)
binaryMatrix |
matrix; binary matrix |
xmin |
integer; minimum x coordinate of the coordinate system |
xmax |
integer; maximum x coordinate of the coordinate system |
ymin |
integer; minimum y coordinate of the coordinate system |
ymax |
integer; maximum y coordinate of the coordinate system |
sf object
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) plot(poly_R)
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) plot(poly_R)
Estimate reconstruction parameters from a set of images
estimateReconstructionParametersSPE( spe, marks, image_col, mark_select = NULL, nimages = NULL, fun = "bw.diggle", dim = 500, ncores = 1, plot_hist = TRUE )
estimateReconstructionParametersSPE( spe, marks, image_col, mark_select = NULL, nimages = NULL, fun = "bw.diggle", dim = 500, ncores = 1, plot_hist = TRUE )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
image_col |
character; name of a column in |
mark_select |
character; name of mark that is to be selected for the reconstruction |
nimages |
integer; number of images for the estimation. Will be randomly sampled |
fun |
character; function to estimate the kernel density. Default bw.diggle. |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speed up computation but lead to less exact reconstruction. Default = 500 |
ncores |
numeric; number of cores for parallel processing using |
plot_hist |
logical; if histogram of estimated densities and thresholds should be plotted. Default = TRUE |
tibble; tibble with estimated intensities and thresholds
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) spe_sel <- spe[, spe[["image_name"]] %in% c("E02", "E03", "E04")] estimateReconstructionParametersSPE(spe_sel, marks = "cell_category", image_col = "image_name", mark_select = "islet", plot_hist = TRUE )
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) spe_sel <- spe[, spe[["image_name"]] %in% c("E02", "E03", "E04")] estimateReconstructionParametersSPE(spe_sel, marks = "cell_category", image_col = "image_name", mark_select = "islet", plot_hist = TRUE )
Estimate the intensity threshold for the reconstruction of spatial structures
findIntensityThreshold(ppp, mark_select = NULL, bndw = NULL, dim, steps = 250)
findIntensityThreshold(ppp, mark_select = NULL, bndw = NULL, dim, steps = 250)
ppp |
point pattern object of class |
mark_select |
character; name of mark that is to be selected for the reconstruction |
bndw |
numeric; bandwith of the sigma parameter in the density estimation,
if no value is given the bandwith is estimated using cross validation with
the |
dim |
numeric; x dimension of the final reconstruction. |
steps |
numeric; value used to filter the density estimates, where only
densities greater than the maximum value divided by |
numeric; estimated intensity threshold
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) ppp <- SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04") findIntensityThreshold(ppp, mark_select = "islet", dim = 250)
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) ppp <- SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04") findIntensityThreshold(ppp, mark_select = "islet", dim = 250)
Function to get the dimension based on dim of y axis
getDimXY(ppp, ydim)
getDimXY(ppp, ydim)
ppp |
point pattern object of class |
ydim |
dimension of y axis |
vector; vector with x and y dimension
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) ppp <- SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04" ) getDimXY(ppp, 500)
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) ppp <- SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04" ) getDimXY(ppp, 500)
Calculate mean shape metrics of a set of polygons
meanShapeMetrics(totalShapeMetricMatrix)
meanShapeMetrics(totalShapeMetricMatrix)
totalShapeMetricMatrix |
matrix of shape metrics |
matrix; matrix of mean shape metrics
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) islet_poly <- reconstructShapeDensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet", dim = 500 ) shape_metrics <- totalShapeMetrics(islet_poly) meanShapeMetrics(shape_metrics)
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) islet_poly <- reconstructShapeDensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet", dim = 500 ) shape_metrics <- totalShapeMetrics(islet_poly) meanShapeMetrics(shape_metrics)
Function to normalize coordinates between zero and one while keep scaling
normalizeCoordinates(coords)
normalizeCoordinates(coords)
coords |
matrix; matrix with coordinates |
matrix; coordinates scaled between 0 and 1
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) coords <- xyCoordinates(matrix_R) normalizeCoordinates(coords)
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) coords <- xyCoordinates(matrix_R) normalizeCoordinates(coords)
This function estimates the density of a spatial point pattern (ppp
),
thresholds the density to create a binary image, and then converts it
to a valid sf
object (polygons).
reconstructShapeDensity( ppp, mark_select = NULL, bndw = NULL, thres = NULL, dim )
reconstructShapeDensity( ppp, mark_select = NULL, bndw = NULL, thres = NULL, dim )
ppp |
point pattern object of class |
mark_select |
character; name of mark that is to be selected for the reconstruction |
bndw |
bandwidth of kernel density estimator |
thres |
intensity threshold for the reconstruction |
dim |
numeric; x dimension of the final reconstruction. |
sf object of class POLYGON
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) ppp <- SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04") thres <- findIntensityThreshold(ppp, mark_select = "islet", dim = 500) islet_poly <- reconstructShapeDensity(ppp, mark_select = "islet", thres = thres, dim = 500) plot(islet_poly)
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) ppp <- SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04") thres <- findIntensityThreshold(ppp, mark_select = "islet", dim = 500) islet_poly <- reconstructShapeDensity(ppp, mark_select = "islet", thres = thres, dim = 500) plot(islet_poly)
Reconstruct structure from spe object with given image id
reconstructShapeDensityImage( spe, marks, image_col, image_id, mark_select, dim = 500, bndw = NULL, thres = NULL )
reconstructShapeDensityImage( spe, marks, image_col, image_id, mark_select, dim = 500, bndw = NULL, thres = NULL )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
image_col |
character; name of a column in |
image_id |
character; image id, must be present in image_col |
mark_select |
character; name of mark that is to be selected for the reconstruction |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speed up computation but lead to less exact reconstruction. Default = 500 |
bndw |
numeric; bandwith of the sigma parameter in the density estimation,
if no value is given the bandwith is estimated using cross validation with
the |
thres |
numeric; intensity threshold for the reconstruction |
sf object of class POLYGON
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) islet_poly <- reconstructShapeDensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet", dim = 500 ) plot(islet_poly)
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) islet_poly <- reconstructShapeDensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet", dim = 500 ) plot(islet_poly)
Reconstruct structure from spatial experiment object per image id
reconstructShapeDensitySPE( spe, marks, image_col, mark_select, dim = 500, bndw = NULL, thres, ncores = 1 )
reconstructShapeDensitySPE( spe, marks, image_col, mark_select, dim = 500, bndw = NULL, thres, ncores = 1 )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
image_col |
character; name of a column in |
mark_select |
character; name of mark that is to be selected for the reconstruction |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speed up computation but lead to less exact reconstruction. Default = 500 |
bndw |
numeric; bandwith of the sigma parameter in the density estimation,
if no value is given the bandwith is estimated using cross validation with
the |
thres |
numeric; intensity threshold for the reconstruction |
ncores |
numeric; number of cores for parallel processing using
|
simple feature collection
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) spe_sel <- spe[, spe[["image_name"]] %in% c("E02", "E03", "E04")] all_islets <- reconstructShapeDensitySPE(spe_sel, marks = "cell_category", image_col = "image_name", mark_select = "islet", bndw = sigma, thres = 0.0025 ) all_islets
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) spe_sel <- spe[, spe[["image_name"]] %in% c("E02", "E03", "E04")] all_islets <- reconstructShapeDensitySPE(spe_sel, marks = "cell_category", image_col = "image_name", mark_select = "islet", bndw = sigma, thres = 0.0025 ) all_islets
This function plots the intensity of a point pattern image and displays a histogram of the intensity values. Note that intensities less than largest intensity value divided by 250 are not displayed in the histogram.
shapeIntensityImage( spe, marks, image_col, image_id, mark_select, bndw = NULL, dim = 500 )
shapeIntensityImage( spe, marks, image_col, image_id, mark_select, bndw = NULL, dim = 500 )
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
image_col |
character; name of a column in |
image_id |
character; image id, must be present in image_col |
mark_select |
character; name of mark that is to be selected for the reconstruction |
bndw |
numeric; bandwith of the sigma parameter in the density estimation,
if no value is given the bandwith is estimated using cross validation with
the |
dim |
numeric; x dimension of the final reconstruction. A lower resolution speeds up computation but lead to less exact reconstruction. Default = 500 |
ggplot object with intensity image and histogram
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) shapeIntensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet" )
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) shapeIntensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet" )
Calculate a set of shape metrics of a polygon
shapeMetrics(sfPoly)
shapeMetrics(sfPoly)
sfPoly |
POLYGON of class sfc |
list; list of shape metrics
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) shapeMetrics(poly_R)
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) shapeMetrics(poly_R)
SpatialExperiment
object to a ppp
objectFunction to convert spatial coordinates of a SpatialExperiment
object to a ppp
object
SPE2ppp(spe, marks, image_col = NULL, image_id = NULL)
SPE2ppp(spe, marks, image_col = NULL, image_id = NULL)
spe |
SpatialExperiment; a object of class |
marks |
character; name of column in |
image_col |
character; name of a column in |
image_id |
character; image id, must be present in image_col |
ppp; object of type ppp
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04")
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) SPE2ppp(spe, marks = "cell_category", image_col = "image_name", image_id = "E04")
Title
st_calculateCurvature(sfPoly, smoothness = 5)
st_calculateCurvature(sfPoly, smoothness = 5)
sfPoly |
|
smoothness |
list; curvature measures |
list; list of curvatures values
https://stackoverflow.com/questions/62250151/calculate-curvature-of-a-closed-object-in-r
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) st_calculateCurvature(poly_R)
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) st_calculateCurvature(poly_R)
Calculate curl of a polygon
st_calculateShapeCurl(sfPoly)
st_calculateShapeCurl(sfPoly)
sfPoly |
|
numeric; the curl of the polygon
matrix_R <- matrix(c( 1, 1, 1, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0, 1, 1, 0, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1 ), nrow = 7, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) st_calculateShapeCurl(poly_R)
matrix_R <- matrix(c( 1, 1, 1, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1, 1, 1, 1, 1, 1, 0, 1, 1, 0, 1, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 0, 0, 1, 1 ), nrow = 7, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) st_calculateShapeCurl(poly_R)
Calculate the length of feature axes of an sf polygon
st_feature_axes(sfPoly)
st_feature_axes(sfPoly)
sfPoly |
|
list; list containing the major and minor axis lengths
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) st_feature_axes(poly_R)
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) poly_R <- binaryImageToSF(matrix_R, xmin = 0, xmax = 1, ymin = 0, ymax = 1) st_feature_axes(poly_R)
Calculate a set of shape metrics of a set of polygons
totalShapeMetrics(sfInput)
totalShapeMetrics(sfInput)
sfInput |
|
Calculate a set of shape metrics of a set of polygons.
The function calculates all metrics that are implemented in the function
shapeMetrics()
matrix; matrix of shape metrics
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) islet_poly <- reconstructShapeDensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet", dim = 500 ) totalShapeMetrics(islet_poly)
spe <- imcdatasets::Damond_2019_Pancreas("spe", full_dataset = FALSE) islet_poly <- reconstructShapeDensityImage(spe, marks = "cell_category", image_col = "image_name", image_id = "E04", mark_select = "islet", dim = 500 ) totalShapeMetrics(islet_poly)
Function to extract x y coordinates from binary image
xyCoordinates(inputMatrix)
xyCoordinates(inputMatrix)
inputMatrix |
a binary matrix |
matrix; matrix with x,y coordinates of the cell of the input matrix
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) xyCoordinates(matrix_R)
matrix_R <- matrix(c( 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 0, 1, 1, 0, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 ), nrow = 9, byrow = TRUE) xyCoordinates(matrix_R)