00001
00009 #include "party.h"
00010
00011
00018 double C_quadformConditionalPvalue(const double tstat, const double df) {
00019 return(pchisq(tstat, df, 0, 0));
00020 }
00021
00022
00029 SEXP R_quadformConditionalPvalue(SEXP tstat, SEXP df) {
00030
00031 SEXP ans;
00032
00033 PROTECT(ans = allocVector(REALSXP, 1));
00034 REAL(ans)[0] = C_quadformConditionalPvalue(REAL(tstat)[0], REAL(df)[0]);
00035 UNPROTECT(1);
00036 return(ans);
00037 }
00038
00039
00052 double C_maxabsConditionalPvalue(const double tstat, const double *Sigma,
00053 const int pq, int *maxpts, double *releps, double *abseps, double *tol) {
00054
00055 int *n, *nu, *inform, i, j, *infin, sub, *index, nonzero, iz, jz;
00056 double *lower, *upper, *delta, *corr, *sd, *myerror,
00057 *prob, ans;
00058
00059
00060 if (pq == 1)
00061 return(2*pnorm(fabs(tstat)*-1.0, 0.0, 1.0, 1, 0));
00062
00063 n = Calloc(1, int);
00064 nu = Calloc(1, int);
00065 myerror = Calloc(1, double);
00066 prob = Calloc(1, double);
00067 nu[0] = 0;
00068 inform = Calloc(1, int);
00069 n[0] = pq;
00070
00071 if (n[0] == 2)
00072 corr = Calloc(1, double);
00073 else
00074 corr = Calloc(n[0] + ((n[0] - 2) * (n[0] - 1))/2, double);
00075
00076 sd = Calloc(n[0], double);
00077 lower = Calloc(n[0], double);
00078 upper = Calloc(n[0], double);
00079 infin = Calloc(n[0], int);
00080 delta = Calloc(n[0], double);
00081 index = Calloc(n[0], int);
00082
00083
00084
00085 nonzero = 0;
00086 for (i = 0; i < n[0]; i++) {
00087 if (Sigma[i*n[0] + i] > tol[0]) {
00088 index[nonzero] = i;
00089 nonzero++;
00090 }
00091 }
00092
00093
00094
00095
00096
00097 for (iz = 0; iz < nonzero; iz++) {
00098
00099
00100 i = index[iz];
00101
00102
00103 sd[i] = sqrt(Sigma[i*n[0] + i]);
00104
00105
00106 lower[iz] = fabs(tstat) * -1.0;
00107 upper[iz] = fabs(tstat);
00108 infin[iz] = 2;
00109 delta[iz] = 0.0;
00110
00111
00112
00113 for (jz = 0; jz < iz; jz++) {
00114 j = index[jz];
00115 sub = (int) (jz + 1) + (double) ((iz - 1) * iz) / 2 - 1;
00116 if (sd[i] == 0.0 || sd[j] == 0.0)
00117 corr[sub] = 0.0;
00118 else
00119 corr[sub] = Sigma[i*n[0] + j] / (sd[i] * sd[j]);
00120 }
00121 }
00122 n[0] = nonzero;
00123
00124
00125 F77_CALL(mvtdst)(n, nu, lower, upper, infin, corr, delta,
00126 maxpts, abseps, releps, myerror, prob, inform);
00127
00128
00129 switch (inform[0]) {
00130 case 0: break;
00131 case 1: warning("cmvnorm: completion with ERROR > EPS"); break;
00132 case 2: warning("cmvnorm: N > 1000 or N < 1");
00133 prob[0] = 0.0;
00134 break;
00135 case 3: warning("cmvnorm: correlation matrix not positive semi-definite");
00136 prob[0] = 0.0;
00137 break;
00138 default: warning("cmvnorm: unknown problem in MVTDST");
00139 prob[0] = 0.0;
00140 }
00141 ans = prob[0];
00142 Free(corr); Free(sd); Free(lower); Free(upper);
00143 Free(infin); Free(delta); Free(myerror); Free(prob);
00144 Free(n); Free(nu); Free(inform);
00145 return(1 - ans);
00146 }
00147
00148
00159 SEXP R_maxabsConditionalPvalue(SEXP tstat, SEXP Sigma, SEXP maxpts,
00160 SEXP releps, SEXP abseps, SEXP tol) {
00161
00162 SEXP ans;
00163 int pq;
00164
00165 pq = nrow(Sigma);
00166
00167 PROTECT(ans = allocVector(REALSXP, 1));
00168 REAL(ans)[0] = C_maxabsConditionalPvalue(REAL(tstat)[0], REAL(Sigma), pq,
00169 INTEGER(maxpts), REAL(releps), REAL(abseps), REAL(tol));
00170 UNPROTECT(1);
00171 return(ans);
00172 }
00173
00174
00186 void C_MonteCarlo(double *criterion, SEXP learnsample, SEXP weights,
00187 SEXP fitmem, SEXP varctrl, SEXP gtctrl, double *ans_pvalues) {
00188
00189 int ninputs, nobs, j, i, k;
00190 SEXP responses, inputs, y, x, xmem, expcovinf;
00191 double sweights, *stats, tmp = 0.0, smax, *dweights;
00192 int m, *counts, b, B, *dummy, *permindex, *index, *permute;
00193
00194 ninputs = get_ninputs(learnsample);
00195 nobs = get_nobs(learnsample);
00196 responses = GET_SLOT(learnsample, PL2_responsesSym);
00197 inputs = GET_SLOT(learnsample, PL2_inputsSym);
00198 dweights = REAL(weights);
00199
00200
00201 B = get_nresample(gtctrl);
00202
00203
00204 y = get_test_trafo(responses);
00205
00206 expcovinf = GET_SLOT(fitmem, PL2_expcovinfSym);
00207
00208 sweights = REAL(GET_SLOT(expcovinf, PL2_sumweightsSym))[0];
00209 m = (int) sweights;
00210
00211 stats = Calloc(ninputs, double);
00212 counts = Calloc(ninputs, int);
00213
00214 dummy = Calloc(m, int);
00215 permute = Calloc(m, int);
00216 index = Calloc(m, int);
00217 permindex = Calloc(m, int);
00218
00219
00220
00221 j = 0;
00222 for (i = 0; i < nobs; i++) {
00223 for (k = 0; k < dweights[i]; k++) {
00224 index[j] = i;
00225 j++;
00226 }
00227 }
00228
00229 for (b = 0; b < B; b++) {
00230
00231
00232 C_SampleNoReplace(dummy, m, m, permute);
00233 for (k = 0; k < m; k++) permindex[k] = index[permute[k]];
00234
00235
00236 for (j = 1; j <= ninputs; j++) {
00237 x = get_transformation(inputs, j);
00238
00239
00240 xmem = get_varmemory(fitmem, j);
00241 if (!has_missings(inputs, j)) {
00242 C_PermutedLinearStatistic(REAL(x), ncol(x), REAL(y), ncol(y),
00243 nobs, m, index, permindex,
00244 REAL(GET_SLOT(xmem, PL2_linearstatisticSym)));
00245 } else {
00246 error("cannot resample with missing values");
00247 }
00248
00249
00250 C_TeststatCriterion(xmem, varctrl, &tmp, &stats[j - 1]);
00251 }
00252
00253
00254 smax = C_max(stats, ninputs);
00255
00256
00257 for (j = 0; j < ninputs; j++) {
00258 if (smax > criterion[j]) counts[j]++;
00259 }
00260 }
00261
00262
00263 for (j = 0; j < ninputs; j++)
00264 ans_pvalues[j] = (double) counts[j] / B;
00265
00266
00267
00268
00269
00270 for (j = 1; j <= ninputs; j++) {
00271 x = get_transformation(inputs, j);
00272
00273 xmem = get_varmemory(fitmem, j);
00274 C_LinearStatistic(REAL(x), ncol(x), REAL(y), ncol(y),
00275 dweights, nobs,
00276 REAL(GET_SLOT(xmem, PL2_linearstatisticSym)));
00277 }
00278
00279
00280 Free(stats); Free(counts); Free(dummy); Free(permute);
00281 Free(index); Free(permindex);
00282 }
00283
00284
00295 SEXP R_MonteCarlo(SEXP criterion, SEXP learnsample, SEXP weights,
00296 SEXP fitmem, SEXP varctrl, SEXP gtctrl) {
00297
00298 SEXP ans;
00299
00300 GetRNGstate();
00301
00302 PROTECT(ans = allocVector(REALSXP, get_ninputs(learnsample)));
00303 C_MonteCarlo(REAL(criterion), learnsample, weights, fitmem, varctrl,
00304 gtctrl, REAL(ans));
00305
00306 PutRNGstate();
00307
00308 UNPROTECT(1);
00309 return(ans);
00310 }