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