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