Actual source code: snesregi.c

  1: #include <petsc/private/snesimpl.h>

  3: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES);
  4: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES);
  5: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES);
  6: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES);
  7: PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES);
  8: PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES);
  9: PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONRSLS(SNES);
 10: PETSC_EXTERN PetscErrorCode SNESCreate_VINEWTONSSLS(SNES);
 11: PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES);
 12: PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES);
 13: PETSC_EXTERN PetscErrorCode SNESCreate_Shell(SNES);
 14: PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES);
 15: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES);
 16: PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES);
 17: PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES);
 18: PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES);
 19: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES);
 20: PETSC_EXTERN PetscErrorCode SNESCreate_ASPIN(SNES);
 21: PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES);
 22: PETSC_EXTERN PetscErrorCode SNESCreate_Patch(SNES);

 24: const char *SNESConvergedReasons_Shifted[] = {" ", "DIVERGED_TR_DELTA", "DIVERGED_JACOBIAN_DOMAIN", "DIVERGED_DTOL", "DIVERGED_LOCAL_MIN", "DIVERGED_INNER", "DIVERGED_LINE_SEARCH", "DIVERGED_MAX_IT", "DIVERGED_FNORM_NAN", "DIVERGED_LINEAR_SOLVE", "DIVERGED_FUNCTION_COUNT", "DIVERGED_FUNCTION_DOMAIN", "CONVERGED_ITERATING", " ", "CONVERGED_FNORM_ABS", "CONVERGED_FNORM_RELATIVE", "CONVERGED_SNORM_RELATIVE", "CONVERGED_ITS", " ", "SNESConvergedReason", "", NULL};
 25: const char *const *SNESConvergedReasons = SNESConvergedReasons_Shifted + 12;

 27: const char              *SNESNormSchedules_Shifted[] = {"DEFAULT", "NONE", "ALWAYS", "INITIALONLY", "FINALONLY", "INITIALFINALONLY", "SNESNormSchedule", "SNES_NORM_", NULL};
 28: const char *const *const SNESNormSchedules           = SNESNormSchedules_Shifted + 1;

 30: const char              *SNESFunctionTypes_Shifted[] = {"DEFAULT", "UNPRECONDITIONED", "PRECONDITIONED", "SNESFunctionType", "SNES_FUNCTION_", NULL};
 31: const char *const *const SNESFunctionTypes           = SNESFunctionTypes_Shifted + 1;

 33: /*@C
 34:   SNESRegisterAll - Registers all of the nonlinear solver methods in the `SNES` package.

 36:   Not Collective

 38:   Level: advanced

 40: .seealso: [](ch_snes), `SNES`, `SNESRegisterDestroy()`
 41: @*/
 42: PetscErrorCode SNESRegisterAll(void)
 43: {
 44:   PetscFunctionBegin;
 45:   if (SNESRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
 46:   SNESRegisterAllCalled = PETSC_TRUE;

 48:   PetscCall(SNESRegister(SNESNEWTONLS, SNESCreate_NEWTONLS));
 49:   PetscCall(SNESRegister(SNESNEWTONTR, SNESCreate_NEWTONTR));
 50:   PetscCall(SNESRegister(SNESNEWTONTRDC, SNESCreate_NEWTONTRDC));
 51:   PetscCall(SNESRegister(SNESNRICHARDSON, SNESCreate_NRichardson));
 52:   PetscCall(SNESRegister(SNESKSPONLY, SNESCreate_KSPONLY));
 53:   PetscCall(SNESRegister(SNESKSPTRANSPOSEONLY, SNESCreate_KSPTRANSPOSEONLY));
 54:   PetscCall(SNESRegister(SNESVINEWTONRSLS, SNESCreate_VINEWTONRSLS));
 55:   PetscCall(SNESRegister(SNESVINEWTONSSLS, SNESCreate_VINEWTONSSLS));
 56:   PetscCall(SNESRegister(SNESNGMRES, SNESCreate_NGMRES));
 57:   PetscCall(SNESRegister(SNESQN, SNESCreate_QN));
 58:   PetscCall(SNESRegister(SNESSHELL, SNESCreate_Shell));
 59:   PetscCall(SNESRegister(SNESNGS, SNESCreate_NGS));
 60:   PetscCall(SNESRegister(SNESNCG, SNESCreate_NCG));
 61:   PetscCall(SNESRegister(SNESFAS, SNESCreate_FAS));
 62:   PetscCall(SNESRegister(SNESMS, SNESCreate_MS));
 63:   PetscCall(SNESRegister(SNESNASM, SNESCreate_NASM));
 64:   PetscCall(SNESRegister(SNESANDERSON, SNESCreate_Anderson));
 65:   PetscCall(SNESRegister(SNESASPIN, SNESCreate_ASPIN));
 66:   PetscCall(SNESRegister(SNESCOMPOSITE, SNESCreate_Composite));
 67:   PetscCall(SNESRegister(SNESPATCH, SNESCreate_Patch));

 69:   PetscCall(KSPMonitorRegister("snes_preconditioned_residual", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, KSPMonitorSNESResidual, NULL, NULL));
 70:   PetscCall(KSPMonitorRegister("snes_preconditioned_residual", PETSCVIEWERDRAW, PETSC_VIEWER_DRAW_LG, KSPMonitorSNESResidualDrawLG, KSPMonitorSNESResidualDrawLGCreate, NULL));
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }