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