Inputparametters.R + SelecModel
[valse.git] / ProcLassoMLE / constructionModelesLassoMLE.c
1 #include "EMGLLF.h"
2 #include "constructionModelesLassoMLE.h"
3 #include <gsl/gsl_linalg.h>
4 #include <omp.h>
5 #include "omp_num_threads.h"
6
7 // TODO: comment on constructionModelesLassoMLE purpose
8 void constructionModelesLassoMLE(
9 // IN parameters
10 const Real* phiInit, // parametre initial de moyenne renormalisé
11 const Real* rhoInit, // parametre initial de variance renormalisé
12 const Real* piInit, // parametre initial des proportions
13 const Real* gamInit, // paramètre initial des probabilités a posteriori de chaque échantillon
14 Int mini, // nombre minimal d'itérations dans l'algorithme EM
15 Int maxi, // nombre maximal d'itérations dans l'algorithme EM
16 Real gamma, // valeur de gamma : puissance des proportions dans la pénalisation pour un Lasso adaptatif
17 const Real* glambda, // valeur des paramètres de régularisation du Lasso
18 const Real* X, // régresseurs
19 const Real* Y, // réponse
20 Real seuil, // seuil pour prendre en compte une variable
21 Real tau, // seuil pour accepter la convergence
22 const Int* A1, // matrice des coefficients des parametres selectionnes
23 const Int* A2, // matrice des coefficients des parametres non selectionnes
24 // OUT parameters
25 Real* phi, // estimateur ainsi calculé par le Lasso
26 Real* rho, // estimateur ainsi calculé par le Lasso
27 Real* pi, // estimateur ainsi calculé par le Lasso
28 Real* lvraisemblance, // estimateur ainsi calculé par le Lasso
29 // additional size parameters
30 mwSize n, // taille de l'echantillon
31 mwSize p, // nombre de covariables
32 mwSize m, // taille de Y (multivarié)
33 mwSize k, // nombre de composantes
34 mwSize L) // taille de glambda
35 {
36 //preparation: phi = 0
37 for (mwSize u=0; u<p*m*k*L; u++)
38 phi[u] = 0.0;
39
40 //initiate parallel section
41 mwSize lambdaIndex;
42 omp_set_num_threads(OMP_NUM_THREADS);
43 #pragma omp parallel default(shared) private(lambdaIndex)
44 {
45 #pragma omp for schedule(dynamic,CHUNK_SIZE) nowait
46 for (lambdaIndex=0; lambdaIndex<L; lambdaIndex++)
47 {
48 //~ a = A1(:,1,lambdaIndex);
49 //~ a(a==0) = [];
50 Int* a = (Int*)malloc(p*sizeof(Int));
51 mwSize lengthA = 0;
52 for (mwSize j=0; j<p; j++)
53 {
54 if (A1[j*(m+1)*L+0*L+lambdaIndex] != 0)
55 a[lengthA++] = A1[j*(m+1)*L+0*L+lambdaIndex] - 1;
56 }
57 if (lengthA == 0)
58 continue;
59
60 //Xa = X(:,a)
61 Real* Xa = (Real*)malloc(n*lengthA*sizeof(Real));
62 for (mwSize i=0; i<n; i++)
63 {
64 for (mwSize j=0; j<lengthA; j++)
65 Xa[i*lengthA+j] = X[i*p+a[j]];
66 }
67
68 //phia = phiInit(a,:,:)
69 Real* phia = (Real*)malloc(lengthA*m*k*sizeof(Real));
70 for (mwSize j=0; j<lengthA; j++)
71 {
72 for (mwSize mm=0; mm<m; mm++)
73 {
74 for (mwSize r=0; r<k; r++)
75 phia[j*m*k+mm*k+r] = phiInit[a[j]*m*k+mm*k+r];
76 }
77 }
78
79 //[phiLambda,rhoLambda,piLambda,~,~] = EMGLLF(...
80 // phiInit(a,:,:),rhoInit,piInit,gamInit,mini,maxi,gamma,0,X(:,a),Y,tau);
81 Real* phiLambda = (Real*)malloc(lengthA*m*k*sizeof(Real));
82 Real* rhoLambda = (Real*)malloc(m*m*k*sizeof(Real));
83 Real* piLambda = (Real*)malloc(k*sizeof(Real));
84 Real* LLF = (Real*)malloc((maxi+1)*sizeof(Real));
85 Real* S = (Real*)malloc(lengthA*m*k*sizeof(Real));
86 EMGLLF(phia,rhoInit,piInit,gamInit,mini,maxi,gamma,0.0,Xa,Y,tau,
87 phiLambda,rhoLambda,piLambda,LLF,S,
88 n,lengthA,m,k);
89 free(Xa);
90 free(phia);
91 free(LLF);
92 free(S);
93
94 //~ for j=1:length(a)
95 //~ phi(a(j),:,:,lambdaIndex) = phiLambda(j,:,:);
96 //~ end
97 for (mwSize j=0; j<lengthA; j++)
98 {
99 for (mwSize mm=0; mm<m; mm++)
100 {
101 for (mwSize r=0; r<k; r++)
102 phi[a[j]*m*k*L+mm*k*L+r*L+lambdaIndex] = phiLambda[j*m*k+mm*k+r];
103 }
104 }
105 free(phiLambda);
106 //~ rho(:,:,:,lambdaIndex) = rhoLambda;
107 for (mwSize u=0; u<m; u++)
108 {
109 for (mwSize v=0; v<m; v++)
110 {
111 for (mwSize r=0; r<k; r++)
112 rho[u*m*k*L+v*k*L+r*L+lambdaIndex] = rhoLambda[u*m*k+v*k+r];
113 }
114 }
115 free(rhoLambda);
116 //~ pi(:,lambdaIndex) = piLambda;
117 for (mwSize r=0; r<k; r++)
118 pi[r*L+lambdaIndex] = piLambda[r];
119 free(piLambda);
120
121 mwSize dimension = 0;
122 Int* b = (Int*)malloc(m*sizeof(Int));
123 for (mwSize j=0; j<p; j++)
124 {
125 //~ b = A2(j,2:end,lambdaIndex);
126 //~ b(b==0) = [];
127 mwSize lengthB = 0;
128 for (mwSize mm=0; mm<m; mm++)
129 {
130 if (A2[j*(m+1)*L+(mm+1)*L+lambdaIndex] != 0)
131 b[lengthB++] = A2[j*(m+1)*L+(mm+1)*L+lambdaIndex] - 1;
132 }
133 //~ if length(b) > 0
134 //~ phi(A2(j,1,lambdaIndex),b,:,lambdaIndex) = 0.0;
135 //~ end
136 if (lengthB > 0)
137 {
138 for (mwSize mm=0; mm<lengthB; mm++)
139 {
140 for (mwSize r=0; r<k; r++)
141 phi[(A2[j*(m+1)*L+0*L+lambdaIndex]-1)*m*k*L + b[mm]*k*L + r*L + lambdaIndex] = 0.0;
142 }
143 }
144
145 //~ c = A1(j,2:end,lambdaIndex);
146 //~ c(c==0) = [];
147 //~ dimension = dimension + length(c);
148 for (mwSize mm=0; mm<m; mm++)
149 {
150 if (A1[j*(m+1)*L+(mm+1)*L+lambdaIndex] != 0)
151 dimension++;
152 }
153 }
154 free(b);
155
156 int signum;
157 Real* densite = (Real*)calloc(L*n,sizeof(Real));
158 Real sumLogDensit = 0.0;
159 gsl_matrix* matrix = gsl_matrix_alloc(m, m);
160 gsl_permutation* permutation = gsl_permutation_alloc(m);
161 Real* YiRhoR = (Real*)malloc(m*sizeof(Real));
162 Real* XiPhiR = (Real*)malloc(m*sizeof(Real));
163 for (mwSize i=0; i<n; i++)
164 {
165 //~ for r=1:k
166 //~ delta = Y(i,:)*rho(:,:,r,lambdaIndex) - (X(i,a)*(phi(a,:,r,lambdaIndex)));
167 //~ densite(i,lambdaIndex) = densite(i,lambdaIndex) +...
168 //~ pi(r,lambdaIndex)*det(rho(:,:,r,lambdaIndex))/(sqrt(2*PI))^m*exp(-dot(delta,delta)/2.0);
169 //~ end
170 for (mwSize r=0; r<k; r++)
171 {
172 //compute det(rho(:,:,r,lambdaIndex)) [TODO: avoid re-computations]
173 for (mwSize u=0; u<m; u++)
174 {
175 for (mwSize v=0; v<m; v++)
176 matrix->data[u*m+v] = rho[u*m*k*L+v*k*L+r*L+lambdaIndex];
177 }
178 gsl_linalg_LU_decomp(matrix, permutation, &signum);
179 Real detRhoR = gsl_linalg_LU_det(matrix, signum);
180
181 //compute Y(i,:)*rho(:,:,r,lambdaIndex)
182 for (mwSize u=0; u<m; u++)
183 {
184 YiRhoR[u] = 0.0;
185 for (mwSize v=0; v<m; v++)
186 YiRhoR[u] += Y[i*m+v] * rho[v*m*k*L+u*k*L+r*L+lambdaIndex];
187 }
188
189 //compute X(i,a)*phi(a,:,r,lambdaIndex)
190 for (mwSize u=0; u<m; u++)
191 {
192 XiPhiR[u] = 0.0;
193 for (mwSize v=0; v<lengthA; v++)
194 XiPhiR[u] += X[i*p+a[v]] * phi[a[v]*m*k*L+u*k*L+r*L+lambdaIndex];
195 }
196 // On peut remplacer X par Xa dans ce dernier calcul, mais je ne sais pas si c'est intéressant ...
197
198 // compute dotProduct < delta . delta >
199 Real dotProduct = 0.0;
200 for (mwSize u=0; u<m; u++)
201 dotProduct += (YiRhoR[u]-XiPhiR[u]) * (YiRhoR[u]-XiPhiR[u]);
202
203 densite[lambdaIndex*n+i] += (pi[r*L+lambdaIndex]*detRhoR/pow(sqrt(2.0*M_PI),m))*exp(-dotProduct/2.0);
204 }
205 sumLogDensit += log(densite[lambdaIndex*n+i]);
206 }
207 lvraisemblance[lambdaIndex*2+0] = sumLogDensit;
208 lvraisemblance[lambdaIndex*2+1] = (dimension+m+1)*k-1;
209
210 free(a);
211 free(YiRhoR);
212 free(XiPhiR);
213 free(densite);
214 gsl_matrix_free(matrix);
215 gsl_permutation_free(permutation);
216 }
217 }
218 }