1a2fc1e05SToby Isaac /*
2a2fc1e05SToby Isaac Defines a direct QR factorization preconditioner for any Mat implementation
3a2fc1e05SToby Isaac Note: this need not be considered a preconditioner since it supplies
4a2fc1e05SToby Isaac a direct solver.
5a2fc1e05SToby Isaac */
6a2fc1e05SToby Isaac #include <../src/ksp/pc/impls/factor/qr/qr.h> /*I "petscpc.h" I*/
7a2fc1e05SToby Isaac
PCSetUp_QR(PC pc)8d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_QR(PC pc)
9d71ae5a4SJacob Faibussowitsch {
10a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data;
11a2fc1e05SToby Isaac MatSolverType stype;
12a2fc1e05SToby Isaac MatFactorError err;
13f023e1d5SPierre Jolivet const char *prefix;
14a2fc1e05SToby Isaac
15a2fc1e05SToby Isaac PetscFunctionBegin;
16f023e1d5SPierre Jolivet PetscCall(PCGetOptionsPrefix(pc, &prefix));
17f023e1d5SPierre Jolivet PetscCall(MatSetOptionsPrefix(pc->pmat, prefix));
18a2fc1e05SToby Isaac pc->failedreason = PC_NOERROR;
19a2fc1e05SToby Isaac if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
20a2fc1e05SToby Isaac
219566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
22a2fc1e05SToby Isaac if (dir->hdr.inplace) {
23a2fc1e05SToby Isaac MatFactorType ftype;
24a2fc1e05SToby Isaac
259566063dSJacob Faibussowitsch PetscCall(MatGetFactorType(pc->pmat, &ftype));
26a2fc1e05SToby Isaac if (ftype == MAT_FACTOR_NONE) {
279566063dSJacob Faibussowitsch PetscCall(MatQRFactor(pc->pmat, dir->col, &((PC_Factor *)dir)->info));
289566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err));
29a2fc1e05SToby Isaac if (err) { /* Factor() fails */
30a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err;
313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
32a2fc1e05SToby Isaac }
33a2fc1e05SToby Isaac }
34a2fc1e05SToby Isaac ((PC_Factor *)dir)->fact = pc->pmat;
35a2fc1e05SToby Isaac } else {
36a2fc1e05SToby Isaac MatInfo info;
37a2fc1e05SToby Isaac
38a2fc1e05SToby Isaac if (!pc->setupcalled) {
39*3a7d0413SPierre Jolivet if (!((PC_Factor *)dir)->fact) PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_QR, &((PC_Factor *)dir)->fact));
409566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->col, &((PC_Factor *)dir)->info));
419566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
42a2fc1e05SToby Isaac dir->hdr.actualfill = info.fill_ratio_needed;
43a2fc1e05SToby Isaac } else if (pc->flag != SAME_NONZERO_PATTERN) {
449566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->col, &((PC_Factor *)dir)->info));
459566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
46a2fc1e05SToby Isaac dir->hdr.actualfill = info.fill_ratio_needed;
47a2fc1e05SToby Isaac } else {
489566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
49a2fc1e05SToby Isaac }
509566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
51a2fc1e05SToby Isaac if (err) { /* FactorSymbolic() fails */
52a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err;
533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
54a2fc1e05SToby Isaac }
55a2fc1e05SToby Isaac
569566063dSJacob Faibussowitsch PetscCall(MatQRFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
579566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
58a2fc1e05SToby Isaac if (err) { /* FactorNumeric() fails */
59a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err;
60a2fc1e05SToby Isaac }
61a2fc1e05SToby Isaac }
62a2fc1e05SToby Isaac
639566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype));
64a2fc1e05SToby Isaac if (!stype) {
65a2fc1e05SToby Isaac MatSolverType solverpackage;
669566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
679566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
68a2fc1e05SToby Isaac }
693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
70a2fc1e05SToby Isaac }
71a2fc1e05SToby Isaac
PCReset_QR(PC pc)72d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_QR(PC pc)
73d71ae5a4SJacob Faibussowitsch {
74a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data;
75a2fc1e05SToby Isaac
76a2fc1e05SToby Isaac PetscFunctionBegin;
779566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
80a2fc1e05SToby Isaac }
81a2fc1e05SToby Isaac
PCDestroy_QR(PC pc)82d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_QR(PC pc)
83d71ae5a4SJacob Faibussowitsch {
84a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data;
85a2fc1e05SToby Isaac
86a2fc1e05SToby Isaac PetscFunctionBegin;
879566063dSJacob Faibussowitsch PetscCall(PCReset_QR(pc));
889566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
899566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
902e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc));
919566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
93a2fc1e05SToby Isaac }
94a2fc1e05SToby Isaac
PCApply_QR(PC pc,Vec x,Vec y)95d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_QR(PC pc, Vec x, Vec y)
96d71ae5a4SJacob Faibussowitsch {
97a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data;
98a2fc1e05SToby Isaac Mat fact;
99a2fc1e05SToby Isaac
100a2fc1e05SToby Isaac PetscFunctionBegin;
101a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1029566063dSJacob Faibussowitsch PetscCall(MatSolve(fact, x, y));
1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
104a2fc1e05SToby Isaac }
105a2fc1e05SToby Isaac
PCMatApply_QR(PC pc,Mat X,Mat Y)106d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_QR(PC pc, Mat X, Mat Y)
107d71ae5a4SJacob Faibussowitsch {
108a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data;
109a2fc1e05SToby Isaac Mat fact;
110a2fc1e05SToby Isaac
111a2fc1e05SToby Isaac PetscFunctionBegin;
112a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1139566063dSJacob Faibussowitsch PetscCall(MatMatSolve(fact, X, Y));
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
115a2fc1e05SToby Isaac }
116a2fc1e05SToby Isaac
PCApplyTranspose_QR(PC pc,Vec x,Vec y)117d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_QR(PC pc, Vec x, Vec y)
118d71ae5a4SJacob Faibussowitsch {
119a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data;
120a2fc1e05SToby Isaac Mat fact;
121a2fc1e05SToby Isaac
122a2fc1e05SToby Isaac PetscFunctionBegin;
123a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1249566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(fact, x, y));
1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
126a2fc1e05SToby Isaac }
127a2fc1e05SToby Isaac
128a2fc1e05SToby Isaac /*MC
129a2fc1e05SToby Isaac PCQR - Uses a direct solver, based on QR factorization, as a preconditioner
130a2fc1e05SToby Isaac
131a2fc1e05SToby Isaac Level: beginner
132a2fc1e05SToby Isaac
133789736e1SBarry Smith Notes:
134789736e1SBarry Smith This is usually used as a preconditioner for `KSPLSQR`
135789736e1SBarry Smith
136789736e1SBarry Smith If the matrix used to construct the preconditioner is a `MATNORMAL` or `MATNORMALHERMITIAN` matrix, then the QR is performed on the
137789736e1SBarry Smith inner matrix to provide the preconditioner.
138a2fc1e05SToby Isaac
139562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSVD`,
140db781477SPatrick Sanan `PCILU`, `PCLU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
141db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
142c2e3fba1SPatrick Sanan `PCFactorSetPivotingInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
143789736e1SBarry Smith `PCFactorReorderForNonzeroDiagonal()`, `KSPLSQR`
144a2fc1e05SToby Isaac M*/
145a2fc1e05SToby Isaac
PCCreate_QR(PC pc)146d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_QR(PC pc)
147d71ae5a4SJacob Faibussowitsch {
148a2fc1e05SToby Isaac PC_QR *dir;
149a2fc1e05SToby Isaac
150a2fc1e05SToby Isaac PetscFunctionBegin;
1514dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir));
152a2fc1e05SToby Isaac pc->data = (void *)dir;
1539566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_QR));
154a2fc1e05SToby Isaac
155a2fc1e05SToby Isaac dir->col = NULL;
156a2fc1e05SToby Isaac pc->ops->reset = PCReset_QR;
157a2fc1e05SToby Isaac pc->ops->destroy = PCDestroy_QR;
158a2fc1e05SToby Isaac pc->ops->apply = PCApply_QR;
159a2fc1e05SToby Isaac pc->ops->matapply = PCMatApply_QR;
160a2fc1e05SToby Isaac pc->ops->applytranspose = PCApplyTranspose_QR;
161a2fc1e05SToby Isaac pc->ops->setup = PCSetUp_QR;
162a2fc1e05SToby Isaac pc->ops->setfromoptions = PCSetFromOptions_Factor;
163a2fc1e05SToby Isaac pc->ops->view = PCView_Factor;
164a2fc1e05SToby Isaac pc->ops->applyrichardson = NULL;
1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
166a2fc1e05SToby Isaac }
167