a45f71ac4f64769f7da037324efef533ffa58ba0
[morpheus.git] / pkg / R / optimParams.R
1 #' Wrapper function for OptimParams class
2 #'
3 #' @param K Number of populations.
4 #' @param link The link type, 'logit' or 'probit'.
5 #' @param X Data matrix of covariables
6 #' @param Y Output as a binary vector
7 #'
8 #' @return An object 'op' of class OptimParams, initialized so that \code{op$run(x0)}
9 #' outputs the list of optimized parameters
10 #' \itemize{
11 #' \item p: proportions, size K
12 #' \item β: regression matrix, size dxK
13 #' \item b: intercepts, size K
14 #' }
15 #' θ0 is a vector containing respectively the K-1 first elements of p, then β by
16 #' columns, and finally b: \code{θ0 = c(p[1:(K-1)],as.double(β),b)}.
17 #'
18 #' @seealso \code{multiRun} to estimate statistics based on β, and
19 #' \code{generateSampleIO} for I/O random generation.
20 #'
21 #' @examples
22 #' # Optimize parameters from estimated μ
23 #' io = generateSampleIO(10000, 1/2, matrix(c(1,-2,3,1),ncol=2), c(0,0), "logit")
24 #' μ = computeMu(io$X, io$Y, list(K=2))
25 #' o <- optimParams(io$X, io$Y, 2, "logit")
26 #' θ0 <- list(p=1/2, β=μ, b=c(0,0))
27 #' par0 <- o$run(θ0)
28 #' # Compare with another starting point
29 #' θ1 <- list(p=1/2, β=2*μ, b=c(0,0))
30 #' par1 <- o$run(θ1)
31 #' o$f( o$linArgs(par0) )
32 #' o$f( o$linArgs(par1) )
33 #' @export
34 optimParams <- function(X, Y, K, link=c("logit","probit"))
35 {
36 # Check arguments
37 if (!is.matrix(X) || any(is.na(X)))
38 stop("X: numeric matrix, no NAs")
39 if (!is.numeric(Y) || any(is.na(Y)) || any(Y!=0 & Y!=1))
40 stop("Y: binary vector with 0 and 1 only")
41 link <- match.arg(link)
42 if (!is.numeric(K) || K!=floor(K) || K < 2)
43 stop("K: integer >= 2")
44
45 # Build and return optimization algorithm object
46 methods::new("OptimParams", "li"=link, "X"=X,
47 "Y"=as.integer(Y), "K"=as.integer(K))
48 }
49
50 #' Encapsulated optimization for p (proportions), β and b (regression parameters)
51 #'
52 #' Optimize the parameters of a mixture of logistic regressions model, possibly using
53 #' \code{mu <- computeMu(...)} as a partial starting point.
54 #'
55 #' @field li Link function, 'logit' or 'probit'
56 #' @field X Data matrix of covariables
57 #' @field Y Output as a binary vector
58 #' @field K Number of populations
59 #' @field d Number of dimensions
60 #' @field W Weights matrix (iteratively refined)
61 #'
62 setRefClass(
63 Class = "OptimParams",
64
65 fields = list(
66 # Inputs
67 li = "character", #link function
68 X = "matrix",
69 Y = "numeric",
70 Mhat = "numeric", #vector of empirical moments
71 # Dimensions
72 K = "integer",
73 n = "integer",
74 d = "integer",
75 # Weights matrix (generalized least square)
76 W = "matrix"
77 ),
78
79 methods = list(
80 initialize = function(...)
81 {
82 "Check args and initialize K, d, W"
83
84 callSuper(...)
85 if (!hasArg("X") || !hasArg("Y") || !hasArg("K") || !hasArg("li"))
86 stop("Missing arguments")
87
88 # Precompute empirical moments
89 M <- computeMoments(X, Y)
90 M1 <- as.double(M[[1]])
91 M2 <- as.double(M[[2]])
92 M3 <- as.double(M[[3]])
93 Mhat <<- c(M1, M2, M3)
94
95 n <<- nrow(X)
96 d <<- length(M1)
97 W <<- diag(d+d^2+d^3) #initialize at W = Identity
98 },
99
100 expArgs = function(v)
101 {
102 "Expand individual arguments from vector v into a list"
103
104 list(
105 # p: dimension K-1, need to be completed
106 "p" = c(v[1:(K-1)], 1-sum(v[1:(K-1)])),
107 "β" = matrix(v[K:(K+d*K-1)], ncol=K),
108 "b" = v[(K+d*K):(K+(d+1)*K-1)])
109 },
110
111 linArgs = function(L)
112 {
113 "Linearize vectors+matrices from list L into a vector"
114
115 c(L$p[1:(K-1)], as.double(L$β), L$b)
116 },
117
118 computeW = function(θ)
119 {
120 #require(MASS)
121 dd <- d + d^2 + d^3
122 W <<- MASS::ginv( matrix( .C("Compute_Omega",
123 X=as.double(X), Y=Y, M=Moments(θ), pn=as.integer(n), pd=as.integer(d),
124 W=as.double(W), PACKAGE="morpheus")$W, nrow=dd, ncol=dd ) )
125 NULL #avoid returning W
126 },
127
128 Moments = function(θ)
129 {
130 "Vector of moments, of size d+d^2+d^3"
131
132 p <- θ$p
133 β <- θ$β
134 λ <- sqrt(colSums(β^2))
135 b <- θ$b
136
137 # Tensorial products β^2 = β2 and β^3 = β3 must be computed from current β1
138 β2 <- apply(β, 2, function(col) col %o% col)
139 β3 <- apply(β, 2, function(col) col %o% col %o% col)
140
141 c(
142 β %*% (p * .G(li,1,λ,b)),
143 β2 %*% (p * .G(li,2,λ,b)),
144 β3 %*% (p * .G(li,3,λ,b)))
145 },
146
147 f = function(θ)
148 {
149 "Product t(hat_Mi - Mi) W (hat_Mi - Mi) with Mi(theta)"
150
151 L <- expArgs(θ)
152 A <- as.matrix(Mhat - Moments(L))
153 t(A) %*% W %*% A
154 },
155
156 grad_f = function(θ)
157 {
158 "Gradient of f, dimension (K-1) + d*K + K = (d+2)*K - 1"
159
160 L <- expArgs(θ)
161 -2 * t(grad_M(L)) %*% W %*% as.matrix((Mhat - Moments(L)))
162 },
163
164 grad_M = function(θ)
165 {
166 "Gradient of the vector of moments, size (dim=)d+d^2+d^3 x K-1+K+d*K"
167
168 p <- θ$p
169 β <- θ$β
170 λ <- sqrt(colSums(β^2))
171 μ <- sweep(β, 2, λ, '/')
172 b <- θ$b
173
174 res <- matrix(nrow=nrow(W), ncol=0)
175
176 # Tensorial products β^2 = β2 and β^3 = β3 must be computed from current β1
177 β2 <- apply(β, 2, function(col) col %o% col)
178 β3 <- apply(β, 2, function(col) col %o% col %o% col)
179
180 # Some precomputations
181 G1 = .G(li,1,λ,b)
182 G2 = .G(li,2,λ,b)
183 G3 = .G(li,3,λ,b)
184 G4 = .G(li,4,λ,b)
185 G5 = .G(li,5,λ,b)
186
187 # Gradient on p: K-1 columns, dim rows
188 km1 = 1:(K-1)
189 res <- cbind(res, rbind(
190 sweep(as.matrix(β [,km1]), 2, G1[km1], '*') - G1[K] * β [,K],
191 sweep(as.matrix(β2[,km1]), 2, G2[km1], '*') - G2[K] * β2[,K],
192 sweep(as.matrix(β3[,km1]), 2, G3[km1], '*') - G3[K] * β3[,K] ))
193
194 for (i in 1:d)
195 {
196 # i determines the derivated matrix dβ[2,3]
197
198 dβ_left <- sweep(β, 2, p * G3 * β[i,], '*')
199 dβ_right <- matrix(0, nrow=d, ncol=K)
200 block <- i
201 dβ_right[block,] <- dβ_right[block,] + 1
202 dβ <- dβ_left + sweep(dβ_right, 2, p * G1, '*')
203
204 dβ2_left <- sweep(β2, 2, p * G4 * β[i,], '*')
205 dβ2_right <- do.call( rbind, lapply(1:d, function(j) {
206 sweep(dβ_right, 2, β[j,], '*')
207 }) )
208 block <- ((i-1)*d+1):(i*d)
209 dβ2_right[block,] <- dβ2_right[block,] + β
210 dβ2 <- dβ2_left + sweep(dβ2_right, 2, p * G2, '*')
211
212 dβ3_left <- sweep(β3, 2, p * G5 * β[i,], '*')
213 dβ3_right <- do.call( rbind, lapply(1:d, function(j) {
214 sweep(dβ2_right, 2, β[j,], '*')
215 }) )
216 block <- ((i-1)*d*d+1):(i*d*d)
217 dβ3_right[block,] <- dβ3_right[block,] + β2
218 dβ3 <- dβ3_left + sweep(dβ3_right, 2, p * G3, '*')
219
220 res <- cbind(res, rbind(dβ, dβ2, dβ3))
221 }
222
223 # Gradient on b
224 res <- cbind(res, rbind(
225 sweep(β, 2, p * G2, '*'),
226 sweep(β2, 2, p * G3, '*'),
227 sweep(β3, 2, p * G4, '*') ))
228
229 res
230 },
231
232 run = function(θ0)
233 {
234 "Run optimization from θ0 with solver..."
235
236 if (!is.list(θ0))
237 stop("θ0: list")
238 if (is.null(θ0$β))
239 stop("At least θ0$β must be provided")
240 if (!is.matrix(θ0$β) || any(is.na(θ0$β)) || ncol(θ0$β) != K)
241 stop("θ0$β: matrix, no NA, ncol == K")
242 if (is.null(θ0$p))
243 θ0$p = rep(1/K, K-1)
244 else if (length(θ0$p) != K-1 || sum(θ0$p) > 1)
245 stop("θ0$p should contain positive integers and sum to < 1")
246 # Next test = heuristic to detect missing b (when matrix is called "beta")
247 if (is.null(θ0$b) || all(θ0$b == θ0$β))
248 θ0$b = rep(0, K)
249 else if (any(is.na(θ0$b)))
250 stop("θ0$b cannot have missing values")
251
252 # TODO: stopping condition? N iterations? Delta <= epsilon ?
253 for (loop in 1:10)
254 {
255 op_res = constrOptim( linArgs(θ0), .self$f, .self$grad_f,
256 ui=cbind(
257 rbind( rep(-1,K-1), diag(K-1) ),
258 matrix(0, nrow=K, ncol=(d+1)*K) ),
259 ci=c(-1,rep(0,K-1)) )
260
261 computeW(expArgs(op_res$par))
262 # debug:
263 #print(W)
264 print(op_res$value)
265 print(expArgs(op_res$par))
266 }
267
268 expArgs(op_res$par)
269 }
270 )
271 )
272
273 # Compute vectorial E[g^{(order)}(<β,x> + b)] with x~N(0,Id) (integral in R^d)
274 # = E[g^{(order)}(z)] with z~N(b,diag(λ))
275 # by numerically evaluating the integral.
276 #
277 # @param link Link, 'logit' or 'probit'
278 # @param order Order of derivative
279 # @param λ Norm of columns of β
280 # @param b Intercept
281 #
282 .G <- function(link, order, λ, b)
283 {
284 # NOTE: weird "integral divergent" error on inputs:
285 # link="probit"; order=2; λ=c(531.8099,586.8893,523.5816); b=c(-118.512674,-3.488020,2.109969)
286 # Switch to pracma package for that (but it seems slow...)
287 sapply( seq_along(λ), function(k) {
288 res <- NULL
289 tryCatch({
290 # Fast code, may fail:
291 res <- stats::integrate(
292 function(z) .deriv[[link]][[order]](λ[k]*z+b[k]) * exp(-z^2/2) / sqrt(2*pi),
293 lower=-Inf, upper=Inf )$value
294 }, error = function(e) {
295 # Robust slow code, no fails observed:
296 sink("/dev/null") #pracma package has some useless printed outputs...
297 res <- pracma::integral(
298 function(z) .deriv[[link]][[order]](λ[k]*z+b[k]) * exp(-z^2/2) / sqrt(2*pi),
299 xmin=-Inf, xmax=Inf, method="Kronrod")
300 sink()
301 })
302 res
303 })
304 }
305
306 # Derivatives list: g^(k)(x) for links 'logit' and 'probit'
307 #
308 .deriv <- list(
309 "probit"=list(
310 # 'probit' derivatives list;
311 # NOTE: exact values for the integral E[g^(k)(λz+b)] could be computed
312 function(x) exp(-x^2/2)/(sqrt(2*pi)), #g'
313 function(x) exp(-x^2/2)/(sqrt(2*pi)) * -x, #g''
314 function(x) exp(-x^2/2)/(sqrt(2*pi)) * ( x^2 - 1), #g^(3)
315 function(x) exp(-x^2/2)/(sqrt(2*pi)) * (-x^3 + 3*x), #g^(4)
316 function(x) exp(-x^2/2)/(sqrt(2*pi)) * ( x^4 - 6*x^2 + 3) #g^(5)
317 ),
318 "logit"=list(
319 # Sigmoid derivatives list, obtained with http://www.derivative-calculator.net/
320 # @seealso http://www.ece.uc.edu/~aminai/papers/minai_sigmoids_NN93.pdf
321 function(x) {e=exp(x); .zin(e /(e+1)^2)}, #g'
322 function(x) {e=exp(x); .zin(e*(-e + 1) /(e+1)^3)}, #g''
323 function(x) {e=exp(x); .zin(e*( e^2 - 4*e + 1) /(e+1)^4)}, #g^(3)
324 function(x) {e=exp(x); .zin(e*(-e^3 + 11*e^2 - 11*e + 1) /(e+1)^5)}, #g^(4)
325 function(x) {e=exp(x); .zin(e*( e^4 - 26*e^3 + 66*e^2 - 26*e + 1)/(e+1)^6)} #g^(5)
326 )
327 )
328
329 # Utility for integration: "[return] zero if [argument is] NaN" (Inf / Inf divs)
330 #
331 # @param x Ratio of polynoms of exponentials, as in .S[[i]]
332 #
333 .zin <- function(x)
334 {
335 x[is.nan(x)] <- 0.
336 x
337 }