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