xref: /petsc/src/ksp/pc/impls/factor/qr/qr.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1a2fc1e05SToby Isaac 
2a2fc1e05SToby Isaac /*
3a2fc1e05SToby Isaac    Defines a direct QR factorization preconditioner for any Mat implementation
4a2fc1e05SToby Isaac    Note: this need not be considered a preconditioner since it supplies
5a2fc1e05SToby Isaac          a direct solver.
6a2fc1e05SToby Isaac */
7a2fc1e05SToby Isaac #include <../src/ksp/pc/impls/factor/qr/qr.h> /*I "petscpc.h" I*/
8a2fc1e05SToby Isaac 
99371c9d4SSatish Balay static PetscErrorCode PCSetUp_QR(PC pc) {
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;
31a2fc1e05SToby Isaac         PetscFunctionReturn(0);
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*4dfa11a4SJacob Faibussowitsch       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;
53a2fc1e05SToby Isaac       PetscFunctionReturn(0);
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   }
69a2fc1e05SToby Isaac   PetscFunctionReturn(0);
70a2fc1e05SToby Isaac }
71a2fc1e05SToby Isaac 
729371c9d4SSatish Balay static PetscErrorCode PCReset_QR(PC pc) {
73a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
74a2fc1e05SToby Isaac 
75a2fc1e05SToby Isaac   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
78a2fc1e05SToby Isaac   PetscFunctionReturn(0);
79a2fc1e05SToby Isaac }
80a2fc1e05SToby Isaac 
819371c9d4SSatish Balay static PetscErrorCode PCDestroy_QR(PC pc) {
82a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
83a2fc1e05SToby Isaac 
84a2fc1e05SToby Isaac   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(PCReset_QR(pc));
869566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
879566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
882e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
899566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
90a2fc1e05SToby Isaac   PetscFunctionReturn(0);
91a2fc1e05SToby Isaac }
92a2fc1e05SToby Isaac 
939371c9d4SSatish Balay static PetscErrorCode PCApply_QR(PC pc, Vec x, Vec y) {
94a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
95a2fc1e05SToby Isaac   Mat    fact;
96a2fc1e05SToby Isaac 
97a2fc1e05SToby Isaac   PetscFunctionBegin;
98a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
999566063dSJacob Faibussowitsch   PetscCall(MatSolve(fact, x, y));
100a2fc1e05SToby Isaac   PetscFunctionReturn(0);
101a2fc1e05SToby Isaac }
102a2fc1e05SToby Isaac 
1039371c9d4SSatish Balay static PetscErrorCode PCMatApply_QR(PC pc, Mat X, Mat Y) {
104a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
105a2fc1e05SToby Isaac   Mat    fact;
106a2fc1e05SToby Isaac 
107a2fc1e05SToby Isaac   PetscFunctionBegin;
108a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1099566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(fact, X, Y));
110a2fc1e05SToby Isaac   PetscFunctionReturn(0);
111a2fc1e05SToby Isaac }
112a2fc1e05SToby Isaac 
1139371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_QR(PC pc, Vec x, Vec y) {
114a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
115a2fc1e05SToby Isaac   Mat    fact;
116a2fc1e05SToby Isaac 
117a2fc1e05SToby Isaac   PetscFunctionBegin;
118a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1199566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(fact, x, y));
120a2fc1e05SToby Isaac   PetscFunctionReturn(0);
121a2fc1e05SToby Isaac }
122a2fc1e05SToby Isaac 
123a2fc1e05SToby Isaac /*MC
124a2fc1e05SToby Isaac    PCQR - Uses a direct solver, based on QR factorization, as a preconditioner
125a2fc1e05SToby Isaac 
126a2fc1e05SToby Isaac    Level: beginner
127a2fc1e05SToby Isaac 
128f1580f4eSBarry Smith    Note:
129a2fc1e05SToby Isaac    Usually this will compute an "exact" solution in one iteration and does
130a2fc1e05SToby Isaac    not need a Krylov method (i.e. you can use -ksp_type preonly, or
131f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
132a2fc1e05SToby Isaac 
133f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSVD`,
134db781477SPatrick Sanan           `PCILU`, `PCLU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
135db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
136c2e3fba1SPatrick Sanan           `PCFactorSetPivotingInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
137db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
138a2fc1e05SToby Isaac M*/
139a2fc1e05SToby Isaac 
1409371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_QR(PC pc) {
141a2fc1e05SToby Isaac   PC_QR *dir;
142a2fc1e05SToby Isaac 
143a2fc1e05SToby Isaac   PetscFunctionBegin;
144*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dir));
145a2fc1e05SToby Isaac   pc->data = (void *)dir;
1469566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_QR));
147a2fc1e05SToby Isaac 
148a2fc1e05SToby Isaac   dir->col                 = NULL;
149a2fc1e05SToby Isaac   pc->ops->reset           = PCReset_QR;
150a2fc1e05SToby Isaac   pc->ops->destroy         = PCDestroy_QR;
151a2fc1e05SToby Isaac   pc->ops->apply           = PCApply_QR;
152a2fc1e05SToby Isaac   pc->ops->matapply        = PCMatApply_QR;
153a2fc1e05SToby Isaac   pc->ops->applytranspose  = PCApplyTranspose_QR;
154a2fc1e05SToby Isaac   pc->ops->setup           = PCSetUp_QR;
155a2fc1e05SToby Isaac   pc->ops->setfromoptions  = PCSetFromOptions_Factor;
156a2fc1e05SToby Isaac   pc->ops->view            = PCView_Factor;
157a2fc1e05SToby Isaac   pc->ops->applyrichardson = NULL;
158a2fc1e05SToby Isaac   PetscFunctionReturn(0);
159a2fc1e05SToby Isaac }
160