Actual source code: tscreate.c

petsc-master 2017-07-18
Report Typos and Errors

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

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

 17: /*@C
 18:   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
 19:        type of solver can then be set with TSSetType().

 21:   Collective on MPI_Comm

 23:   Input Parameter:
 24: . comm - The communicator

 26:   Output Parameter:
 27: . ts   - The TS

 29:   Level: beginner

 31:   Developer Notes:  TS essentially always creates a SNES object even though explicit methods do not use it. This is
 32:                     unfortunate and should be fixed at some point. The flag snes->usessnes indicates if the
 33:                     particular method does use SNES and regulates if the information about the SNES is printed
 34:                     in TSView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
 35:                     by help messages about meaningless SNES options.

 37: .keywords: TS, create
 38: .seealso: TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
 39: @*/
 40: PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
 41: {
 42:   TS             t;

 47:   *ts = NULL;
 48:   TSInitializePackage();

 50:   PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);

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

 79:   t->atol             = 1e-4;
 80:   t->rtol             = 1e-4;
 81:   t->cfltime          = PETSC_MAX_REAL;
 82:   t->cfltime_local    = PETSC_MAX_REAL;
 83:   t->exact_final_time = TS_EXACTFINALTIME_UNSPECIFIED;
 84:   t->vec_costintegral = NULL;
 85:   t->trajectory       = NULL;

 87:   /* All methods that do adaptivity should specify
 88:    * its preferred adapt type in their constructor */
 89:   t->default_adapt_type = TSADAPTNONE;

 91:   *ts = t;
 92:   return(0);
 93: }