Actual source code: gmres2.c

  1: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>

  3: /*@C
  4:   KSPGMRESSetOrthogonalization - Sets the orthogonalization routine used by `KSPGMRES` and `KSPFGMRES`.

  6:   Logically Collective

  8:   Input Parameters:
  9: + ksp - iterative context obtained from `KSPCreate()`
 10: - fcn - orthogonalization function

 12:   Calling sequence of `fcn`:
 13: + ksp - the solver context
 14: - it  - the current iteration

 16:   Options Database Keys:
 17: + -ksp_gmres_classicalgramschmidt - Activates KSPGMRESClassicalGramSchmidtOrthogonalization() (default)
 18: - -ksp_gmres_modifiedgramschmidt  - Activates KSPGMRESModifiedGramSchmidtOrthogonalization()

 20:   Level: intermediate

 22:   Notes:
 23:   Two orthogonalization routines are predefined, `KSPGMRESModifiedGramSchmidtOrthogonalization()` and the default
 24:   `KSPGMRESClassicalGramSchmidtOrthogonalization()`.

 26:   Use `KSPGMRESSetCGSRefinementType()` to determine if iterative refinement is used to increase stability.

 28: .seealso: [](ch_ksp), `KSPGMRESSetRestart()`, `KSPGMRESSetPreAllocateVectors()`,
 29: `KSPGMRESSetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
 30: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESGetCGSRefinementType()`
 31: @*/
 32: PetscErrorCode KSPGMRESSetOrthogonalization(KSP ksp, PetscErrorCode (*fcn)(KSP ksp, PetscInt it))
 33: {
 34:   PetscFunctionBegin;
 36:   PetscTryMethod(ksp, "KSPGMRESSetOrthogonalization_C", (KSP, PetscErrorCode (*)(KSP, PetscInt)), (ksp, fcn));
 37:   PetscFunctionReturn(PETSC_SUCCESS);
 38: }

 40: /*@C
 41:   KSPGMRESGetOrthogonalization - Gets the orthogonalization routine used by `KSPGMRES` and `KSPFGMRES`.

 43:   Not Collective

 45:   Input Parameter:
 46: . ksp - iterative context obtained from `KSPCreate()`

 48:   Output Parameter:
 49: . fcn - orthogonalization function

 51:   Calling sequence of `fcn`:
 52: + ksp - the solver context
 53: - it  - the current iteration

 55:   Level: intermediate

 57:   Notes:
 58:   Two orthogonalization routines are predefined,  `KSPGMRESModifiedGramSchmidtOrthogonalization()`, and the default
 59:   `KSPGMRESClassicalGramSchmidtOrthogonalization()`

 61:   Use `KSPGMRESSetCGSRefinementType()` to determine if iterative refinement is used to increase stability.

 63: .seealso: [](ch_ksp), `KSPGMRESSetRestart()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESSetOrthogonalization()`,
 64:           `KSPGMRESModifiedGramSchmidtOrthogonalization()`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESGetCGSRefinementType()`
 65: @*/
 66: PetscErrorCode KSPGMRESGetOrthogonalization(KSP ksp, PetscErrorCode (**fcn)(KSP ksp, PetscInt it))
 67: {
 68:   PetscFunctionBegin;
 70:   PetscUseMethod(ksp, "KSPGMRESGetOrthogonalization_C", (KSP, PetscErrorCode (**)(KSP, PetscInt)), (ksp, fcn));
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }