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