Actual source code: tsreg.c

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

  3: PetscFunctionList TSList              = NULL;
  4: PetscBool         TSRegisterAllCalled = PETSC_FALSE;

  6: /*@
  7:   TSSetType - Sets the algorithm/method to be used for integrating the ODE with the given `TS`.

  9:   Collective

 11:   Input Parameters:
 12: + ts   - The `TS` context
 13: - type - A known method

 15:   Options Database Key:
 16: . -ts_type <type> - Sets the method; use -help for a list of available methods (for instance, euler)

 18:   Level: intermediate

 20:   Notes:
 21:   See `TSType` for available methods (for instance)
 22: +  TSEULER - Euler
 23: .  TSBEULER - Backward Euler
 24: -  TSPSEUDO - Pseudo-timestepping

 26:   Normally, it is best to use the `TSSetFromOptions()` command and
 27:   then set the `TS` type from the options database rather than by using
 28:   this routine.  Using the options database provides the user with
 29:   maximum flexibility in evaluating the many different solvers.
 30:   The TSSetType() routine is provided for those situations where it
 31:   is necessary to set the timestepping solver independently of the
 32:   command line or options database.  This might be the case, for example,
 33:   when the choice of solver changes during the execution of the
 34:   program, and the user's application is taking responsibility for
 35:   choosing the appropriate method.  In other words, this routine is
 36:   not for beginners.

 38: .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSCreate()`, `TSSetFromOptions()`, `TSDestroy()`, `TSType`
 39: @*/
 40: PetscErrorCode TSSetType(TS ts, TSType type)
 41: {
 42:   PetscErrorCode (*r)(TS);
 43:   PetscBool match;

 45:   PetscFunctionBegin;
 47:   PetscAssertPointer(type, 2);
 48:   PetscCall(PetscObjectTypeCompare((PetscObject)ts, type, &match));
 49:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

 51:   PetscCall(PetscFunctionListFind(TSList, type, &r));
 52:   PetscCheck(r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TS type: %s", type);
 53:   PetscTryTypeMethod(ts, destroy);
 54:   PetscCall(PetscMemzero(ts->ops, sizeof(*ts->ops)));
 55:   ts->usessnes           = PETSC_FALSE;
 56:   ts->default_adapt_type = TSADAPTNONE;

 58:   ts->setupcalled = PETSC_FALSE;

 60:   PetscCall(PetscObjectChangeTypeName((PetscObject)ts, type));
 61:   PetscCall((*r)(ts));
 62:   PetscFunctionReturn(PETSC_SUCCESS);
 63: }

 65: /*@
 66:   TSGetType - Gets the `TS` method type (as a string) that is being used to solve the ODE with the given `TS`

 68:   Not Collective

 70:   Input Parameter:
 71: . ts - The `TS`

 73:   Output Parameter:
 74: . type - The `TSType`

 76:   Level: intermediate

 78: .seealso: [](ch_ts), `TS`, `TSType`, `TSSetType()`
 79: @*/
 80: PetscErrorCode TSGetType(TS ts, TSType *type)
 81: {
 82:   PetscFunctionBegin;
 84:   PetscAssertPointer(type, 2);
 85:   *type = ((PetscObject)ts)->type_name;
 86:   PetscFunctionReturn(PETSC_SUCCESS);
 87: }

 89: /*@C
 90:   TSRegister - Adds a creation method to the `TS` package.

 92:   Not Collective, No Fortran Support

 94:   Input Parameters:
 95: + sname    - The name of a new user-defined creation routine
 96: - function - The creation routine itself

 98:   Calling sequence of `function`:
 99: . ts - the `TS` being setup for the new `TSType` being registered

101:   Level: advanced

103:   Notes:
104:   `TSRegister()` may be called multiple times to add several user-defined tses.

106:   Example Usage:
107: .vb
108:   TSRegister("my_ts",  MyTSCreate);
109: .ve

111:   Then, your ts type can be chosen with the procedural interface via
112: .vb
113:     TS ts;
114:     TSCreate(MPI_Comm, &ts);
115:     TSSetType(ts, "my_ts")
116: .ve
117:   or at runtime via the option
118: .vb
119:     -ts_type my_ts
120: .ve

122: .seealso: [](ch_ts), `TSSetType()`, `TSType`, `TSRegisterAll()`, `TSRegisterDestroy()`
123: @*/
124: PetscErrorCode TSRegister(const char sname[], PetscErrorCode (*function)(TS ts))
125: {
126:   PetscFunctionBegin;
127:   PetscCall(TSInitializePackage());
128:   PetscCall(PetscFunctionListAdd(&TSList, sname, function));
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }