#' Find (if possible) an interior point of a polytope by solving a linear program #' #' Solves the LP: max over {w,eps} {eps | SAw - eps >= Sb, 0 < eps <= epsbound} #' #' where S is diag(s). If at the solution eps > 0, then w is a valid interior point #' otherwise the LP fails to find an interior point, and another s must be tried. #' #' @param A is a n by m matrix of hyperplane slope coefficients #' @param b is an n vector of hyperplane intercept coefficients #' @param s is an n vector of signs #' @param epsbound is a scalar tolerance controlling how close the witness point #' can be to an edge of the polytope #' @param epstol is a scalar tolerance for the LP convergence #' @param presolve controls whether Mosek should presolve the LP #' @param verb controls verbosity of Mosek solution #' #' @return List with components: #' itemize{ #' item w proposed interior point at solution #' item fail indicator of whether w is a valid interior point #' } #' @importFrom Matrix Rmasek #' @export witness <- function(A,b, s,epsbound = 1,epstol = 1e-07, presolve = 1, verb = 0){
if(!is.list(s)) s <- list(s) m = ncol(A) n = nrow(A) if (is.null(m)){ m = length(A) n = 1 } # Interior defines a Linear Programming Problem Interior <- list() Interior$sense <- "max" Interior$c <- c(rep(0, m), 1) blx <- c(rep(-Inf,m),0) # lower bound of eps set as 0 bux <- c(rep(Inf,m),epsbound) Interior$bx <- rbind(blx,bux) lapply(s, FUN=function(item){ S <- if (length(item) > 1) diag(item) else item A = cbind(S %*% A, matrix(-1, nrow = n,1)) Interior$A <- Matrix::Matrix(A, sparse = TRUE) buc <- c(rep(Inf, n)) blc <- c(S %*% b) Interior$bc <- rbind(blc,buc) Interior$iparam$MSK_IPAR_PRESOLVE_USE <- presolve r <- Rmosek::mosek(Interior,opts = list(verbose = verb)) eps <- r$sol$bas$xx[3] fail <- eps <= epstol list(w = r$sol$bas$xx[1:2], fail = fail, s = item, eps = eps) })
} #' New Incremental Cell Enumeration (in) R #' #' Modified version of the algorithm of Rada and Cerny (2018). #' The main modifications include preprocessing as hyperplanes are added #' to determine which new cells are created, thereby reducing the number of #' calls to the witness function to solve LPs, and treatment of degenerate #' configurations as well as those in “general position.” When the hyperplanes #' are in general position the number of cells (polytopes) is determined by the #' elegant formula of Zazlavsky (1975) deqn{m = {n choose d} - n + 1}. In #' degenerate cases, i.e. when hyperplanes are not in general position, the #' number of cells is more complicated as considered by Alexanderson and Wetzel (1981). #' The function code{polycount} is provided to check agreement with their results #' in an effort to aid in the selection of tolerances for the code{witness} function. #' Current version is intended for use with eqn{d = 2}, but the algorithm is adaptable to #' eqn{d > 2}. #' #' @param A is a n by m matrix of hyperplane slope coefficients #' @param b is an n vector of hyperplane intercept coefficients #' @param initial origin for the interior point vectors code{w} #' @param epsbound is a scalar tolerance controlling how close the witness point #' can be to an edge of the polytope #' @param epstol is a scalar tolerance for the LP convergence #' @param verb controls verbosity of Mosek solution #' #' @return A list with components: #' itemize{ #' item SignVector a n by m matrix of signs determining position of cell relative #' to each hyperplane. #' item w a d by m matrix of interior points for the m cells #' } #' @references #' Alexanderson, G.L and J.E. Wetzel, (1981) Arrangements of planes in space, #' Discrete Math, 34, 219–240. #' Rada, M. and M. Cerny (2018) A new algorithm for the enumeration of cells #' of hyperplane arrangements and a comparison with Avis and Fukada's reverse #' search, SIAM J. of Discrete Math, 32, 455-473. #' Zaslavsky, T. (1975) Facing up to arrangements: Face-Count Formulas for #' Partitions of Space by Hyperplanes, Memoirs of the AMS, Number 154. #' @export #' NICER <- function(A, b, initial = c(0,0), verb = TRUE, epsbound = 1, epstol = 1e-07){
if (sum(abs(A %*% initial - b)) < epstol) stop("Error: initial point falls on some of the hyperplanes") n <- nrow(A) maxdim <- choose(n,2) + n + 1 # maximum column dimension of SignVector SignVector <- matrix(NA, n, maxdim) w <- matrix(0, 2, maxdim) eps <- matrix(0,1,maxdim) mode(SignVector) <- "integer" i <- 0 if (t(as.vector(A[(1:(i+1)),]))%*%as.vector(initial)>b[(1:(i+1))]) SignVector[(i+1),1] <- 1L else SignVector[(i+1),1] <- -1L w[,1] <- initial eps[,1] <- epsbound s <- SignVector[1:(i+1),1] * c(rep(1, i),-1) test <- witness(A[(1:(i+1)),], b[(1:(i+1))], s = s, epsbound = epsbound, epstol = epstol) if (!test[[1]]$fail){ SignVector[(i+1),2] <- test[[1]]$s w[,2] <- c(test[[1]]$w) eps[,2] <- c(test[[1]]$eps) } # for each i columns of SignVector and w are filled with values # test if the corresponding w is on the upper side of the hyperplane for (i in 1:(n-1)){ # adding hyperplane 2 to n tempdim = sum(!is.na(SignVector[1,]) ) SignVector[i+1,1:tempdim] <- sign(as.vector(t(w[,1:tempdim])%*%A[i+1,] - b[i+1])) # the following lines avoid the situation that the newly added line crosses # any of the existing interior points. If so, find a new interior point # N.B. seems using the same epsbound does not help, so I've tentatively # set it at 0.1, this needs further testing) if(any(SignVector[i+1,1:tempdim]==0)){ for(k in which(SignVector[i+1,1:tempdim]==0)){ testk <- witness(A[1:(i+1),], b[1:(i+1)], s= c(SignVector[1:i, k],1), epsbound=epsbound, epstol=epstol) if(!testk[[1]]$fail){ w[,k] <- c(testk[[1]]$w) eps[,k] <- c(testk[[1]]$eps) signs <- sign(t(w[,k])%*%A[i+1,] - b[i+1]) mode(signs) <- "integer" SignVector[i+1,k] <- signs } } } # if the newly added line does not coincide with previous existing lines # pre-process to see which polygons creates new cells with the addition the hyperplane # otherwise, no new polygons should be added and move on. dup <- duplicated(cbind(A[1:(i+1),]/A[1:(i+1),1], b[1:(i+1)]/A[1:(i+1),1]), MARGIN=1) if (dup[i+1] == 0){ if (i > 1){ vertex <- matrix(NA, nrow = 2, ncol = i) for (k in 1:i){ if (A[k,1]/A[k,2] != A[(i+1),1]/A[(i+1),2]) # new line not parallel to existing lines vertex[,k] <- solve(rbind(A[k,], A[(i+1),]), rbind(b[k], b[(i+1)])) } if (sum(!is.na(vertex))==0) colset = c(1:tempdim) else{ # feas records sign constraints of each existing hyperplane wrt to each vertex feas <- matrix(NA_integer_, nrow = i, ncol = i) Atemp = A[1:i,] btemp = b[1:i] for (j in 1:i){ if (!is.na(vertex[1,j])) feas[j, -j] <- sign(c(Atemp[-j,] %*% vertex[,j]) - btemp[-j]) else feas[j,] <- NA_integer_ } rownames(feas) <- 1:i feas <- feas[which(rowSums(is.na(feas))!=i),,drop=FALSE] feas[which(feas==0)] <- NA_integer_ feas <- unique(feas,MARGIN=1) tempSignVector <- SignVector[1:i, 1:tempdim] # if all rows have only one NA (i.e. general position case) # enumerate the NA by 1 and -1 exhausts all possible polygons that may have been crossed # otherwise, distinguish rows that have one NA (set1na) and more than one NA (set2na) # set2na only occurs when the new added hyperplane crossed an existing vertex, # colset collects all polygons that are crossed by the newly added hyperplane if (sum(is.na(feas))==nrow(feas)){ fea1 = feas fea1[is.na(fea1)] = 1L fea2 = feas fea2[is.na(fea2)] = -1L feaset = rbind(fea1,fea2) feaset = unique(t(feaset),MARGIN = 2) colset <-which(duplicated(cbind(tempSignVector,feaset),fromLast=TRUE,MARGIN=2)) } else{ set1na <- rowSums(is.na(feas))==1 set2na <- rowSums(is.na(feas))>1 fea1na1 <- feas[set1na,] fea1na1[is.na(fea1na1)] = 1L fea1na2 <- feas[set1na,] fea1na2[is.na(fea1na2)] = -1L fea1naset = rbind(fea1na1,fea1na2) fea1naset = unique(t(fea1naset),MARGIN=2) colset1na <- which(duplicated(cbind(tempSignVector,fea1naset), fromLast = TRUE, MARGIN = 2)) fea2na <- feas[set2na,,drop=FALSE] constraintmat <- abs(fea2na)==1 constraintmat <- constraintmat & !is.na(constraintmat) # picking out location of constraints in feas step <- 0.1*c(1,-1)/A[i+1,] colsettemp <- lapply(1:nrow(constraintmat), FUN=function(j){ constraint <- constraintmat[j,] realj <- as.integer(rownames(constraintmat)[j]) clines <- which(!constraint) clines <- clines[clines!=realj] linecols <- which(colSums(abs(tempSignVector[constraint,,drop=FALSE]-fea2na[j,constraint]))==0) remcols <- logical(length(linecols)) for(line in clines){ vert <- solve(A[c(i+1,line),],b[c(i+1,line)]) vert_m <- vert - step vert_p <- vert + step goodsigns <- sign(A[c(realj,line),]%*%cbind(vert_m,vert_p) - b[c(realj,line)]) csum1 <- colSums(abs(tempSignVector[c(realj,line),linecols, drop = FALSE]-goodsigns[,1])) != 0 csum2 <- colSums(abs(tempSignVector[c(realj,line),linecols, drop = FALSE]-goodsigns[,2])) != 0 remcols <- remcols | csum1 & csum2 } linecols <- linecols[!remcols] }) colset2na <- unique(unlist(colsettemp)) colset <- unique(c(colset1na, colset2na)) } } } else colset = c(1:2) # Find new interior point by flipping the sign of the last entry in relevant SignVector s = lapply(colset, FUN=function(x) SignVector[1:(i+1),x] * c(rep(1,i),-1)) testlist <- witness(A[(1:(i+1)),], b[(1:(i+1))], s = s, epsbound = epsbound, epstol = epstol) testlist.keep = which(!unlist(lapply(testlist,"[[",2))) # keep only those that has fail = FALSE if (length(testlist.keep)>0){ for(l in 1:length(testlist.keep)){ SignVector[1:(i+1),tempdim+l] <- testlist[[testlist.keep[l]]]$s w[,tempdim+l] <- c(testlist[[testlist.keep[l]]]$w) eps[,tempdim+l] <- c(testlist[[testlist.keep[l]]]$eps) } } } if ((verb != 0) && round((i+1)/10)==(i+1)/10) print(paste0(i+1, " :: Polygons: ", sum(!is.na(SignVector[1,]) )), quote=FALSE) } w = w[,colSums(is.na(SignVector))<nrow(SignVector)] SignVector = SignVector[,colSums(is.na(SignVector))<nrow(SignVector)] eps = eps[,1:ncol(w)] if (dim(SignVector)[2] != polycount(A,b)) print("warning: polycount does not match, consider changing epsbound.") list(w = w, SignVector = SignVector, eps = eps)
} #' Check Neighbouring Cell Counts #' #' Compare cell counts for each cell with its neighbours and return indices #' of the locally maximal cells. #' #' @param SignVector n by m matrix of signs produced by NICER #' #' @return Column indices of the cells that are locally maximal, i.e. those #' whose neighbours have strictly fewer cell counts. The corresponding #' interior points of these cells can be used as potential mass points #' for the NPMLE function code{rcbr.fit.KW}. #' @export neighbours <- function(SignVector){
SVcolsums <- colSums(SignVector==1L) SVorder <- order(SVcolsums,decreasing=TRUE) SVcolsums_od <- SVcolsums[SVorder] SVordered <- SignVector[,SVorder] elim <- logical(nc <- ncol(SignVector)) for(i in 1:nc){ if(!elim[i]){ count.i <- SVcolsums_od[i] tem_big <- SVordered[,i] - SVordered[,SVcolsums_od-count.i==1,drop=FALSE] tem_small <- SVordered[,i] - SVordered[,SVcolsums_od-count.i==-1,drop=FALSE] nbr_big <- which(colSums(abs(tem_big))==2L) if(!is.na(ind_b <- match(count.i+1,SVcolsums_od))) nbr_big <- nbr_big + ind_b - 1 nbr_small <- which(colSums(abs(tem_small))==2L) + match(count.i-1,SVcolsums_od) - 1 if(i > 1 & length(nbr_big)) elim[i] <- count.i < max(SVcolsums_od[nbr_big]) if(i < nc & length(nbr_small)) elim[nbr_small[SVcolsums_od[nbr_small] < count.i]] <- TRUE } } elim[SVorder] <- elim !elim
} #' Check Cell Count for degenerate hyperplane arrangements #' #' When the hyperplane arrangement is degenerate, i.e. not in general position, #' the number of distinct cells can be checked against the formula of #' Alexanderson and Wetzel (1981). #' #' @param A is a n by m matrix of hyperplane slope coefficients #' @param b is an n vector of hyperplane intercept coefficients #' @param maxints is maximum number of lines allowed to cross at the same vertex #' #' @return number of distinct cells #' @references #' Alexanderson, G.L and J.E. Wetzel, (1981) Arrangements of planes in space, #' Discrete Math, 34, 219–240. #' @export #' polycount <- function(A, b, maxints=10){
Ab <- cbind(A,b) dup <- which(duplicated(Ab, MARGIN=1)) if (length(dup)>0){ A = A[-dup,,drop=FALSE] b = b[-dup,] } Ab = cbind(A,b) Ab2 = cbind(-A,-b) dup2 <- which(duplicated(rbind(Ab,Ab2),MARGIN=1))-nrow(A) if (length(dup2)>0){ duploc <- (1:(length(dup2)/2))*2 Ab = Ab[-dup2[duploc],,drop=FALSE] } A <- Ab[,1:2,drop=FALSE] b <- Ab[,3] combs <- combn(nrow(A),2) # A here must have redundant rows removed to work vertices <- vapply(1:ncol(combs), FUN=function(j){ if(A[combs[1,j],1]/A[combs[1,j],2]!= A[combs[2,j],1]/A[combs[2,j],2]) solve(A[combs[,j],],b[combs[,j]]) else rep(NA,2)}, FUN.VALUE=numeric(2)) vertices.orig <- vertices vertices <- vertices[,!is.na(colSums(vertices))] vertices <- round(vertices,digits=12) vertices <- vertices[,order(vertices[1,],vertices[2,])] whichdup <- which(duplicated(vertices,MARGIN=2)) diffwhich <- diff(whichdup); diffwhich[diffwhich!=1] <- 0 reps <- vapply(strsplit(paste0(diffwhich,collapse=""),"0")[[1]], FUN=nchar, FUN.VALUE=numeric(1))+2 if(!all(as.numeric(names(table(reps))) %in% choose(3:maxints,2))) stop("Error: inappropriate number of repetitions") counts <- c(NA, NA, sapply(3:maxints, FUN=function(i) sum(reps==choose(i,2)))) counts[2] <- if(maxints==2) ncol(vertices) else ncol(vertices) - sum(choose(3:maxints,2)*counts[3:maxints]) 1 + nrow(A) + sum((1:(maxints-1))*counts[2:maxints])
}