static void updateRegressionForNewCommonScale(SEXP regression, MERCache* cache) { // the update is only required if the total sum of squares term depends on the common scale, // itself being a question of whether or not the unmodeled coefficient prior is SEXP unmodeledCoefPrior = GET_SLOT(regression, blme_unmodeledCoefficientPriorSym); if (PRIOR_TYPE_SLOT(unmodeledCoefPrior) != PRIOR_TYPE_DIRECT || PRIOR_FAMILIES_SLOT(unmodeledCoefPrior)[0] != PRIOR_FAMILY_GAUSSIAN || getCommonScaleBit(PRIOR_SCALES_SLOT(unmodeledCoefPrior)[0]) != PRIOR_COMMON_SCALE_FALSE) return; const int* dims = DIMS_SLOT(regression); double* deviances = DEV_SLOT(regression); int numUnmodeledCoefs = dims[p_POS]; // we need to refactor (X'X - Rzx'Rzx + sigma^2 / sigma_beta^2 * I) double* lowerRightBlockRightFactorization = RX_SLOT(regression); // recover the cached version of X'X - Rzx'Rzx Memcpy(lowerRightBlockRightFactorization, (const double*) cache->downdatedDenseCrossproduct, numUnmodeledCoefs * numUnmodeledCoefs); addGaussianContributionToDenseBlock(regression, lowerRightBlockRightFactorization, deviances[dims[isREML_POS] ? sigmaREML_POS : sigmaML_POS]); int choleskyResult = getDenseCholeskyDecomposition(lowerRightBlockRightFactorization, numUnmodeledCoefs, TRIANGLE_TYPE_UPPER); if (choleskyResult > 0) error("Leading minor %d of downdated X'X is not positive definite.", choleskyResult); if (choleskyResult < 0) error("Illegal argument %d to cholesky decomposition (dpotrf).", -choleskyResult); deviances[ldRX2_POS] = 0.0; for (int j = 0; j < numUnmodeledCoefs; ++j) { deviances[ldRX2_POS] += 2.0 * log(lowerRightBlockRightFactorization[j * (numUnmodeledCoefs + 1)]); } // at this point, we have the correct Rx stored. now we need to // compute the new half projection and the new sum of squares double* unmodeledCoefProjection = cache->unmodeledCoefProjection; // copy in (X'Y - Rzx' theta half projection); beta half projection is Rx^-1 times that Memcpy(unmodeledCoefProjection, (const double*) cache->downdatedDenseResponseRotation, numUnmodeledCoefs); int i_one = 1; // solve A'x = b for A an Upper triangular, Tranposed, Non-unit matrix F77_CALL(dtrsv)("U", "T", "N", &numUnmodeledCoefs, lowerRightBlockRightFactorization, &numUnmodeledCoefs, unmodeledCoefProjection, &i_one); // now update the sums of squares double newSumOfSquares = getSumOfSquares(unmodeledCoefProjection, numUnmodeledCoefs); cache->totalSumOfSquares -= newSumOfSquares - cache->unmodeledCoefProjectionSumOfSquares; cache->unmodeledCoefProjectionSumOfSquares = newSumOfSquares; }
/** * Update the fixed effects and the orthogonal random effects in an MCMC sample * from an mer object. * * @param x an mer object * @param sigma current standard deviation of the per-observation * noise terms. * @param fvals pointer to memory in which to store the updated beta * @param rvals pointer to memory in which to store the updated b (may * be (double*)NULL) */ static void MCMC_beta_u(SEXP x, double sigma, double *fvals, double *rvals) { int *dims = DIMS_SLOT(x); int i1 = 1, p = dims[p_POS], q = dims[q_POS]; double *V = V_SLOT(x), *fixef = FIXEF_SLOT(x), *muEta = MUETA_SLOT(x), *u = U_SLOT(x), mone[] = {-1,0}, one[] = {1,0}; CHM_FR L = L_SLOT(x); double *del1 = Calloc(q, double), *del2 = Alloca(p, double); CHM_DN sol, rhs = N_AS_CHM_DN(del1, q, 1); R_CheckStack(); if (V || muEta) { error(_("Update not yet written")); } else { /* Linear mixed model */ update_L(x); update_RX(x); lmm_update_fixef_u(x); /* Update beta */ for (int j = 0; j < p; j++) del2[j] = sigma * norm_rand(); F77_CALL(dtrsv)("U", "N", "N", &p, RX_SLOT(x), &p, del2, &i1); for (int j = 0; j < p; j++) fixef[j] += del2[j]; /* Update u */ for (int j = 0; j < q; j++) del1[j] = sigma * norm_rand(); F77_CALL(dgemv)("N", &q, &p, mone, RZX_SLOT(x), &q, del2, &i1, one, del1, &i1); sol = M_cholmod_solve(CHOLMOD_Lt, L, rhs, &c); for (int j = 0; j < q; j++) u[j] += ((double*)(sol->x))[j]; M_cholmod_free_dense(&sol, &c); update_mu(x); /* and parts of the deviance slot */ } Memcpy(fvals, fixef, p); if (rvals) { update_ranef(x); Memcpy(rvals, RANEF_SLOT(x), q); } Free(del1); }
static void getDerivativesOfSumOfSquares(SEXP regression, MERCache* cache, double* firstDerivative, double* secondDerivative) { const int* dims = DIMS_SLOT(regression); double* deviances = DEV_SLOT(regression); int i_one = 1; double d_one = 1.0; int numUnmodeledCoefs = dims[p_POS]; SEXP unmodeledCoefPrior = GET_SLOT(regression, blme_unmodeledCoefficientPriorSym); const double* hyperparameters = PRIOR_HYPERPARAMETERS_SLOT(unmodeledCoefPrior) + 1; // skip over the log det of the covar, not needed here unsigned int numHyperparameters = LENGTH(GET_SLOT(unmodeledCoefPrior, blme_prior_hyperparametersSym)) - 1; // take Rx and get Rx^-1 const double* lowerRightFactor = RX_SLOT(regression); double rightFactorInverse[numUnmodeledCoefs * numUnmodeledCoefs]; // Rx^-1 invertUpperTriangularMatrix(lowerRightFactor, numUnmodeledCoefs, rightFactorInverse); // calculate Lbeta^-1 * Rx^-1 int factorIsTriangular = TRUE; if (numHyperparameters == 1) { // multiply by a scalar // printMatrix(lowerRightFactor, numUnmodeledCoefs, numUnmodeledCoefs); for (int col = 0; col < numUnmodeledCoefs; ++col) { int offset = col * numUnmodeledCoefs; for (int row = 0; row <= col; ++row) { rightFactorInverse[offset++] *= hyperparameters[0]; } } // printMatrix(rightFactorInverse, numUnmodeledCoefs, numUnmodeledCoefs); } else if (numHyperparameters == numUnmodeledCoefs) { // left multiply by a diagonal matrix const double *diagonal = hyperparameters; for (int col = 0; col < numUnmodeledCoefs; ++col) { int offset = col * numUnmodeledCoefs; for (int row = 0; row <= col; ++row) { rightFactorInverse[offset++] *= diagonal[row]; } } } else { const double* priorLeftFactorInverse = hyperparameters; // want L * R // Left multiply, Lower triangluar matrix, No-transpose, Non-unit F77_CALL(dtrmm)("L", "L", "N", "N", &numUnmodeledCoefs, &numUnmodeledCoefs, &d_one, (double*) priorLeftFactorInverse, &numUnmodeledCoefs, rightFactorInverse, &numUnmodeledCoefs); factorIsTriangular = FALSE; } double projectionRotation[numUnmodeledCoefs]; Memcpy(projectionRotation, (const double *) cache->unmodeledCoefProjection, numUnmodeledCoefs); // this step corresponds to Rx^-1 * unmodeled coef projection if (factorIsTriangular) { // X := A x, A triangular F77_CALL(dtrmv)("Upper triangular", "Non transposed", "Non unit diagonal", &numUnmodeledCoefs, rightFactorInverse, &numUnmodeledCoefs, projectionRotation, &i_one); } else { applyMatrixToVector(rightFactorInverse, numUnmodeledCoefs, numUnmodeledCoefs, FALSE, projectionRotation, projectionRotation); } double firstRotationSumOfSquares = getSumOfSquares(projectionRotation, numUnmodeledCoefs); // now for Rx^-T Rx^-1 * modeled coef projection if (factorIsTriangular) { // X: = A' x, A triangular F77_CALL(dtrmv)("Upper triangular", "Transposed", "Non unit diagonal", &numUnmodeledCoefs, rightFactorInverse, &numUnmodeledCoefs, projectionRotation, &i_one); } else { applyMatrixToVector(rightFactorInverse, numUnmodeledCoefs, numUnmodeledCoefs, TRUE, projectionRotation, projectionRotation); } double secondRotationSumOfSquares = getSumOfSquares(projectionRotation, numUnmodeledCoefs); // in general, DoF depends on unmodeled coefficient prior scale, as we can get back those // lost DoF. However, in that case we can't get here, where optimization is required. double sigma = deviances[dims[isREML_POS] ? sigmaREML_POS : sigmaML_POS]; double sigma_sq = sigma * sigma; *firstDerivative -= firstRotationSumOfSquares / sigma; *secondDerivative += 3.0 * firstRotationSumOfSquares / sigma_sq + 4.0 * secondRotationSumOfSquares; // From here, done unless REML. REML involves taking the derivative of // the log determinant of LxLx' (with some unmodeled covariance terms), // which is just the trace of the product. The second derivative is // the trace of the "square" of that product. if (dims[isREML_POS]) { int covarianceMatrixLength = numUnmodeledCoefs * numUnmodeledCoefs; double crossproduct[covarianceMatrixLength]; // we square the left factor Lx^-T * Lx^-1. the trace of this is immediately // useful, but we also need the trace of its square. Fortunately, the trace // of AA' is simply the sum of the squares of all of the elements. if (factorIsTriangular) { // want UU' singleTriangularMatrixCrossproduct(rightFactorInverse, numUnmodeledCoefs, TRUE, TRIANGLE_TYPE_UPPER, crossproduct); } else { singleMatrixCrossproduct(rightFactorInverse, numUnmodeledCoefs, numUnmodeledCoefs, crossproduct, TRUE, TRIANGLE_TYPE_UPPER); } double firstOrderTrace = 0.0; double secondOrderTrace = 0.0; int offset; // as the cross product is symmetric, we only have to use its upper // triangle and the diagonal for (int col = 0; col < numUnmodeledCoefs; ++col) { offset = col * numUnmodeledCoefs; for (int row = 0; row < col; ++row) { secondOrderTrace += 2.0 * crossproduct[offset] * crossproduct[offset]; ++offset; } firstOrderTrace += crossproduct[offset]; secondOrderTrace += crossproduct[offset] * crossproduct[offset]; } *firstDerivative -= sigma * firstOrderTrace; *secondDerivative += -firstOrderTrace + 2.0 * sigma_sq * secondOrderTrace; } }