xref: /petsc/src/ts/interface/tscreate.c (revision 69f65d41fd6b4b8a63ef21431c8781fe8a588e09)
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 /*@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().
20 
21   Collective on MPI_Comm
22 
23   Input Parameter:
24 . comm - The communicator
25 
26   Output Parameter:
27 . ts   - The TS
28 
29   Level: beginner
30 
31 .keywords: TS, create
32 .seealso: TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
33 @*/
34 PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
35 {
36   TS             t;
37   PetscErrorCode ierr;
38 
39   PetscFunctionBegin;
40   PetscValidPointer(ts,1);
41   *ts = NULL;
42   ierr = TSInitializePackage();CHKERRQ(ierr);
43 
44   ierr = PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);CHKERRQ(ierr);
45 
46   /* General TS description */
47   t->problem_type      = TS_NONLINEAR;
48   t->vec_sol           = NULL;
49   t->numbermonitors    = 0;
50   t->snes              = NULL;
51   t->setupcalled       = 0;
52   t->data              = NULL;
53   t->user              = NULL;
54   t->ptime             = 0.0;
55   t->time_step         = 0.1;
56   t->max_time          = 5.0;
57   t->steprollback      = PETSC_FALSE;
58   t->steprestart       = PETSC_FALSE;
59   t->steps             = 0;
60   t->max_steps         = 5000;
61   t->ksp_its           = 0;
62   t->snes_its          = 0;
63   t->work              = NULL;
64   t->nwork             = 0;
65   t->max_snes_failures = 1;
66   t->max_reject        = 10;
67   t->errorifstepfailed = PETSC_TRUE;
68   t->rhsjacobian.time  = -1e20;
69   t->rhsjacobian.scale = 1.;
70   t->ijacobian.shift   = 1.;
71   t->equation_type     = TS_EQ_UNSPECIFIED;
72 
73   t->atol             = 1e-4;
74   t->rtol             = 1e-4;
75   t->cfltime          = PETSC_MAX_REAL;
76   t->cfltime_local    = PETSC_MAX_REAL;
77   t->exact_final_time = TS_EXACTFINALTIME_UNSPECIFIED;
78   t->vec_costintegral = NULL;
79   t->trajectory       = NULL;
80 
81   /* All methods that do adaptivity should specify
82    * its preferred adapt type in their constructor */
83   t->default_adapt_type = TSADAPTNONE;
84 
85   *ts = t;
86   PetscFunctionReturn(0);
87 }
88