1*9b54502bSHong Zhang /* 2*9b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 3*9b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 4*9b54502bSHong Zhang a direct solver. 5*9b54502bSHong Zhang */ 6*9b54502bSHong Zhang #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 7*9b54502bSHong Zhang 8*9b54502bSHong Zhang typedef struct { 9*9b54502bSHong Zhang Mat fact; /* factored matrix */ 10*9b54502bSHong Zhang PetscReal actualfill; /* actual fill in factor */ 11*9b54502bSHong Zhang PetscTruth inplace; /* flag indicating in-place factorization */ 12*9b54502bSHong Zhang IS row,col; /* index sets used for reordering */ 13*9b54502bSHong Zhang MatOrderingType ordering; /* matrix ordering */ 14*9b54502bSHong Zhang PetscTruth reuseordering; /* reuses previous reordering computed */ 15*9b54502bSHong Zhang PetscTruth reusefill; /* reuse fill from previous LU */ 16*9b54502bSHong Zhang MatFactorInfo info; 17*9b54502bSHong Zhang PetscTruth nonzerosalongdiagonal; 18*9b54502bSHong Zhang PetscReal nonzerosalongdiagonaltol; 19*9b54502bSHong Zhang } PC_LU; 20*9b54502bSHong Zhang 21*9b54502bSHong Zhang 22*9b54502bSHong Zhang EXTERN_C_BEGIN 23*9b54502bSHong Zhang #undef __FUNCT__ 24*9b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal_LU" 25*9b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 26*9b54502bSHong Zhang { 27*9b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 28*9b54502bSHong Zhang 29*9b54502bSHong Zhang PetscFunctionBegin; 30*9b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 31*9b54502bSHong Zhang if (z == PETSC_DECIDE) { 32*9b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 33*9b54502bSHong Zhang } else { 34*9b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 35*9b54502bSHong Zhang } 36*9b54502bSHong Zhang PetscFunctionReturn(0); 37*9b54502bSHong Zhang } 38*9b54502bSHong Zhang EXTERN_C_END 39*9b54502bSHong Zhang 40*9b54502bSHong Zhang EXTERN_C_BEGIN 41*9b54502bSHong Zhang #undef __FUNCT__ 42*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetSetZeroPivot" 43*9b54502bSHong Zhang PetscErrorCode PCLUSetZeroPivot_LU(PC pc,PetscReal z) 44*9b54502bSHong Zhang { 45*9b54502bSHong Zhang PC_LU *lu; 46*9b54502bSHong Zhang 47*9b54502bSHong Zhang PetscFunctionBegin; 48*9b54502bSHong Zhang lu = (PC_LU*)pc->data; 49*9b54502bSHong Zhang lu->info.zeropivot = z; 50*9b54502bSHong Zhang PetscFunctionReturn(0); 51*9b54502bSHong Zhang } 52*9b54502bSHong Zhang EXTERN_C_END 53*9b54502bSHong Zhang 54*9b54502bSHong Zhang EXTERN_C_BEGIN 55*9b54502bSHong Zhang #undef __FUNCT__ 56*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering_LU" 57*9b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag) 58*9b54502bSHong Zhang { 59*9b54502bSHong Zhang PC_LU *lu; 60*9b54502bSHong Zhang 61*9b54502bSHong Zhang PetscFunctionBegin; 62*9b54502bSHong Zhang lu = (PC_LU*)pc->data; 63*9b54502bSHong Zhang lu->reuseordering = flag; 64*9b54502bSHong Zhang PetscFunctionReturn(0); 65*9b54502bSHong Zhang } 66*9b54502bSHong Zhang EXTERN_C_END 67*9b54502bSHong Zhang 68*9b54502bSHong Zhang EXTERN_C_BEGIN 69*9b54502bSHong Zhang #undef __FUNCT__ 70*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill_LU" 71*9b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill_LU(PC pc,PetscTruth flag) 72*9b54502bSHong Zhang { 73*9b54502bSHong Zhang PC_LU *lu; 74*9b54502bSHong Zhang 75*9b54502bSHong Zhang PetscFunctionBegin; 76*9b54502bSHong Zhang lu = (PC_LU*)pc->data; 77*9b54502bSHong Zhang lu->reusefill = flag; 78*9b54502bSHong Zhang PetscFunctionReturn(0); 79*9b54502bSHong Zhang } 80*9b54502bSHong Zhang EXTERN_C_END 81*9b54502bSHong Zhang 82*9b54502bSHong Zhang EXTERN_C_BEGIN 83*9b54502bSHong Zhang #undef __FUNCT__ 84*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetShift_LU" 85*9b54502bSHong Zhang PetscErrorCode PCLUSetShift_LU(PC pc,PetscTruth shift) 86*9b54502bSHong Zhang { 87*9b54502bSHong Zhang PC_LU *dir; 88*9b54502bSHong Zhang 89*9b54502bSHong Zhang PetscFunctionBegin; 90*9b54502bSHong Zhang dir = (PC_LU*)pc->data; 91*9b54502bSHong Zhang dir->info.shift = shift; 92*9b54502bSHong Zhang if (shift) dir->info.shift_fraction = 0.0; 93*9b54502bSHong Zhang PetscFunctionReturn(0); 94*9b54502bSHong Zhang } 95*9b54502bSHong Zhang EXTERN_C_END 96*9b54502bSHong Zhang 97*9b54502bSHong Zhang #undef __FUNCT__ 98*9b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 99*9b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 100*9b54502bSHong Zhang { 101*9b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 102*9b54502bSHong Zhang PetscErrorCode ierr; 103*9b54502bSHong Zhang PetscTruth flg,set; 104*9b54502bSHong Zhang char tname[256]; 105*9b54502bSHong Zhang PetscFList ordlist; 106*9b54502bSHong Zhang PetscReal tol; 107*9b54502bSHong Zhang 108*9b54502bSHong Zhang PetscFunctionBegin; 109*9b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 110*9b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 111*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);CHKERRQ(ierr); 112*9b54502bSHong Zhang if (flg) { 113*9b54502bSHong Zhang ierr = PCLUSetUseInPlace(pc);CHKERRQ(ierr); 114*9b54502bSHong Zhang } 115*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 116*9b54502bSHong Zhang 117*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",&flg);CHKERRQ(ierr); 118*9b54502bSHong Zhang if (flg) { 119*9b54502bSHong Zhang ierr = PCLUSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 120*9b54502bSHong Zhang } 121*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);CHKERRQ(ierr); 122*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_shift","Manteuffel shift applied to diagonal","PCLUSetShift",&flg);CHKERRQ(ierr); 123*9b54502bSHong Zhang if (flg) { 124*9b54502bSHong Zhang ierr = PCLUSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 125*9b54502bSHong Zhang } 126*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 127*9b54502bSHong Zhang 128*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);CHKERRQ(ierr); 129*9b54502bSHong Zhang if (flg) { 130*9b54502bSHong Zhang ierr = PCLUSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 131*9b54502bSHong Zhang } 132*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);CHKERRQ(ierr); 133*9b54502bSHong Zhang if (flg) { 134*9b54502bSHong Zhang ierr = PCLUSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 135*9b54502bSHong Zhang } 136*9b54502bSHong Zhang 137*9b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 138*9b54502bSHong Zhang ierr = PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 139*9b54502bSHong Zhang if (flg) { 140*9b54502bSHong Zhang ierr = PCLUSetMatOrdering(pc,tname);CHKERRQ(ierr); 141*9b54502bSHong Zhang } 142*9b54502bSHong Zhang 143*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 144*9b54502bSHong Zhang if (flg) { 145*9b54502bSHong Zhang tol = PETSC_DECIDE; 146*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 147*9b54502bSHong Zhang ierr = PCLUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 148*9b54502bSHong Zhang } 149*9b54502bSHong Zhang 150*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 151*9b54502bSHong Zhang 152*9b54502bSHong Zhang flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 153*9b54502bSHong Zhang ierr = PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 154*9b54502bSHong Zhang if (set) { 155*9b54502bSHong Zhang ierr = PCLUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 156*9b54502bSHong Zhang } 157*9b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 158*9b54502bSHong Zhang PetscFunctionReturn(0); 159*9b54502bSHong Zhang } 160*9b54502bSHong Zhang 161*9b54502bSHong Zhang #undef __FUNCT__ 162*9b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 163*9b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 164*9b54502bSHong Zhang { 165*9b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 166*9b54502bSHong Zhang PetscErrorCode ierr; 167*9b54502bSHong Zhang PetscTruth iascii,isstring; 168*9b54502bSHong Zhang 169*9b54502bSHong Zhang PetscFunctionBegin; 170*9b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 171*9b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 172*9b54502bSHong Zhang if (iascii) { 173*9b54502bSHong Zhang MatInfo info; 174*9b54502bSHong Zhang 175*9b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 176*9b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 177*9b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 178*9b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %g\n",lu->info.zeropivot);CHKERRQ(ierr); 179*9b54502bSHong Zhang if (lu->info.shift) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 180*9b54502bSHong Zhang if (lu->fact) { 181*9b54502bSHong Zhang ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 182*9b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU nonzeros %g\n",info.nz_used);CHKERRQ(ierr); 183*9b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr); 184*9b54502bSHong Zhang ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 185*9b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 186*9b54502bSHong Zhang } 187*9b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 188*9b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 189*9b54502bSHong Zhang } else if (isstring) { 190*9b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 191*9b54502bSHong Zhang } else { 192*9b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 193*9b54502bSHong Zhang } 194*9b54502bSHong Zhang PetscFunctionReturn(0); 195*9b54502bSHong Zhang } 196*9b54502bSHong Zhang 197*9b54502bSHong Zhang #undef __FUNCT__ 198*9b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_LU" 199*9b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat) 200*9b54502bSHong Zhang { 201*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 202*9b54502bSHong Zhang 203*9b54502bSHong Zhang PetscFunctionBegin; 204*9b54502bSHong Zhang if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 205*9b54502bSHong Zhang *mat = dir->fact; 206*9b54502bSHong Zhang PetscFunctionReturn(0); 207*9b54502bSHong Zhang } 208*9b54502bSHong Zhang 209*9b54502bSHong Zhang #undef __FUNCT__ 210*9b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 211*9b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 212*9b54502bSHong Zhang { 213*9b54502bSHong Zhang PetscErrorCode ierr; 214*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 215*9b54502bSHong Zhang 216*9b54502bSHong Zhang PetscFunctionBegin; 217*9b54502bSHong Zhang if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 218*9b54502bSHong Zhang 219*9b54502bSHong Zhang if (dir->inplace) { 220*9b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 221*9b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 222*9b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 223*9b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 224*9b54502bSHong Zhang ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 225*9b54502bSHong Zhang dir->fact = pc->pmat; 226*9b54502bSHong Zhang } else { 227*9b54502bSHong Zhang MatInfo info; 228*9b54502bSHong Zhang if (!pc->setupcalled) { 229*9b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 230*9b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 231*9b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 232*9b54502bSHong Zhang } 233*9b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 234*9b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 235*9b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 236*9b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 237*9b54502bSHong Zhang PetscLogObjectParent(pc,dir->fact); 238*9b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 239*9b54502bSHong Zhang if (!dir->reuseordering) { 240*9b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 241*9b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 242*9b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 243*9b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 244*9b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 245*9b54502bSHong Zhang } 246*9b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 247*9b54502bSHong Zhang } 248*9b54502bSHong Zhang ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 249*9b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 250*9b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 251*9b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 252*9b54502bSHong Zhang PetscLogObjectParent(pc,dir->fact); 253*9b54502bSHong Zhang } 254*9b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 255*9b54502bSHong Zhang } 256*9b54502bSHong Zhang PetscFunctionReturn(0); 257*9b54502bSHong Zhang } 258*9b54502bSHong Zhang 259*9b54502bSHong Zhang #undef __FUNCT__ 260*9b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 261*9b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 262*9b54502bSHong Zhang { 263*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 264*9b54502bSHong Zhang PetscErrorCode ierr; 265*9b54502bSHong Zhang 266*9b54502bSHong Zhang PetscFunctionBegin; 267*9b54502bSHong Zhang if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 268*9b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 269*9b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 270*9b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 271*9b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 272*9b54502bSHong Zhang PetscFunctionReturn(0); 273*9b54502bSHong Zhang } 274*9b54502bSHong Zhang 275*9b54502bSHong Zhang #undef __FUNCT__ 276*9b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 277*9b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 278*9b54502bSHong Zhang { 279*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 280*9b54502bSHong Zhang PetscErrorCode ierr; 281*9b54502bSHong Zhang 282*9b54502bSHong Zhang PetscFunctionBegin; 283*9b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 284*9b54502bSHong Zhang else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 285*9b54502bSHong Zhang PetscFunctionReturn(0); 286*9b54502bSHong Zhang } 287*9b54502bSHong Zhang 288*9b54502bSHong Zhang #undef __FUNCT__ 289*9b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 290*9b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 291*9b54502bSHong Zhang { 292*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 293*9b54502bSHong Zhang PetscErrorCode ierr; 294*9b54502bSHong Zhang 295*9b54502bSHong Zhang PetscFunctionBegin; 296*9b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 297*9b54502bSHong Zhang else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 298*9b54502bSHong Zhang PetscFunctionReturn(0); 299*9b54502bSHong Zhang } 300*9b54502bSHong Zhang 301*9b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 302*9b54502bSHong Zhang 303*9b54502bSHong Zhang EXTERN_C_BEGIN 304*9b54502bSHong Zhang #undef __FUNCT__ 305*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill_LU" 306*9b54502bSHong Zhang PetscErrorCode PCLUSetFill_LU(PC pc,PetscReal fill) 307*9b54502bSHong Zhang { 308*9b54502bSHong Zhang PC_LU *dir; 309*9b54502bSHong Zhang 310*9b54502bSHong Zhang PetscFunctionBegin; 311*9b54502bSHong Zhang dir = (PC_LU*)pc->data; 312*9b54502bSHong Zhang dir->info.fill = fill; 313*9b54502bSHong Zhang PetscFunctionReturn(0); 314*9b54502bSHong Zhang } 315*9b54502bSHong Zhang EXTERN_C_END 316*9b54502bSHong Zhang 317*9b54502bSHong Zhang EXTERN_C_BEGIN 318*9b54502bSHong Zhang #undef __FUNCT__ 319*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetDamping_LU" 320*9b54502bSHong Zhang PetscErrorCode PCLUSetDamping_LU(PC pc,PetscReal damping) 321*9b54502bSHong Zhang { 322*9b54502bSHong Zhang PC_LU *dir; 323*9b54502bSHong Zhang 324*9b54502bSHong Zhang PetscFunctionBegin; 325*9b54502bSHong Zhang dir = (PC_LU*)pc->data; 326*9b54502bSHong Zhang if (damping == (PetscReal) PETSC_DECIDE) { 327*9b54502bSHong Zhang dir->info.damping = 1.e-12; 328*9b54502bSHong Zhang } else { 329*9b54502bSHong Zhang dir->info.damping = damping; 330*9b54502bSHong Zhang } 331*9b54502bSHong Zhang PetscFunctionReturn(0); 332*9b54502bSHong Zhang } 333*9b54502bSHong Zhang EXTERN_C_END 334*9b54502bSHong Zhang 335*9b54502bSHong Zhang EXTERN_C_BEGIN 336*9b54502bSHong Zhang #undef __FUNCT__ 337*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace_LU" 338*9b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace_LU(PC pc) 339*9b54502bSHong Zhang { 340*9b54502bSHong Zhang PC_LU *dir; 341*9b54502bSHong Zhang 342*9b54502bSHong Zhang PetscFunctionBegin; 343*9b54502bSHong Zhang dir = (PC_LU*)pc->data; 344*9b54502bSHong Zhang dir->inplace = PETSC_TRUE; 345*9b54502bSHong Zhang PetscFunctionReturn(0); 346*9b54502bSHong Zhang } 347*9b54502bSHong Zhang EXTERN_C_END 348*9b54502bSHong Zhang 349*9b54502bSHong Zhang EXTERN_C_BEGIN 350*9b54502bSHong Zhang #undef __FUNCT__ 351*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering_LU" 352*9b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering) 353*9b54502bSHong Zhang { 354*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 355*9b54502bSHong Zhang PetscErrorCode ierr; 356*9b54502bSHong Zhang 357*9b54502bSHong Zhang PetscFunctionBegin; 358*9b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 359*9b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 360*9b54502bSHong Zhang PetscFunctionReturn(0); 361*9b54502bSHong Zhang } 362*9b54502bSHong Zhang EXTERN_C_END 363*9b54502bSHong Zhang 364*9b54502bSHong Zhang EXTERN_C_BEGIN 365*9b54502bSHong Zhang #undef __FUNCT__ 366*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting_LU" 367*9b54502bSHong Zhang PetscErrorCode PCLUSetPivoting_LU(PC pc,PetscReal dtcol) 368*9b54502bSHong Zhang { 369*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 370*9b54502bSHong Zhang 371*9b54502bSHong Zhang PetscFunctionBegin; 372*9b54502bSHong Zhang if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",dtcol); 373*9b54502bSHong Zhang dir->info.dtcol = dtcol; 374*9b54502bSHong Zhang PetscFunctionReturn(0); 375*9b54502bSHong Zhang } 376*9b54502bSHong Zhang EXTERN_C_END 377*9b54502bSHong Zhang 378*9b54502bSHong Zhang EXTERN_C_BEGIN 379*9b54502bSHong Zhang #undef __FUNCT__ 380*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks_LU" 381*9b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 382*9b54502bSHong Zhang { 383*9b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 384*9b54502bSHong Zhang 385*9b54502bSHong Zhang PetscFunctionBegin; 386*9b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 387*9b54502bSHong Zhang PetscFunctionReturn(0); 388*9b54502bSHong Zhang } 389*9b54502bSHong Zhang EXTERN_C_END 390*9b54502bSHong Zhang 391*9b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 392*9b54502bSHong Zhang 393*9b54502bSHong Zhang #undef __FUNCT__ 394*9b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal" 395*9b54502bSHong Zhang /*@ 396*9b54502bSHong Zhang PCLUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 397*9b54502bSHong Zhang 398*9b54502bSHong Zhang Collective on PC 399*9b54502bSHong Zhang 400*9b54502bSHong Zhang Input Parameters: 401*9b54502bSHong Zhang + pc - the preconditioner context 402*9b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 403*9b54502bSHong Zhang 404*9b54502bSHong Zhang Options Database Key: 405*9b54502bSHong Zhang . -pc_lu_nonzeros_along_diagonal 406*9b54502bSHong Zhang 407*9b54502bSHong Zhang Level: intermediate 408*9b54502bSHong Zhang 409*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 410*9b54502bSHong Zhang 411*9b54502bSHong Zhang .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot(), MatReorderForNonzeroDiagonal() 412*9b54502bSHong Zhang @*/ 413*9b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 414*9b54502bSHong Zhang { 415*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 416*9b54502bSHong Zhang 417*9b54502bSHong Zhang PetscFunctionBegin; 418*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 419*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 420*9b54502bSHong Zhang if (f) { 421*9b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 422*9b54502bSHong Zhang } 423*9b54502bSHong Zhang PetscFunctionReturn(0); 424*9b54502bSHong Zhang } 425*9b54502bSHong Zhang 426*9b54502bSHong Zhang #undef __FUNCT__ 427*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetZeroPivot" 428*9b54502bSHong Zhang /*@ 429*9b54502bSHong Zhang PCLUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 430*9b54502bSHong Zhang 431*9b54502bSHong Zhang Collective on PC 432*9b54502bSHong Zhang 433*9b54502bSHong Zhang Input Parameters: 434*9b54502bSHong Zhang + pc - the preconditioner context 435*9b54502bSHong Zhang - zero - all pivots smaller than this will be considered zero 436*9b54502bSHong Zhang 437*9b54502bSHong Zhang Options Database Key: 438*9b54502bSHong Zhang . -pc_lu_zeropivot <zero> - Sets the zero pivot size 439*9b54502bSHong Zhang 440*9b54502bSHong Zhang Level: intermediate 441*9b54502bSHong Zhang 442*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 443*9b54502bSHong Zhang 444*9b54502bSHong Zhang .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot() 445*9b54502bSHong Zhang @*/ 446*9b54502bSHong Zhang PetscErrorCode PCLUSetZeroPivot(PC pc,PetscReal zero) 447*9b54502bSHong Zhang { 448*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 449*9b54502bSHong Zhang 450*9b54502bSHong Zhang PetscFunctionBegin; 451*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 452*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 453*9b54502bSHong Zhang if (f) { 454*9b54502bSHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 455*9b54502bSHong Zhang } 456*9b54502bSHong Zhang PetscFunctionReturn(0); 457*9b54502bSHong Zhang } 458*9b54502bSHong Zhang 459*9b54502bSHong Zhang #undef __FUNCT__ 460*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetShift" 461*9b54502bSHong Zhang /*@ 462*9b54502bSHong Zhang PCLUSetShift - specify whether to use Manteuffel shifting of LU. 463*9b54502bSHong Zhang If an LU factorisation breaks down because of nonpositive pivots, 464*9b54502bSHong Zhang adding sufficient identity to the diagonal will remedy this. 465*9b54502bSHong Zhang Setting this causes a bisection method to find the minimum shift that 466*9b54502bSHong Zhang will lead to a well-defined LU. 467*9b54502bSHong Zhang 468*9b54502bSHong Zhang Input parameters: 469*9b54502bSHong Zhang + pc - the preconditioner context 470*9b54502bSHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 471*9b54502bSHong Zhang 472*9b54502bSHong Zhang Options Database Key: 473*9b54502bSHong Zhang . -pc_lu_shift - Activate PCLUSetShift() 474*9b54502bSHong Zhang 475*9b54502bSHong Zhang Level: intermediate 476*9b54502bSHong Zhang 477*9b54502bSHong Zhang .keywords: PC, indefinite, factorization 478*9b54502bSHong Zhang 479*9b54502bSHong Zhang .seealso: PCLUSetDamping(), PCILUSetShift() 480*9b54502bSHong Zhang @*/ 481*9b54502bSHong Zhang PetscErrorCode PCLUSetShift(PC pc,PetscTruth shifting) 482*9b54502bSHong Zhang { 483*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 484*9b54502bSHong Zhang 485*9b54502bSHong Zhang PetscFunctionBegin; 486*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 487*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 488*9b54502bSHong Zhang if (f) { 489*9b54502bSHong Zhang ierr = (*f)(pc,shifting);CHKERRQ(ierr); 490*9b54502bSHong Zhang } 491*9b54502bSHong Zhang PetscFunctionReturn(0); 492*9b54502bSHong Zhang } 493*9b54502bSHong Zhang 494*9b54502bSHong Zhang 495*9b54502bSHong Zhang #undef __FUNCT__ 496*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering" 497*9b54502bSHong Zhang /*@ 498*9b54502bSHong Zhang PCLUSetReuseOrdering - When similar matrices are factored, this 499*9b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 500*9b54502bSHong Zhang following factors; applies to both fill and drop tolerance LUs. 501*9b54502bSHong Zhang 502*9b54502bSHong Zhang Collective on PC 503*9b54502bSHong Zhang 504*9b54502bSHong Zhang Input Parameters: 505*9b54502bSHong Zhang + pc - the preconditioner context 506*9b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 507*9b54502bSHong Zhang 508*9b54502bSHong Zhang Options Database Key: 509*9b54502bSHong Zhang . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 510*9b54502bSHong Zhang 511*9b54502bSHong Zhang Level: intermediate 512*9b54502bSHong Zhang 513*9b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 514*9b54502bSHong Zhang 515*9b54502bSHong Zhang .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill() 516*9b54502bSHong Zhang @*/ 517*9b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering(PC pc,PetscTruth flag) 518*9b54502bSHong Zhang { 519*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 520*9b54502bSHong Zhang 521*9b54502bSHong Zhang PetscFunctionBegin; 522*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 523*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 524*9b54502bSHong Zhang if (f) { 525*9b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 526*9b54502bSHong Zhang } 527*9b54502bSHong Zhang PetscFunctionReturn(0); 528*9b54502bSHong Zhang } 529*9b54502bSHong Zhang 530*9b54502bSHong Zhang #undef __FUNCT__ 531*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill" 532*9b54502bSHong Zhang /*@ 533*9b54502bSHong Zhang PCLUSetReuseFill - When matrices with same nonzero structure are LU factored, 534*9b54502bSHong Zhang this causes later ones to use the fill computed in the initial factorization. 535*9b54502bSHong Zhang 536*9b54502bSHong Zhang Collective on PC 537*9b54502bSHong Zhang 538*9b54502bSHong Zhang Input Parameters: 539*9b54502bSHong Zhang + pc - the preconditioner context 540*9b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 541*9b54502bSHong Zhang 542*9b54502bSHong Zhang Options Database Key: 543*9b54502bSHong Zhang . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 544*9b54502bSHong Zhang 545*9b54502bSHong Zhang Level: intermediate 546*9b54502bSHong Zhang 547*9b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 548*9b54502bSHong Zhang 549*9b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill() 550*9b54502bSHong Zhang @*/ 551*9b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill(PC pc,PetscTruth flag) 552*9b54502bSHong Zhang { 553*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 554*9b54502bSHong Zhang 555*9b54502bSHong Zhang PetscFunctionBegin; 556*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 557*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 558*9b54502bSHong Zhang if (f) { 559*9b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 560*9b54502bSHong Zhang } 561*9b54502bSHong Zhang PetscFunctionReturn(0); 562*9b54502bSHong Zhang } 563*9b54502bSHong Zhang 564*9b54502bSHong Zhang #undef __FUNCT__ 565*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill" 566*9b54502bSHong Zhang /*@ 567*9b54502bSHong Zhang PCLUSetFill - Indicate the amount of fill you expect in the factored matrix, 568*9b54502bSHong Zhang fill = number nonzeros in factor/number nonzeros in original matrix. 569*9b54502bSHong Zhang 570*9b54502bSHong Zhang Collective on PC 571*9b54502bSHong Zhang 572*9b54502bSHong Zhang Input Parameters: 573*9b54502bSHong Zhang + pc - the preconditioner context 574*9b54502bSHong Zhang - fill - amount of expected fill 575*9b54502bSHong Zhang 576*9b54502bSHong Zhang Options Database Key: 577*9b54502bSHong Zhang . -pc_lu_fill <fill> - Sets fill amount 578*9b54502bSHong Zhang 579*9b54502bSHong Zhang Level: intermediate 580*9b54502bSHong Zhang 581*9b54502bSHong Zhang Note: 582*9b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 583*9b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 584*9b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 585*9b54502bSHong Zhang future runs. Default PETSc uses a value of 5.0 586*9b54502bSHong Zhang 587*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 588*9b54502bSHong Zhang 589*9b54502bSHong Zhang .seealso: PCILUSetFill() 590*9b54502bSHong Zhang @*/ 591*9b54502bSHong Zhang PetscErrorCode PCLUSetFill(PC pc,PetscReal fill) 592*9b54502bSHong Zhang { 593*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 594*9b54502bSHong Zhang 595*9b54502bSHong Zhang PetscFunctionBegin; 596*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 597*9b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 598*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 599*9b54502bSHong Zhang if (f) { 600*9b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 601*9b54502bSHong Zhang } 602*9b54502bSHong Zhang PetscFunctionReturn(0); 603*9b54502bSHong Zhang } 604*9b54502bSHong Zhang 605*9b54502bSHong Zhang #undef __FUNCT__ 606*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetDamping" 607*9b54502bSHong Zhang /*@ 608*9b54502bSHong Zhang PCLUSetDamping - adds this quantity to the diagonal of the matrix during the 609*9b54502bSHong Zhang LU numerical factorization 610*9b54502bSHong Zhang 611*9b54502bSHong Zhang Collective on PC 612*9b54502bSHong Zhang 613*9b54502bSHong Zhang Input Parameters: 614*9b54502bSHong Zhang + pc - the preconditioner context 615*9b54502bSHong Zhang - damping - amount of damping (use PETSC_DECIDE for default of 1.e-12) 616*9b54502bSHong Zhang 617*9b54502bSHong Zhang Options Database Key: 618*9b54502bSHong Zhang . -pc_lu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 619*9b54502bSHong Zhang 620*9b54502bSHong Zhang Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 621*9b54502bSHong Zhang pivot, then the damping is doubled until this is alleviated. 622*9b54502bSHong Zhang 623*9b54502bSHong Zhang Level: intermediate 624*9b54502bSHong Zhang 625*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 626*9b54502bSHong Zhang .seealso: PCILUSetFill(), PCILUSetDamp() 627*9b54502bSHong Zhang @*/ 628*9b54502bSHong Zhang PetscErrorCode PCLUSetDamping(PC pc,PetscReal damping) 629*9b54502bSHong Zhang { 630*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 631*9b54502bSHong Zhang 632*9b54502bSHong Zhang PetscFunctionBegin; 633*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 634*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 635*9b54502bSHong Zhang if (f) { 636*9b54502bSHong Zhang ierr = (*f)(pc,damping);CHKERRQ(ierr); 637*9b54502bSHong Zhang } 638*9b54502bSHong Zhang PetscFunctionReturn(0); 639*9b54502bSHong Zhang } 640*9b54502bSHong Zhang 641*9b54502bSHong Zhang #undef __FUNCT__ 642*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace" 643*9b54502bSHong Zhang /*@ 644*9b54502bSHong Zhang PCLUSetUseInPlace - Tells the system to do an in-place factorization. 645*9b54502bSHong Zhang For dense matrices, this enables the solution of much larger problems. 646*9b54502bSHong Zhang For sparse matrices the factorization cannot be done truly in-place 647*9b54502bSHong Zhang so this does not save memory during the factorization, but after the matrix 648*9b54502bSHong Zhang is factored, the original unfactored matrix is freed, thus recovering that 649*9b54502bSHong Zhang space. 650*9b54502bSHong Zhang 651*9b54502bSHong Zhang Collective on PC 652*9b54502bSHong Zhang 653*9b54502bSHong Zhang Input Parameters: 654*9b54502bSHong Zhang . pc - the preconditioner context 655*9b54502bSHong Zhang 656*9b54502bSHong Zhang Options Database Key: 657*9b54502bSHong Zhang . -pc_lu_in_place - Activates in-place factorization 658*9b54502bSHong Zhang 659*9b54502bSHong Zhang Notes: 660*9b54502bSHong Zhang PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when 661*9b54502bSHong Zhang a different matrix is provided for the multiply and the preconditioner in 662*9b54502bSHong Zhang a call to KSPSetOperators(). 663*9b54502bSHong Zhang This is because the Krylov space methods require an application of the 664*9b54502bSHong Zhang matrix multiplication, which is not possible here because the matrix has 665*9b54502bSHong Zhang been factored in-place, replacing the original matrix. 666*9b54502bSHong Zhang 667*9b54502bSHong Zhang Level: intermediate 668*9b54502bSHong Zhang 669*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU 670*9b54502bSHong Zhang 671*9b54502bSHong Zhang .seealso: PCILUSetUseInPlace() 672*9b54502bSHong Zhang @*/ 673*9b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace(PC pc) 674*9b54502bSHong Zhang { 675*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 676*9b54502bSHong Zhang 677*9b54502bSHong Zhang PetscFunctionBegin; 678*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 679*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 680*9b54502bSHong Zhang if (f) { 681*9b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 682*9b54502bSHong Zhang } 683*9b54502bSHong Zhang PetscFunctionReturn(0); 684*9b54502bSHong Zhang } 685*9b54502bSHong Zhang 686*9b54502bSHong Zhang #undef __FUNCT__ 687*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering" 688*9b54502bSHong Zhang /*@C 689*9b54502bSHong Zhang PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to 690*9b54502bSHong Zhang be used in the LU factorization. 691*9b54502bSHong Zhang 692*9b54502bSHong Zhang Collective on PC 693*9b54502bSHong Zhang 694*9b54502bSHong Zhang Input Parameters: 695*9b54502bSHong Zhang + pc - the preconditioner context 696*9b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 697*9b54502bSHong Zhang 698*9b54502bSHong Zhang Options Database Key: 699*9b54502bSHong Zhang . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 700*9b54502bSHong Zhang 701*9b54502bSHong Zhang Level: intermediate 702*9b54502bSHong Zhang 703*9b54502bSHong Zhang Notes: nested dissection is used by default 704*9b54502bSHong Zhang 705*9b54502bSHong Zhang .seealso: PCILUSetMatOrdering() 706*9b54502bSHong Zhang @*/ 707*9b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering(PC pc,MatOrderingType ordering) 708*9b54502bSHong Zhang { 709*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 710*9b54502bSHong Zhang 711*9b54502bSHong Zhang PetscFunctionBegin; 712*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 713*9b54502bSHong Zhang if (f) { 714*9b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 715*9b54502bSHong Zhang } 716*9b54502bSHong Zhang PetscFunctionReturn(0); 717*9b54502bSHong Zhang } 718*9b54502bSHong Zhang 719*9b54502bSHong Zhang #undef __FUNCT__ 720*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting" 721*9b54502bSHong Zhang /*@ 722*9b54502bSHong Zhang PCLUSetPivoting - Determines when pivoting is done during LU. 723*9b54502bSHong Zhang For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 724*9b54502bSHong Zhang it is never done. For the Matlab and SuperLU factorization this is used. 725*9b54502bSHong Zhang 726*9b54502bSHong Zhang Collective on PC 727*9b54502bSHong Zhang 728*9b54502bSHong Zhang Input Parameters: 729*9b54502bSHong Zhang + pc - the preconditioner context 730*9b54502bSHong Zhang - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 731*9b54502bSHong Zhang 732*9b54502bSHong Zhang Options Database Key: 733*9b54502bSHong Zhang . -pc_lu_pivoting <dtcol> 734*9b54502bSHong Zhang 735*9b54502bSHong Zhang Level: intermediate 736*9b54502bSHong Zhang 737*9b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks() 738*9b54502bSHong Zhang @*/ 739*9b54502bSHong Zhang PetscErrorCode PCLUSetPivoting(PC pc,PetscReal dtcol) 740*9b54502bSHong Zhang { 741*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 742*9b54502bSHong Zhang 743*9b54502bSHong Zhang PetscFunctionBegin; 744*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 745*9b54502bSHong Zhang if (f) { 746*9b54502bSHong Zhang ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 747*9b54502bSHong Zhang } 748*9b54502bSHong Zhang PetscFunctionReturn(0); 749*9b54502bSHong Zhang } 750*9b54502bSHong Zhang 751*9b54502bSHong Zhang #undef __FUNCT__ 752*9b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks" 753*9b54502bSHong Zhang /*@ 754*9b54502bSHong Zhang PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block 755*9b54502bSHong Zhang with BAIJ or SBAIJ matrices 756*9b54502bSHong Zhang 757*9b54502bSHong Zhang Collective on PC 758*9b54502bSHong Zhang 759*9b54502bSHong Zhang Input Parameters: 760*9b54502bSHong Zhang + pc - the preconditioner context 761*9b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 762*9b54502bSHong Zhang 763*9b54502bSHong Zhang Options Database Key: 764*9b54502bSHong Zhang . -pc_lu_pivot_in_blocks <true,false> 765*9b54502bSHong Zhang 766*9b54502bSHong Zhang Level: intermediate 767*9b54502bSHong Zhang 768*9b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivoting() 769*9b54502bSHong Zhang @*/ 770*9b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks(PC pc,PetscTruth pivot) 771*9b54502bSHong Zhang { 772*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 773*9b54502bSHong Zhang 774*9b54502bSHong Zhang PetscFunctionBegin; 775*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 776*9b54502bSHong Zhang if (f) { 777*9b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 778*9b54502bSHong Zhang } 779*9b54502bSHong Zhang PetscFunctionReturn(0); 780*9b54502bSHong Zhang } 781*9b54502bSHong Zhang 782*9b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 783*9b54502bSHong Zhang 784*9b54502bSHong Zhang /*MC 785*9b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 786*9b54502bSHong Zhang 787*9b54502bSHong Zhang Options Database Keys: 788*9b54502bSHong Zhang + -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 789*9b54502bSHong Zhang . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 790*9b54502bSHong Zhang . -pc_lu_fill <fill> - Sets fill amount 791*9b54502bSHong Zhang . -pc_lu_damping <damping> - Sets damping amount 792*9b54502bSHong Zhang . -pc_lu_in_place - Activates in-place factorization 793*9b54502bSHong Zhang . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 794*9b54502bSHong Zhang . -pc_lu_pivoting <dtcol> 795*9b54502bSHong Zhang - -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 796*9b54502bSHong Zhang stability of factorization. 797*9b54502bSHong Zhang 798*9b54502bSHong Zhang Notes: Not all options work for all matrix formats 799*9b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 800*9b54502bSHong Zhang algorithms 801*9b54502bSHong Zhang 802*9b54502bSHong Zhang Level: beginner 803*9b54502bSHong Zhang 804*9b54502bSHong Zhang Concepts: LU factorization, direct solver 805*9b54502bSHong Zhang 806*9b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 807*9b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 808*9b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 809*9b54502bSHong Zhang 810*9b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 811*9b54502bSHong Zhang PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(), 812*9b54502bSHong Zhang PCLUSetFill(), PCLUSetDamping(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCLUSetPivoting(), 813*9b54502bSHong Zhang PCLUSetPivotingInBlocks() 814*9b54502bSHong Zhang M*/ 815*9b54502bSHong Zhang 816*9b54502bSHong Zhang EXTERN_C_BEGIN 817*9b54502bSHong Zhang #undef __FUNCT__ 818*9b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 819*9b54502bSHong Zhang PetscErrorCode PCCreate_LU(PC pc) 820*9b54502bSHong Zhang { 821*9b54502bSHong Zhang PetscErrorCode ierr; 822*9b54502bSHong Zhang PetscMPIInt size; 823*9b54502bSHong Zhang PC_LU *dir; 824*9b54502bSHong Zhang 825*9b54502bSHong Zhang PetscFunctionBegin; 826*9b54502bSHong Zhang ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr); 827*9b54502bSHong Zhang PetscLogObjectMemory(pc,sizeof(PC_LU)); 828*9b54502bSHong Zhang 829*9b54502bSHong Zhang ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 830*9b54502bSHong Zhang dir->fact = 0; 831*9b54502bSHong Zhang dir->inplace = PETSC_FALSE; 832*9b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 833*9b54502bSHong Zhang 834*9b54502bSHong Zhang dir->info.fill = 5.0; 835*9b54502bSHong Zhang dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 836*9b54502bSHong Zhang dir->info.damping = 0.0; 837*9b54502bSHong Zhang dir->info.zeropivot = 1.e-12; 838*9b54502bSHong Zhang dir->info.pivotinblocks = 1.0; 839*9b54502bSHong Zhang dir->info.shift = PETSC_FALSE; 840*9b54502bSHong Zhang dir->info.shift_fraction = 0.0; 841*9b54502bSHong Zhang dir->col = 0; 842*9b54502bSHong Zhang dir->row = 0; 843*9b54502bSHong Zhang ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 844*9b54502bSHong Zhang if (size == 1) { 845*9b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 846*9b54502bSHong Zhang } else { 847*9b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 848*9b54502bSHong Zhang } 849*9b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 850*9b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 851*9b54502bSHong Zhang pc->data = (void*)dir; 852*9b54502bSHong Zhang 853*9b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 854*9b54502bSHong Zhang pc->ops->apply = PCApply_LU; 855*9b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 856*9b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 857*9b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 858*9b54502bSHong Zhang pc->ops->view = PCView_LU; 859*9b54502bSHong Zhang pc->ops->applyrichardson = 0; 860*9b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU; 861*9b54502bSHong Zhang 862*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU", 863*9b54502bSHong Zhang PCLUSetFill_LU);CHKERRQ(ierr); 864*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU", 865*9b54502bSHong Zhang PCLUSetDamping_LU);CHKERRQ(ierr); 866*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetShift_C","PCLUSetShift_LU", 867*9b54502bSHong Zhang PCLUSetShift_LU);CHKERRQ(ierr); 868*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU", 869*9b54502bSHong Zhang PCLUSetUseInPlace_LU);CHKERRQ(ierr); 870*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU", 871*9b54502bSHong Zhang PCLUSetMatOrdering_LU);CHKERRQ(ierr); 872*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU", 873*9b54502bSHong Zhang PCLUSetReuseOrdering_LU);CHKERRQ(ierr); 874*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU", 875*9b54502bSHong Zhang PCLUSetReuseFill_LU);CHKERRQ(ierr); 876*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU", 877*9b54502bSHong Zhang PCLUSetPivoting_LU);CHKERRQ(ierr); 878*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU", 879*9b54502bSHong Zhang PCLUSetPivotInBlocks_LU);CHKERRQ(ierr); 880*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetZeroPivot_C","PCLUSetZeroPivot_LU", 881*9b54502bSHong Zhang PCLUSetZeroPivot_LU);CHKERRQ(ierr); 882*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C","PCLUReorderForNonzeroDiagonal_LU", 883*9b54502bSHong Zhang PCLUReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 884*9b54502bSHong Zhang PetscFunctionReturn(0); 885*9b54502bSHong Zhang } 886*9b54502bSHong Zhang EXTERN_C_END 887