59be2f734c64697ae045c79ce15af23ebaf16d7f
[valse.git] / src / sources / constructionModelesLassoRank.c
1 #include <stdlib.h>
2 #include <omp.h>
3 #include <gsl/gsl_linalg.h>
4 #include "EMGrank.h"
5 #include "utils.h"
6
7 // TODO: comment on constructionModelesLassoRank purpose
8 void constructionModelesLassoRank_core(
9 // IN parameters
10 const Real* pi,// parametre initial des proportions
11 const Real* rho, // parametre initial de variance renormalisé
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 const Real* X,// régresseurs
15 const Real* Y,// réponse
16 Real tau, // seuil pour accepter la convergence
17 const int* A1, // matrice des coefficients des parametres selectionnes
18 int rangmin, //rang minimum autorisé
19 int rangmax, //rang maximum autorisé
20 // OUT parameters (all pointers, to be modified)
21 Real* phi,// estimateur ainsi calculé par le Lasso
22 Real* llh,// estimateur ainsi calculé par le Lasso
23 // additional size parameters
24 int n,// taille de l'echantillon
25 int p,// nombre de covariables
26 int m,// taille de Y (multivarié)
27 int k,// nombre de composantes
28 int L)// taille de glambda
29 {
30 //On cherche les rangs possiblement intéressants
31 int deltaRank = rangmax-rangmin+1;
32 int Size = (int)pow(deltaRank,k);
33 int* Rank = (int*)malloc(Size*k*sizeof(int));
34 for (int r=0; r<k; r++)
35 {
36 // On veut le tableau de toutes les combinaisons de rangs possibles
37 // Dans la première colonne : on répète (rangmax-rangmin)^(k-1) chaque chiffre :
38 // ça remplit la colonne
39 // Dans la deuxieme : on répète (rangmax-rangmin)^(k-2) chaque chiffre,
40 // et on fait ça (rangmax-rangmin)^2 fois
41 // ...
42 // Dans la dernière, on répète chaque chiffre une fois,
43 // et on fait ça (rangmin-rangmax)^(k-1) fois.
44 int indexInRank = 0;
45 int value = 0;
46 while (indexInRank < Size)
47 {
48 int upperBoundU = (int)pow(deltaRank,k-r-1);
49 for (int u=0; u<upperBoundU; u++)
50 Rank[mi(indexInRank++,r,Size,k)] = rangmin + value;
51 value = (value+1) % deltaRank;
52 }
53 }
54
55 //Initialize phi to zero, because unactive variables won't be assigned
56 for (int i=0; i<p*m*k*L*Size; i++)
57 phi[i] = 0.0;
58
59 //initiate parallel section
60 int lambdaIndex;
61 omp_set_num_threads(OMP_NUM_THREADS);
62 #pragma omp parallel default(shared) private(lambdaIndex)
63 {
64 #pragma omp for schedule(dynamic,CHUNK_SIZE) nowait
65 for (lambdaIndex=0; lambdaIndex<L; lambdaIndex++)
66 {
67 //On ne garde que les colonnes actives : active sera l'ensemble des variables informatives
68 int* active = (int*)malloc(p*sizeof(int));
69 int longueurActive = 0;
70 for (int j=0; j<p; j++)
71 {
72 if (A1[mi(j,lambdaIndex,p,L)] != 0)
73 active[longueurActive++] = A1[mi(j,lambdaIndex,p,L)] - 1;
74 }
75 if (longueurActive == 0)
76 continue;
77
78 //from now on, longueurActive > 0
79 Real* phiLambda = (Real*)malloc(longueurActive*m*k*sizeof(Real));
80 Real LLF;
81 for (int j=0; j<Size; j++)
82 {
83 int* rank = (int*)malloc(k*sizeof(int));
84 for (int r=0; r<k; r++)
85 rank[r] = Rank[mi(j,r,Size,k)];
86 Real* Xactive = (Real*)malloc(n*longueurActive*sizeof(Real));
87 for (int i=0; i<n; i++)
88 {
89 for (int jj=0; jj<longueurActive; jj++)
90 Xactive[mi(i,jj,n,longueurActive)] = X[mi(i,active[jj],n,p)];
91 }
92 Real* piLambda = (Real*)malloc(k*sizeof(Real));
93 for (int r=0; r<k; r++)
94 piLambda[r] = pi[mi(r,lambdaIndex,k,L)];
95 Real* rhoLambda = (Real*)malloc(m*m*k*sizeof(Real));
96 for (int u=0; u<m; u++)
97 {
98 for (int v=0; v<m; v++)
99 {
100 for (int r=0; r<k; r++)
101 rhoLambda[ai(u,v,r,m,m,k)] = rho[ai4(u,v,r,lambdaIndex,m,m,k,L)];
102 }
103 }
104 EMGrank_core(piLambda,rhoLambda,mini,maxi,Xactive,Y,tau,rank,
105 phiLambda,&LLF,
106 n,longueurActive,m,k);
107 free(rank);
108 free(Xactive);
109 free(piLambda);
110 free(rhoLambda);
111 //llh[(lambdaIndex-1)*Size+j,] = c(LLF, ...)
112 llh[mi(lambdaIndex*Size+j,0,L*Size,2)] = LLF;
113 //sum(Rank[j,] * (length(active)- Rank[j,] + m)) )
114 Real dotProduct = 0.0;
115 for (int r=0; r<k; r++)
116 dotProduct += Rank[mi(j,r,Size,k)] * (longueurActive-Rank[mi(j,r,Size,k)]+m);
117 llh[mi(lambdaIndex*Size+j,1,Size*L,2)] = dotProduct;
118 //phi[active,,,(lambdaIndex-1)*Size+j] = res$phi
119 for (int jj=0; jj<longueurActive; jj++)
120 {
121 for (int mm=0; mm<m; mm++)
122 {
123 for (int r=0; r<k; r++)
124 {
125 phi[ai4(active[jj],mm,r,lambdaIndex*Size+j,p,m,k,L*Size)] =
126 phiLambda[ai(jj,mm,r,longueurActive,m,k)];
127 }
128 }
129 }
130 }
131 free(active);
132 free(phiLambda);
133 }
134 }
135 free(Rank);
136 }