void solveSparseCholeskySystem(int operationType, const_CHM_FR leftHandSide, const double *rightHandSide, int numColsRightHandSide, double *target) { int commonDim = leftHandSide->n; CHM_DN chol_rightHandSide = N_AS_CHM_DN((double *) rightHandSide, commonDim, numColsRightHandSide); CHM_DN solution; R_CheckStack(); solution = M_cholmod_solve(operationType, leftHandSide, chol_rightHandSide, &cholmodCommon); if (!solution) error("cholmod_solve (%d) failed", operationType); // copy result into the dense, non-cholmod matrix Memcpy(target, (const double *) (solution->x), commonDim * numColsRightHandSide); M_cholmod_free_dense(&solution, &cholmodCommon); }
/** * 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); }
/** * Update the Xb, Zu, eta and mu slots in cpglmm according to x * * @param x pointer to the vector of values for beta or u * @param is_beta indicates whether x contains the values for beta or u. * 1: x contains values of beta, and Zu is not updated. If x is null, the fixef slot is used; * 0: x contains values of u, and Xb is not updated. If x is null, the u slot is used; * -1: x is ignored, and the fixef and u slots are used. * @param da an SEXP object * */ static void cpglmm_fitted(double *x, int is_beta, SEXP da){ int *dm = DIMS_SLOT(da) ; int nO = dm[nO_POS], nB = dm[nB_POS], nU = dm[nU_POS]; double *X = X_SLOT(da), *eta = ETA_SLOT(da), *mu = MU_SLOT(da), *beta = FIXEF_SLOT(da), *u = U_SLOT(da), *offset= OFFSET_SLOT(da), *Xb = XB_SLOT(da), *Zu = ZU_SLOT(da), lp = LKP_SLOT(da)[0], one[] = {1, 0}, zero[] = {0, 0}; if (is_beta == -1) x = NULL ; /* update from the fixef and u slots */ /* update Xb */ if (is_beta == 1 || is_beta == -1){ /* beta is updated */ if (x) beta = x ; /* point beta to x if x is not NULL */ mult_mv("N", nO, nB, X, beta, Xb) ; /* Xb = x * beta */ } /* update Zu */ if (is_beta == 0 || is_beta == -1){ /* u is updated */ SEXP tmp; /* create an SEXP object to be coerced to CHM_DN */ PROTECT(tmp = allocVector(REALSXP, nU)); if (x) Memcpy(REAL(tmp), x, nU); else Memcpy(REAL(tmp), u, nU); CHM_DN ceta, us = AS_CHM_DN(tmp); CHM_SP Zt = Zt_SLOT(da); R_CheckStack(); ceta = N_AS_CHM_DN(Zu, nO, 1); /* update Zu */ R_CheckStack(); /* Y = alpha * A * X + beta * Y */ if (!M_cholmod_sdmult(Zt, 1 , one, zero, us, ceta, &c)) error(_("cholmod_sdmult error returned")); UNPROTECT(1) ; } for (int i = 0; i < nO; i++){ /* update mu */ eta[i] = Xb[i] + Zu[i] + offset[i]; /* eta = Xb + Z * u + offset*/ mu[i] = link_inv(eta[i], lp); } }