Actual source code: tsreg.c

petsc-3.3-p7 2013-05-11
  1: #include <petsc-private/tsimpl.h>      /*I "petscts.h"  I*/

  3: PetscFList TSList                       = PETSC_NULL;
  4: PetscBool  TSRegisterAllCalled          = PETSC_FALSE;

  8: /*@C
  9:   TSSetType - Sets the method for the timestepping solver.

 11:   Collective on TS

 13:   Input Parameters:
 14: + ts   - The TS context
 15: - type - A known method

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

 20:    Notes:
 21:    See "petsc/include/petscts.h" for available methods (for instance)
 22: +  TSEULER - Euler
 23: .  TSSUNDIALS - SUNDIALS interface
 24: .  TSBEULER - Backward Euler
 25: -  TSPSEUDO - Pseudo-timestepping

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

 39:    Level: intermediate

 41: .keywords: TS, set, type

 43: @*/
 44: PetscErrorCode  TSSetType(TS ts,const TSType type)
 45: {
 46:   PetscErrorCode (*r)(TS);
 47:   PetscBool      match;

 52:   PetscObjectTypeCompare((PetscObject) ts, type, &match);
 53:   if (match) return(0);

 55:   PetscFListFind( TSList,((PetscObject)ts)->comm, type,PETSC_TRUE, (void (**)(void)) &r);
 56:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TS type: %s", type);
 57:   if (ts->ops->destroy) {
 58:     (*(ts)->ops->destroy)(ts);
 59:     ts->ops->destroy = PETSC_NULL;
 60:   }
 61:   PetscMemzero(ts->ops,sizeof(*ts->ops));
 62:   ts->setupcalled = PETSC_FALSE;
 63:   PetscObjectChangeTypeName((PetscObject)ts, type);
 64:   (*r)(ts);
 65: #if defined(PETSC_HAVE_AMS)
 66:   if (PetscAMSPublishAll) {
 67:     PetscObjectAMSPublish((PetscObject)ts);
 68:   }
 69: #endif
 70:   return(0);
 71: }

 75: /*@C
 76:   TSGetType - Gets the TS method type (as a string).

 78:   Not Collective

 80:   Input Parameter:
 81: . ts - The TS

 83:   Output Parameter:
 84: . type - The name of TS method

 86:   Level: intermediate

 88: .keywords: TS, timestepper, get, type, name
 89: .seealso TSSetType()
 90: @*/
 91: PetscErrorCode  TSGetType(TS ts, const TSType *type)
 92: {
 96:   *type = ((PetscObject)ts)->type_name;
 97:   return(0);
 98: }

100: /*--------------------------------------------------------------------------------------------------------------------*/

104: /*@C
105:   TSRegister - See TSRegisterDynamic()

107:   Level: advanced
108: @*/
109: PetscErrorCode  TSRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(TS))
110: {
111:   char           fullname[PETSC_MAX_PATH_LEN];

115:   PetscStrcpy(fullname, path);
116:   PetscStrcat(fullname, ":");
117:   PetscStrcat(fullname, name);
118:   PetscFListAdd(&TSList, sname, fullname, (void (*)(void)) function);
119:   return(0);
120: }

122: /*-------------------------------------------------------------------------------------------------------------------*/
125: /*@C
126:    TSRegisterDestroy - Frees the list of timestepping routines that were registered by TSRegister()/TSRegisterDynamic().

128:    Not Collective

130:    Level: advanced

132: .keywords: TS, timestepper, register, destroy
133: .seealso: TSRegister(), TSRegisterAll(), TSRegisterDynamic()
134: @*/
135: PetscErrorCode  TSRegisterDestroy(void)
136: {

140:   PetscFListDestroy(&TSList);
141:   TSRegisterAllCalled = PETSC_FALSE;
142:   return(0);
143: }