xref: /petsc/src/ts/interface/tscreate.c (revision feed9e9da6c2d07546ffbbe4cafa0424dc138dae)
1bdad233fSMatthew Knepley 
2af0996ceSBarry Smith #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",
103118ae5eSBarry Smith   "CONVERGED_USER",
113118ae5eSBarry Smith   "CONVERGED_EVENT",
123118ae5eSBarry Smith   "CONVERGED_PSEUDO_FATOL",
133118ae5eSBarry Smith   "CONVERGED_PSEUDO_FATOL",
14193ac0bcSJed Brown   "TSConvergedReason","TS_",0};
15193ac0bcSJed Brown const char *const*TSConvergedReasons = TSConvergedReasons_Shifted + 2;
16193ac0bcSJed Brown 
17bdad233fSMatthew Knepley #undef  __FUNCT__
18bdad233fSMatthew Knepley #define __FUNCT__ "TSCreate"
1960893bc3SSatish Balay /*@C
20bd6a702fSBarry Smith   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
21bd6a702fSBarry Smith        type of solver can then be set with TSSetType().
22bdad233fSMatthew Knepley 
23bdad233fSMatthew Knepley   Collective on MPI_Comm
24bdad233fSMatthew Knepley 
25bdad233fSMatthew Knepley   Input Parameter:
26bdad233fSMatthew Knepley . comm - The communicator
27bdad233fSMatthew Knepley 
28bdad233fSMatthew Knepley   Output Parameter:
29bdad233fSMatthew Knepley . ts   - The TS
30bdad233fSMatthew Knepley 
31bdad233fSMatthew Knepley   Level: beginner
32bdad233fSMatthew Knepley 
33bdad233fSMatthew Knepley .keywords: TS, create
3404ceabe5SMatthew G. Knepley .seealso: TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
35bdad233fSMatthew Knepley @*/
360adebc6cSBarry Smith PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
370adebc6cSBarry Smith {
38bdad233fSMatthew Knepley   TS             t;
39dfbe8321SBarry Smith   PetscErrorCode ierr;
40bdad233fSMatthew Knepley 
41bdad233fSMatthew Knepley   PetscFunctionBegin;
424482741eSBarry Smith   PetscValidPointer(ts,1);
430298fd71SBarry Smith   *ts = NULL;
44607a6623SBarry Smith   ierr = TSInitializePackage();CHKERRQ(ierr);
45bdad233fSMatthew Knepley 
4673107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);CHKERRQ(ierr);
47bdad233fSMatthew Knepley 
48bdad233fSMatthew Knepley   /* General TS description */
49089b2837SJed Brown   t->problem_type      = TS_NONLINEAR;
500298fd71SBarry Smith   t->vec_sol           = NULL;
51bdad233fSMatthew Knepley   t->numbermonitors    = 0;
520298fd71SBarry Smith   t->snes              = NULL;
53bdad233fSMatthew Knepley   t->setupcalled       = 0;
540298fd71SBarry Smith   t->data              = NULL;
550298fd71SBarry Smith   t->user              = NULL;
56cdbf8f93SLisandro Dalcin   t->ptime             = 0.0;
57cdbf8f93SLisandro Dalcin   t->time_step         = 0.1;
5831748224SBarry Smith   t->time_step_orig    = 0.1;
59bdad233fSMatthew Knepley   t->max_time          = 5.0;
608392e04aSShri Abhyankar   t->steprollback      = PETSC_FALSE;
61bdad233fSMatthew Knepley   t->steps             = 0;
62cdbf8f93SLisandro Dalcin   t->max_steps         = 5000;
635ef26d82SJed Brown   t->ksp_its           = 0;
645ef26d82SJed Brown   t->snes_its          = 0;
650298fd71SBarry Smith   t->work              = NULL;
66bdad233fSMatthew Knepley   t->nwork             = 0;
67193ac0bcSJed Brown   t->max_snes_failures = 1;
68193ac0bcSJed Brown   t->max_reject        = 10;
69193ac0bcSJed Brown   t->errorifstepfailed = PETSC_TRUE;
700e4ef248SJed Brown   t->rhsjacobian.time  = -1e20;
71e1244c69SJed Brown   t->rhsjacobian.scale = 1.;
72b41af12eSJed Brown   t->ijacobian.shift   = 1.;
73e817cc15SEmil Constantinescu   t->equation_type     = TS_EQ_UNSPECIFIED;
74bdad233fSMatthew Knepley 
751c3436cfSJed Brown   t->atol             = 1e-4;
761c3436cfSJed Brown   t->rtol             = 1e-4;
778d59e960SJed Brown   t->cfltime          = PETSC_MAX_REAL;
788d59e960SJed Brown   t->cfltime_local    = PETSC_MAX_REAL;
79*feed9e9dSBarry Smith   t->exact_final_time = TS_EXACTFINALTIME_UNSPECIFIED;
80814b06c0SHong Zhang   t->vec_costintegral = NULL;
8168bece0bSHong Zhang   t->trajectory       = NULL;
8268bece0bSHong Zhang 
83bdad233fSMatthew Knepley   *ts = t;
84bdad233fSMatthew Knepley   PetscFunctionReturn(0);
85bdad233fSMatthew Knepley }
86