1: ! 2: ! 3: ! Include file for Fortran use of the KSP package in PETSc 4: ! 5: #if !defined (__PETSCKSP_H) 8: #define KSP PetscFortranAddr 9: #define KSPType character*(80) 10: #define KSPCGType PetscEnum 11: #define KSPConvergedReason PetscEnum 12: #define KSPNormType PetscEnum 13: #define KSPGMRESCGSRefinementType PetscEnum 14: ! 15: ! Various Krylov subspace methods 16: ! 17: #define KSPRICHARDSON 'richardson' 18: #define KSPCHEBYCHEV 'chebychev' 19: #define KSPCG 'cg' 20: #define KSPCGNE 'cgne' 21: #define KSPSTCG 'stcg' 22: #define KSPGLTR 'gltr' 23: #define KSPGMRES 'gmres' 24: #define KSPFGMRES 'fgmres' 25: #define KSPLGMRES 'lgmres' 26: #define KSPTCQMR 'tcqmr' 27: #define KSPBCGS 'bcgs' 28: #define KSPBCGSL 'bcgsl' 29: #define KSPCGS 'cgs' 30: #define KSPTFQMR 'tfqmr' 31: #define KSPCR 'cr' 32: #define KSPLSQR 'lsqr' 33: #define KSPPREONLY 'preonly' 34: #define KSPQCG 'qcg' 35: #define KSPBICG 'bicg' 36: #define KSPMINRES 'minres' 37: #define KSPSYMMLQ 'symmlq' 38: #define KSPLCD 'lcd' 39: #endif 42: #if !defined (PETSC_AVOID_DECLARATIONS) 44: ! 45: ! CG Types 46: ! 47: PetscEnum KSP_CG_SYMMETRIC 48: PetscEnum KSP_CG_HERMITIAN 50: parameter (KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1) 52: PetscEnum KSP_CONVERGED_RTOL 53: PetscEnum KSP_CONVERGED_ATOL 54: PetscEnum KSP_CONVERGED_ITS 55: PetscEnum KSP_DIVERGED_NULL 56: PetscEnum KSP_DIVERGED_ITS 57: PetscEnum KSP_DIVERGED_DTOL 58: PetscEnum KSP_DIVERGED_BREAKDOWN 59: PetscEnum KSP_CONVERGED_ITERATING 60: PetscEnum KSP_CONVERGED_CG_NEG_CURVE 61: PetscEnum KSP_CONVERGED_CG_CONSTRAINED 62: PetscEnum KSP_CONVERGED_STEP_LENGTH 63: PetscEnum KSP_CONVERGED_HAPPY_BREAKDOWN 64: PetscEnum KSP_DIVERGED_BREAKDOWN_BICG 65: PetscEnum KSP_DIVERGED_NONSYMMETRIC 66: PetscEnum KSP_DIVERGED_INDEFINITE_PC 67: PetscEnum KSP_DIVERGED_NAN 68: PetscEnum KSP_DIVERGED_INDEFINITE_MAT 70: parameter (KSP_CONVERGED_RTOL = 2) 71: parameter (KSP_CONVERGED_ATOL = 3) 72: parameter (KSP_CONVERGED_ITS = 4) 73: parameter (KSP_CONVERGED_CG_NEG_CURVE = 5) 74: parameter (KSP_CONVERGED_CG_CONSTRAINED = 6) 75: parameter (KSP_CONVERGED_STEP_LENGTH = 7) 76: parameter (KSP_CONVERGED_HAPPY_BREAKDOWN = 8) 78: parameter (KSP_DIVERGED_NULL = -2) 79: parameter (KSP_DIVERGED_ITS = -3) 80: parameter (KSP_DIVERGED_DTOL = -4) 81: parameter (KSP_DIVERGED_BREAKDOWN = -5) 82: parameter (KSP_DIVERGED_BREAKDOWN_BICG = -6) 83: parameter (KSP_DIVERGED_NONSYMMETRIC = -7) 84: parameter (KSP_DIVERGED_INDEFINITE_PC = -8) 85: parameter (KSP_DIVERGED_NAN = -9) 86: parameter (KSP_DIVERGED_INDEFINITE_MAT = -10) 88: parameter (KSP_CONVERGED_ITERATING = 0) 89: ! 90: ! Possible arguments to KSPSetNormType() 91: ! 92: PetscEnum KSP_NORM_NO 93: PetscEnum KSP_NORM_PRECONDITIONED 94: PetscEnum KSP_NORM_UNPRECONDITIONED 95: PetscEnum KSP_NORM_NATURAL 96: 97: parameter (KSP_NORM_NO=0) 98: parameter (KSP_NORM_PRECONDITIONED=1) 99: parameter (KSP_NORM_UNPRECONDITIONED=2) 100: parameter (KSP_NORM_NATURAL=3) 101: ! 102: ! Possible arguments to KSPMonitorSet() 103: ! 104: external KSPDEFAULTCONVERGED 105: external KSPMONITORDEFAULT 106: external KSPMONITORTRUERESIDUALNORM 107: external KSPMONITORLG 108: external KSPMONITORLGTRUERESIDUALNORM 109: external KSPMONITORSOLUTION 110: external KSPMONITORSINGULARVALUE 111: external KSPGMRESMONITORKRYLOV 112: ! 113: ! Possible arguments to KSPGMRESSetRefinementType() 114: ! 115: PetscEnum KSP_GMRES_CGS_REFINE_NEVER 116: PetscEnum KSP_GMRES_CGS_REFINE_IFNEEDED 117: PetscEnum KSP_GMRES_CGS_REFINE_ALWAYS 118: ! 119: parameter (KSP_GMRES_CGS_REFINE_NEVER = 0) 120: parameter (KSP_GMRES_CGS_REFINE_IFNEEDED = 1) 121: parameter (KSP_GMRES_CGS_REFINE_ALWAYS = 2) 122: ! 123: !PETSC_DEC_ATTRIBUTES(KSPDEFAULTCONVERGED,'_KSPDEFAULTCONVERGED') 124: !PETSC_DEC_ATTRIBUTES(KSPMONITORDEFAULT,'_KSPMONITORDEFAULT') 125: !PETSC_DEC_ATTRIBUTES(KSPMONITORTRUERESIDUALNORM,'_KSPMONITORTRUERESIDUALNORM') 126: !PETSC_DEC_ATTRIBUTES(KSPMONITORLG,'_KSPMONITORLG') 127: !PETSC_DEC_ATTRIBUTES(KSPMONITORLGTRUERESIDUALNORM,'_KSPMONITORLGTRUERESIDUALNORM') 128: !PETSC_DEC_ATTRIBUTES(KSPMONITORSOLUTION,'_KSPMONITORSOLUTION') 129: !PETSC_DEC_ATTRIBUTES(KSPMONITORSINGULARVALUE,'_KSPMONITORSINGULARVALUE') 130: !PETSC_DEC_ATTRIBUTES(KSPGMRESMONITORKRYLOV,'_KSPGMRESMONITORKRYLOV') 132: ! 133: ! End of Fortran include file for the KSP package in PETSc 134: ! 136: #endif