00001
00009 #include "party.h"
00010
00011
00023 void C_LinearStatistic (const double *x, const int p,
00024 const double *y, const int q,
00025 const double *weights, const int n,
00026 double *ans) {
00027
00028 int i, j, k, kp, kn;
00029 double tmp;
00030
00031 for (k = 0; k < q; k++) {
00032
00033 kn = k * n;
00034 kp = k * p;
00035 for (j = 0; j < p; j++) ans[kp + j] = 0.0;
00036
00037 for (i = 0; i < n; i++) {
00038
00039
00040 if (weights[i] == 0.0) continue;
00041
00042 tmp = y[kn + i] * weights[i];
00043
00044 for (j = 0; j < p; j++)
00045 ans[kp + j] += x[j*n + i] * tmp;
00046 }
00047 }
00048 }
00049
00050
00058 SEXP R_LinearStatistic(SEXP x, SEXP y, SEXP weights) {
00059
00060
00061 SEXP ans;
00062
00063
00064 int n, p, q;
00065
00066
00067
00068
00069
00070
00071 if (!isReal(x) || !isReal(y) || !isReal(weights))
00072 error("LinStat: arguments are not of type REALSXP");
00073
00074 n = nrow(y);
00075 if (nrow(x) != n || LENGTH(weights) != n)
00076 error("LinStat: dimensions don't match");
00077
00078 q = ncol(y);
00079 p = ncol(x);
00080
00081 PROTECT(ans = allocVector(REALSXP, p*q));
00082
00083 C_LinearStatistic(REAL(x), p, REAL(y), q, REAL(weights), n,
00084 REAL(ans));
00085
00086 UNPROTECT(1);
00087 return(ans);
00088 }
00089
00090
00100 void C_ExpectCovarInfluence(const double* y, const int q,
00101 const double* weights, const int n,
00102 SEXP ans) {
00103
00104 int i, j, k, jq;
00105
00106
00107 double *dExp_y, *dCov_y, *dsweights, tmp;
00108
00109
00110 dExp_y = REAL(GET_SLOT(ans, PL2_expectationSym));
00111 for (j = 0; j < q; j++) dExp_y[j] = 0.0;
00112
00113 dCov_y = REAL(GET_SLOT(ans, PL2_covarianceSym));
00114 for (j = 0; j < q*q; j++) dCov_y[j] = 0.0;
00115
00116 dsweights = REAL(GET_SLOT(ans, PL2_sumweightsSym));
00117
00118
00119
00120 dsweights[0] = 0;
00121 for (i = 0; i < n; i++) dsweights[0] += weights[i];
00122 if (dsweights[0] <= 1)
00123 error("C_ExpectCovarInfluence: sum of weights is less than one");
00124
00125
00126
00127
00128
00129 for (i = 0; i < n; i++) {
00130
00131
00132
00133 if (weights[i] == 0.0) continue;
00134
00135 for (j = 0; j < q; j++)
00136 dExp_y[j] += weights[i] * y[j * n + i];
00137 }
00138
00139 for (j = 0; j < q; j++)
00140 dExp_y[j] = dExp_y[j] / dsweights[0];
00141
00142
00143
00144
00145
00146
00147 for (i = 0; i < n; i++) {
00148
00149 if (weights[i] == 0.0) continue;
00150
00151 for (j = 0; j < q; j++) {
00152 tmp = weights[i] * (y[j * n + i] - dExp_y[j]);
00153 jq = j * q;
00154 for (k = 0; k < q; k++)
00155 dCov_y[jq + k] += tmp * (y[k * n + i] - dExp_y[k]);
00156 }
00157 }
00158
00159 for (j = 0; j < q*q; j++)
00160 dCov_y[j] = dCov_y[j] / dsweights[0];
00161 }
00162
00163
00170 SEXP R_ExpectCovarInfluence(SEXP y, SEXP weights) {
00171
00172 SEXP ans;
00173 int q, n;
00174
00175 if (!isReal(y) || !isReal(weights))
00176 error("R_ExpectCovarInfluence: arguments are not of type REALSXP");
00177
00178 n = nrow(y);
00179 q = ncol(y);
00180
00181 if (LENGTH(weights) != n)
00182 error("R_ExpectCovarInfluence: vector of case weights does not have %d elements", n);
00183
00184
00185 PROTECT(ans = NEW_OBJECT(MAKE_CLASS("ExpectCovarInfluence")));
00186 SET_SLOT(ans, PL2_expectationSym,
00187 PROTECT(allocVector(REALSXP, q)));
00188 SET_SLOT(ans, PL2_covarianceSym,
00189 PROTECT(allocMatrix(REALSXP, q, q)));
00190 SET_SLOT(ans, PL2_sumweightsSym,
00191 PROTECT(allocVector(REALSXP, 1)));
00192
00193 C_ExpectCovarInfluence(REAL(y), q, REAL(weights), n, ans);
00194
00195 UNPROTECT(4);
00196 return(ans);
00197 }
00198
00199
00212 void C_ExpectCovarLinearStatistic(const double* x, const int p,
00213 const double* y, const int q,
00214 const double* weights, const int n,
00215 const SEXP expcovinf, SEXP ans) {
00216
00217 int i, j, k, pq;
00218 double sweights = 0.0, f1, f2, tmp;
00219 double *swx, *CT1, *CT2, *Covy_x_swx,
00220 *dExp_y, *dCov_y, *dExp_T, *dCov_T;
00221
00222 pq = p * q;
00223
00224
00225 dExp_y = REAL(GET_SLOT(expcovinf, PL2_expectationSym));
00226 dCov_y = REAL(GET_SLOT(expcovinf, PL2_covarianceSym));
00227 sweights = REAL(GET_SLOT(expcovinf, PL2_sumweightsSym))[0];
00228
00229 if (sweights <= 1.0)
00230 error("C_ExpectCovarLinearStatistic: sum of weights is less than one");
00231
00232
00233 dExp_T = REAL(GET_SLOT(ans, PL2_expectationSym));
00234 dCov_T = REAL(GET_SLOT(ans, PL2_covarianceSym));
00235
00236
00237 swx = Calloc(p, double);
00238 CT1 = Calloc(p * p, double);
00239
00240 for (i = 0; i < n; i++) {
00241
00242
00243 if (weights[i] == 0.0) continue;
00244
00245 for (k = 0; k < p; k++) {
00246 tmp = weights[i] * x[k * n + i];
00247 swx[k] += tmp;
00248
00249
00250 for (j = 0; j < p; j++) {
00251 CT1[j * p + k] += tmp * x[j * n + i];
00252 }
00253 }
00254 }
00255
00256
00257
00258
00259
00260 for (k = 0; k < p; k++) {
00261 for (j = 0; j < q; j++)
00262 dExp_T[j * p + k] = swx[k] * dExp_y[j];
00263 }
00264
00265
00266
00267
00268
00269 f1 = sweights/(sweights - 1);
00270 f2 = (1/(sweights - 1));
00271
00272 if (pq == 1) {
00273 dCov_T[0] = f1 * dCov_y[0] * CT1[0];
00274 dCov_T[0] -= f2 * dCov_y[0] * swx[0] * swx[0];
00275 } else {
00276
00277 CT2 = Calloc(pq * pq, double);
00278 Covy_x_swx = Calloc(pq * q, double);
00279
00280 C_kronecker(dCov_y, q, q, CT1, p, p, dCov_T);
00281 C_kronecker(dCov_y, q, q, swx, p, 1, Covy_x_swx);
00282 C_kronecker(Covy_x_swx, pq, q, swx, 1, p, CT2);
00283
00284 for (k = 0; k < (pq * pq); k++)
00285 dCov_T[k] = f1 * dCov_T[k] - f2 * CT2[k];
00286
00287
00288 Free(CT2); Free(Covy_x_swx);
00289 }
00290
00291
00292 Free(swx); Free(CT1);
00293 }
00294
00295
00304 SEXP R_ExpectCovarLinearStatistic(SEXP x, SEXP y, SEXP weights,
00305 SEXP expcovinf) {
00306
00307 SEXP ans;
00308 int n, p, q, pq;
00309
00310
00311
00312 n = nrow(x);
00313 p = ncol(x);
00314 q = ncol(y);
00315 pq = p * q;
00316
00317 if (nrow(y) != n)
00318 error("y does not have %d rows", n);
00319 if (LENGTH(weights) != n)
00320 error("vector of case weights does not have %d elements", n);
00321
00322 PROTECT(ans = NEW_OBJECT(MAKE_CLASS("ExpectCovar")));
00323 SET_SLOT(ans, PL2_expectationSym,
00324 PROTECT(allocVector(REALSXP, pq)));
00325 SET_SLOT(ans, PL2_covarianceSym,
00326 PROTECT(allocMatrix(REALSXP, pq, pq)));
00327
00328 C_ExpectCovarLinearStatistic(REAL(x), p, REAL(y), q,
00329 REAL(weights), n, expcovinf, ans);
00330
00331 UNPROTECT(3);
00332 return(ans);
00333 }
00334
00335
00349 void C_PermutedLinearStatistic(const double *x, const int p,
00350 const double *y, const int q,
00351 const int n, const int nperm,
00352 const int *indx, const int *perm,
00353 double *ans) {
00354
00355 int i, j, k, kp, kn, knpi;
00356
00357 for (k = 0; k < q; k++) {
00358
00359 kn = k * n;
00360 kp = k * p;
00361 for (j = 0; j < p; j++) ans[kp + j] = 0.0;
00362
00363 for (i = 0; i < nperm; i++) {
00364
00365 knpi = kn + perm[i];
00366
00367 for (j = 0; j < p; j++)
00368 ans[kp + j] += x[j*n + indx[i]] * y[knpi];
00369 }
00370 }
00371 }
00372
00373
00382 SEXP R_PermutedLinearStatistic(SEXP x, SEXP y, SEXP indx, SEXP perm) {
00383
00384 SEXP ans;
00385 int n, nperm, p, q, i, *iperm, *iindx;
00386
00387
00388
00389
00390
00391 if (!isReal(x) || !isReal(y))
00392 error("R_PermutedLinearStatistic: arguments are not of type REALSXP");
00393
00394 if (!isInteger(perm))
00395 error("R_PermutedLinearStatistic: perm is not of type INTSXP");
00396 if (!isInteger(indx))
00397 error("R_PermutedLinearStatistic: indx is not of type INTSXP");
00398
00399 n = nrow(y);
00400 nperm = LENGTH(perm);
00401 iperm = INTEGER(perm);
00402 if (LENGTH(indx) != nperm)
00403 error("R_PermutedLinearStatistic: dimensions don't match");
00404 iindx = INTEGER(indx);
00405
00406 if (nrow(x) != n)
00407 error("R_PermutedLinearStatistic: dimensions don't match");
00408
00409 for (i = 0; i < nperm; i++) {
00410 if (iperm[i] < 0 || iperm[i] > (n - 1) )
00411 error("R_PermutedLinearStatistic: perm is not between 1 and nobs");
00412 if (iindx[i] < 0 || iindx[i] > (n - 1) )
00413 error("R_PermutedLinearStatistic: indx is not between 1 and nobs");
00414 }
00415
00416 q = ncol(y);
00417 p = ncol(x);
00418
00419 PROTECT(ans = allocVector(REALSXP, p*q));
00420
00421 C_PermutedLinearStatistic(REAL(x), p, REAL(y), q, n, nperm,
00422 iindx, iperm, REAL(ans));
00423
00424 UNPROTECT(1);
00425 return(ans);
00426 }