Fix errors in EMGLLF.c :)
[valse.git] / pkg / src / sources / EMGLLF.c
1 #include "utils.h"
2 #include <stdlib.h>
3 #include <gsl/gsl_linalg.h>
4
5 // TODO: don't recompute indexes ai(...) and mi(...) when possible
6 void EMGLLF_core(
7 // IN parameters
8 const Real* phiInit, // parametre initial de moyenne renormalisé
9 const Real* rhoInit, // parametre initial de variance renormalisé
10 const Real* piInit, // parametre initial des proportions
11 const Real* gamInit, // paramètre initial des probabilités a posteriori de chaque échantillon
12 int mini, // nombre minimal d'itérations dans l'algorithme EM
13 int maxi, // nombre maximal d'itérations dans l'algorithme EM
14 Real gamma, // puissance des proportions dans la pénalisation pour un Lasso adaptatif
15 Real lambda, // valeur du paramètre de régularisation du Lasso
16 const Real* X, // régresseurs
17 const Real* Y, // réponse
18 Real tau, // seuil pour accepter la convergence
19 // OUT parameters (all pointers, to be modified)
20 Real* phi, // parametre de moyenne renormalisé, calculé par l'EM
21 Real* rho, // parametre de variance renormalisé, calculé par l'EM
22 Real* pi, // parametre des proportions renormalisé, calculé par l'EM
23 Real* LLF, // log vraisemblance associée à cet échantillon, pour les valeurs estimées des paramètres
24 Real* S,
25 // additional size parameters
26 int n, // nombre d'echantillons
27 int p, // nombre de covariables
28 int m, // taille de Y (multivarié)
29 int k) // nombre de composantes dans le mélange
30 {
31 //Initialize outputs
32 copyArray(phiInit, phi, p*m*k);
33 copyArray(rhoInit, rho, m*m*k);
34 copyArray(piInit, pi, k);
35 zeroArray(LLF, maxi);
36 //S is already allocated, and doesn't need to be 'zeroed'
37
38 //Other local variables: same as in R
39 Real* gam = (Real*)malloc(n*k*sizeof(Real));
40 copyArray(gamInit, gam, n*k);
41 Real* Gram2 = (Real*)malloc(p*p*k*sizeof(Real));
42 Real* ps2 = (Real*)malloc(p*m*k*sizeof(Real));
43 Real* b = (Real*)malloc(k*sizeof(Real));
44 Real* X2 = (Real*)malloc(n*p*k*sizeof(Real));
45 Real* Y2 = (Real*)malloc(n*m*k*sizeof(Real));
46 Real dist = 0.;
47 Real dist2 = 0.;
48 int ite = 0;
49 Real* pi2 = (Real*)malloc(k*sizeof(Real));
50 Real* ps = (Real*)malloc(m*k*sizeof(Real));
51 Real* nY2 = (Real*)malloc(m*k*sizeof(Real));
52 Real* ps1 = (Real*)malloc(n*m*k*sizeof(Real));
53 Real* Gam = (Real*)malloc(n*k*sizeof(Real));
54 const Real EPS = 1e-15;
55 // Additional (not at this place, in R file)
56 Real* gam2 = (Real*)malloc(k*sizeof(Real));
57 Real* sqNorm2 = (Real*)malloc(k*sizeof(Real));
58 gsl_matrix* matrix = gsl_matrix_alloc(m, m);
59 gsl_permutation* permutation = gsl_permutation_alloc(m);
60 Real* YiRhoR = (Real*)malloc(m*sizeof(Real));
61 Real* XiPhiR = (Real*)malloc(m*sizeof(Real));
62 const Real gaussConstM = pow(2.*M_PI,m/2.);
63 Real* Phi = (Real*)malloc(p*m*k*sizeof(Real));
64 Real* Rho = (Real*)malloc(m*m*k*sizeof(Real));
65 Real* Pi = (Real*)malloc(k*sizeof(Real));
66
67 while (ite < mini || (ite < maxi && (dist >= tau || dist2 >= sqrt(tau))))
68 {
69 copyArray(phi, Phi, p*m*k);
70 copyArray(rho, Rho, m*m*k);
71 copyArray(pi, Pi, k);
72
73 // Calculs associés a Y et X
74 for (int r=0; r<k; r++)
75 {
76 for (int mm=0; mm<m; mm++)
77 {
78 //Y2[,mm,r] = sqrt(gam[,r]) * Y[,mm]
79 for (int u=0; u<n; u++)
80 Y2[ai(u,mm,r,n,m,k)] = sqrt(gam[mi(u,r,n,k)]) * Y[mi(u,mm,n,m)];
81 }
82 for (int i=0; i<n; i++)
83 {
84 //X2[i,,r] = sqrt(gam[i,r]) * X[i,]
85 for (int u=0; u<p; u++)
86 X2[ai(i,u,r,n,p,k)] = sqrt(gam[mi(i,r,n,k)]) * X[mi(i,u,n,p)];
87 }
88 for (int mm=0; mm<m; mm++)
89 {
90 //ps2[,mm,r] = crossprod(X2[,,r],Y2[,mm,r])
91 for (int u=0; u<p; u++)
92 {
93 Real dotProduct = 0.;
94 for (int v=0; v<n; v++)
95 dotProduct += X2[ai(v,u,r,n,p,k)] * Y2[ai(v,mm,r,n,m,k)];
96 ps2[ai(u,mm,r,p,m,k)] = dotProduct;
97 }
98 }
99 for (int j=0; j<p; j++)
100 {
101 for (int s=0; s<p; s++)
102 {
103 //Gram2[j,s,r] = crossprod(X2[,j,r], X2[,s,r])
104 Real dotProduct = 0.;
105 for (int u=0; u<n; u++)
106 dotProduct += X2[ai(u,j,r,n,p,k)] * X2[ai(u,s,r,n,p,k)];
107 Gram2[ai(j,s,r,p,p,k)] = dotProduct;
108 }
109 }
110 }
111
112 /////////////
113 // Etape M //
114 /////////////
115
116 // Pour pi
117 for (int r=0; r<k; r++)
118 {
119 //b[r] = sum(abs(phi[,,r]))
120 Real sumAbsPhi = 0.;
121 for (int u=0; u<p; u++)
122 for (int v=0; v<m; v++)
123 sumAbsPhi += fabs(phi[ai(u,v,r,p,m,k)]);
124 b[r] = sumAbsPhi;
125 }
126 //gam2 = colSums(gam)
127 for (int u=0; u<k; u++)
128 {
129 Real sumOnColumn = 0.;
130 for (int v=0; v<n; v++)
131 sumOnColumn += gam[mi(v,u,n,k)];
132 gam2[u] = sumOnColumn;
133 }
134 //a = sum(gam %*% log(pi))
135 Real a = 0.;
136 for (int u=0; u<n; u++)
137 {
138 Real dotProduct = 0.;
139 for (int v=0; v<k; v++)
140 dotProduct += gam[mi(u,v,n,k)] * log(pi[v]);
141 a += dotProduct;
142 }
143
144 //tant que les proportions sont negatives
145 int kk = 0;
146 int pi2AllPositive = 0;
147 Real invN = 1./n;
148 while (!pi2AllPositive)
149 {
150 //pi2 = pi + 0.1^kk * ((1/n)*gam2 - pi)
151 Real pow_01_kk = pow(0.1,kk);
152 for (int r=0; r<k; r++)
153 pi2[r] = pi[r] + pow_01_kk * (invN*gam2[r] - pi[r]);
154 //pi2AllPositive = all(pi2 >= 0)
155 pi2AllPositive = 1;
156 for (int r=0; r<k; r++)
157 {
158 if (pi2[r] < 0)
159 {
160 pi2AllPositive = 0;
161 break;
162 }
163 }
164 kk++;
165 }
166
167 //sum(pi^gamma * b)
168 Real piPowGammaDotB = 0.;
169 for (int v=0; v<k; v++)
170 piPowGammaDotB += pow(pi[v],gamma) * b[v];
171 //sum(pi2^gamma * b)
172 Real pi2PowGammaDotB = 0.;
173 for (int v=0; v<k; v++)
174 pi2PowGammaDotB += pow(pi2[v],gamma) * b[v];
175 //sum(gam2 * log(pi2))
176 Real gam2DotLogPi2 = 0.;
177 for (int v=0; v<k; v++)
178 gam2DotLogPi2 += gam2[v] * log(pi2[v]);
179
180 //t(m) la plus grande valeur dans la grille O.1^k tel que ce soit décroissante ou constante
181 while (-invN*a + lambda*piPowGammaDotB < -invN*gam2DotLogPi2 + lambda*pi2PowGammaDotB
182 && kk<1000)
183 {
184 Real pow_01_kk = pow(0.1,kk);
185 //pi2 = pi + 0.1^kk * (1/n*gam2 - pi)
186 for (int v=0; v<k; v++)
187 pi2[v] = pi[v] + pow_01_kk * (invN*gam2[v] - pi[v]);
188 //pi2 was updated, so we recompute pi2PowGammaDotB and gam2DotLogPi2
189 pi2PowGammaDotB = 0.;
190 for (int v=0; v<k; v++)
191 pi2PowGammaDotB += pow(pi2[v],gamma) * b[v];
192 gam2DotLogPi2 = 0.;
193 for (int v=0; v<k; v++)
194 gam2DotLogPi2 += gam2[v] * log(pi2[v]);
195 kk++;
196 }
197 Real t = pow(0.1,kk);
198 //sum(pi + t*(pi2-pi))
199 Real sumPiPlusTbyDiff = 0.;
200 for (int v=0; v<k; v++)
201 sumPiPlusTbyDiff += (pi[v] + t*(pi2[v] - pi[v]));
202 //pi = (pi + t*(pi2-pi)) / sum(pi + t*(pi2-pi))
203 for (int v=0; v<k; v++)
204 pi[v] = (pi[v] + t*(pi2[v] - pi[v])) / sumPiPlusTbyDiff;
205
206 //Pour phi et rho
207 for (int r=0; r<k; r++)
208 {
209 for (int mm=0; mm<m; mm++)
210 {
211 for (int i=0; i<n; i++)
212 {
213 //< X2[i,,r] , phi[,mm,r] >
214 Real dotProduct = 0.;
215 for (int u=0; u<p; u++)
216 dotProduct += X2[ai(i,u,r,n,p,k)] * phi[ai(u,mm,r,p,m,k)];
217 //ps1[i,mm,r] = Y2[i,mm,r] * sum(X2[i,,r] * phi[,mm,r])
218 ps1[ai(i,mm,r,n,m,k)] = Y2[ai(i,mm,r,n,m,k)] * dotProduct;
219 }
220 //ps[mm,r] = sum(ps1[,mm,r])
221 Real sumPs1 = 0.;
222 for (int u=0; u<n; u++)
223 sumPs1 += ps1[ai(u,mm,r,n,m,k)];
224 ps[mi(mm,r,m,k)] = sumPs1;
225 //nY2[mm,r] = sum(Y2[,mm,r]^2)
226 Real sumY2 = 0.;
227 for (int u=0; u<n; u++)
228 sumY2 += Y2[ai(u,mm,r,n,m,k)] * Y2[ai(u,mm,r,n,m,k)];
229 nY2[mi(mm,r,m,k)] = sumY2;
230 //rho[mm,mm,r] = (ps[mm,r]+sqrt(ps[mm,r]^2+4*nY2[mm,r]*(gam2[r]))) / (2*nY2[mm,r])
231 rho[ai(mm,mm,r,m,m,k)] = ( ps[mi(mm,r,m,k)] + sqrt( ps[mi(mm,r,m,k)]*ps[mi(mm,r,m,k)]
232 + 4*nY2[mi(mm,r,m,k)] * gam2[r] ) ) / (2*nY2[mi(mm,r,m,k)]);
233 }
234 }
235
236 for (int r=0; r<k; r++)
237 {
238 for (int j=0; j<p; j++)
239 {
240 for (int mm=0; mm<m; mm++)
241 {
242 //sum(phi[-j,mm,r] * Gram2[j, setdiff(1:p,j),r])
243 Real phiDotGram2 = 0.;
244 for (int u=0; u<p; u++)
245 {
246 if (u != j)
247 phiDotGram2 += phi[ai(u,mm,r,p,m,k)] * Gram2[ai(j,u,r,p,p,k)];
248 }
249 //S[j,mm,r] = -rho[mm,mm,r]*ps2[j,mm,r] + sum(phi[-j,mm,r] * Gram2[j,-j,r])
250 S[ai(j,mm,r,p,m,k)] = -rho[ai(mm,mm,r,m,m,k)] * ps2[ai(j,mm,r,p,m,k)] + phiDotGram2;
251 Real pirPowGamma = pow(pi[r],gamma);
252 if (fabs(S[ai(j,mm,r,p,m,k)]) <= n*lambda*pirPowGamma)
253 phi[ai(j,mm,r,p,m,k)] = 0.;
254 else if (S[ai(j,mm,r,p,m,k)] > n*lambda*pirPowGamma)
255 {
256 phi[ai(j,mm,r,p,m,k)] = (n*lambda*pirPowGamma - S[ai(j,mm,r,p,m,k)])
257 / Gram2[ai(j,j,r,p,p,k)];
258 }
259 else
260 {
261 phi[ai(j,mm,r,p,m,k)] = -(n*lambda*pirPowGamma + S[ai(j,mm,r,p,m,k)])
262 / Gram2[ai(j,j,r,p,p,k)];
263 }
264 }
265 }
266 }
267
268 /////////////
269 // Etape E //
270 /////////////
271
272 int signum;
273 Real sumLogLLF2 = 0.;
274 for (int i=0; i<n; i++)
275 {
276 for (int r=0; r<k; r++)
277 {
278 //compute Y[i,]%*%rho[,,r]
279 for (int u=0; u<m; u++)
280 {
281 YiRhoR[u] = 0.;
282 for (int v=0; v<m; v++)
283 YiRhoR[u] += Y[mi(i,v,n,m)] * rho[ai(v,u,r,m,m,k)];
284 }
285
286 //compute X[i,]%*%phi[,,r]
287 for (int u=0; u<m; u++)
288 {
289 XiPhiR[u] = 0.;
290 for (int v=0; v<p; v++)
291 XiPhiR[u] += X[mi(i,v,n,p)] * phi[ai(v,u,r,p,m,k)];
292 }
293
294 //compute sq norm || Y(:,i)*rho(:,:,r)-X(i,:)*phi(:,:,r) ||_2^2
295 sqNorm2[r] = 0.;
296 for (int u=0; u<m; u++)
297 sqNorm2[r] += (YiRhoR[u]-XiPhiR[u]) * (YiRhoR[u]-XiPhiR[u]);
298 }
299
300 Real sumLLF1 = 0.;
301 Real sumGamI = 0.;
302 for (int r=0; r<k; r++)
303 {
304 //compute det(rho[,,r]) [TODO: avoid re-computations]
305 for (int u=0; u<m; u++)
306 {
307 for (int v=0; v<m; v++)
308 matrix->data[u*m+v] = rho[ai(u,v,r,m,m,k)];
309 }
310 gsl_linalg_LU_decomp(matrix, permutation, &signum);
311 Real detRhoR = gsl_linalg_LU_det(matrix, signum);
312 Gam[mi(i,r,n,k)] = pi[r] * exp(-.5*sqNorm2[r]) * detRhoR;
313 sumLLF1 += Gam[mi(i,r,n,k)] / gaussConstM;
314 sumGamI += Gam[mi(i,r,n,k)];
315 }
316
317 sumLogLLF2 += log(sumLLF1);
318 for (int r=0; r<k; r++)
319 {
320 //gam[i,] = Gam[i,] / sumGamI
321 gam[mi(i,r,n,k)] = sumGamI > EPS ? Gam[mi(i,r,n,k)] / sumGamI : 0.;
322 }
323 }
324
325 //sumPen = sum(pi^gamma * b)
326 Real sumPen = 0.;
327 for (int r=0; r<k; r++)
328 sumPen += pow(pi[r],gamma) * b[r];
329 //LLF[ite] = -sumLogLLF2/n + lambda*sumPen
330 LLF[ite] = -invN * sumLogLLF2 + lambda * sumPen;
331 dist = ite==0 ? LLF[ite] : (LLF[ite] - LLF[ite-1]) / (1. + fabs(LLF[ite]));
332
333 //Dist1 = max( abs(phi-Phi) / (1+abs(phi)) )
334 Real Dist1 = 0.;
335 for (int u=0; u<p; u++)
336 {
337 for (int v=0; v<m; v++)
338 {
339 for (int w=0; w<k; w++)
340 {
341 Real tmpDist = fabs(phi[ai(u,v,w,p,m,k)]-Phi[ai(u,v,w,p,m,k)])
342 / (1.+fabs(phi[ai(u,v,w,p,m,k)]));
343 if (tmpDist > Dist1)
344 Dist1 = tmpDist;
345 }
346 }
347 }
348 //Dist2 = max( (abs(rho-Rho)) / (1+abs(rho)) )
349 Real Dist2 = 0.;
350 for (int u=0; u<m; u++)
351 {
352 for (int v=0; v<m; v++)
353 {
354 for (int w=0; w<k; w++)
355 {
356 Real tmpDist = fabs(rho[ai(u,v,w,m,m,k)]-Rho[ai(u,v,w,m,m,k)])
357 / (1.+fabs(rho[ai(u,v,w,m,m,k)]));
358 if (tmpDist > Dist2)
359 Dist2 = tmpDist;
360 }
361 }
362 }
363 //Dist3 = max( (abs(pi-Pi)) / (1+abs(Pi)))
364 Real Dist3 = 0.;
365 for (int u=0; u<n; u++)
366 {
367 for (int v=0; v<k; v++)
368 {
369 Real tmpDist = fabs(pi[v]-Pi[v]) / (1.+fabs(pi[v]));
370 if (tmpDist > Dist3)
371 Dist3 = tmpDist;
372 }
373 }
374 //dist2=max([max(Dist1),max(Dist2),max(Dist3)]);
375 dist2 = Dist1;
376 if (Dist2 > dist2)
377 dist2 = Dist2;
378 if (Dist3 > dist2)
379 dist2 = Dist3;
380
381 ite++;
382 }
383
384 //TODO: affec = apply(gam, 1,which.max) à traduire, et retourner affec aussi.
385
386 //free memory
387 free(b);
388 free(gam);
389 free(Gam);
390 free(Phi);
391 free(Rho);
392 free(Pi);
393 free(ps);
394 free(nY2);
395 free(ps1);
396 free(Gram2);
397 free(ps2);
398 gsl_matrix_free(matrix);
399 gsl_permutation_free(permutation);
400 free(XiPhiR);
401 free(YiRhoR);
402 free(gam2);
403 free(pi2);
404 free(X2);
405 free(Y2);
406 free(sqNorm2);
407 }