1a2fc1e05SToby Isaac #include <petscsys.h> 2a2fc1e05SToby Isaac #include <../src/mat/impls/aij/seq/aij.h> 3a2fc1e05SToby Isaac #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 4a2fc1e05SToby Isaac 5a2fc1e05SToby Isaac EXTERN_C_BEGIN 6a2fc1e05SToby Isaac #include <SuiteSparseQR_C.h> 7a2fc1e05SToby Isaac EXTERN_C_END 8a2fc1e05SToby Isaac 9*789736e1SBarry Smith /* 10*789736e1SBarry Smith If A is a MATNORMALHERMITIAN or MATNORMAL format that it pulls the underlying matrix out 11*789736e1SBarry Smith and performs the factorization on that matrix 12*789736e1SBarry Smith */ 13d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_SPQR_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) 14d71ae5a4SJacob Faibussowitsch { 15a2fc1e05SToby Isaac Mat_SeqAIJ *aij; 16e2b46ddfSPierre Jolivet Mat AT, B = NULL; 17b9c875b8SPierre Jolivet Vec left, right; 18e2b46ddfSPierre Jolivet const PetscScalar *aa, *L = NULL; 19b9c875b8SPierre Jolivet PetscScalar *ca, scale; 20a2fc1e05SToby Isaac const PetscInt *ai, *aj; 21a2fc1e05SToby Isaac PetscInt n = A->cmap->n, i, j, k, nz; 22a2fc1e05SToby Isaac SuiteSparse_long *ci, *cj; /* SuiteSparse_long is the only choice for SPQR */ 2302ef7dfaSPierre Jolivet PetscBool vain = PETSC_FALSE, flg; 24a2fc1e05SToby Isaac 25a2fc1e05SToby Isaac PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMALHERMITIAN, &flg)); 2765f45395SPierre Jolivet if (flg) { 28e2b46ddfSPierre Jolivet PetscCall(MatNormalHermitianGetMat(A, &AT)); 2965f45395SPierre Jolivet } else if (!PetscDefined(USE_COMPLEX)) { 309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMAL, &flg)); 31e2b46ddfSPierre Jolivet if (flg) PetscCall(MatNormalGetMat(A, &AT)); 32e2b46ddfSPierre Jolivet } 33e2b46ddfSPierre Jolivet if (flg) { 34e2b46ddfSPierre Jolivet B = A; 35e2b46ddfSPierre Jolivet A = AT; 36b9c875b8SPierre Jolivet PetscCall(MatShellGetScalingShifts(B, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, &left, &right, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 37b9c875b8SPierre Jolivet PetscCheck(left == right, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Cannot call SuiteSparseQR if MatDiagonalScale() has been called on the input Mat with L != R"); 38b9c875b8SPierre Jolivet if (values && left) { 39b9c875b8SPierre Jolivet PetscCall(VecGetArrayRead(left, &L)); 40e2b46ddfSPierre Jolivet #if PetscDefined(USE_COMPLEX) 41e2b46ddfSPierre Jolivet for (j = 0; j < n; j++) 42e2b46ddfSPierre Jolivet PetscCheck(PetscAbsReal(PetscImaginaryPart(L[j])) < PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Cannot call SuiteSparseQR if MatDiagonalScale() has been called on the input Mat with a complex Vec"); 43e2b46ddfSPierre Jolivet #endif 44e2b46ddfSPierre Jolivet } 4565f45395SPierre Jolivet } 46a2fc1e05SToby Isaac /* cholmod_sparse is compressed sparse column */ 47b94d7dedSBarry Smith PetscCall(MatIsSymmetric(A, 0.0, &flg)); 4802ef7dfaSPierre Jolivet if (flg) { 499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 50a2fc1e05SToby Isaac AT = A; 51a2fc1e05SToby Isaac } else { 529566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &AT)); 53a2fc1e05SToby Isaac } 54a2fc1e05SToby Isaac aij = (Mat_SeqAIJ *)AT->data; 55a2fc1e05SToby Isaac ai = aij->j; 56a2fc1e05SToby Isaac aj = aij->i; 57a2fc1e05SToby Isaac for (j = 0, nz = 0; j < n; j++) nz += aj[j + 1] - aj[j]; 589566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n + 1, &cj, nz, &ci)); 59a2fc1e05SToby Isaac if (values) { 60a2fc1e05SToby Isaac vain = PETSC_TRUE; 619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz, &ca)); 629566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(AT, &aa)); 63a2fc1e05SToby Isaac } 64a2fc1e05SToby Isaac for (j = 0, k = 0; j < n; j++) { 65a2fc1e05SToby Isaac cj[j] = k; 66a2fc1e05SToby Isaac for (i = aj[j]; i < aj[j + 1]; i++, k++) { 67a2fc1e05SToby Isaac ci[k] = ai[i]; 68e2b46ddfSPierre Jolivet if (values) { 69e2b46ddfSPierre Jolivet ca[k] = aa[i]; 70e2b46ddfSPierre Jolivet if (L) ca[k] *= L[j]; 71e2b46ddfSPierre Jolivet } 72a2fc1e05SToby Isaac } 73a2fc1e05SToby Isaac } 74a2fc1e05SToby Isaac cj[j] = k; 75a2fc1e05SToby Isaac *aijalloc = PETSC_TRUE; 76a2fc1e05SToby Isaac *valloc = vain; 77e2b46ddfSPierre Jolivet if (values) { 78e2b46ddfSPierre Jolivet PetscCall(MatSeqAIJRestoreArrayRead(AT, &aa)); 79b9c875b8SPierre Jolivet if (L) PetscCall(VecRestoreArrayRead(left, &L)); 80e2b46ddfSPierre Jolivet } 81a2fc1e05SToby Isaac 829566063dSJacob Faibussowitsch PetscCall(PetscMemzero(C, sizeof(*C))); 83a2fc1e05SToby Isaac 8402ef7dfaSPierre Jolivet C->nrow = (size_t)AT->cmap->n; 8502ef7dfaSPierre Jolivet C->ncol = (size_t)AT->rmap->n; 86a2fc1e05SToby Isaac C->nzmax = (size_t)nz; 87a2fc1e05SToby Isaac C->p = cj; 88a2fc1e05SToby Isaac C->i = ci; 89a2fc1e05SToby Isaac C->x = values ? ca : 0; 90a2fc1e05SToby Isaac C->stype = 0; 91a2fc1e05SToby Isaac C->itype = CHOLMOD_LONG; 92a2fc1e05SToby Isaac C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 93a2fc1e05SToby Isaac C->dtype = CHOLMOD_DOUBLE; 94a2fc1e05SToby Isaac C->sorted = 1; 95a2fc1e05SToby Isaac C->packed = 1; 96a2fc1e05SToby Isaac 979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT)); 983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 99a2fc1e05SToby Isaac } 100a2fc1e05SToby Isaac 101d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_SPQR(Mat A, MatSolverType *type) 102d71ae5a4SJacob Faibussowitsch { 103a2fc1e05SToby Isaac PetscFunctionBegin; 104a2fc1e05SToby Isaac *type = MATSOLVERSPQR; 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106a2fc1e05SToby Isaac } 107a2fc1e05SToby Isaac 108a2fc1e05SToby Isaac #define GET_ARRAY_READ 0 109a2fc1e05SToby Isaac #define GET_ARRAY_WRITE 1 110a2fc1e05SToby Isaac 111d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_SPQR_Internal(Mat F, cholmod_dense *cholB, cholmod_dense **_Y_handle) 112d71ae5a4SJacob Faibussowitsch { 113a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 11402ef7dfaSPierre Jolivet cholmod_dense *Y_handle = NULL, *QTB_handle = NULL, *Z_handle = NULL; 115a2fc1e05SToby Isaac 116a2fc1e05SToby Isaac PetscFunctionBegin; 11702ef7dfaSPierre Jolivet if (!chol->normal) { 118a2fc1e05SToby Isaac QTB_handle = SuiteSparseQR_C_qmult(SPQR_QTX, chol->spqrfact, cholB, chol->common); 11928b400f6SJacob Faibussowitsch PetscCheck(QTB_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_qmult failed"); 120a2fc1e05SToby Isaac Y_handle = SuiteSparseQR_C_solve(SPQR_RETX_EQUALS_B, chol->spqrfact, QTB_handle, chol->common); 12128b400f6SJacob Faibussowitsch PetscCheck(Y_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed"); 12202ef7dfaSPierre Jolivet } else { 12302ef7dfaSPierre Jolivet Z_handle = SuiteSparseQR_C_solve(SPQR_RTX_EQUALS_ETB, chol->spqrfact, cholB, chol->common); 12428b400f6SJacob Faibussowitsch PetscCheck(Z_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed"); 12502ef7dfaSPierre Jolivet Y_handle = SuiteSparseQR_C_solve(SPQR_RETX_EQUALS_B, chol->spqrfact, Z_handle, chol->common); 12628b400f6SJacob Faibussowitsch PetscCheck(Y_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed"); 1273ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &Z_handle, chol->common); 12802ef7dfaSPierre Jolivet } 129a2fc1e05SToby Isaac *_Y_handle = Y_handle; 1303ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &QTB_handle, chol->common); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132a2fc1e05SToby Isaac } 133a2fc1e05SToby Isaac 134d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_SPQR(Mat F, Vec B, Vec X) 135d71ae5a4SJacob Faibussowitsch { 136a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 137a2fc1e05SToby Isaac cholmod_dense cholB, *Y_handle = NULL; 138a2fc1e05SToby Isaac PetscInt n; 139a2fc1e05SToby Isaac PetscScalar *v; 140a2fc1e05SToby Isaac 141a2fc1e05SToby Isaac PetscFunctionBegin; 1429566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB)); 1439566063dSJacob Faibussowitsch PetscCall(MatSolve_SPQR_Internal(F, &cholB, &Y_handle)); 1449566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 1459566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &v)); 146f4f49eeaSPierre Jolivet PetscCall(PetscArraycpy(v, (PetscScalar *)Y_handle->x, n)); 1479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, &v)); 1483ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common); 1499566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 150e2b46ddfSPierre Jolivet PetscCall(VecScale(X, chol->scale)); 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152a2fc1e05SToby Isaac } 153a2fc1e05SToby Isaac 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_SPQR(Mat F, Mat B, Mat X) 155d71ae5a4SJacob Faibussowitsch { 156a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 157a2fc1e05SToby Isaac cholmod_dense cholB, *Y_handle = NULL; 158a2fc1e05SToby Isaac PetscScalar *v; 159a2fc1e05SToby Isaac PetscInt lda; 160a2fc1e05SToby Isaac 161a2fc1e05SToby Isaac PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB)); 1639566063dSJacob Faibussowitsch PetscCall(MatSolve_SPQR_Internal(F, &cholB, &Y_handle)); 1649566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(X, &v)); 1659566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X, &lda)); 166a2fc1e05SToby Isaac if ((size_t)lda == Y_handle->d) { 167f4f49eeaSPierre Jolivet PetscCall(PetscArraycpy(v, (PetscScalar *)Y_handle->x, lda * Y_handle->ncol)); 168a2fc1e05SToby Isaac } else { 16948a46eb9SPierre Jolivet for (size_t j = 0; j < Y_handle->ncol; j++) PetscCall(PetscArraycpy(&v[j * lda], &(((PetscScalar *)Y_handle->x)[j * Y_handle->d]), Y_handle->nrow)); 170a2fc1e05SToby Isaac } 1719566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(X, &v)); 1723ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common); 1739566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 174e2b46ddfSPierre Jolivet PetscCall(MatScale(X, chol->scale)); 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176a2fc1e05SToby Isaac } 177a2fc1e05SToby Isaac 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_SPQR_Internal(Mat F, cholmod_dense *cholB, cholmod_dense **_Y_handle) 179d71ae5a4SJacob Faibussowitsch { 180a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 181a2fc1e05SToby Isaac cholmod_dense *Y_handle = NULL, *RTB_handle = NULL; 182a2fc1e05SToby Isaac 183a2fc1e05SToby Isaac PetscFunctionBegin; 184a2fc1e05SToby Isaac RTB_handle = SuiteSparseQR_C_solve(SPQR_RTX_EQUALS_ETB, chol->spqrfact, cholB, chol->common); 18528b400f6SJacob Faibussowitsch PetscCheck(RTB_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed"); 186a2fc1e05SToby Isaac Y_handle = SuiteSparseQR_C_qmult(SPQR_QX, chol->spqrfact, RTB_handle, chol->common); 18728b400f6SJacob Faibussowitsch PetscCheck(Y_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_qmult failed"); 188a2fc1e05SToby Isaac *_Y_handle = Y_handle; 1893ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &RTB_handle, chol->common); 1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 191a2fc1e05SToby Isaac } 192a2fc1e05SToby Isaac 193d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_SPQR(Mat F, Vec B, Vec X) 194d71ae5a4SJacob Faibussowitsch { 195a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 196a2fc1e05SToby Isaac cholmod_dense cholB, *Y_handle = NULL; 197a2fc1e05SToby Isaac PetscInt n; 198a2fc1e05SToby Isaac PetscScalar *v; 199a2fc1e05SToby Isaac 200a2fc1e05SToby Isaac PetscFunctionBegin; 2019566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB)); 2029566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose_SPQR_Internal(F, &cholB, &Y_handle)); 2039566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 2049566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &v)); 2059566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(v, (PetscScalar *)Y_handle->x, n)); 2069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, &v)); 2073ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common); 2089566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210a2fc1e05SToby Isaac } 211a2fc1e05SToby Isaac 212d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolveTranspose_SPQR(Mat F, Mat B, Mat X) 213d71ae5a4SJacob Faibussowitsch { 214a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 215a2fc1e05SToby Isaac cholmod_dense cholB, *Y_handle = NULL; 216a2fc1e05SToby Isaac PetscScalar *v; 217a2fc1e05SToby Isaac PetscInt lda; 218a2fc1e05SToby Isaac 219a2fc1e05SToby Isaac PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB)); 2219566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose_SPQR_Internal(F, &cholB, &Y_handle)); 2229566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(X, &v)); 2239566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X, &lda)); 224a2fc1e05SToby Isaac if ((size_t)lda == Y_handle->d) { 2259566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(v, (PetscScalar *)Y_handle->x, lda * Y_handle->ncol)); 226a2fc1e05SToby Isaac } else { 22748a46eb9SPierre Jolivet for (size_t j = 0; j < Y_handle->ncol; j++) PetscCall(PetscArraycpy(&v[j * lda], &(((PetscScalar *)Y_handle->x)[j * Y_handle->d]), Y_handle->nrow)); 228a2fc1e05SToby Isaac } 2299566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(X, &v)); 2303ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common); 2319566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233a2fc1e05SToby Isaac } 234a2fc1e05SToby Isaac 235d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatQRFactorNumeric_SPQR(Mat F, Mat A, const MatFactorInfo *info) 236d71ae5a4SJacob Faibussowitsch { 237a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 238a2fc1e05SToby Isaac cholmod_sparse cholA; 23965f45395SPierre Jolivet PetscBool aijalloc, valloc; 240d0609cedSBarry Smith int err; 241a2fc1e05SToby Isaac 242a2fc1e05SToby Isaac PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMALHERMITIAN, &chol->normal)); 24448a46eb9SPierre Jolivet if (!chol->normal && !PetscDefined(USE_COMPLEX)) PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMAL, &chol->normal)); 2459566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc)); 246d0609cedSBarry Smith err = !SuiteSparseQR_C_numeric(PETSC_SMALL, &cholA, chol->spqrfact, chol->common); 247d0609cedSBarry Smith PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SPQR factorization failed with status %d", chol->common->status); 248a2fc1e05SToby Isaac 2499566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 2509566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 251a2fc1e05SToby Isaac 252a2fc1e05SToby Isaac F->ops->solve = MatSolve_SPQR; 253a2fc1e05SToby Isaac F->ops->matsolve = MatMatSolve_SPQR; 25402ef7dfaSPierre Jolivet if (chol->normal) { 255b9c875b8SPierre Jolivet PetscScalar scale; 256b9c875b8SPierre Jolivet 257b9c875b8SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, NULL, NULL, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 25802ef7dfaSPierre Jolivet F->ops->solvetranspose = MatSolve_SPQR; 25902ef7dfaSPierre Jolivet F->ops->matsolvetranspose = MatMatSolve_SPQR; 260b9c875b8SPierre Jolivet chol->scale = 1.0 / scale; 26102ef7dfaSPierre Jolivet } else if (A->cmap->n == A->rmap->n) { 262a2fc1e05SToby Isaac F->ops->solvetranspose = MatSolveTranspose_SPQR; 263a2fc1e05SToby Isaac F->ops->matsolvetranspose = MatMatSolveTranspose_SPQR; 264e2b46ddfSPierre Jolivet chol->scale = 1.0; 265a2fc1e05SToby Isaac } 2663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 267a2fc1e05SToby Isaac } 268a2fc1e05SToby Isaac 269d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatQRFactorSymbolic_SPQR(Mat F, Mat A, IS perm, const MatFactorInfo *info) 270d71ae5a4SJacob Faibussowitsch { 271a2fc1e05SToby Isaac Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 272a2fc1e05SToby Isaac cholmod_sparse cholA; 27365f45395SPierre Jolivet PetscBool aijalloc, valloc; 274a2fc1e05SToby Isaac 275a2fc1e05SToby Isaac PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMALHERMITIAN, &chol->normal)); 27748a46eb9SPierre Jolivet if (!chol->normal && !PetscDefined(USE_COMPLEX)) PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMAL, &chol->normal)); 2789566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc)); 2793ba16761SJacob Faibussowitsch if (PetscDefined(USE_DEBUG)) PetscCallExternal(!cholmod_l_check_sparse, &cholA, chol->common); 2803ba16761SJacob Faibussowitsch if (chol->spqrfact) PetscCallExternal(!SuiteSparseQR_C_free, &chol->spqrfact, chol->common); 281d0f82a0bSPierre Jolivet chol->spqrfact = SuiteSparseQR_C_symbolic(SPQR_ORDERING_DEFAULT, 1, &cholA, chol->common); 28228b400f6SJacob Faibussowitsch PetscCheck(chol->spqrfact, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status); 283a2fc1e05SToby Isaac 2849566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 2859566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 286a2fc1e05SToby Isaac 2879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", MatQRFactorNumeric_SPQR)); 2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 289a2fc1e05SToby Isaac } 290a2fc1e05SToby Isaac 291a2fc1e05SToby Isaac /*MC 292a2fc1e05SToby Isaac MATSOLVERSPQR 293a2fc1e05SToby Isaac 294*789736e1SBarry Smith A matrix solver type providing direct solvers, QR factorizations, for sequential `MATSEQAIJ` and `MATNORMAL` matrices (that contain a `MATSEQAIJ`). 295a2fc1e05SToby Isaac via the external package SPQR. 296a2fc1e05SToby Isaac 2972ef1f0ffSBarry Smith Use `./configure --download-suitesparse` to install PETSc to use SPQR 298a2fc1e05SToby Isaac 299a2fc1e05SToby Isaac Level: beginner 300a2fc1e05SToby Isaac 30111a5261eSBarry Smith Note: 3021d27aa22SBarry Smith SPQR is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html> 303a2fc1e05SToby Isaac 3041cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCQR`, `PCFactorSetMatSolverType()`, `MatSolverType` 305a2fc1e05SToby Isaac M*/ 306a2fc1e05SToby Isaac 307d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat A, MatFactorType ftype, Mat *F) 308d71ae5a4SJacob Faibussowitsch { 309a2fc1e05SToby Isaac Mat B; 310a2fc1e05SToby Isaac Mat_CHOLMOD *chol; 311a2fc1e05SToby Isaac PetscInt m = A->rmap->n, n = A->cmap->n; 312a2fc1e05SToby Isaac const char *prefix; 313a2fc1e05SToby Isaac 314a2fc1e05SToby Isaac PetscFunctionBegin; 315a2fc1e05SToby Isaac /* Create the factorization matrix F */ 3169566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 3179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 3189566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("spqr", &((PetscObject)B)->type_name)); 3199566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A, &prefix)); 3209566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, prefix)); 3219566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 3224dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&chol)); 323a2fc1e05SToby Isaac 324a2fc1e05SToby Isaac chol->Wrap = MatWrapCholmod_SPQR_seqaij; 325a2fc1e05SToby Isaac B->data = chol; 326a2fc1e05SToby Isaac 327a2fc1e05SToby Isaac B->ops->getinfo = MatGetInfo_CHOLMOD; 328a2fc1e05SToby Isaac B->ops->view = MatView_CHOLMOD; 329a2fc1e05SToby Isaac B->ops->destroy = MatDestroy_CHOLMOD; 330a2fc1e05SToby Isaac 3319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_SPQR)); 3329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatQRFactorSymbolic_C", MatQRFactorSymbolic_SPQR)); 333a2fc1e05SToby Isaac 334a2fc1e05SToby Isaac B->factortype = MAT_FACTOR_QR; 335a2fc1e05SToby Isaac B->assembled = PETSC_TRUE; 336a2fc1e05SToby Isaac B->preallocated = PETSC_TRUE; 337a2fc1e05SToby Isaac 3389566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 3399566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype)); 340a2fc1e05SToby Isaac B->canuseordering = PETSC_FALSE; 3419566063dSJacob Faibussowitsch PetscCall(CholmodStart(B)); 342a2fc1e05SToby Isaac chol->common->itype = CHOLMOD_LONG; 343a2fc1e05SToby Isaac *F = B; 3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 345a2fc1e05SToby Isaac } 346