1bdad233fSMatthew Knepley 2c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 3bdad233fSMatthew Knepley 4193ac0bcSJed Brown const char *const TSConvergedReasons_Shifted[] = { 5193ac0bcSJed Brown "DIVERGED_STEP_REJECTED", 6193ac0bcSJed Brown "DIVERGED_NONLINEAR_SOLVE", 7193ac0bcSJed Brown "CONVERGED_ITERATING", 8193ac0bcSJed Brown "CONVERGED_TIME", 9193ac0bcSJed Brown "CONVERGED_ITS", 10193ac0bcSJed Brown "TSConvergedReason","TS_",0}; 11193ac0bcSJed Brown const char *const*TSConvergedReasons = TSConvergedReasons_Shifted + 2; 12193ac0bcSJed Brown 137adad957SLisandro Dalcin #if 0 14bdad233fSMatthew Knepley #undef __FUNCT__ 15bdad233fSMatthew Knepley #define __FUNCT__ "TSPublish_Petsc" 166849ba73SBarry Smith static PetscErrorCode TSPublish_Petsc(PetscObject obj) 17bdad233fSMatthew Knepley { 18bdad233fSMatthew Knepley PetscFunctionBegin; 19bdad233fSMatthew Knepley PetscFunctionReturn(0); 20bdad233fSMatthew Knepley } 217adad957SLisandro Dalcin #endif 22bdad233fSMatthew Knepley 23bdad233fSMatthew Knepley #undef __FUNCT__ 24bdad233fSMatthew Knepley #define __FUNCT__ "TSCreate" 2560893bc3SSatish Balay /*@C 26bd6a702fSBarry Smith TSCreate - This function creates an empty timestepper. The problem type can then be set with TSSetProblemType() and the 27bd6a702fSBarry Smith type of solver can then be set with TSSetType(). 28bdad233fSMatthew Knepley 29bdad233fSMatthew Knepley Collective on MPI_Comm 30bdad233fSMatthew Knepley 31bdad233fSMatthew Knepley Input Parameter: 32bdad233fSMatthew Knepley . comm - The communicator 33bdad233fSMatthew Knepley 34bdad233fSMatthew Knepley Output Parameter: 35bdad233fSMatthew Knepley . ts - The TS 36bdad233fSMatthew Knepley 37bdad233fSMatthew Knepley Level: beginner 38bdad233fSMatthew Knepley 39bdad233fSMatthew Knepley .keywords: TS, create 40bd6a702fSBarry Smith .seealso: TSSetType(), TSSetUp(), TSDestroy(), MeshCreate(), TSSetProblemType() 41bdad233fSMatthew Knepley @*/ 427087cfbeSBarry Smith PetscErrorCode TSCreate(MPI_Comm comm, TS *ts) { 43bdad233fSMatthew Knepley TS t; 44dfbe8321SBarry Smith PetscErrorCode ierr; 45bdad233fSMatthew Knepley 46bdad233fSMatthew Knepley PetscFunctionBegin; 474482741eSBarry Smith PetscValidPointer(ts,1); 48bdad233fSMatthew Knepley *ts = PETSC_NULL; 49bdad233fSMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 50bdad233fSMatthew Knepley ierr = TSInitializePackage(PETSC_NULL);CHKERRQ(ierr); 51bdad233fSMatthew Knepley #endif 52bdad233fSMatthew Knepley 530700a824SBarry Smith ierr = PetscHeaderCreate(t, _p_TS, struct _TSOps, TS_CLASSID, -1, "TS", comm, TSDestroy, TSView);CHKERRQ(ierr); 54bdad233fSMatthew Knepley ierr = PetscMemzero(t->ops, sizeof(struct _TSOps));CHKERRQ(ierr); 55bdad233fSMatthew Knepley 56bdad233fSMatthew Knepley /* General TS description */ 57089b2837SJed Brown t->problem_type = TS_NONLINEAR; 58bdad233fSMatthew Knepley t->vec_sol = PETSC_NULL; 59bdad233fSMatthew Knepley t->numbermonitors = 0; 60bdad233fSMatthew Knepley t->snes = PETSC_NULL; 61bdad233fSMatthew Knepley t->funP = PETSC_NULL; 62bdad233fSMatthew Knepley t->jacP = PETSC_NULL; 63bdad233fSMatthew Knepley t->setupcalled = 0; 64bdad233fSMatthew Knepley t->data = PETSC_NULL; 65bdad233fSMatthew Knepley t->user = PETSC_NULL; 66bdad233fSMatthew Knepley t->max_steps = 5000; 67bdad233fSMatthew Knepley t->max_time = 5.0; 68bdad233fSMatthew Knepley t->steps = 0; 69bdad233fSMatthew Knepley t->linear_its = 0; 70bdad233fSMatthew Knepley t->nonlinear_its = 0; 71bdad233fSMatthew Knepley t->work = PETSC_NULL; 72bdad233fSMatthew Knepley t->nwork = 0; 73193ac0bcSJed Brown t->max_snes_failures = 1; 74193ac0bcSJed Brown t->max_reject = 10; 75193ac0bcSJed Brown t->errorifstepfailed = PETSC_TRUE; 76*0e4ef248SJed Brown t->rhsjacobian.time = -1e20; 77bdad233fSMatthew Knepley 78e144a568SJed Brown ierr = TSSetInitialTimeStep(t,0.,0.1);CHKERRQ(ierr); 79e144a568SJed Brown 80bdad233fSMatthew Knepley *ts = t; 81bdad233fSMatthew Knepley PetscFunctionReturn(0); 82bdad233fSMatthew Knepley } 83