Can execute any R function at every voxel for a set of MINC files

mincApply(filenames, function.string, mask = NULL, maskval = NULL,
  reduce = FALSE)

Arguments

filenames

The MINC files over which to apply the function. Have to be the same sampling.

function.string

The function which to apply. Can only take a single argument, which has to be 'x'. See details and examples.

mask

The filename of a mask - function will only be evaluated inside the mask.

maskval

An integer specifying the value inside the mask where to apply the function. If left blank (the default) then anything above 0.5 will be considered inside the mask. This argument only works for mincApply, not pMincApply.

reduce

Whether to filter masked voxels from the resulting object

Value

out The output is a matrix with the same number of rows a the file sizes and the same number of columns as output by the function that was applied. Cast into one of the MINC classes to make writing it out with mincWriteVolume easier.

Details

mincApply allows one to execute any R function at every voxel of a set of files. There are two variants: mincApply, which works inside the current R session. Unless the function to be applied takes a single argument a wrapper function has to be written. This is illustrated in the examples; briefly, the wrapper function takes a single argument, called 'x', which will take on the voxel values at every voxel. The function has to be turned into a string; the quote function can be very handy. The output of this function also has to be a simple vector or scalar. Note that interpreted R can be very slow. Mindnumbingly slow. It therefore pays to write optimal functions or, whenever available, use the optimized equivalents. In other words and to give one example, use mincLm rather than applying lm, and if lm has to really be applied, try to use lm.fit rather than plain lm. When using the pbs method, one can also set the options --> TMPDIR,MAX_NODES and WORKDIR
If you encounter memory issues, it could be due to minc file caching. Consider trying with the environment variable MINC_FILE_CACHE_MB set to a small value like 1.

See also

mincWriteVolume,mincMean,mincSd,mincLm,mincAnova

Examples

# NOT RUN {
getRMINCTestData()
gf <- read.csv("/tmp/rminctestdata/test_data_set.csv")
ma <- mincApply(gf$jacobians_fixed_2,quote(mean(x)));
mincWriteVolume(ma, "means.mnc")
### run the whole thing in parallel on the local machine
testFunc <- function (x) { return(c(1,2))}
pout <- pMincApply(gf$jacobians_fixed_2,
                   quote(testFunc(x)),
                   workers = 4,
                   global = c('gf','testFunc'))
mincWriteVolume(pout, "pmincOut.mnc")

### pbs example 1 (hpf)
getRMINCTestData()
gf <- read.csv("/tmp/rminctestdata/test_data_set.csv")
testFunc <- function (x) { return(c(1,2))}
options(TMPDIR="/localhd/$PBS_ID")
pout <- pMincApply(gf$jacobians_fixed_2,
                   quote(testFunc(x)),
                   workers = 4,
                   method="pbs",
                   global = c('gf','testFunc'))

### pbs example 2 (scinet)
getRMINCTestData()
gf <- read.csv("/tmp/rminctestdata/test_data_set.csv")
testFunc <- function (x) { return(c(1,2))}
options(WORKDIR="$SCRATCH")
options(MAX_NODES=8)
options(TMP_DIR="/tmp")
pout <- pMincApply(gf$jacobians_fixed_2,
                   quote(testFunc(x)),
                   modules=c("intel","openmpi","R/3.1.1"),
                   workers = 4,
                   method="pbs",
                   global = c('gf','testFunc'))

#NOTE 1: On SCINET jobs are limited to 32*48 hours of cpu time
#NOTE 2: On SCINET the Rmpi library must be compiled for use with R 3.1.1
# }