xref: /petsc/src/ksp/pc/impls/factor/qr/qr.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
9*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_QR(PC pc)
10*d71ae5a4SJacob Faibussowitsch {
11a2fc1e05SToby Isaac   PC_QR         *dir = (PC_QR *)pc->data;
12a2fc1e05SToby Isaac   MatSolverType  stype;
13a2fc1e05SToby Isaac   MatFactorError err;
14f023e1d5SPierre Jolivet   const char    *prefix;
15a2fc1e05SToby Isaac 
16a2fc1e05SToby Isaac   PetscFunctionBegin;
17f023e1d5SPierre Jolivet   PetscCall(PCGetOptionsPrefix(pc, &prefix));
18f023e1d5SPierre Jolivet   PetscCall(MatSetOptionsPrefix(pc->pmat, prefix));
19a2fc1e05SToby Isaac   pc->failedreason = PC_NOERROR;
20a2fc1e05SToby Isaac   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
21a2fc1e05SToby Isaac 
229566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
23a2fc1e05SToby Isaac   if (dir->hdr.inplace) {
24a2fc1e05SToby Isaac     MatFactorType ftype;
25a2fc1e05SToby Isaac 
269566063dSJacob Faibussowitsch     PetscCall(MatGetFactorType(pc->pmat, &ftype));
27a2fc1e05SToby Isaac     if (ftype == MAT_FACTOR_NONE) {
289566063dSJacob Faibussowitsch       PetscCall(MatQRFactor(pc->pmat, dir->col, &((PC_Factor *)dir)->info));
299566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat, &err));
30a2fc1e05SToby Isaac       if (err) { /* Factor() fails */
31a2fc1e05SToby Isaac         pc->failedreason = (PCFailedReason)err;
32a2fc1e05SToby Isaac         PetscFunctionReturn(0);
33a2fc1e05SToby Isaac       }
34a2fc1e05SToby Isaac     }
35a2fc1e05SToby Isaac     ((PC_Factor *)dir)->fact = pc->pmat;
36a2fc1e05SToby Isaac   } else {
37a2fc1e05SToby Isaac     MatInfo info;
38a2fc1e05SToby Isaac 
39a2fc1e05SToby Isaac     if (!pc->setupcalled) {
404dfa11a4SJacob Faibussowitsch       if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_QR, &((PC_Factor *)dir)->fact)); }
419566063dSJacob Faibussowitsch       PetscCall(MatQRFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->col, &((PC_Factor *)dir)->info));
429566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
43a2fc1e05SToby Isaac       dir->hdr.actualfill = info.fill_ratio_needed;
44a2fc1e05SToby Isaac     } else if (pc->flag != SAME_NONZERO_PATTERN) {
459566063dSJacob Faibussowitsch       PetscCall(MatQRFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->col, &((PC_Factor *)dir)->info));
469566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
47a2fc1e05SToby Isaac       dir->hdr.actualfill = info.fill_ratio_needed;
48a2fc1e05SToby Isaac     } else {
499566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
50a2fc1e05SToby Isaac     }
519566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
52a2fc1e05SToby Isaac     if (err) { /* FactorSymbolic() fails */
53a2fc1e05SToby Isaac       pc->failedreason = (PCFailedReason)err;
54a2fc1e05SToby Isaac       PetscFunctionReturn(0);
55a2fc1e05SToby Isaac     }
56a2fc1e05SToby Isaac 
579566063dSJacob Faibussowitsch     PetscCall(MatQRFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
589566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
59a2fc1e05SToby Isaac     if (err) { /* FactorNumeric() fails */
60a2fc1e05SToby Isaac       pc->failedreason = (PCFailedReason)err;
61a2fc1e05SToby Isaac     }
62a2fc1e05SToby Isaac   }
63a2fc1e05SToby Isaac 
649566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
65a2fc1e05SToby Isaac   if (!stype) {
66a2fc1e05SToby Isaac     MatSolverType solverpackage;
679566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
689566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
69a2fc1e05SToby Isaac   }
70a2fc1e05SToby Isaac   PetscFunctionReturn(0);
71a2fc1e05SToby Isaac }
72a2fc1e05SToby Isaac 
73*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_QR(PC pc)
74*d71ae5a4SJacob Faibussowitsch {
75a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
76a2fc1e05SToby Isaac 
77a2fc1e05SToby Isaac   PetscFunctionBegin;
789566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
799566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
80a2fc1e05SToby Isaac   PetscFunctionReturn(0);
81a2fc1e05SToby Isaac }
82a2fc1e05SToby Isaac 
83*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_QR(PC pc)
84*d71ae5a4SJacob Faibussowitsch {
85a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
86a2fc1e05SToby Isaac 
87a2fc1e05SToby Isaac   PetscFunctionBegin;
889566063dSJacob Faibussowitsch   PetscCall(PCReset_QR(pc));
899566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
909566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
912e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
929566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
93a2fc1e05SToby Isaac   PetscFunctionReturn(0);
94a2fc1e05SToby Isaac }
95a2fc1e05SToby Isaac 
96*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_QR(PC pc, Vec x, Vec y)
97*d71ae5a4SJacob Faibussowitsch {
98a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
99a2fc1e05SToby Isaac   Mat    fact;
100a2fc1e05SToby Isaac 
101a2fc1e05SToby Isaac   PetscFunctionBegin;
102a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1039566063dSJacob Faibussowitsch   PetscCall(MatSolve(fact, x, y));
104a2fc1e05SToby Isaac   PetscFunctionReturn(0);
105a2fc1e05SToby Isaac }
106a2fc1e05SToby Isaac 
107*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_QR(PC pc, Mat X, Mat Y)
108*d71ae5a4SJacob Faibussowitsch {
109a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
110a2fc1e05SToby Isaac   Mat    fact;
111a2fc1e05SToby Isaac 
112a2fc1e05SToby Isaac   PetscFunctionBegin;
113a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1149566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(fact, X, Y));
115a2fc1e05SToby Isaac   PetscFunctionReturn(0);
116a2fc1e05SToby Isaac }
117a2fc1e05SToby Isaac 
118*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_QR(PC pc, Vec x, Vec y)
119*d71ae5a4SJacob Faibussowitsch {
120a2fc1e05SToby Isaac   PC_QR *dir = (PC_QR *)pc->data;
121a2fc1e05SToby Isaac   Mat    fact;
122a2fc1e05SToby Isaac 
123a2fc1e05SToby Isaac   PetscFunctionBegin;
124a2fc1e05SToby Isaac   fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact;
1259566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(fact, x, y));
126a2fc1e05SToby Isaac   PetscFunctionReturn(0);
127a2fc1e05SToby Isaac }
128a2fc1e05SToby Isaac 
129a2fc1e05SToby Isaac /*MC
130a2fc1e05SToby Isaac    PCQR - Uses a direct solver, based on QR factorization, as a preconditioner
131a2fc1e05SToby Isaac 
132a2fc1e05SToby Isaac    Level: beginner
133a2fc1e05SToby Isaac 
134f1580f4eSBarry Smith    Note:
135a2fc1e05SToby Isaac    Usually this will compute an "exact" solution in one iteration and does
136a2fc1e05SToby Isaac    not need a Krylov method (i.e. you can use -ksp_type preonly, or
137f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
138a2fc1e05SToby Isaac 
139f1580f4eSBarry Smith .seealso: `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()`
143db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
144a2fc1e05SToby Isaac M*/
145a2fc1e05SToby Isaac 
146*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_QR(PC pc)
147*d71ae5a4SJacob 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;
165a2fc1e05SToby Isaac   PetscFunctionReturn(0);
166a2fc1e05SToby Isaac }
167