xref: /petsc/src/ts/interface/tscreate.c (revision 67c2884ef72a42a4f3a1babc25bcfa28f347172e)
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 
137adad957SLisandro Dalcin #if 0
14bdad233fSMatthew Knepley #undef __FUNCT__
15bdad233fSMatthew Knepley #define __FUNCT__ "TSPublish_Petsc"
166849ba73SBarry Smith static PetscErrorCode TSPublish_Petsc(PetscObject obj)
17bdad233fSMatthew Knepley {
18bdad233fSMatthew Knepley   PetscFunctionBegin;
19bdad233fSMatthew Knepley   PetscFunctionReturn(0);
20bdad233fSMatthew Knepley }
217adad957SLisandro Dalcin #endif
22bdad233fSMatthew Knepley 
23bdad233fSMatthew Knepley #undef  __FUNCT__
24bdad233fSMatthew Knepley #define __FUNCT__ "TSCreate"
2560893bc3SSatish Balay /*@C
26bd6a702fSBarry Smith   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
27bd6a702fSBarry Smith        type of solver can then be set with TSSetType().
28bdad233fSMatthew Knepley 
29bdad233fSMatthew Knepley   Collective on MPI_Comm
30bdad233fSMatthew Knepley 
31bdad233fSMatthew Knepley   Input Parameter:
32bdad233fSMatthew Knepley . comm - The communicator
33bdad233fSMatthew Knepley 
34bdad233fSMatthew Knepley   Output Parameter:
35bdad233fSMatthew Knepley . ts   - The TS
36bdad233fSMatthew Knepley 
37bdad233fSMatthew Knepley   Level: beginner
38bdad233fSMatthew Knepley 
39bdad233fSMatthew Knepley .keywords: TS, create
40bd6a702fSBarry Smith .seealso: TSSetType(), TSSetUp(), TSDestroy(), MeshCreate(), TSSetProblemType()
41bdad233fSMatthew Knepley @*/
420adebc6cSBarry Smith PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
430adebc6cSBarry Smith {
44bdad233fSMatthew Knepley   TS             t;
45dfbe8321SBarry Smith   PetscErrorCode ierr;
46bdad233fSMatthew Knepley 
47bdad233fSMatthew Knepley   PetscFunctionBegin;
484482741eSBarry Smith   PetscValidPointer(ts,1);
490298fd71SBarry Smith   *ts = NULL;
50519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
510298fd71SBarry Smith   ierr = TSInitializePackage(NULL);CHKERRQ(ierr);
52bdad233fSMatthew Knepley #endif
53bdad233fSMatthew Knepley 
54*67c2884eSBarry Smith   ierr = PetscHeaderCreate(t, _p_TS, struct _TSOps, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);CHKERRQ(ierr);
55bdad233fSMatthew Knepley   ierr = PetscMemzero(t->ops, sizeof(struct _TSOps));CHKERRQ(ierr);
56bdad233fSMatthew Knepley 
57bdad233fSMatthew Knepley   /* General TS description */
58089b2837SJed Brown   t->problem_type      = TS_NONLINEAR;
590298fd71SBarry Smith   t->vec_sol           = NULL;
60bdad233fSMatthew Knepley   t->numbermonitors    = 0;
610298fd71SBarry Smith   t->snes              = NULL;
62bdad233fSMatthew Knepley   t->setupcalled       = 0;
630298fd71SBarry Smith   t->data              = NULL;
640298fd71SBarry Smith   t->user              = NULL;
65cdbf8f93SLisandro Dalcin   t->ptime             = 0.0;
66cdbf8f93SLisandro Dalcin   t->time_step         = 0.1;
6731748224SBarry Smith   t->time_step_orig    = 0.1;
68bdad233fSMatthew Knepley   t->max_time          = 5.0;
69bdad233fSMatthew Knepley   t->steps             = 0;
70cdbf8f93SLisandro Dalcin   t->max_steps         = 5000;
715ef26d82SJed Brown   t->ksp_its           = 0;
725ef26d82SJed Brown   t->snes_its          = 0;
730298fd71SBarry Smith   t->work              = NULL;
74bdad233fSMatthew Knepley   t->nwork             = 0;
75193ac0bcSJed Brown   t->max_snes_failures = 1;
76193ac0bcSJed Brown   t->max_reject        = 10;
77193ac0bcSJed Brown   t->errorifstepfailed = PETSC_TRUE;
780e4ef248SJed Brown   t->rhsjacobian.time  = -1e20;
790026cea9SSean Farley   t->ijacobian.time    = -1e20;
80e817cc15SEmil Constantinescu   t->equation_type     = TS_EQ_UNSPECIFIED;
81bdad233fSMatthew Knepley 
821c3436cfSJed Brown   t->atol             = 1e-4;
831c3436cfSJed Brown   t->rtol             = 1e-4;
848d59e960SJed Brown   t->cfltime          = PETSC_MAX_REAL;
858d59e960SJed Brown   t->cfltime_local    = PETSC_MAX_REAL;
8649354f04SShri Abhyankar   t->exact_final_time = TS_EXACTFINALTIME_STEPOVER;
87e144a568SJed Brown 
88bdad233fSMatthew Knepley   *ts = t;
89bdad233fSMatthew Knepley   PetscFunctionReturn(0);
90bdad233fSMatthew Knepley }
91