1dba47a55SKris Buschelman 29b54502bSHong Zhang /* 39b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 49b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 59b54502bSHong Zhang a direct solver. 69b54502bSHong Zhang */ 7ee45ca4aSHong Zhang 8c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/ 99b54502bSHong Zhang 10680c5173SHong Zhang 117087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 129b54502bSHong Zhang { 139b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 149b54502bSHong Zhang 159b54502bSHong Zhang PetscFunctionBegin; 169b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 172fa5cd67SKarl Rupp if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 182fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 199b54502bSHong Zhang PetscFunctionReturn(0); 209b54502bSHong Zhang } 219b54502bSHong Zhang 224416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc) 239b54502bSHong Zhang { 249b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 259b54502bSHong Zhang PetscErrorCode ierr; 26ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 279b54502bSHong Zhang PetscReal tol; 289b54502bSHong Zhang 299b54502bSHong Zhang PetscFunctionBegin; 30e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr); 31e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 325c9eb25fSBarry Smith 332401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 349b54502bSHong Zhang if (flg) { 359b54502bSHong Zhang tol = PETSC_DECIDE; 362401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 372401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 389b54502bSHong Zhang } 399b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 409b54502bSHong Zhang PetscFunctionReturn(0); 419b54502bSHong Zhang } 429b54502bSHong Zhang 439b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 449b54502bSHong Zhang { 459b54502bSHong Zhang PetscErrorCode ierr; 469b54502bSHong Zhang 479b54502bSHong Zhang PetscFunctionBegin; 48914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 499b54502bSHong Zhang PetscFunctionReturn(0); 509b54502bSHong Zhang } 519b54502bSHong Zhang 529b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 539b54502bSHong Zhang { 549b54502bSHong Zhang PetscErrorCode ierr; 559b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 56ea799195SBarry Smith MatSolverType stype; 5700e125f8SBarry Smith MatFactorError err; 583d1c1ea0SBarry Smith 599b54502bSHong Zhang PetscFunctionBegin; 60c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 613d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill; 629b54502bSHong Zhang 6384d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 643d1c1ea0SBarry Smith if (dir->hdr.inplace) { 65*9d76b4d0SMatthew G. Knepley MatFactorType ftype; 66*9d76b4d0SMatthew G. Knepley 67*9d76b4d0SMatthew G. Knepley ierr = MatGetFactorType(pc->pmat, &ftype);CHKERRQ(ierr); 68*9d76b4d0SMatthew G. Knepley if (ftype == MAT_FACTOR_NONE) { 69fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 70fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 71075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 7203c60df9SBarry Smith if (dir->row) { 733bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 743bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 7503c60df9SBarry Smith } 76075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 7700e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 7800e125f8SBarry Smith if (err) { /* Factor() fails */ 7900e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 806baea169SHong Zhang PetscFunctionReturn(0); 816baea169SHong Zhang } 82*9d76b4d0SMatthew G. Knepley } 83075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 849b54502bSHong Zhang } else { 859b54502bSHong Zhang MatInfo info; 8600e125f8SBarry Smith 879b54502bSHong Zhang if (!pc->setupcalled) { 88075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 899b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 909b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 919b54502bSHong Zhang } 9203c60df9SBarry Smith if (dir->row) { 933bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 943bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 9503c60df9SBarry Smith } 96d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 97075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 98a1f19f5aSHong Zhang } 99075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 100075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1013d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 1023bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1039b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1043d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 105fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 106fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 107075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1089b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1099b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1109b54502bSHong Zhang } 11103c60df9SBarry Smith if (dir->row) { 1123bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 1133bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 11403c60df9SBarry Smith } 1159b54502bSHong Zhang } 1166bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 117075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 118075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 119075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1203d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 1213bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 12204545d6dSBarry Smith } else { 123b8b68cfdSBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 124160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 125b8b68cfdSBarry Smith ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 126b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 12704545d6dSBarry Smith } 1289b54502bSHong Zhang } 12900e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 13000e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 13100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1328c1cd74cSHong Zhang PetscFunctionReturn(0); 1338c1cd74cSHong Zhang } 1348c1cd74cSHong Zhang 135075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 13600e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 13700e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 13800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1398c1cd74cSHong Zhang } 140680c5173SHong Zhang 1419b54502bSHong Zhang } 14200c67f3bSHong Zhang 1433ca39a21SBarry Smith ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 14400c67f3bSHong Zhang if (!stype) { 145ea799195SBarry Smith MatSolverType solverpackage; 1463ca39a21SBarry Smith ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 1473ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 14800c67f3bSHong Zhang } 1499b54502bSHong Zhang PetscFunctionReturn(0); 1509b54502bSHong Zhang } 1519b54502bSHong Zhang 152574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1539b54502bSHong Zhang { 1549b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1559b54502bSHong Zhang PetscErrorCode ierr; 1569b54502bSHong Zhang 1579b54502bSHong Zhang PetscFunctionBegin; 1583d1c1ea0SBarry Smith if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 159fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 160fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 161574deadeSBarry Smith PetscFunctionReturn(0); 162574deadeSBarry Smith } 163574deadeSBarry Smith 164574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 165574deadeSBarry Smith { 166574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 167574deadeSBarry Smith PetscErrorCode ierr; 168574deadeSBarry Smith 169574deadeSBarry Smith PetscFunctionBegin; 170574deadeSBarry Smith ierr = PCReset_LU(pc);CHKERRQ(ierr); 171503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 172503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 173c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1749b54502bSHong Zhang PetscFunctionReturn(0); 1759b54502bSHong Zhang } 1769b54502bSHong Zhang 1779b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 1789b54502bSHong Zhang { 1799b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1809b54502bSHong Zhang PetscErrorCode ierr; 1819b54502bSHong Zhang 1829b54502bSHong Zhang PetscFunctionBegin; 1833d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1842fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 1852fa5cd67SKarl Rupp } else { 1862fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1872fa5cd67SKarl Rupp } 1889b54502bSHong Zhang PetscFunctionReturn(0); 1899b54502bSHong Zhang } 1909b54502bSHong Zhang 1919b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 1929b54502bSHong Zhang { 1939b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1949b54502bSHong Zhang PetscErrorCode ierr; 1959b54502bSHong Zhang 1969b54502bSHong Zhang PetscFunctionBegin; 1973d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1982fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 1992fa5cd67SKarl Rupp } else { 2002fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2012fa5cd67SKarl Rupp } 2029b54502bSHong Zhang PetscFunctionReturn(0); 2039b54502bSHong Zhang } 2049b54502bSHong Zhang 2059b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2069b54502bSHong Zhang 2079b54502bSHong Zhang /*MC 2089b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2099b54502bSHong Zhang 2109b54502bSHong Zhang Options Database Keys: 2112401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2123ca39a21SBarry Smith . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 2132401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 21455ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2152401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2162401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2172401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2189b54502bSHong Zhang stability of factorization. 219145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 220145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 221e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 222e22d95b2SBarry Smith diagonal. 2239b54502bSHong Zhang 22495452b02SPatrick Sanan Notes: 22595452b02SPatrick Sanan Not all options work for all matrix formats 2269b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2279b54502bSHong Zhang algorithms 2289b54502bSHong Zhang 2299b54502bSHong Zhang Level: beginner 2309b54502bSHong Zhang 2319b54502bSHong Zhang Concepts: LU factorization, direct solver 2329b54502bSHong Zhang 23395452b02SPatrick Sanan Notes: 23495452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2359b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2369b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2379b54502bSHong Zhang 2389b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 239a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2408ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 241145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2428ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2439b54502bSHong Zhang M*/ 2449b54502bSHong Zhang 2458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 2469b54502bSHong Zhang { 2479b54502bSHong Zhang PetscErrorCode ierr; 2489b54502bSHong Zhang PetscMPIInt size; 2499b54502bSHong Zhang PC_LU *dir; 2509b54502bSHong Zhang 2519b54502bSHong Zhang PetscFunctionBegin; 252b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 2533d1c1ea0SBarry Smith pc->data = (void*)dir; 2543d1c1ea0SBarry Smith ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 255d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 2569b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2579b54502bSHong Zhang 258075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 259075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 260f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 2619b54502bSHong Zhang dir->col = 0; 2629b54502bSHong Zhang dir->row = 0; 2635c9eb25fSBarry Smith 264ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 2659b54502bSHong Zhang if (size == 1) { 26619fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 2679b54502bSHong Zhang } else { 26819fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 2699b54502bSHong Zhang } 2709b54502bSHong Zhang 271574deadeSBarry Smith pc->ops->reset = PCReset_LU; 2729b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 2739b54502bSHong Zhang pc->ops->apply = PCApply_LU; 2749b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 2759b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 2769b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 2779b54502bSHong Zhang pc->ops->view = PCView_LU; 2789b54502bSHong Zhang pc->ops->applyrichardson = 0; 279bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 2809b54502bSHong Zhang PetscFunctionReturn(0); 2819b54502bSHong Zhang } 282