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) { 39a2fc1e05SToby Isaac if (!((PC_Factor *)dir)->fact) { 409566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_QR, &((PC_Factor *)dir)->fact)); 419566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)dir)->fact)); 42a2fc1e05SToby Isaac } 439566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->col, &((PC_Factor *)dir)->info)); 449566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 45a2fc1e05SToby Isaac dir->hdr.actualfill = info.fill_ratio_needed; 46a2fc1e05SToby Isaac } else if (pc->flag != SAME_NONZERO_PATTERN) { 479566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->col, &((PC_Factor *)dir)->info)); 489566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 49a2fc1e05SToby Isaac dir->hdr.actualfill = info.fill_ratio_needed; 50a2fc1e05SToby Isaac } else { 519566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 52a2fc1e05SToby Isaac } 539566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 54a2fc1e05SToby Isaac if (err) { /* FactorSymbolic() fails */ 55a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err; 56a2fc1e05SToby Isaac PetscFunctionReturn(0); 57a2fc1e05SToby Isaac } 58a2fc1e05SToby Isaac 599566063dSJacob Faibussowitsch PetscCall(MatQRFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 609566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 61a2fc1e05SToby Isaac if (err) { /* FactorNumeric() fails */ 62a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err; 63a2fc1e05SToby Isaac } 64a2fc1e05SToby Isaac } 65a2fc1e05SToby Isaac 669566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 67a2fc1e05SToby Isaac if (!stype) { 68a2fc1e05SToby Isaac MatSolverType solverpackage; 699566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 709566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 71a2fc1e05SToby Isaac } 72a2fc1e05SToby Isaac PetscFunctionReturn(0); 73a2fc1e05SToby Isaac } 74a2fc1e05SToby Isaac 759371c9d4SSatish Balay static PetscErrorCode PCReset_QR(PC pc) { 76a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data; 77a2fc1e05SToby Isaac 78a2fc1e05SToby Isaac PetscFunctionBegin; 799566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 81a2fc1e05SToby Isaac PetscFunctionReturn(0); 82a2fc1e05SToby Isaac } 83a2fc1e05SToby Isaac 849371c9d4SSatish Balay static PetscErrorCode PCDestroy_QR(PC pc) { 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 969371c9d4SSatish Balay static PetscErrorCode PCApply_QR(PC pc, Vec x, Vec y) { 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)); 103a2fc1e05SToby Isaac PetscFunctionReturn(0); 104a2fc1e05SToby Isaac } 105a2fc1e05SToby Isaac 1069371c9d4SSatish Balay static PetscErrorCode PCMatApply_QR(PC pc, Mat X, Mat Y) { 107a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data; 108a2fc1e05SToby Isaac Mat fact; 109a2fc1e05SToby Isaac 110a2fc1e05SToby Isaac PetscFunctionBegin; 111a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact; 1129566063dSJacob Faibussowitsch PetscCall(MatMatSolve(fact, X, Y)); 113a2fc1e05SToby Isaac PetscFunctionReturn(0); 114a2fc1e05SToby Isaac } 115a2fc1e05SToby Isaac 1169371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_QR(PC pc, Vec x, Vec y) { 117a2fc1e05SToby Isaac PC_QR *dir = (PC_QR *)pc->data; 118a2fc1e05SToby Isaac Mat fact; 119a2fc1e05SToby Isaac 120a2fc1e05SToby Isaac PetscFunctionBegin; 121a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor *)dir)->fact; 1229566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(fact, x, y)); 123a2fc1e05SToby Isaac PetscFunctionReturn(0); 124a2fc1e05SToby Isaac } 125a2fc1e05SToby Isaac 126a2fc1e05SToby Isaac /*MC 127a2fc1e05SToby Isaac PCQR - Uses a direct solver, based on QR factorization, as a preconditioner 128a2fc1e05SToby Isaac 129a2fc1e05SToby Isaac Level: beginner 130a2fc1e05SToby Isaac 131*f1580f4eSBarry Smith Note: 132a2fc1e05SToby Isaac Usually this will compute an "exact" solution in one iteration and does 133a2fc1e05SToby Isaac not need a Krylov method (i.e. you can use -ksp_type preonly, or 134*f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 135a2fc1e05SToby Isaac 136*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSVD`, 137db781477SPatrick Sanan `PCILU`, `PCLU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 138db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`, 139c2e3fba1SPatrick Sanan `PCFactorSetPivotingInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 140db781477SPatrick Sanan `PCFactorReorderForNonzeroDiagonal()` 141a2fc1e05SToby Isaac M*/ 142a2fc1e05SToby Isaac 1439371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_QR(PC pc) { 144a2fc1e05SToby Isaac PC_QR *dir; 145a2fc1e05SToby Isaac 146a2fc1e05SToby Isaac PetscFunctionBegin; 1479566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &dir)); 148a2fc1e05SToby Isaac pc->data = (void *)dir; 1499566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_QR)); 150a2fc1e05SToby Isaac 151a2fc1e05SToby Isaac dir->col = NULL; 152a2fc1e05SToby Isaac pc->ops->reset = PCReset_QR; 153a2fc1e05SToby Isaac pc->ops->destroy = PCDestroy_QR; 154a2fc1e05SToby Isaac pc->ops->apply = PCApply_QR; 155a2fc1e05SToby Isaac pc->ops->matapply = PCMatApply_QR; 156a2fc1e05SToby Isaac pc->ops->applytranspose = PCApplyTranspose_QR; 157a2fc1e05SToby Isaac pc->ops->setup = PCSetUp_QR; 158a2fc1e05SToby Isaac pc->ops->setfromoptions = PCSetFromOptions_Factor; 159a2fc1e05SToby Isaac pc->ops->view = PCView_Factor; 160a2fc1e05SToby Isaac pc->ops->applyrichardson = NULL; 161a2fc1e05SToby Isaac PetscFunctionReturn(0); 162a2fc1e05SToby Isaac } 163