xref: /petsc/src/mat/impls/aij/seq/spqr/aijspqr.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
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 
9d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_SPQR_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc)
10d71ae5a4SJacob Faibussowitsch {
11a2fc1e05SToby Isaac   Mat_SeqAIJ        *aij;
12a2fc1e05SToby Isaac   Mat                AT;
13a2fc1e05SToby Isaac   const PetscScalar *aa;
14a2fc1e05SToby Isaac   PetscScalar       *ca;
15a2fc1e05SToby Isaac   const PetscInt    *ai, *aj;
16a2fc1e05SToby Isaac   PetscInt           n = A->cmap->n, i, j, k, nz;
17a2fc1e05SToby Isaac   SuiteSparse_long  *ci, *cj; /* SuiteSparse_long is the only choice for SPQR */
1802ef7dfaSPierre Jolivet   PetscBool          vain = PETSC_FALSE, flg;
19a2fc1e05SToby Isaac 
20a2fc1e05SToby Isaac   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMALHERMITIAN, &flg));
2265f45395SPierre Jolivet   if (flg) {
239566063dSJacob Faibussowitsch     PetscCall(MatNormalHermitianGetMat(A, &A));
2465f45395SPierre Jolivet   } else if (!PetscDefined(USE_COMPLEX)) {
259566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMAL, &flg));
2648a46eb9SPierre Jolivet     if (flg) PetscCall(MatNormalGetMat(A, &A));
2765f45395SPierre Jolivet   }
28a2fc1e05SToby Isaac   /* cholmod_sparse is compressed sparse column */
29b94d7dedSBarry Smith   PetscCall(MatIsSymmetric(A, 0.0, &flg));
3002ef7dfaSPierre Jolivet   if (flg) {
319566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
32a2fc1e05SToby Isaac     AT = A;
33a2fc1e05SToby Isaac   } else {
349566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &AT));
35a2fc1e05SToby Isaac   }
36a2fc1e05SToby Isaac   aij = (Mat_SeqAIJ *)AT->data;
37a2fc1e05SToby Isaac   ai  = aij->j;
38a2fc1e05SToby Isaac   aj  = aij->i;
39a2fc1e05SToby Isaac   for (j = 0, nz = 0; j < n; j++) nz += aj[j + 1] - aj[j];
409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n + 1, &cj, nz, &ci));
41a2fc1e05SToby Isaac   if (values) {
42a2fc1e05SToby Isaac     vain = PETSC_TRUE;
439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &ca));
449566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArrayRead(AT, &aa));
45a2fc1e05SToby Isaac   }
46a2fc1e05SToby Isaac   for (j = 0, k = 0; j < n; j++) {
47a2fc1e05SToby Isaac     cj[j] = k;
48a2fc1e05SToby Isaac     for (i = aj[j]; i < aj[j + 1]; i++, k++) {
49a2fc1e05SToby Isaac       ci[k] = ai[i];
50a2fc1e05SToby Isaac       if (values) ca[k] = aa[i];
51a2fc1e05SToby Isaac     }
52a2fc1e05SToby Isaac   }
53a2fc1e05SToby Isaac   cj[j]     = k;
54a2fc1e05SToby Isaac   *aijalloc = PETSC_TRUE;
55a2fc1e05SToby Isaac   *valloc   = vain;
5648a46eb9SPierre Jolivet   if (values) PetscCall(MatSeqAIJRestoreArrayRead(AT, &aa));
57a2fc1e05SToby Isaac 
589566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(C, sizeof(*C)));
59a2fc1e05SToby Isaac 
6002ef7dfaSPierre Jolivet   C->nrow   = (size_t)AT->cmap->n;
6102ef7dfaSPierre Jolivet   C->ncol   = (size_t)AT->rmap->n;
62a2fc1e05SToby Isaac   C->nzmax  = (size_t)nz;
63a2fc1e05SToby Isaac   C->p      = cj;
64a2fc1e05SToby Isaac   C->i      = ci;
65a2fc1e05SToby Isaac   C->x      = values ? ca : 0;
66a2fc1e05SToby Isaac   C->stype  = 0;
67a2fc1e05SToby Isaac   C->itype  = CHOLMOD_LONG;
68a2fc1e05SToby Isaac   C->xtype  = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
69a2fc1e05SToby Isaac   C->dtype  = CHOLMOD_DOUBLE;
70a2fc1e05SToby Isaac   C->sorted = 1;
71a2fc1e05SToby Isaac   C->packed = 1;
72a2fc1e05SToby Isaac 
739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AT));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75a2fc1e05SToby Isaac }
76a2fc1e05SToby Isaac 
77d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_SPQR(Mat A, MatSolverType *type)
78d71ae5a4SJacob Faibussowitsch {
79a2fc1e05SToby Isaac   PetscFunctionBegin;
80a2fc1e05SToby Isaac   *type = MATSOLVERSPQR;
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82a2fc1e05SToby Isaac }
83a2fc1e05SToby Isaac 
84a2fc1e05SToby Isaac #define GET_ARRAY_READ  0
85a2fc1e05SToby Isaac #define GET_ARRAY_WRITE 1
86a2fc1e05SToby Isaac 
87d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_SPQR_Internal(Mat F, cholmod_dense *cholB, cholmod_dense **_Y_handle)
88d71ae5a4SJacob Faibussowitsch {
89a2fc1e05SToby Isaac   Mat_CHOLMOD   *chol     = (Mat_CHOLMOD *)F->data;
9002ef7dfaSPierre Jolivet   cholmod_dense *Y_handle = NULL, *QTB_handle = NULL, *Z_handle = NULL;
91a2fc1e05SToby Isaac 
92a2fc1e05SToby Isaac   PetscFunctionBegin;
9302ef7dfaSPierre Jolivet   if (!chol->normal) {
94a2fc1e05SToby Isaac     QTB_handle = SuiteSparseQR_C_qmult(SPQR_QTX, chol->spqrfact, cholB, chol->common);
9528b400f6SJacob Faibussowitsch     PetscCheck(QTB_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_qmult failed");
96a2fc1e05SToby Isaac     Y_handle = SuiteSparseQR_C_solve(SPQR_RETX_EQUALS_B, chol->spqrfact, QTB_handle, chol->common);
9728b400f6SJacob Faibussowitsch     PetscCheck(Y_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed");
9802ef7dfaSPierre Jolivet   } else {
9902ef7dfaSPierre Jolivet     Z_handle = SuiteSparseQR_C_solve(SPQR_RTX_EQUALS_ETB, chol->spqrfact, cholB, chol->common);
10028b400f6SJacob Faibussowitsch     PetscCheck(Z_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed");
10102ef7dfaSPierre Jolivet     Y_handle = SuiteSparseQR_C_solve(SPQR_RETX_EQUALS_B, chol->spqrfact, Z_handle, chol->common);
10228b400f6SJacob Faibussowitsch     PetscCheck(Y_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed");
1033ba16761SJacob Faibussowitsch     PetscCallExternal(!cholmod_l_free_dense, &Z_handle, chol->common);
10402ef7dfaSPierre Jolivet   }
105a2fc1e05SToby Isaac   *_Y_handle = Y_handle;
1063ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_l_free_dense, &QTB_handle, chol->common);
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108a2fc1e05SToby Isaac }
109a2fc1e05SToby Isaac 
110d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_SPQR(Mat F, Vec B, Vec X)
111d71ae5a4SJacob Faibussowitsch {
112a2fc1e05SToby Isaac   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
113a2fc1e05SToby Isaac   cholmod_dense cholB, *Y_handle = NULL;
114a2fc1e05SToby Isaac   PetscInt      n;
115a2fc1e05SToby Isaac   PetscScalar  *v;
116a2fc1e05SToby Isaac 
117a2fc1e05SToby Isaac   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB));
1199566063dSJacob Faibussowitsch   PetscCall(MatSolve_SPQR_Internal(F, &cholB, &Y_handle));
1209566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
1219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(X, &v));
1229566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(v, (PetscScalar *)(Y_handle->x), n));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(X, &v));
1243ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common);
1259566063dSJacob Faibussowitsch   PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
127a2fc1e05SToby Isaac }
128a2fc1e05SToby Isaac 
129d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_SPQR(Mat F, Mat B, Mat X)
130d71ae5a4SJacob Faibussowitsch {
131a2fc1e05SToby Isaac   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
132a2fc1e05SToby Isaac   cholmod_dense cholB, *Y_handle = NULL;
133a2fc1e05SToby Isaac   PetscScalar  *v;
134a2fc1e05SToby Isaac   PetscInt      lda;
135a2fc1e05SToby Isaac 
136a2fc1e05SToby Isaac   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB));
1389566063dSJacob Faibussowitsch   PetscCall(MatSolve_SPQR_Internal(F, &cholB, &Y_handle));
1399566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(X, &v));
1409566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
141a2fc1e05SToby Isaac   if ((size_t)lda == Y_handle->d) {
1429566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(v, (PetscScalar *)(Y_handle->x), lda * Y_handle->ncol));
143a2fc1e05SToby Isaac   } else {
14448a46eb9SPierre 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));
145a2fc1e05SToby Isaac   }
1469566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(X, &v));
1473ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common);
1489566063dSJacob Faibussowitsch   PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150a2fc1e05SToby Isaac }
151a2fc1e05SToby Isaac 
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_SPQR_Internal(Mat F, cholmod_dense *cholB, cholmod_dense **_Y_handle)
153d71ae5a4SJacob Faibussowitsch {
154a2fc1e05SToby Isaac   Mat_CHOLMOD   *chol     = (Mat_CHOLMOD *)F->data;
155a2fc1e05SToby Isaac   cholmod_dense *Y_handle = NULL, *RTB_handle = NULL;
156a2fc1e05SToby Isaac 
157a2fc1e05SToby Isaac   PetscFunctionBegin;
158a2fc1e05SToby Isaac   RTB_handle = SuiteSparseQR_C_solve(SPQR_RTX_EQUALS_ETB, chol->spqrfact, cholB, chol->common);
15928b400f6SJacob Faibussowitsch   PetscCheck(RTB_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_solve failed");
160a2fc1e05SToby Isaac   Y_handle = SuiteSparseQR_C_qmult(SPQR_QX, chol->spqrfact, RTB_handle, chol->common);
16128b400f6SJacob Faibussowitsch   PetscCheck(Y_handle, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SuiteSparseQR_C_qmult failed");
162a2fc1e05SToby Isaac   *_Y_handle = Y_handle;
1633ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_l_free_dense, &RTB_handle, chol->common);
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165a2fc1e05SToby Isaac }
166a2fc1e05SToby Isaac 
167d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_SPQR(Mat F, Vec B, Vec X)
168d71ae5a4SJacob Faibussowitsch {
169a2fc1e05SToby Isaac   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
170a2fc1e05SToby Isaac   cholmod_dense cholB, *Y_handle = NULL;
171a2fc1e05SToby Isaac   PetscInt      n;
172a2fc1e05SToby Isaac   PetscScalar  *v;
173a2fc1e05SToby Isaac 
174a2fc1e05SToby Isaac   PetscFunctionBegin;
1759566063dSJacob Faibussowitsch   PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB));
1769566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose_SPQR_Internal(F, &cholB, &Y_handle));
1779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(X, &v));
1799566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(v, (PetscScalar *)Y_handle->x, n));
1809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(X, &v));
1813ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common);
1829566063dSJacob Faibussowitsch   PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184a2fc1e05SToby Isaac }
185a2fc1e05SToby Isaac 
186d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolveTranspose_SPQR(Mat F, Mat B, Mat X)
187d71ae5a4SJacob Faibussowitsch {
188a2fc1e05SToby Isaac   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
189a2fc1e05SToby Isaac   cholmod_dense cholB, *Y_handle = NULL;
190a2fc1e05SToby Isaac   PetscScalar  *v;
191a2fc1e05SToby Isaac   PetscInt      lda;
192a2fc1e05SToby Isaac 
193a2fc1e05SToby Isaac   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB));
1959566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose_SPQR_Internal(F, &cholB, &Y_handle));
1969566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(X, &v));
1979566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
198a2fc1e05SToby Isaac   if ((size_t)lda == Y_handle->d) {
1999566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(v, (PetscScalar *)Y_handle->x, lda * Y_handle->ncol));
200a2fc1e05SToby Isaac   } else {
20148a46eb9SPierre 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));
202a2fc1e05SToby Isaac   }
2039566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(X, &v));
2043ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_l_free_dense, &Y_handle, chol->common);
2059566063dSJacob Faibussowitsch   PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207a2fc1e05SToby Isaac }
208a2fc1e05SToby Isaac 
209d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatQRFactorNumeric_SPQR(Mat F, Mat A, const MatFactorInfo *info)
210d71ae5a4SJacob Faibussowitsch {
211a2fc1e05SToby Isaac   Mat_CHOLMOD   *chol = (Mat_CHOLMOD *)F->data;
212a2fc1e05SToby Isaac   cholmod_sparse cholA;
21365f45395SPierre Jolivet   PetscBool      aijalloc, valloc;
214d0609cedSBarry Smith   int            err;
215a2fc1e05SToby Isaac 
216a2fc1e05SToby Isaac   PetscFunctionBegin;
2179566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMALHERMITIAN, &chol->normal));
21848a46eb9SPierre Jolivet   if (!chol->normal && !PetscDefined(USE_COMPLEX)) PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMAL, &chol->normal));
2199566063dSJacob Faibussowitsch   PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc));
220d0609cedSBarry Smith   err = !SuiteSparseQR_C_numeric(PETSC_SMALL, &cholA, chol->spqrfact, chol->common);
221d0609cedSBarry Smith   PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "SPQR factorization failed with status %d", chol->common->status);
222a2fc1e05SToby Isaac 
2239566063dSJacob Faibussowitsch   if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
2249566063dSJacob Faibussowitsch   if (valloc) PetscCall(PetscFree(cholA.x));
225a2fc1e05SToby Isaac 
226a2fc1e05SToby Isaac   F->ops->solve    = MatSolve_SPQR;
227a2fc1e05SToby Isaac   F->ops->matsolve = MatMatSolve_SPQR;
22802ef7dfaSPierre Jolivet   if (chol->normal) {
22902ef7dfaSPierre Jolivet     F->ops->solvetranspose    = MatSolve_SPQR;
23002ef7dfaSPierre Jolivet     F->ops->matsolvetranspose = MatMatSolve_SPQR;
23102ef7dfaSPierre Jolivet   } else if (A->cmap->n == A->rmap->n) {
232a2fc1e05SToby Isaac     F->ops->solvetranspose    = MatSolveTranspose_SPQR;
233a2fc1e05SToby Isaac     F->ops->matsolvetranspose = MatMatSolveTranspose_SPQR;
234a2fc1e05SToby Isaac   }
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
236a2fc1e05SToby Isaac }
237a2fc1e05SToby Isaac 
238d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatQRFactorSymbolic_SPQR(Mat F, Mat A, IS perm, const MatFactorInfo *info)
239d71ae5a4SJacob Faibussowitsch {
240a2fc1e05SToby Isaac   Mat_CHOLMOD   *chol = (Mat_CHOLMOD *)F->data;
241a2fc1e05SToby Isaac   cholmod_sparse cholA;
24265f45395SPierre Jolivet   PetscBool      aijalloc, valloc;
243a2fc1e05SToby Isaac 
244a2fc1e05SToby Isaac   PetscFunctionBegin;
2459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMALHERMITIAN, &chol->normal));
24648a46eb9SPierre Jolivet   if (!chol->normal && !PetscDefined(USE_COMPLEX)) PetscCall(PetscObjectTypeCompare((PetscObject)A, MATNORMAL, &chol->normal));
2479566063dSJacob Faibussowitsch   PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc));
2483ba16761SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) PetscCallExternal(!cholmod_l_check_sparse, &cholA, chol->common);
2493ba16761SJacob Faibussowitsch   if (chol->spqrfact) PetscCallExternal(!SuiteSparseQR_C_free, &chol->spqrfact, chol->common);
250d0f82a0bSPierre Jolivet   chol->spqrfact = SuiteSparseQR_C_symbolic(SPQR_ORDERING_DEFAULT, 1, &cholA, chol->common);
25128b400f6SJacob Faibussowitsch   PetscCheck(chol->spqrfact, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status);
252a2fc1e05SToby Isaac 
2539566063dSJacob Faibussowitsch   if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
2549566063dSJacob Faibussowitsch   if (valloc) PetscCall(PetscFree(cholA.x));
255a2fc1e05SToby Isaac 
2569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", MatQRFactorNumeric_SPQR));
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258a2fc1e05SToby Isaac }
259a2fc1e05SToby Isaac 
260a2fc1e05SToby Isaac /*MC
261a2fc1e05SToby Isaac   MATSOLVERSPQR
262a2fc1e05SToby Isaac 
26311a5261eSBarry Smith   A matrix type providing direct solvers, QR factorizations, for sequential matrices
264a2fc1e05SToby Isaac   via the external package SPQR.
265a2fc1e05SToby Isaac 
2662ef1f0ffSBarry Smith   Use `./configure --download-suitesparse` to install PETSc to use SPQR
267a2fc1e05SToby Isaac 
26811a5261eSBarry Smith   Consult SPQR documentation for more information about the common parameters
269a2fc1e05SToby Isaac   which correspond to the options database keys below.
270a2fc1e05SToby Isaac 
271a2fc1e05SToby Isaac    Level: beginner
272a2fc1e05SToby Isaac 
27311a5261eSBarry Smith    Note:
274*1d27aa22SBarry Smith    SPQR is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html>
275a2fc1e05SToby Isaac 
2761cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCQR`, `PCFactorSetMatSolverType()`, `MatSolverType`
277a2fc1e05SToby Isaac M*/
278a2fc1e05SToby Isaac 
279d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat A, MatFactorType ftype, Mat *F)
280d71ae5a4SJacob Faibussowitsch {
281a2fc1e05SToby Isaac   Mat          B;
282a2fc1e05SToby Isaac   Mat_CHOLMOD *chol;
283a2fc1e05SToby Isaac   PetscInt     m = A->rmap->n, n = A->cmap->n;
284a2fc1e05SToby Isaac   const char  *prefix;
285a2fc1e05SToby Isaac 
286a2fc1e05SToby Isaac   PetscFunctionBegin;
287a2fc1e05SToby Isaac   /* Create the factorization matrix F */
2889566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2899566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
2909566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("spqr", &((PetscObject)B)->type_name));
2919566063dSJacob Faibussowitsch   PetscCall(MatGetOptionsPrefix(A, &prefix));
2929566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(B, prefix));
2939566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
2944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&chol));
295a2fc1e05SToby Isaac 
296a2fc1e05SToby Isaac   chol->Wrap = MatWrapCholmod_SPQR_seqaij;
297a2fc1e05SToby Isaac   B->data    = chol;
298a2fc1e05SToby Isaac 
299a2fc1e05SToby Isaac   B->ops->getinfo = MatGetInfo_CHOLMOD;
300a2fc1e05SToby Isaac   B->ops->view    = MatView_CHOLMOD;
301a2fc1e05SToby Isaac   B->ops->destroy = MatDestroy_CHOLMOD;
302a2fc1e05SToby Isaac 
3039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_SPQR));
3049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatQRFactorSymbolic_C", MatQRFactorSymbolic_SPQR));
305a2fc1e05SToby Isaac 
306a2fc1e05SToby Isaac   B->factortype   = MAT_FACTOR_QR;
307a2fc1e05SToby Isaac   B->assembled    = PETSC_TRUE;
308a2fc1e05SToby Isaac   B->preallocated = PETSC_TRUE;
309a2fc1e05SToby Isaac 
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
3119566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype));
312a2fc1e05SToby Isaac   B->canuseordering = PETSC_FALSE;
3139566063dSJacob Faibussowitsch   PetscCall(CholmodStart(B));
314a2fc1e05SToby Isaac   chol->common->itype = CHOLMOD_LONG;
315a2fc1e05SToby Isaac   *F                  = B;
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
317a2fc1e05SToby Isaac }
318