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