xref: /petsc/src/ts/interface/tscreate.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 #include <petsc/private/tsimpl.h>      /*I "petscts.h"  I*/
2 
3 const char *const TSConvergedReasons_Shifted[] = {
4   "ADJOINT_DIVERGED_LINEAR_SOLVE",
5   "FORWARD_DIVERGED_LINEAR_SOLVE",
6   "DIVERGED_STEP_REJECTED",
7   "DIVERGED_NONLINEAR_SOLVE",
8   "CONVERGED_ITERATING",
9   "CONVERGED_TIME",
10   "CONVERGED_ITS",
11   "CONVERGED_USER",
12   "CONVERGED_EVENT",
13   "CONVERGED_PSEUDO_FATOL",
14   "CONVERGED_PSEUDO_FATOL",
15   "TSConvergedReason","TS_",NULL};
16 const char *const*TSConvergedReasons = TSConvergedReasons_Shifted + 4;
17 
18 /*@C
19   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
20        type of solver can then be set with TSSetType().
21 
22   Collective
23 
24   Input Parameter:
25 . comm - The communicator
26 
27   Output Parameter:
28 . ts   - The TS
29 
30   Level: beginner
31 
32   Developer Notes:
33     TS essentially always creates a SNES object even though explicit methods do not use it. This is
34                     unfortunate and should be fixed at some point. The flag snes->usessnes indicates if the
35                     particular method does use SNES and regulates if the information about the SNES is printed
36                     in TSView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
37                     by help messages about meaningless SNES options.
38 
39 .seealso: TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
40 @*/
41 PetscErrorCode  TSCreate(MPI_Comm comm, TS *ts)
42 {
43   TS             t;
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   PetscValidPointer(ts,2);
48   *ts = NULL;
49   ierr = TSInitializePackage();CHKERRQ(ierr);
50 
51   ierr = PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);CHKERRQ(ierr);
52 
53   /* General TS description */
54   t->problem_type      = TS_NONLINEAR;
55   t->equation_type     = TS_EQ_UNSPECIFIED;
56 
57   t->ptime             = 0.0;
58   t->time_step         = 0.1;
59   t->max_time          = PETSC_MAX_REAL;
60   t->exact_final_time  = TS_EXACTFINALTIME_UNSPECIFIED;
61   t->steps             = 0;
62   t->max_steps         = PETSC_MAX_INT;
63   t->steprestart       = PETSC_TRUE;
64 
65   t->max_snes_failures = 1;
66   t->max_reject        = 10;
67   t->errorifstepfailed = PETSC_TRUE;
68 
69   t->rhsjacobian.time  = PETSC_MIN_REAL;
70   t->rhsjacobian.scale = 1.0;
71   t->ijacobian.shift   = 1.0;
72 
73   /* All methods that do adaptivity should specify
74    * its preferred adapt type in their constructor */
75   t->default_adapt_type = TSADAPTNONE;
76   t->atol               = 1e-4;
77   t->rtol               = 1e-4;
78   t->cfltime            = PETSC_MAX_REAL;
79   t->cfltime_local      = PETSC_MAX_REAL;
80 
81   t->num_rhs_splits     = 0;
82 
83   t->axpy_pattern       = UNKNOWN_NONZERO_PATTERN;
84   *ts = t;
85   PetscFunctionReturn(0);
86 }
87