xref: /petsc/src/ts/interface/tscreate.c (revision eb2eea041fe0c038f64a45e6c11b98f844c78178)
1 
2 #include <petsc/private/tsimpl.h>      /*I "petscts.h"  I*/
3 
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;
16 
17 #undef  __FUNCT__
18 #define __FUNCT__ "TSCreate"
19 /*@C
20   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
21        type of solver can then be set with TSSetType().
22 
23   Collective on MPI_Comm
24 
25   Input Parameter:
26 . comm - The communicator
27 
28   Output Parameter:
29 . ts   - The TS
30 
31   Level: beginner
32 
33 .keywords: TS, create
34 .seealso: TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
35 @*/
36 PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
37 {
38   TS             t;
39   PetscErrorCode ierr;
40 
41   PetscFunctionBegin;
42   PetscValidPointer(ts,1);
43   *ts = NULL;
44   ierr = TSInitializePackage();CHKERRQ(ierr);
45 
46   ierr = PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);CHKERRQ(ierr);
47 
48   /* General TS description */
49   t->problem_type      = TS_NONLINEAR;
50   t->vec_sol           = NULL;
51   t->numbermonitors    = 0;
52   t->snes              = NULL;
53   t->setupcalled       = 0;
54   t->data              = NULL;
55   t->user              = NULL;
56   t->ptime             = 0.0;
57   t->time_step         = 0.1;
58   t->time_step_orig    = 0.1;
59   t->max_time          = 5.0;
60   t->steprollback      = PETSC_FALSE;
61   t->steps             = 0;
62   t->max_steps         = 5000;
63   t->ksp_its           = 0;
64   t->snes_its          = 0;
65   t->work              = NULL;
66   t->nwork             = 0;
67   t->max_snes_failures = 1;
68   t->max_reject        = 10;
69   t->errorifstepfailed = PETSC_TRUE;
70   t->rhsjacobian.time  = -1e20;
71   t->rhsjacobian.scale = 1.;
72   t->ijacobian.shift   = 1.;
73   t->equation_type     = TS_EQ_UNSPECIFIED;
74 
75   t->atol             = 1e-4;
76   t->rtol             = 1e-4;
77   t->cfltime          = PETSC_MAX_REAL;
78   t->cfltime_local    = PETSC_MAX_REAL;
79   t->exact_final_time = TS_EXACTFINALTIME_STEPOVER;
80 
81   t->trajectory       = NULL;
82 
83   *ts = t;
84   PetscFunctionReturn(0);
85 }
86