1 #define PETSCKSP_DLL 2 3 #include "../src/ksp/pc/impls/factor/factor.h" /*I "petscpc.h" I*/ 4 5 /* ------------------------------------------------------------------------------------------*/ 6 7 EXTERN_C_BEGIN 8 #undef __FUNCT__ 9 #define __FUNCT__ "PCFactorSetZeroPivot_Factor" 10 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Factor(PC pc,PetscReal z) 11 { 12 PC_Factor *ilu = (PC_Factor*)pc->data; 13 14 PetscFunctionBegin; 15 ilu->info.zeropivot = z; 16 PetscFunctionReturn(0); 17 } 18 EXTERN_C_END 19 20 EXTERN_C_BEGIN 21 #undef __FUNCT__ 22 #define __FUNCT__ "PCFactorSetShiftType_Factor" 23 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype) 24 { 25 PC_Factor *dir = (PC_Factor*)pc->data; 26 27 PetscFunctionBegin; 28 if (shifttype == (MatFactorShiftType)PETSC_DECIDE){ 29 dir->info.shifttype = MAT_SHIFT_NONE; 30 } else { 31 dir->info.shifttype = shifttype; 32 } 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "PCFactorSetShiftAmount_Factor" 38 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount) 39 { 40 PC_Factor *dir = (PC_Factor*)pc->data; 41 42 PetscFunctionBegin; 43 if (shiftamount == (PetscReal) PETSC_DECIDE){ 44 dir->info.shiftamount = 1.e-12; 45 } else { 46 dir->info.shiftamount = shiftamount; 47 } 48 PetscFunctionReturn(0); 49 } 50 EXTERN_C_END 51 52 EXTERN_C_BEGIN 53 #undef __FUNCT__ 54 #define __FUNCT__ "PCFactorSetDropTolerance_Factor" 55 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 56 { 57 PC_Factor *ilu = (PC_Factor*)pc->data; 58 59 PetscFunctionBegin; 60 if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 61 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use"); 62 } 63 ilu->info.usedt = PETSC_TRUE; 64 ilu->info.dt = dt; 65 ilu->info.dtcol = dtcol; 66 ilu->info.dtcount = dtcount; 67 ilu->info.fill = PETSC_DEFAULT; 68 PetscFunctionReturn(0); 69 } 70 EXTERN_C_END 71 72 EXTERN_C_BEGIN 73 #undef __FUNCT__ 74 #define __FUNCT__ "PCFactorSetFill_Factor" 75 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Factor(PC pc,PetscReal fill) 76 { 77 PC_Factor *dir = (PC_Factor*)pc->data; 78 79 PetscFunctionBegin; 80 dir->info.fill = fill; 81 PetscFunctionReturn(0); 82 } 83 EXTERN_C_END 84 85 EXTERN_C_BEGIN 86 #undef __FUNCT__ 87 #define __FUNCT__ "PCFactorSetMatOrderingType_Factor" 88 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering) 89 { 90 PC_Factor *dir = (PC_Factor*)pc->data; 91 PetscErrorCode ierr; 92 PetscTruth flg; 93 94 PetscFunctionBegin; 95 if (!pc->setupcalled) { 96 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 97 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 98 } else { 99 ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 100 if (!flg) { 101 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use"); 102 } 103 } 104 PetscFunctionReturn(0); 105 } 106 EXTERN_C_END 107 108 EXTERN_C_BEGIN 109 #undef __FUNCT__ 110 #define __FUNCT__ "PCFactorSetLevels_Factor" 111 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_Factor(PC pc,PetscInt levels) 112 { 113 PC_Factor *ilu = (PC_Factor*)pc->data; 114 115 PetscFunctionBegin; 116 if (!pc->setupcalled) { 117 ilu->info.levels = levels; 118 } else if (ilu->info.usedt || ilu->info.levels != levels) { 119 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use"); 120 } 121 PetscFunctionReturn(0); 122 } 123 EXTERN_C_END 124 125 EXTERN_C_BEGIN 126 #undef __FUNCT__ 127 #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor" 128 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_Factor(PC pc) 129 { 130 PC_Factor *dir = (PC_Factor*)pc->data; 131 132 PetscFunctionBegin; 133 dir->info.diagonal_fill = 1; 134 PetscFunctionReturn(0); 135 } 136 EXTERN_C_END 137 138 139 /* ------------------------------------------------------------------------------------------*/ 140 141 EXTERN_C_BEGIN 142 #undef __FUNCT__ 143 #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor" 144 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_Factor(PC pc,PetscTruth pivot) 145 { 146 PC_Factor *dir = (PC_Factor*)pc->data; 147 148 PetscFunctionBegin; 149 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 150 PetscFunctionReturn(0); 151 } 152 EXTERN_C_END 153 154 EXTERN_C_BEGIN 155 #undef __FUNCT__ 156 #define __FUNCT__ "PCFactorGetMatrix_Factor" 157 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix_Factor(PC pc,Mat *mat) 158 { 159 PC_Factor *ilu = (PC_Factor*)pc->data; 160 161 PetscFunctionBegin; 162 if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 163 *mat = ilu->fact; 164 PetscFunctionReturn(0); 165 } 166 EXTERN_C_END 167 168 EXTERN_C_BEGIN 169 #undef __FUNCT__ 170 #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor" 171 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype) 172 { 173 PetscErrorCode ierr; 174 PC_Factor *lu = (PC_Factor*)pc->data; 175 176 PetscFunctionBegin; 177 if (lu->fact) { 178 const MatSolverPackage ltype; 179 PetscTruth flg; 180 ierr = MatFactorGetSolverPackage(lu->fact,<ype);CHKERRQ(ierr); 181 ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr); 182 if (!flg) { 183 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used"); 184 } 185 } 186 ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr); 187 ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr); 188 PetscFunctionReturn(0); 189 } 190 EXTERN_C_END 191 192 EXTERN_C_BEGIN 193 #undef __FUNCT__ 194 #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor" 195 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype) 196 { 197 PC_Factor *lu = (PC_Factor*)pc->data; 198 199 PetscFunctionBegin; 200 *stype = lu->solvertype; 201 PetscFunctionReturn(0); 202 } 203 EXTERN_C_END 204 205 EXTERN_C_BEGIN 206 #undef __FUNCT__ 207 #define __FUNCT__ "PCFactorSetColumnPivot_Factor" 208 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol) 209 { 210 PC_Factor *dir = (PC_Factor*)pc->data; 211 212 PetscFunctionBegin; 213 if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol); 214 dir->info.dtcol = dtcol; 215 PetscFunctionReturn(0); 216 } 217 EXTERN_C_END 218 219 EXTERN_C_BEGIN 220 #undef __FUNCT__ 221 #define __FUNCT__ "PCSetFromOptions_Factor" 222 PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions_Factor(PC pc) 223 { 224 PC_Factor *factor = (PC_Factor*)pc->data; 225 PetscErrorCode ierr; 226 PetscTruth flg = PETSC_FALSE,set; 227 char tname[256], solvertype[64]; 228 PetscFList ordlist; 229 230 PetscFunctionBegin; 231 if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 232 ierr = PetscOptionsTruth("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 233 if (flg) { 234 ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 235 } 236 ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);CHKERRQ(ierr); 237 238 ierr = PetscOptionsEnum("-pc_factor_shift_type","Shift added to diagonal","PCFactorSetShiftType", 239 MatFactorShiftTypes,(PetscEnum)((PC_Factor*)factor)->info.shifttype,(PetscEnum*)&((PC_Factor*)factor)->info.shifttype,&flg);CHKERRQ(ierr); 240 241 ierr = PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);CHKERRQ(ierr); 242 243 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);CHKERRQ(ierr); 244 ierr = PetscOptionsReal("-pc_factor_column_pivot","Column pivot tolerance (used only for some factorization)","PCFactorSetColumnPivot",((PC_Factor*)factor)->info.dtcol,&((PC_Factor*)factor)->info.dtcol,&flg);CHKERRQ(ierr); 245 246 flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 247 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 248 if (set) { 249 ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 250 } 251 252 flg = PETSC_FALSE; 253 ierr = PetscOptionsTruth("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 254 if (flg) { 255 ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 256 } 257 flg = PETSC_FALSE; 258 ierr = PetscOptionsTruth("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 259 if (flg) { 260 ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 261 } 262 263 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 264 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr); 265 if (flg) { 266 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 267 } 268 269 /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 270 ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 271 if (flg) { 272 ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 273 } 274 PetscFunctionReturn(0); 275 } 276 EXTERN_C_END 277