## ----install-BiocManager, eval = FALSE----------------------------------------
# if (!requireNamespace("BiocManager", quietly = TRUE)) {
#   install.packages("BiocManager")
# }
# BiocManager::install("ImageArray")

## ----load-libs, message=FALSE, warning=FALSE----------------------------------
library(ImageArray)
library(BiocFileCache)
library(RBioFormats)
library(EBImage)
library(magick)
library(ggplot2)
library(shiny)

## ----read---------------------------------------------------------------------
# make random EBImage image
f <- system.file("images", "sample.png", package = "EBImage")

# read image with EBImage
img <- readImage(f)

## ----write--------------------------------------------------------------------
output_h5 <- tempfile(fileext = ".h5")
imgarray <- writeImageArray(img, output = output_h5, n.levels = 2)
imgarray

## ----visualize_read-----------------------------------------------------------
bfa.raster <- as.raster(imgarray, level = 2)
plot(bfa.raster)

## ----visualize_read2----------------------------------------------------------
# visualize
bfa.raster <- as.raster(imgarray, max.pixel.size = 400)
dim(bfa.raster)

## ----rotate-------------------------------------------------------------------
imgarray_rotated <- rotate(imgarray, angle = 90)
imgarray_rotated

# visualize
bfa.raster <- as.raster(imgarray_rotated, max.pixel.size = 400)
plot(bfa.raster)

## ----crop---------------------------------------------------------------------
imgarray_cropped <- imgarray[300:500, 200:300]
imgarray_cropped

# visualize
bfa.raster <- as.raster(imgarray_cropped, max.pixel.size = 120)
plot(bfa.raster)

# cropping with indices
imgarray_cropped <- crop(imgarray, index = list(300:500, 200:300))
imgarray_cropped

## ----crop_named---------------------------------------------------------------
# cropping with named indices
imgarray_cropped <- crop(imgarray, index = list(x=300:500, y=200:300))
imgarray_cropped

## ----crop_named2--------------------------------------------------------------
# index in yx order, results in same slice as xy order
imgarray_cropped <- crop(imgarray, index = list(y=200:300, x=300:500))
imgarray_cropped

# subsetting with `[` provides a different slice
imgarray_cropped <- imgarray[200:300, 300:500]
imgarray_cropped

## ----read3--------------------------------------------------------------------
# make random EBImage image
f <- system.file("images", "sample.png", package = "EBImage")

# read image with magick
img <- image_read(f)

# create image array
output_h5 <- tempfile(fileext = ".h5")
imgarray <- writeImageArray(img, output = output_h5)
imgarray

# plot
bfa.raster <- as.raster(imgarray, level = 2)
plot(bfa.raster)

## ----ometiff------------------------------------------------------------------
ome.tiff.file <- system.file("extdata", "xy_12bit__plant.ome.tiff",
                             package = "ImageArray")
read.metadata(ome.tiff.file)

# define ImageArray object
imgarray <- createImageArray(ome.tiff.file, series = 1, resolution = 1:2)
imgarray

## ----ometiffvis2, out.width="50%"---------------------------------------------
bfa.raster <- as.raster(imgarray, max.pixel.size = 300)
plot(bfa.raster)
dim(bfa.raster)

bfa.raster <- as.raster(imgarray, min.pixel.size = 300)
plot(bfa.raster)
dim(bfa.raster)

bfa.raster <- as.raster(imgarray, level = 1)
plot(bfa.raster)
dim(bfa.raster)

## ----he-----------------------------------------------------------------------
image_file <- paste(
  "https://cf.10xgenomics.com/samples/xenium/1.0.1",
  "Xenium_FFPE_Human_Breast_Cancer_Rep1",
  "Xenium_FFPE_Human_Breast_Cancer_Rep1_he_image.ome.tif",
  sep = "/")
library(BiocFileCache)
bfc <- BiocFileCache()
image_file <- bfcrpath(bfc, image_file)

## ----he_imgarray--------------------------------------------------------------
imgarray <- createImageArray(image_file,
                             series = 1,
                             resolution = 1:6)
imgarray

## ----speedy_vis, out.width="70%"----------------------------------------------
# crop image
imgarray_sub <- crop(imgarray, index = list(16000:19000, 7000:10000, NULL))

# convert to raster
img_raster <- as.raster(imgarray_sub, max.pixel.size = 800)
dim(img_raster)

# plot with ggplot
ggplot(data.frame(x = 0, y = 0), aes(x, y)) +
  coord_fixed(expand = FALSE,
              xlim = c(0, dim(img_raster)[2]),
              ylim = c(0, dim(img_raster)[1])) +
  annotation_raster(img_raster, 
                    0, dim(img_raster)[2],
                    dim(img_raster)[1], 0, interpolate = FALSE)

## ----speedy_vis_shiny---------------------------------------------------------
if (interactive()) {

  # variables
  dimimg <- dim(imgarray)
  max.pixel.size <- 800

  # Define UI
  ui <- fluidPage(
    sliderInput("x_slider", label = "X", min = 1,
                max = dimimg[1], value = c(16000, 19000)),
    sliderInput("y_slider", label = "Y", min = 1,
                max = dimimg[2], value = c(7000, 10000)),
    mainPanel(
      plotOutput(outputId = "scatterPlot")
    )
  )

  # Define server logic
  server <- function(input, output) {

    output$scatterPlot <- renderPlot({

      # crop image
      indx <- seq(input$x_slider[1], input$x_slider[2])
      indy <- seq(input$y_slider[1], input$y_slider[2])
      imgarray_sub <- crop(imgarray, index = list(indx, indy, NULL))

      # convert to raster
      img_raster <- as.raster(imgarray_sub, max.pixel.size = max.pixel.size)

      # plot with ggplot
      imgggplot <- ggplot(data.frame(x = 0, y = 0), aes(x, y)) +
        coord_fixed(expand = FALSE,
                    xlim = c(0, dim(img_raster)[2]),
                    ylim = c(0, dim(img_raster)[1])) +
        annotation_raster(img_raster, 
                          0, dim(img_raster)[2],
                          dim(img_raster)[1], 0, interpolate = FALSE)
      imgggplot
    })
  }

  # Run the app
  shinyApp(ui = ui, server = server)
}

## ----sessionInfo, echo=FALSE--------------------------------------------------
sessionInfo()

