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