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