Actual source code: tscreate.c
petsc-3.7.3 2016-08-01
2: #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
4: const char *const TSConvergedReasons_Shifted[] = {
5: "DIVERGED_STEP_REJECTED",
6: "DIVERGED_NONLINEAR_SOLVE",
7: "CONVERGED_ITERATING",
8: "CONVERGED_TIME",
9: "CONVERGED_ITS",
10: "CONVERGED_USER",
11: "CONVERGED_EVENT",
12: "CONVERGED_PSEUDO_FATOL",
13: "CONVERGED_PSEUDO_FATOL",
14: "TSConvergedReason","TS_",0};
15: const char *const*TSConvergedReasons = TSConvergedReasons_Shifted + 2;
17: #undef __FUNCT__
19: /*@C
20: TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
21: type of solver can then be set with TSSetType().
23: Collective on MPI_Comm
25: Input Parameter:
26: . comm - The communicator
28: Output Parameter:
29: . ts - The TS
31: Level: beginner
33: .keywords: TS, create
34: .seealso: TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
35: @*/
36: PetscErrorCode TSCreate(MPI_Comm comm, TS *ts)
37: {
38: TS t;
43: *ts = NULL;
44: TSInitializePackage();
46: PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView);
48: /* General TS description */
49: t->problem_type = TS_NONLINEAR;
50: t->vec_sol = NULL;
51: t->numbermonitors = 0;
52: t->snes = NULL;
53: t->setupcalled = 0;
54: t->data = NULL;
55: t->user = NULL;
56: t->ptime = 0.0;
57: t->time_step = 0.1;
58: t->max_time = 5.0;
59: t->steprollback = PETSC_FALSE;
60: t->steps = 0;
61: t->max_steps = 5000;
62: t->ksp_its = 0;
63: t->snes_its = 0;
64: t->work = NULL;
65: t->nwork = 0;
66: t->max_snes_failures = 1;
67: t->max_reject = 10;
68: t->errorifstepfailed = PETSC_TRUE;
69: t->rhsjacobian.time = -1e20;
70: t->rhsjacobian.scale = 1.;
71: t->ijacobian.shift = 1.;
72: t->equation_type = TS_EQ_UNSPECIFIED;
74: t->atol = 1e-4;
75: t->rtol = 1e-4;
76: t->cfltime = PETSC_MAX_REAL;
77: t->cfltime_local = PETSC_MAX_REAL;
78: t->exact_final_time = TS_EXACTFINALTIME_UNSPECIFIED;
79: t->vec_costintegral = NULL;
80: t->trajectory = NULL;
82: /* All methods that do not do adaptivity at all
83: * should delete this object in their constructor */
84: TSGetAdapt(t,&t->adapt);
86: *ts = t;
87: return(0);
88: }