Actual source code: tscreate.c

petsc-3.4.5 2014-06-29
  2: #include <petsc-private/tsimpl.h>      /*I "petscts.h"  I*/

  4: const char *const TSConvergedReasons_Shifted[] = {
  5:   "DIVERGED_STEP_REJECTED",
  6:   "DIVERGED_NONLINEAR_SOLVE",
  7:   "CONVERGED_ITERATING",
  8:   "CONVERGED_TIME",
  9:   "CONVERGED_ITS",
 10:   "TSConvergedReason","TS_",0};
 11: const char *const*TSConvergedReasons = TSConvergedReasons_Shifted + 2;

 13: #undef  __FUNCT__
 15: /*@C
 16:   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
 17:        type of solver can then be set with TSSetType().

 19:   Collective on MPI_Comm

 21:   Input Parameter:
 22: . comm - The communicator

 24:   Output Parameter:
 25: . ts   - The TS

 27:   Level: beginner

 29: .keywords: TS, create
 30: .seealso: TSSetType(), TSSetUp(), TSDestroy(), MeshCreate(), TSSetProblemType()
 31: @*/
 32: PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
 33: {
 34:   TS             t;

 39:   *ts = NULL;
 40: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
 41:   TSInitializePackage();
 42: #endif

 44:   PetscHeaderCreate(t, _p_TS, struct _TSOps, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);
 45:   PetscMemzero(t->ops, sizeof(struct _TSOps));

 47:   /* General TS description */
 48:   t->problem_type      = TS_NONLINEAR;
 49:   t->vec_sol           = NULL;
 50:   t->numbermonitors    = 0;
 51:   t->snes              = NULL;
 52:   t->setupcalled       = 0;
 53:   t->data              = NULL;
 54:   t->user              = NULL;
 55:   t->ptime             = 0.0;
 56:   t->time_step         = 0.1;
 57:   t->time_step_orig    = 0.1;
 58:   t->max_time          = 5.0;
 59:   t->steps             = 0;
 60:   t->max_steps         = 5000;
 61:   t->ksp_its           = 0;
 62:   t->snes_its          = 0;
 63:   t->work              = NULL;
 64:   t->nwork             = 0;
 65:   t->max_snes_failures = 1;
 66:   t->max_reject        = 10;
 67:   t->errorifstepfailed = PETSC_TRUE;
 68:   t->rhsjacobian.time  = -1e20;
 69:   t->rhsjacobian.scale = 1.;
 70:   t->ijacobian.shift   = 1.;
 71:   t->equation_type     = TS_EQ_UNSPECIFIED;

 73:   t->atol             = 1e-4;
 74:   t->rtol             = 1e-4;
 75:   t->cfltime          = PETSC_MAX_REAL;
 76:   t->cfltime_local    = PETSC_MAX_REAL;
 77:   t->exact_final_time = TS_EXACTFINALTIME_STEPOVER;

 79:   *ts = t;
 80:   return(0);
 81: }