Actual source code: snessaws.c

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

  4: typedef struct {
  5:   PetscViewer viewer;
  6: } SNESMonitor_SAWs;

  8: /*@C
  9:   SNESMonitorSAWsCreate - create an SAWs monitor context for `SNES`

 11:   Collective

 13:   Input Parameter:
 14: . snes - `SNES` to monitor

 16:   Output Parameter:
 17: . ctx - context for monitor

 19:   Level: developer

 21: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNES`, `SNESMonitorSAWs()`, `SNESMonitorSAWsDestroy()`
 22: @*/
 23: PetscErrorCode SNESMonitorSAWsCreate(SNES snes, void **ctx)
 24: {
 25:   SNESMonitor_SAWs *mon;

 27:   PetscFunctionBegin;
 28:   PetscCall(PetscNew(&mon));
 29:   mon->viewer = PETSC_VIEWER_SAWS_(PetscObjectComm((PetscObject)snes));
 30:   PetscCheck(mon->viewer, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Cannot create SAWs default viewer");
 31:   *ctx = (void *)mon;
 32:   PetscFunctionReturn(PETSC_SUCCESS);
 33: }

 35: /*@C
 36:   SNESMonitorSAWsDestroy - destroy a monitor context created with `SNESMonitorSAWsCreate()`

 38:   Collective

 40:   Input Parameter:
 41: . ctx - monitor context

 43:   Level: developer

 45: .seealso: [](ch_snes), `SNESMonitorSAWsCreate()`
 46: @*/
 47: PetscErrorCode SNESMonitorSAWsDestroy(void **ctx)
 48: {
 49:   PetscFunctionBegin;
 50:   PetscCall(PetscFree(*ctx));
 51:   PetscFunctionReturn(PETSC_SUCCESS);
 52: }

 54: /*@C
 55:   SNESMonitorSAWs - monitor solution process of `SNES` using SAWs

 57:   Collective

 59:   Input Parameters:
 60: + snes  - iterative context
 61: . n     - iteration number
 62: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
 63: - ctx   - `PetscViewer` of type `PETSCVIEWERSAWS`

 65:   Level: advanced

 67: .seealso: [](ch_snes), `PetscViewerSAWsOpen()`, `SNESMonitorSAWsDestroy()`, `SNESMonitorSAWsCreate()`
 68: @*/
 69: PetscErrorCode SNESMonitorSAWs(SNES snes, PetscInt n, PetscReal rnorm, void *ctx)
 70: {
 71:   PetscMPIInt rank;

 73:   PetscFunctionBegin;

 76:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
 77:   if (rank == 0) {
 78:     PetscCallSAWs(SAWs_Register, ("/PETSc/snes_monitor_saws/its", &snes->iter, 1, SAWs_READ, SAWs_INT));
 79:     PetscCallSAWs(SAWs_Register, ("/PETSc/snes_monitor_saws/rnorm", &snes->norm, 1, SAWs_READ, SAWs_DOUBLE));
 80:     PetscCall(PetscSAWsBlock());
 81:   }
 82:   PetscFunctionReturn(PETSC_SUCCESS);
 83: }