xref: /petsc/src/ts/interface/tscreate.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2bdad233fSMatthew Knepley 
39371c9d4SSatish Balay 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};
422b60d99SBarry Smith const char *const *TSConvergedReasons = TSConvergedReasons_Shifted + 4;
5193ac0bcSJed Brown 
660893bc3SSatish Balay /*@C
7bd6a702fSBarry Smith   TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the
8bd6a702fSBarry Smith        type of solver can then be set with TSSetType().
9bdad233fSMatthew Knepley 
10d083f849SBarry Smith   Collective
11bdad233fSMatthew Knepley 
12bdad233fSMatthew Knepley   Input Parameter:
13bdad233fSMatthew Knepley . comm - The communicator
14bdad233fSMatthew Knepley 
15bdad233fSMatthew Knepley   Output Parameter:
16bdad233fSMatthew Knepley . ts   - The TS
17bdad233fSMatthew Knepley 
18bdad233fSMatthew Knepley   Level: beginner
19bdad233fSMatthew Knepley 
2095452b02SPatrick Sanan   Developer Notes:
2195452b02SPatrick Sanan     TS essentially always creates a SNES object even though explicit methods do not use it. This is
22efd4aadfSBarry Smith                     unfortunate and should be fixed at some point. The flag snes->usessnes indicates if the
23efd4aadfSBarry Smith                     particular method does use SNES and regulates if the information about the SNES is printed
24efd4aadfSBarry Smith                     in TSView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
25efd4aadfSBarry Smith                     by help messages about meaningless SNES options.
26efd4aadfSBarry Smith 
27db781477SPatrick Sanan .seealso: `TSSetType()`, `TSSetUp()`, `TSDestroy()`, `TSSetProblemType()`
28bdad233fSMatthew Knepley @*/
29*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSCreate(MPI_Comm comm, TS *ts)
30*d71ae5a4SJacob Faibussowitsch {
31bdad233fSMatthew Knepley   TS t;
32bdad233fSMatthew Knepley 
33bdad233fSMatthew Knepley   PetscFunctionBegin;
34064a246eSJacob Faibussowitsch   PetscValidPointer(ts, 2);
350298fd71SBarry Smith   *ts = NULL;
369566063dSJacob Faibussowitsch   PetscCall(TSInitializePackage());
37bdad233fSMatthew Knepley 
389566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", comm, TSDestroy, TSView));
39bdad233fSMatthew Knepley 
40bdad233fSMatthew Knepley   /* General TS description */
41089b2837SJed Brown   t->problem_type  = TS_NONLINEAR;
42ef85077eSLisandro Dalcin   t->equation_type = TS_EQ_UNSPECIFIED;
43ef85077eSLisandro Dalcin 
44cdbf8f93SLisandro Dalcin   t->ptime            = 0.0;
45cdbf8f93SLisandro Dalcin   t->time_step        = 0.1;
46ef85077eSLisandro Dalcin   t->max_time         = PETSC_MAX_REAL;
47ef85077eSLisandro Dalcin   t->exact_final_time = TS_EXACTFINALTIME_UNSPECIFIED;
48bdad233fSMatthew Knepley   t->steps            = 0;
49ef85077eSLisandro Dalcin   t->max_steps        = PETSC_MAX_INT;
50ef85077eSLisandro Dalcin   t->steprestart      = PETSC_TRUE;
51ef85077eSLisandro Dalcin 
52193ac0bcSJed Brown   t->max_snes_failures = 1;
53193ac0bcSJed Brown   t->max_reject        = 10;
54193ac0bcSJed Brown   t->errorifstepfailed = PETSC_TRUE;
55bdad233fSMatthew Knepley 
56ef85077eSLisandro Dalcin   t->rhsjacobian.time  = PETSC_MIN_REAL;
57ef85077eSLisandro Dalcin   t->rhsjacobian.scale = 1.0;
58ef85077eSLisandro Dalcin   t->ijacobian.shift   = 1.0;
5968bece0bSHong Zhang 
60b92453a8SLisandro Dalcin   /* All methods that do adaptivity should specify
61b92453a8SLisandro Dalcin    * its preferred adapt type in their constructor */
62b92453a8SLisandro Dalcin   t->default_adapt_type = TSADAPTNONE;
63ef85077eSLisandro Dalcin   t->atol               = 1e-4;
64ef85077eSLisandro Dalcin   t->rtol               = 1e-4;
65ef85077eSLisandro Dalcin   t->cfltime            = PETSC_MAX_REAL;
66ef85077eSLisandro Dalcin   t->cfltime_local      = PETSC_MAX_REAL;
67715f1b00SHong Zhang 
6887f4e208SHong Zhang   t->num_rhs_splits = 0;
6987f4e208SHong Zhang 
70d60b7d5cSBarry Smith   t->axpy_pattern = UNKNOWN_NONZERO_PATTERN;
71bdad233fSMatthew Knepley   *ts             = t;
72bdad233fSMatthew Knepley   PetscFunctionReturn(0);
73bdad233fSMatthew Knepley }
74