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 9a2fc1e05SToby Isaac static PetscErrorCode PCSetUp_QR(PC pc) 10a2fc1e05SToby Isaac { 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) { 40a2fc1e05SToby Isaac if (!((PC_Factor*)dir)->fact) { 419566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_QR,&((PC_Factor*)dir)->fact)); 429566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact)); 43a2fc1e05SToby Isaac } 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 if (pc->flag != SAME_NONZERO_PATTERN) { 489566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->col,&((PC_Factor*)dir)->info)); 499566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info)); 50a2fc1e05SToby Isaac dir->hdr.actualfill = info.fill_ratio_needed; 51a2fc1e05SToby Isaac } else { 529566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err)); 53a2fc1e05SToby Isaac } 549566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err)); 55a2fc1e05SToby Isaac if (err) { /* FactorSymbolic() fails */ 56a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err; 57a2fc1e05SToby Isaac PetscFunctionReturn(0); 58a2fc1e05SToby Isaac } 59a2fc1e05SToby Isaac 609566063dSJacob Faibussowitsch PetscCall(MatQRFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info)); 619566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err)); 62a2fc1e05SToby Isaac if (err) { /* FactorNumeric() fails */ 63a2fc1e05SToby Isaac pc->failedreason = (PCFailedReason)err; 64a2fc1e05SToby Isaac } 65a2fc1e05SToby Isaac } 66a2fc1e05SToby Isaac 679566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,&stype)); 68a2fc1e05SToby Isaac if (!stype) { 69a2fc1e05SToby Isaac MatSolverType solverpackage; 709566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage)); 719566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solverpackage)); 72a2fc1e05SToby Isaac } 73a2fc1e05SToby Isaac PetscFunctionReturn(0); 74a2fc1e05SToby Isaac } 75a2fc1e05SToby Isaac 76a2fc1e05SToby Isaac static PetscErrorCode PCReset_QR(PC pc) 77a2fc1e05SToby Isaac { 78a2fc1e05SToby Isaac PC_QR *dir = (PC_QR*)pc->data; 79a2fc1e05SToby Isaac 80a2fc1e05SToby Isaac PetscFunctionBegin; 819566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) PetscCall(MatDestroy(&((PC_Factor*)dir)->fact)); 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 83a2fc1e05SToby Isaac PetscFunctionReturn(0); 84a2fc1e05SToby Isaac } 85a2fc1e05SToby Isaac 86a2fc1e05SToby Isaac static PetscErrorCode PCDestroy_QR(PC pc) 87a2fc1e05SToby Isaac { 88a2fc1e05SToby Isaac PC_QR *dir = (PC_QR*)pc->data; 89a2fc1e05SToby Isaac 90a2fc1e05SToby Isaac PetscFunctionBegin; 919566063dSJacob Faibussowitsch PetscCall(PCReset_QR(pc)); 929566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)dir)->ordering)); 939566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)dir)->solvertype)); 949566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 95a2fc1e05SToby Isaac PetscFunctionReturn(0); 96a2fc1e05SToby Isaac } 97a2fc1e05SToby Isaac 98a2fc1e05SToby Isaac static PetscErrorCode PCApply_QR(PC pc,Vec x,Vec y) 99a2fc1e05SToby Isaac { 100a2fc1e05SToby Isaac PC_QR *dir = (PC_QR*)pc->data; 101a2fc1e05SToby Isaac Mat fact; 102a2fc1e05SToby Isaac 103a2fc1e05SToby Isaac PetscFunctionBegin; 104a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor*)dir)->fact; 1059566063dSJacob Faibussowitsch PetscCall(MatSolve(fact,x,y)); 106a2fc1e05SToby Isaac PetscFunctionReturn(0); 107a2fc1e05SToby Isaac } 108a2fc1e05SToby Isaac 109a2fc1e05SToby Isaac static PetscErrorCode PCMatApply_QR(PC pc,Mat X,Mat Y) 110a2fc1e05SToby Isaac { 111a2fc1e05SToby Isaac PC_QR *dir = (PC_QR*)pc->data; 112a2fc1e05SToby Isaac Mat fact; 113a2fc1e05SToby Isaac 114a2fc1e05SToby Isaac PetscFunctionBegin; 115a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor*)dir)->fact; 1169566063dSJacob Faibussowitsch PetscCall(MatMatSolve(fact,X,Y)); 117a2fc1e05SToby Isaac PetscFunctionReturn(0); 118a2fc1e05SToby Isaac } 119a2fc1e05SToby Isaac 120a2fc1e05SToby Isaac static PetscErrorCode PCApplyTranspose_QR(PC pc,Vec x,Vec y) 121a2fc1e05SToby Isaac { 122a2fc1e05SToby Isaac PC_QR *dir = (PC_QR*)pc->data; 123a2fc1e05SToby Isaac Mat fact; 124a2fc1e05SToby Isaac 125a2fc1e05SToby Isaac PetscFunctionBegin; 126a2fc1e05SToby Isaac fact = dir->hdr.inplace ? pc->pmat : ((PC_Factor*)dir)->fact; 1279566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(fact,x,y)); 128a2fc1e05SToby Isaac PetscFunctionReturn(0); 129a2fc1e05SToby Isaac } 130a2fc1e05SToby Isaac 131a2fc1e05SToby Isaac /* -----------------------------------------------------------------------------------*/ 132a2fc1e05SToby Isaac 133a2fc1e05SToby Isaac /*MC 134a2fc1e05SToby Isaac PCQR - Uses a direct solver, based on QR factorization, as a preconditioner 135a2fc1e05SToby Isaac 136a2fc1e05SToby Isaac Level: beginner 137a2fc1e05SToby Isaac 138a2fc1e05SToby Isaac Notes: 139a2fc1e05SToby Isaac Usually this will compute an "exact" solution in one iteration and does 140a2fc1e05SToby Isaac not need a Krylov method (i.e. you can use -ksp_type preonly, or 141a2fc1e05SToby Isaac KSPSetType(ksp,KSPPREONLY) for the Krylov method 142a2fc1e05SToby Isaac 143db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 144db781477SPatrick Sanan `PCILU`, `PCLU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 145db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`, 146*c2e3fba1SPatrick Sanan `PCFactorSetPivotingInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 147db781477SPatrick Sanan `PCFactorReorderForNonzeroDiagonal()` 148a2fc1e05SToby Isaac M*/ 149a2fc1e05SToby Isaac 150a2fc1e05SToby Isaac PETSC_EXTERN PetscErrorCode PCCreate_QR(PC pc) 151a2fc1e05SToby Isaac { 152a2fc1e05SToby Isaac PC_QR *dir; 153a2fc1e05SToby Isaac 154a2fc1e05SToby Isaac PetscFunctionBegin; 1559566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&dir)); 156a2fc1e05SToby Isaac pc->data = (void*)dir; 1579566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_QR)); 158a2fc1e05SToby Isaac 159a2fc1e05SToby Isaac dir->col = NULL; 160a2fc1e05SToby Isaac pc->ops->reset = PCReset_QR; 161a2fc1e05SToby Isaac pc->ops->destroy = PCDestroy_QR; 162a2fc1e05SToby Isaac pc->ops->apply = PCApply_QR; 163a2fc1e05SToby Isaac pc->ops->matapply = PCMatApply_QR; 164a2fc1e05SToby Isaac pc->ops->applytranspose = PCApplyTranspose_QR; 165a2fc1e05SToby Isaac pc->ops->setup = PCSetUp_QR; 166a2fc1e05SToby Isaac pc->ops->setfromoptions = PCSetFromOptions_Factor; 167a2fc1e05SToby Isaac pc->ops->view = PCView_Factor; 168a2fc1e05SToby Isaac pc->ops->applyrichardson = NULL; 169a2fc1e05SToby Isaac PetscFunctionReturn(0); 170a2fc1e05SToby Isaac } 171