Actual source code: factor.h

  1: /*
  2:    Private data structure for ILU/ICC/LU/Cholesky preconditioners.
  3: */
  4: #pragma once

  6: #include <petsc/private/pcimpl.h>

  8: typedef struct {
  9:   Mat             fact; /* factored matrix */
 10:   MatFactorInfo   info;
 11:   MatOrderingType ordering; /* matrix reordering */
 12:   char           *solvertype;
 13:   MatFactorType   factortype;
 14:   PetscReal       actualfill;
 15:   PetscBool       inplace;       /* flag indicating in-place factorization */
 16:   PetscBool       reuseordering; /* reuses previous reordering computed */
 17:   PetscBool       reusefill;     /* reuse fill from previous LU */
 18: } PC_Factor;

 20: PETSC_INTERN PetscErrorCode PCFactorInitialize(PC, MatFactorType);
 21: PETSC_INTERN PetscErrorCode PCFactorGetMatrix_Factor(PC, Mat *);

 23: PETSC_INTERN PetscErrorCode PCFactorSetZeroPivot_Factor(PC, PetscReal);
 24: PETSC_INTERN PetscErrorCode PCFactorGetZeroPivot_Factor(PC, PetscReal *);
 25: PETSC_INTERN PetscErrorCode PCFactorSetShiftType_Factor(PC, MatFactorShiftType);
 26: PETSC_INTERN PetscErrorCode PCFactorGetShiftType_Factor(PC, MatFactorShiftType *);
 27: PETSC_INTERN PetscErrorCode PCFactorSetShiftAmount_Factor(PC, PetscReal);
 28: PETSC_INTERN PetscErrorCode PCFactorGetShiftAmount_Factor(PC, PetscReal *);
 29: PETSC_INTERN PetscErrorCode PCFactorSetDropTolerance_Factor(PC, PetscReal, PetscReal, PetscInt);
 30: PETSC_INTERN PetscErrorCode PCFactorSetFill_Factor(PC, PetscReal);
 31: PETSC_INTERN PetscErrorCode PCFactorSetMatOrderingType_Factor(PC, MatOrderingType);
 32: PETSC_INTERN PetscErrorCode PCFactorGetLevels_Factor(PC, PetscInt *);
 33: PETSC_INTERN PetscErrorCode PCFactorSetLevels_Factor(PC, PetscInt);
 34: PETSC_INTERN PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC, PetscBool);
 35: PETSC_INTERN PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC, PetscBool *);
 36: PETSC_INTERN PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC, PetscBool);
 37: PETSC_INTERN PetscErrorCode PCFactorSetMatSolverType_Factor(PC, MatSolverType);
 38: PETSC_INTERN PetscErrorCode PCFactorSetUpMatSolverType_Factor(PC);
 39: PETSC_INTERN PetscErrorCode PCFactorGetMatSolverType_Factor(PC, MatSolverType *);
 40: PETSC_INTERN PetscErrorCode PCFactorSetColumnPivot_Factor(PC, PetscReal);
 41: PETSC_INTERN PetscErrorCode PCSetFromOptions_Factor(PC, PetscOptionItems *PetscOptionsObject);
 42: PETSC_INTERN PetscErrorCode PCView_Factor(PC, PetscViewer);
 43: PETSC_INTERN PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC);
 44: PETSC_INTERN PetscErrorCode PCFactorClearComposedFunctions(PC);