xref: /petsc/src/ts/interface/tscreate.c (revision b41af12e0cf5c8cfcc61fcc96f5ea4cd576f4469)
1bdad233fSMatthew Knepley 
2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>      /*I "petscts.h"  I*/
3bdad233fSMatthew Knepley 
4193ac0bcSJed Brown const char *const TSConvergedReasons_Shifted[] = {
5193ac0bcSJed Brown   "DIVERGED_STEP_REJECTED",
6193ac0bcSJed Brown   "DIVERGED_NONLINEAR_SOLVE",
7193ac0bcSJed Brown   "CONVERGED_ITERATING",
8193ac0bcSJed Brown   "CONVERGED_TIME",
9193ac0bcSJed Brown   "CONVERGED_ITS",
10193ac0bcSJed Brown   "TSConvergedReason","TS_",0};
11193ac0bcSJed Brown const char *const*TSConvergedReasons = TSConvergedReasons_Shifted + 2;
12193ac0bcSJed Brown 
13bdad233fSMatthew Knepley #undef  __FUNCT__
14bdad233fSMatthew Knepley #define __FUNCT__ "TSCreate"
1560893bc3SSatish Balay /*@C
16bd6a702fSBarry Smith   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
17bd6a702fSBarry Smith        type of solver can then be set with TSSetType().
18bdad233fSMatthew Knepley 
19bdad233fSMatthew Knepley   Collective on MPI_Comm
20bdad233fSMatthew Knepley 
21bdad233fSMatthew Knepley   Input Parameter:
22bdad233fSMatthew Knepley . comm - The communicator
23bdad233fSMatthew Knepley 
24bdad233fSMatthew Knepley   Output Parameter:
25bdad233fSMatthew Knepley . ts   - The TS
26bdad233fSMatthew Knepley 
27bdad233fSMatthew Knepley   Level: beginner
28bdad233fSMatthew Knepley 
29bdad233fSMatthew Knepley .keywords: TS, create
30bd6a702fSBarry Smith .seealso: TSSetType(), TSSetUp(), TSDestroy(), MeshCreate(), TSSetProblemType()
31bdad233fSMatthew Knepley @*/
320adebc6cSBarry Smith PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
330adebc6cSBarry Smith {
34bdad233fSMatthew Knepley   TS             t;
35dfbe8321SBarry Smith   PetscErrorCode ierr;
36bdad233fSMatthew Knepley 
37bdad233fSMatthew Knepley   PetscFunctionBegin;
384482741eSBarry Smith   PetscValidPointer(ts,1);
390298fd71SBarry Smith   *ts = NULL;
40519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
41607a6623SBarry Smith   ierr = TSInitializePackage();CHKERRQ(ierr);
42bdad233fSMatthew Knepley #endif
43bdad233fSMatthew Knepley 
4467c2884eSBarry Smith   ierr = PetscHeaderCreate(t, _p_TS, struct _TSOps, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);CHKERRQ(ierr);
45bdad233fSMatthew Knepley   ierr = PetscMemzero(t->ops, sizeof(struct _TSOps));CHKERRQ(ierr);
46bdad233fSMatthew Knepley 
47bdad233fSMatthew Knepley   /* General TS description */
48089b2837SJed Brown   t->problem_type      = TS_NONLINEAR;
490298fd71SBarry Smith   t->vec_sol           = NULL;
50bdad233fSMatthew Knepley   t->numbermonitors    = 0;
510298fd71SBarry Smith   t->snes              = NULL;
52bdad233fSMatthew Knepley   t->setupcalled       = 0;
530298fd71SBarry Smith   t->data              = NULL;
540298fd71SBarry Smith   t->user              = NULL;
55cdbf8f93SLisandro Dalcin   t->ptime             = 0.0;
56cdbf8f93SLisandro Dalcin   t->time_step         = 0.1;
5731748224SBarry Smith   t->time_step_orig    = 0.1;
58bdad233fSMatthew Knepley   t->max_time          = 5.0;
59bdad233fSMatthew Knepley   t->steps             = 0;
60cdbf8f93SLisandro Dalcin   t->max_steps         = 5000;
615ef26d82SJed Brown   t->ksp_its           = 0;
625ef26d82SJed Brown   t->snes_its          = 0;
630298fd71SBarry Smith   t->work              = NULL;
64bdad233fSMatthew Knepley   t->nwork             = 0;
65193ac0bcSJed Brown   t->max_snes_failures = 1;
66193ac0bcSJed Brown   t->max_reject        = 10;
67193ac0bcSJed Brown   t->errorifstepfailed = PETSC_TRUE;
680e4ef248SJed Brown   t->rhsjacobian.time  = -1e20;
69e1244c69SJed Brown   t->rhsjacobian.scale = 1.;
70*b41af12eSJed Brown   t->ijacobian.shift   = 1.;
71e817cc15SEmil Constantinescu   t->equation_type     = TS_EQ_UNSPECIFIED;
72bdad233fSMatthew Knepley 
731c3436cfSJed Brown   t->atol             = 1e-4;
741c3436cfSJed Brown   t->rtol             = 1e-4;
758d59e960SJed Brown   t->cfltime          = PETSC_MAX_REAL;
768d59e960SJed Brown   t->cfltime_local    = PETSC_MAX_REAL;
7749354f04SShri Abhyankar   t->exact_final_time = TS_EXACTFINALTIME_STEPOVER;
78e144a568SJed Brown 
79bdad233fSMatthew Knepley   *ts = t;
80bdad233fSMatthew Knepley   PetscFunctionReturn(0);
81bdad233fSMatthew Knepley }
82