Actual source code: snesshell.c

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

  3: typedef struct {
  4:   PetscErrorCode (*solve)(SNES, Vec);
  5:   void *ctx;
  6: } SNES_Shell;

  8: /*@C
  9:   SNESShellSetSolve - Sets routine to apply as solver to a `SNESSHELL` `SNES` object

 11:   Logically Collective

 13:   Input Parameters:
 14: + snes  - the `SNES` nonlinear solver context
 15: - solve - the application-provided solver routine

 17:   Calling sequence of `apply`:
 18: + snes - the preconditioner, get the application context with `SNESShellGetContext()` provided with `SNESShellSetContext()`
 19: - xout - solution vector

 21:   Level: advanced

 23: .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESShellSetContext()`, `SNESShellGetContext()`
 24: @*/
 25: PetscErrorCode SNESShellSetSolve(SNES snes, PetscErrorCode (*solve)(SNES snes, Vec xout))
 26: {
 27:   PetscFunctionBegin;
 29:   PetscTryMethod(snes, "SNESShellSetSolve_C", (SNES, PetscErrorCode (*)(SNES, Vec)), (snes, solve));
 30:   PetscFunctionReturn(PETSC_SUCCESS);
 31: }

 33: static PetscErrorCode SNESDestroy_Shell(SNES snes)
 34: {
 35:   PetscFunctionBegin;
 36:   PetscCall(PetscFree(snes->data));
 37:   PetscFunctionReturn(PETSC_SUCCESS);
 38: }

 40: /*@
 41:   SNESShellGetContext - Returns the user-provided context associated with a `SNESSHELL`

 43:   Not Collective

 45:   Input Parameter:
 46: . snes - should have been created with `SNESSetType`(snes,`SNESSHELL`);

 48:   Output Parameter:
 49: . ctx - the user provided context

 51:   Level: advanced

 53: .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESCreateShell()`, `SNESShellSetContext()`
 54: @*/
 55: PetscErrorCode SNESShellGetContext(SNES snes, void *ctx)
 56: {
 57:   PetscBool flg;

 59:   PetscFunctionBegin;
 61:   PetscAssertPointer(ctx, 2);
 62:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESSHELL, &flg));
 63:   if (!flg) *(void **)ctx = NULL;
 64:   else *(void **)ctx = ((SNES_Shell *)snes->data)->ctx;
 65:   PetscFunctionReturn(PETSC_SUCCESS);
 66: }

 68: /*@
 69:   SNESShellSetContext - sets the context for a `SNESSHELL`

 71:   Logically Collective

 73:   Input Parameters:
 74: + snes - the `SNESSHELL`
 75: - ctx  - the context

 77:   Level: advanced

 79: .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESCreateShell()`, `SNESShellGetContext()`
 80: @*/
 81: PetscErrorCode SNESShellSetContext(SNES snes, void *ctx)
 82: {
 83:   SNES_Shell *shell = (SNES_Shell *)snes->data;
 84:   PetscBool   flg;

 86:   PetscFunctionBegin;
 88:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESSHELL, &flg));
 89:   if (flg) shell->ctx = ctx;
 90:   PetscFunctionReturn(PETSC_SUCCESS);
 91: }

 93: static PetscErrorCode SNESSolve_Shell(SNES snes)
 94: {
 95:   SNES_Shell *shell = (SNES_Shell *)snes->data;

 97:   PetscFunctionBegin;
 98:   PetscCheck(shell->solve, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESShellSetSolve() first");
 99:   snes->reason = SNES_CONVERGED_ITS;
100:   PetscCall((*shell->solve)(snes, snes->vec_sol));
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: static PetscErrorCode SNESShellSetSolve_Shell(SNES snes, PetscErrorCode (*solve)(SNES, Vec))
105: {
106:   SNES_Shell *shell = (SNES_Shell *)snes->data;

108:   PetscFunctionBegin;
109:   shell->solve = solve;
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: /*MC
114:   SNESSHELL - a user provided nonlinear solver

116:    Level: advanced

118: .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESShellGetContext()`, `SNESShellSetContext()`, `SNESShellSetSolve()`
119: M*/

121: PETSC_EXTERN PetscErrorCode SNESCreate_Shell(SNES snes)
122: {
123:   SNES_Shell *shell;

125:   PetscFunctionBegin;
126:   snes->ops->destroy = SNESDestroy_Shell;
127:   snes->ops->solve   = SNESSolve_Shell;

129:   snes->usesksp = PETSC_FALSE;
130:   snes->usesnpc = PETSC_FALSE;

132:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

134:   PetscCall(SNESParametersInitialize(snes));

136:   PetscCall(PetscNew(&shell));
137:   snes->data = (void *)shell;
138:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESShellSetSolve_C", SNESShellSetSolve_Shell));
139:   PetscFunctionReturn(PETSC_SUCCESS);
140: }