163dd3a1aSKris Buschelman #define PETSCTS_DLL 263dd3a1aSKris Buschelman 37c4f633dSBarry Smith #include "private/tsimpl.h" /*I "petscts.h" I*/ 4d763cef2SBarry Smith 5d5ba7fb7SMatthew Knepley /* Logging support */ 6166c7f25SBarry Smith PetscCookie PETSCTS_DLLEXPORT TS_COOKIE; 7166c7f25SBarry Smith PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 8d405a339SMatthew Knepley 94a2ae208SSatish Balay #undef __FUNCT__ 10bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions" 11bdad233fSMatthew Knepley /* 12bdad233fSMatthew Knepley TSSetTypeFromOptions - Sets the type of ts from user options. 13bdad233fSMatthew Knepley 14bdad233fSMatthew Knepley Collective on TS 15bdad233fSMatthew Knepley 16bdad233fSMatthew Knepley Input Parameter: 17bdad233fSMatthew Knepley . ts - The ts 18bdad233fSMatthew Knepley 19bdad233fSMatthew Knepley Level: intermediate 20bdad233fSMatthew Knepley 21bdad233fSMatthew Knepley .keywords: TS, set, options, database, type 22bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType() 23bdad233fSMatthew Knepley */ 246849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts) 25bdad233fSMatthew Knepley { 26bdad233fSMatthew Knepley PetscTruth opt; 272fc52814SBarry Smith const char *defaultType; 28bdad233fSMatthew Knepley char typeName[256]; 29dfbe8321SBarry Smith PetscErrorCode ierr; 30bdad233fSMatthew Knepley 31bdad233fSMatthew Knepley PetscFunctionBegin; 327adad957SLisandro Dalcin if (((PetscObject)ts)->type_name) { 337adad957SLisandro Dalcin defaultType = ((PetscObject)ts)->type_name; 34bdad233fSMatthew Knepley } else { 35bdad233fSMatthew Knepley defaultType = TS_EULER; 36bdad233fSMatthew Knepley } 37bdad233fSMatthew Knepley 38cce0b1b2SLisandro Dalcin if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 39bdad233fSMatthew Knepley ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 40a7cc72afSBarry Smith if (opt) { 41bdad233fSMatthew Knepley ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 42bdad233fSMatthew Knepley } else { 43bdad233fSMatthew Knepley ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 44bdad233fSMatthew Knepley } 45bdad233fSMatthew Knepley PetscFunctionReturn(0); 46bdad233fSMatthew Knepley } 47bdad233fSMatthew Knepley 48bdad233fSMatthew Knepley #undef __FUNCT__ 49bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions" 50bdad233fSMatthew Knepley /*@ 51bdad233fSMatthew Knepley TSSetFromOptions - Sets various TS parameters from user options. 52bdad233fSMatthew Knepley 53bdad233fSMatthew Knepley Collective on TS 54bdad233fSMatthew Knepley 55bdad233fSMatthew Knepley Input Parameter: 56bdad233fSMatthew Knepley . ts - the TS context obtained from TSCreate() 57bdad233fSMatthew Knepley 58bdad233fSMatthew Knepley Options Database Keys: 590f3b3ca1SHong Zhang + -ts_type <type> - TS_EULER, TS_BEULER, TS_SUNDIALS, TS_PSEUDO, TS_CRANK_NICHOLSON 60bdad233fSMatthew Knepley . -ts_max_steps maxsteps - maximum number of time-steps to take 61bdad233fSMatthew Knepley . -ts_max_time time - maximum time to compute to 62bdad233fSMatthew Knepley . -ts_dt dt - initial time step 63bdad233fSMatthew Knepley . -ts_monitor - print information at each timestep 64a6570f20SBarry Smith - -ts_monitor_draw - plot information at each timestep 65bdad233fSMatthew Knepley 66bdad233fSMatthew Knepley Level: beginner 67bdad233fSMatthew Knepley 68bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database 69bdad233fSMatthew Knepley 70a313700dSBarry Smith .seealso: TSGetType() 71bdad233fSMatthew Knepley @*/ 7263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS ts) 73bdad233fSMatthew Knepley { 74bdad233fSMatthew Knepley PetscReal dt; 75eabae89aSBarry Smith PetscTruth opt,flg; 76dfbe8321SBarry Smith PetscErrorCode ierr; 77a34d58ebSBarry Smith PetscViewerASCIIMonitor monviewer; 78eabae89aSBarry Smith char monfilename[PETSC_MAX_PATH_LEN]; 79bdad233fSMatthew Knepley 80bdad233fSMatthew Knepley PetscFunctionBegin; 814482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 827adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr); 83bdad233fSMatthew Knepley 84bdad233fSMatthew Knepley /* Handle generic TS options */ 85bdad233fSMatthew Knepley ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 86bdad233fSMatthew Knepley ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 87bdad233fSMatthew Knepley ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr); 88bdad233fSMatthew Knepley ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr); 89a7cc72afSBarry Smith if (opt) { 90bdad233fSMatthew Knepley ts->initial_time_step = ts->time_step = dt; 91bdad233fSMatthew Knepley } 92bdad233fSMatthew Knepley 93bdad233fSMatthew Knepley /* Monitor options */ 94a6570f20SBarry Smith ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 95eabae89aSBarry Smith if (flg) { 96050a712dSBarry Smith ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,monfilename,((PetscObject)ts)->tablevel,&monviewer);CHKERRQ(ierr); 97a34d58ebSBarry Smith ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 98bdad233fSMatthew Knepley } 9990d69ab7SBarry Smith opt = PETSC_FALSE; 10090d69ab7SBarry Smith ierr = PetscOptionsTruth("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 101a7cc72afSBarry Smith if (opt) { 102a6570f20SBarry Smith ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 103bdad233fSMatthew Knepley } 10490d69ab7SBarry Smith opt = PETSC_FALSE; 10590d69ab7SBarry Smith ierr = PetscOptionsTruth("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 106a7cc72afSBarry Smith if (opt) { 107a6570f20SBarry Smith ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 108bdad233fSMatthew Knepley } 109bdad233fSMatthew Knepley 110bdad233fSMatthew Knepley /* Handle TS type options */ 111bdad233fSMatthew Knepley ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 112bdad233fSMatthew Knepley 113bdad233fSMatthew Knepley /* Handle specific TS options */ 114abc0a331SBarry Smith if (ts->ops->setfromoptions) { 115bdad233fSMatthew Knepley ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 116bdad233fSMatthew Knepley } 117bdad233fSMatthew Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 118bdad233fSMatthew Knepley 119bdad233fSMatthew Knepley /* Handle subobject options */ 120bdad233fSMatthew Knepley switch(ts->problem_type) { 121156fc9a6SMatthew Knepley /* Should check for implicit/explicit */ 122bdad233fSMatthew Knepley case TS_LINEAR: 123abc0a331SBarry Smith if (ts->ksp) { 1248beb423aSHong Zhang ierr = KSPSetOperators(ts->ksp,ts->Arhs,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 12594b7f48cSBarry Smith ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr); 126156fc9a6SMatthew Knepley } 127bdad233fSMatthew Knepley break; 128bdad233fSMatthew Knepley case TS_NONLINEAR: 129abc0a331SBarry Smith if (ts->snes) { 1307c236d22SBarry Smith /* this is a bit of a hack, but it gets the matrix information into SNES earlier 1317c236d22SBarry Smith so that SNES and KSP have more information to pick reasonable defaults 1327c0b301bSJed Brown before they allow users to set options 1337c0b301bSJed Brown * If ts->A has been set at this point, we are probably using the implicit form 1347c0b301bSJed Brown and Arhs will never be used. */ 1357c0b301bSJed Brown ierr = SNESSetJacobian(ts->snes,ts->A?ts->A:ts->Arhs,ts->B,0,ts);CHKERRQ(ierr); 136bdad233fSMatthew Knepley ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 137156fc9a6SMatthew Knepley } 138bdad233fSMatthew Knepley break; 139bdad233fSMatthew Knepley default: 14077431f27SBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type); 141bdad233fSMatthew Knepley } 142bdad233fSMatthew Knepley 143bdad233fSMatthew Knepley PetscFunctionReturn(0); 144bdad233fSMatthew Knepley } 145bdad233fSMatthew Knepley 146bdad233fSMatthew Knepley #undef __FUNCT__ 147bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions" 148bdad233fSMatthew Knepley /*@ 149bdad233fSMatthew Knepley TSViewFromOptions - This function visualizes the ts based upon user options. 150bdad233fSMatthew Knepley 151bdad233fSMatthew Knepley Collective on TS 152bdad233fSMatthew Knepley 153bdad233fSMatthew Knepley Input Parameter: 154bdad233fSMatthew Knepley . ts - The ts 155bdad233fSMatthew Knepley 156bdad233fSMatthew Knepley Level: intermediate 157bdad233fSMatthew Knepley 158bdad233fSMatthew Knepley .keywords: TS, view, options, database 159bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView() 160bdad233fSMatthew Knepley @*/ 16163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[]) 162bdad233fSMatthew Knepley { 163bdad233fSMatthew Knepley PetscViewer viewer; 164bdad233fSMatthew Knepley PetscDraw draw; 16590d69ab7SBarry Smith PetscTruth opt = PETSC_FALSE; 166e10c49a3SBarry Smith char fileName[PETSC_MAX_PATH_LEN]; 167dfbe8321SBarry Smith PetscErrorCode ierr; 168bdad233fSMatthew Knepley 169bdad233fSMatthew Knepley PetscFunctionBegin; 1707adad957SLisandro Dalcin ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 171eabae89aSBarry Smith if (opt && !PetscPreLoadingOn) { 1727adad957SLisandro Dalcin ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 173bdad233fSMatthew Knepley ierr = TSView(ts, viewer);CHKERRQ(ierr); 174bdad233fSMatthew Knepley ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 175bdad233fSMatthew Knepley } 1768e83347fSKai Germaschewski opt = PETSC_FALSE; 17790d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 178a7cc72afSBarry Smith if (opt) { 1797adad957SLisandro Dalcin ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 180bdad233fSMatthew Knepley ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 181a7cc72afSBarry Smith if (title) { 1821836bdbcSSatish Balay ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 183bdad233fSMatthew Knepley } else { 184bdad233fSMatthew Knepley ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 1857adad957SLisandro Dalcin ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 186bdad233fSMatthew Knepley } 187bdad233fSMatthew Knepley ierr = TSView(ts, viewer);CHKERRQ(ierr); 188bdad233fSMatthew Knepley ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 189bdad233fSMatthew Knepley ierr = PetscDrawPause(draw);CHKERRQ(ierr); 190bdad233fSMatthew Knepley ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 191bdad233fSMatthew Knepley } 192bdad233fSMatthew Knepley PetscFunctionReturn(0); 193bdad233fSMatthew Knepley } 194bdad233fSMatthew Knepley 195bdad233fSMatthew Knepley #undef __FUNCT__ 1964a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian" 197a7a1495cSBarry Smith /*@ 1988c385f81SBarry Smith TSComputeRHSJacobian - Computes the Jacobian matrix that has been 199a7a1495cSBarry Smith set with TSSetRHSJacobian(). 200a7a1495cSBarry Smith 201a7a1495cSBarry Smith Collective on TS and Vec 202a7a1495cSBarry Smith 203a7a1495cSBarry Smith Input Parameters: 204316643e7SJed Brown + ts - the TS context 205a7a1495cSBarry Smith . t - current timestep 206a7a1495cSBarry Smith - x - input vector 207a7a1495cSBarry Smith 208a7a1495cSBarry Smith Output Parameters: 209a7a1495cSBarry Smith + A - Jacobian matrix 210a7a1495cSBarry Smith . B - optional preconditioning matrix 211a7a1495cSBarry Smith - flag - flag indicating matrix structure 212a7a1495cSBarry Smith 213a7a1495cSBarry Smith Notes: 214a7a1495cSBarry Smith Most users should not need to explicitly call this routine, as it 215a7a1495cSBarry Smith is used internally within the nonlinear solvers. 216a7a1495cSBarry Smith 21794b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 218a7a1495cSBarry Smith flag parameter. 219a7a1495cSBarry Smith 220a7a1495cSBarry Smith TSComputeJacobian() is valid only for TS_NONLINEAR 221a7a1495cSBarry Smith 222a7a1495cSBarry Smith Level: developer 223a7a1495cSBarry Smith 224a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix 225a7a1495cSBarry Smith 22694b7f48cSBarry Smith .seealso: TSSetRHSJacobian(), KSPSetOperators() 227a7a1495cSBarry Smith @*/ 22863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 229a7a1495cSBarry Smith { 230dfbe8321SBarry Smith PetscErrorCode ierr; 231a7a1495cSBarry Smith 232a7a1495cSBarry Smith PetscFunctionBegin; 2334482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 2344482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,3); 235c9780b6fSBarry Smith PetscCheckSameComm(ts,1,X,3); 236a7a1495cSBarry Smith if (ts->problem_type != TS_NONLINEAR) { 23729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only"); 238a7a1495cSBarry Smith } 239000e7ae3SMatthew Knepley if (ts->ops->rhsjacobian) { 240d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 241a7a1495cSBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 242a7a1495cSBarry Smith PetscStackPush("TS user Jacobian function"); 243000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 244a7a1495cSBarry Smith PetscStackPop; 245d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 246a7a1495cSBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 2474482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,4); 2484482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,5); 249ef66eb69SBarry Smith } else { 250ef66eb69SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 251ef66eb69SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252ef66eb69SBarry Smith if (*A != *B) { 253ef66eb69SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 254ef66eb69SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255ef66eb69SBarry Smith } 256ef66eb69SBarry Smith } 257a7a1495cSBarry Smith PetscFunctionReturn(0); 258a7a1495cSBarry Smith } 259a7a1495cSBarry Smith 2604a2ae208SSatish Balay #undef __FUNCT__ 2614a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction" 262316643e7SJed Brown /*@ 263d763cef2SBarry Smith TSComputeRHSFunction - Evaluates the right-hand-side function. 264d763cef2SBarry Smith 265316643e7SJed Brown Collective on TS and Vec 266316643e7SJed Brown 267316643e7SJed Brown Input Parameters: 268316643e7SJed Brown + ts - the TS context 269316643e7SJed Brown . t - current time 270316643e7SJed Brown - x - state vector 271316643e7SJed Brown 272316643e7SJed Brown Output Parameter: 273316643e7SJed Brown . y - right hand side 274316643e7SJed Brown 275316643e7SJed Brown Note: 276316643e7SJed Brown Most users should not need to explicitly call this routine, as it 277316643e7SJed Brown is used internally within the nonlinear solvers. 278316643e7SJed Brown 279316643e7SJed Brown If the user did not provide a function but merely a matrix, 280d763cef2SBarry Smith this routine applies the matrix. 281316643e7SJed Brown 282316643e7SJed Brown Level: developer 283316643e7SJed Brown 284316643e7SJed Brown .keywords: TS, compute 285316643e7SJed Brown 286316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction() 287316643e7SJed Brown @*/ 288dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 289d763cef2SBarry Smith { 290dfbe8321SBarry Smith PetscErrorCode ierr; 291d763cef2SBarry Smith 292d763cef2SBarry Smith PetscFunctionBegin; 2934482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 294316643e7SJed Brown PetscValidHeaderSpecific(x,VEC_COOKIE,3); 295316643e7SJed Brown PetscValidHeaderSpecific(y,VEC_COOKIE,4); 296d763cef2SBarry Smith 297d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 298000e7ae3SMatthew Knepley if (ts->ops->rhsfunction) { 299d763cef2SBarry Smith PetscStackPush("TS user right-hand-side function"); 300000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 301d763cef2SBarry Smith PetscStackPop; 3027c922b88SBarry Smith } else { 303000e7ae3SMatthew Knepley if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 304d763cef2SBarry Smith MatStructure flg; 305d763cef2SBarry Smith PetscStackPush("TS user right-hand-side matrix function"); 3068beb423aSHong Zhang ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 307d763cef2SBarry Smith PetscStackPop; 308d763cef2SBarry Smith } 3098beb423aSHong Zhang ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr); 3107c922b88SBarry Smith } 311d763cef2SBarry Smith 312d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 313d763cef2SBarry Smith 314d763cef2SBarry Smith PetscFunctionReturn(0); 315d763cef2SBarry Smith } 316d763cef2SBarry Smith 3174a2ae208SSatish Balay #undef __FUNCT__ 318316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction" 319316643e7SJed Brown /*@ 320316643e7SJed Brown TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 321316643e7SJed Brown 322316643e7SJed Brown Collective on TS and Vec 323316643e7SJed Brown 324316643e7SJed Brown Input Parameters: 325316643e7SJed Brown + ts - the TS context 326316643e7SJed Brown . t - current time 327316643e7SJed Brown . X - state vector 328316643e7SJed Brown - Xdot - time derivative of state vector 329316643e7SJed Brown 330316643e7SJed Brown Output Parameter: 331316643e7SJed Brown . Y - right hand side 332316643e7SJed Brown 333316643e7SJed Brown Note: 334316643e7SJed Brown Most users should not need to explicitly call this routine, as it 335316643e7SJed Brown is used internally within the nonlinear solvers. 336316643e7SJed Brown 337316643e7SJed Brown If the user did did not write their equations in implicit form, this 338316643e7SJed Brown function recasts them in implicit form. 339316643e7SJed Brown 340316643e7SJed Brown Level: developer 341316643e7SJed Brown 342316643e7SJed Brown .keywords: TS, compute 343316643e7SJed Brown 344316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction() 345316643e7SJed Brown @*/ 346316643e7SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y) 347316643e7SJed Brown { 348316643e7SJed Brown PetscErrorCode ierr; 349316643e7SJed Brown 350316643e7SJed Brown PetscFunctionBegin; 351316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 352316643e7SJed Brown PetscValidHeaderSpecific(X,VEC_COOKIE,3); 353316643e7SJed Brown PetscValidHeaderSpecific(Xdot,VEC_COOKIE,4); 354316643e7SJed Brown PetscValidHeaderSpecific(Y,VEC_COOKIE,5); 355316643e7SJed Brown 356316643e7SJed Brown ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 357316643e7SJed Brown if (ts->ops->ifunction) { 358316643e7SJed Brown PetscStackPush("TS user implicit function"); 359316643e7SJed Brown ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 360316643e7SJed Brown PetscStackPop; 361316643e7SJed Brown } else { 362316643e7SJed Brown PetscScalar *y; 363316643e7SJed Brown PetscInt i,n; 364316643e7SJed Brown 365316643e7SJed Brown if (ts->ops->rhsfunction) { 366316643e7SJed Brown PetscStackPush("TS user right-hand-side function"); 367316643e7SJed Brown ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr); 368316643e7SJed Brown PetscStackPop; 369316643e7SJed Brown } else { 370316643e7SJed Brown if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 371316643e7SJed Brown MatStructure flg; 372316643e7SJed Brown PetscStackPush("TS user right-hand-side matrix function"); 373316643e7SJed Brown ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 374316643e7SJed Brown PetscStackPop; 375316643e7SJed Brown } 376316643e7SJed Brown ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr); 377316643e7SJed Brown } 378316643e7SJed Brown 379316643e7SJed Brown /* Convert to implicit form */ 380*ace68cafSJed Brown ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 381316643e7SJed Brown } 382316643e7SJed Brown ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 383316643e7SJed Brown PetscFunctionReturn(0); 384316643e7SJed Brown } 385316643e7SJed Brown 386316643e7SJed Brown #undef __FUNCT__ 387316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian" 388316643e7SJed Brown /*@ 389316643e7SJed Brown TSComputeIJacobian - Evaluates the Jacobian of the DAE 390316643e7SJed Brown 391316643e7SJed Brown Collective on TS and Vec 392316643e7SJed Brown 393316643e7SJed Brown Input 394316643e7SJed Brown Input Parameters: 395316643e7SJed Brown + ts - the TS context 396316643e7SJed Brown . t - current timestep 397316643e7SJed Brown . X - state vector 398316643e7SJed Brown . Xdot - time derivative of state vector 399316643e7SJed Brown - shift - shift to apply, see note below 400316643e7SJed Brown 401316643e7SJed Brown Output Parameters: 402316643e7SJed Brown + A - Jacobian matrix 403316643e7SJed Brown . B - optional preconditioning matrix 404316643e7SJed Brown - flag - flag indicating matrix structure 405316643e7SJed Brown 406316643e7SJed Brown Notes: 407316643e7SJed Brown If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 408316643e7SJed Brown 409316643e7SJed Brown dF/dX + shift*dF/dXdot 410316643e7SJed Brown 411316643e7SJed Brown Most users should not need to explicitly call this routine, as it 412316643e7SJed Brown is used internally within the nonlinear solvers. 413316643e7SJed Brown 414316643e7SJed Brown TSComputeIJacobian() is valid only for TS_NONLINEAR 415316643e7SJed Brown 416316643e7SJed Brown Level: developer 417316643e7SJed Brown 418316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix 419316643e7SJed Brown 420316643e7SJed Brown .seealso: TSSetIJacobian() 42163495f91SJed Brown @*/ 422316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg) 423316643e7SJed Brown { 424316643e7SJed Brown PetscErrorCode ierr; 425316643e7SJed Brown 426316643e7SJed Brown PetscFunctionBegin; 427316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 428316643e7SJed Brown PetscValidHeaderSpecific(X,VEC_COOKIE,3); 429316643e7SJed Brown PetscValidHeaderSpecific(Xdot,VEC_COOKIE,4); 430316643e7SJed Brown PetscValidPointer(A,6); 431316643e7SJed Brown PetscValidHeaderSpecific(*A,MAT_COOKIE,6); 432316643e7SJed Brown PetscValidPointer(B,7); 433316643e7SJed Brown PetscValidHeaderSpecific(*B,MAT_COOKIE,7); 434316643e7SJed Brown PetscValidPointer(flg,8); 435316643e7SJed Brown 436316643e7SJed Brown ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 437316643e7SJed Brown if (ts->ops->ijacobian) { 438316643e7SJed Brown PetscStackPush("TS user implicit Jacobian"); 439316643e7SJed Brown ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 440316643e7SJed Brown PetscStackPop; 441316643e7SJed Brown } else { 442316643e7SJed Brown if (ts->ops->rhsjacobian) { 443316643e7SJed Brown PetscStackPush("TS user right-hand-side Jacobian"); 444316643e7SJed Brown ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 445316643e7SJed Brown PetscStackPop; 446316643e7SJed Brown } else { 447316643e7SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 448316643e7SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 449316643e7SJed Brown if (*A != *B) { 450316643e7SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 451316643e7SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 452316643e7SJed Brown } 453316643e7SJed Brown } 454316643e7SJed Brown 455316643e7SJed Brown /* Convert to implicit form */ 456316643e7SJed Brown /* inefficient because these operations will normally traverse all matrix elements twice */ 457316643e7SJed Brown ierr = MatScale(*A,-1);CHKERRQ(ierr); 458316643e7SJed Brown ierr = MatShift(*A,shift);CHKERRQ(ierr); 459316643e7SJed Brown if (*A != *B) { 460316643e7SJed Brown ierr = MatScale(*B,-1);CHKERRQ(ierr); 461316643e7SJed Brown ierr = MatShift(*B,shift);CHKERRQ(ierr); 462316643e7SJed Brown } 463316643e7SJed Brown } 464316643e7SJed Brown ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 465316643e7SJed Brown PetscFunctionReturn(0); 466316643e7SJed Brown } 467316643e7SJed Brown 468316643e7SJed Brown #undef __FUNCT__ 4694a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction" 470d763cef2SBarry Smith /*@C 471d763cef2SBarry Smith TSSetRHSFunction - Sets the routine for evaluating the function, 472d763cef2SBarry Smith F(t,u), where U_t = F(t,u). 473d763cef2SBarry Smith 474d763cef2SBarry Smith Collective on TS 475d763cef2SBarry Smith 476d763cef2SBarry Smith Input Parameters: 477d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 478d763cef2SBarry Smith . f - routine for evaluating the right-hand-side function 479d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 480d763cef2SBarry Smith function evaluation routine (may be PETSC_NULL) 481d763cef2SBarry Smith 482d763cef2SBarry Smith Calling sequence of func: 48387828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 484d763cef2SBarry Smith 485d763cef2SBarry Smith + t - current timestep 486d763cef2SBarry Smith . u - input vector 487d763cef2SBarry Smith . F - function vector 488d763cef2SBarry Smith - ctx - [optional] user-defined function context 489d763cef2SBarry Smith 490d763cef2SBarry Smith Important: 49176f2fa84SHong Zhang The user MUST call either this routine or TSSetMatrices(). 492d763cef2SBarry Smith 493d763cef2SBarry Smith Level: beginner 494d763cef2SBarry Smith 495d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function 496d763cef2SBarry Smith 49776f2fa84SHong Zhang .seealso: TSSetMatrices() 498d763cef2SBarry Smith @*/ 49963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 500d763cef2SBarry Smith { 501d763cef2SBarry Smith PetscFunctionBegin; 502d763cef2SBarry Smith 5034482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 504d763cef2SBarry Smith if (ts->problem_type == TS_LINEAR) { 50529bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem"); 506d763cef2SBarry Smith } 507000e7ae3SMatthew Knepley ts->ops->rhsfunction = f; 508d763cef2SBarry Smith ts->funP = ctx; 509d763cef2SBarry Smith PetscFunctionReturn(0); 510d763cef2SBarry Smith } 511d763cef2SBarry Smith 5124a2ae208SSatish Balay #undef __FUNCT__ 51395f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices" 51495f0b562SHong Zhang /*@C 51595f0b562SHong Zhang TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs, 51695f0b562SHong Zhang where Alhs(t) U_t = Arhs(t) U. 51795f0b562SHong Zhang 51895f0b562SHong Zhang Collective on TS 51995f0b562SHong Zhang 52095f0b562SHong Zhang Input Parameters: 52195f0b562SHong Zhang + ts - the TS context obtained from TSCreate() 52295f0b562SHong Zhang . Arhs - matrix 52395f0b562SHong Zhang . frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 52495f0b562SHong Zhang if Arhs is not a function of t. 52595f0b562SHong Zhang . Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix. 52695f0b562SHong Zhang . flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 52795f0b562SHong Zhang if Alhs is not a function of t. 52895f0b562SHong Zhang . flag - flag indicating information about the matrix structure of Arhs and Alhs. 52995f0b562SHong Zhang The available options are 53095f0b562SHong Zhang SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs 53195f0b562SHong Zhang DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs 53295f0b562SHong Zhang - ctx - [optional] user-defined context for private data for the 53395f0b562SHong Zhang matrix evaluation routine (may be PETSC_NULL) 53495f0b562SHong Zhang 53595f0b562SHong Zhang Calling sequence of func: 53695f0b562SHong Zhang $ func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx); 53795f0b562SHong Zhang 53895f0b562SHong Zhang + t - current timestep 53995f0b562SHong Zhang . A - matrix A, where U_t = A(t) U 54095f0b562SHong Zhang . B - preconditioner matrix, usually the same as A 54195f0b562SHong Zhang . flag - flag indicating information about the preconditioner matrix 54295f0b562SHong Zhang structure (same as flag in KSPSetOperators()) 54395f0b562SHong Zhang - ctx - [optional] user-defined context for matrix evaluation routine 54495f0b562SHong Zhang 54595f0b562SHong Zhang Notes: 54695f0b562SHong Zhang The routine func() takes Mat* as the matrix arguments rather than Mat. 54795f0b562SHong Zhang This allows the matrix evaluation routine to replace Arhs or Alhs with a 54895f0b562SHong Zhang completely new new matrix structure (not just different matrix elements) 54995f0b562SHong Zhang when appropriate, for instance, if the nonzero structure is changing 55095f0b562SHong Zhang throughout the global iterations. 55195f0b562SHong Zhang 55295f0b562SHong Zhang Important: 55395f0b562SHong Zhang The user MUST call either this routine or TSSetRHSFunction(). 55495f0b562SHong Zhang 55595f0b562SHong Zhang Level: beginner 55695f0b562SHong Zhang 55795f0b562SHong Zhang .keywords: TS, timestep, set, matrix 55895f0b562SHong Zhang 55995f0b562SHong Zhang .seealso: TSSetRHSFunction() 56095f0b562SHong Zhang @*/ 56195f0b562SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSetMatrices(TS ts,Mat Arhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx) 56295f0b562SHong Zhang { 56395f0b562SHong Zhang PetscFunctionBegin; 56495f0b562SHong Zhang PetscValidHeaderSpecific(ts,TS_COOKIE,1); 56592af4f6aSHong Zhang if (Arhs){ 56695f0b562SHong Zhang PetscValidHeaderSpecific(Arhs,MAT_COOKIE,2); 56795f0b562SHong Zhang PetscCheckSameComm(ts,1,Arhs,2); 56895f0b562SHong Zhang ts->Arhs = Arhs; 56992af4f6aSHong Zhang ts->ops->rhsmatrix = frhs; 57092af4f6aSHong Zhang } 57192af4f6aSHong Zhang if (Alhs){ 57292af4f6aSHong Zhang PetscValidHeaderSpecific(Alhs,MAT_COOKIE,4); 57392af4f6aSHong Zhang PetscCheckSameComm(ts,1,Alhs,4); 57495f0b562SHong Zhang ts->Alhs = Alhs; 57592af4f6aSHong Zhang ts->ops->lhsmatrix = flhs; 57692af4f6aSHong Zhang } 57792af4f6aSHong Zhang 57892af4f6aSHong Zhang ts->jacP = ctx; 57995f0b562SHong Zhang ts->matflg = flag; 58095f0b562SHong Zhang PetscFunctionReturn(0); 58195f0b562SHong Zhang } 582d763cef2SBarry Smith 583aa644b49SHong Zhang #undef __FUNCT__ 584cda39b92SHong Zhang #define __FUNCT__ "TSGetMatrices" 585cda39b92SHong Zhang /*@C 586cda39b92SHong Zhang TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep, 587cda39b92SHong Zhang where Alhs(t) U_t = Arhs(t) U. 588cda39b92SHong Zhang 589cda39b92SHong Zhang Not Collective, but parallel objects are returned if TS is parallel 590cda39b92SHong Zhang 591cda39b92SHong Zhang Input Parameter: 592cda39b92SHong Zhang . ts - The TS context obtained from TSCreate() 593cda39b92SHong Zhang 594cda39b92SHong Zhang Output Parameters: 595cda39b92SHong Zhang + Arhs - The right-hand side matrix 596cda39b92SHong Zhang . Alhs - The left-hand side matrix 597cda39b92SHong Zhang - ctx - User-defined context for matrix evaluation routine 598cda39b92SHong Zhang 599cda39b92SHong Zhang Notes: You can pass in PETSC_NULL for any return argument you do not need. 600cda39b92SHong Zhang 601cda39b92SHong Zhang Level: intermediate 602cda39b92SHong Zhang 603cda39b92SHong Zhang .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 604cda39b92SHong Zhang 605cda39b92SHong Zhang .keywords: TS, timestep, get, matrix 606cda39b92SHong Zhang 607cda39b92SHong Zhang @*/ 608cda39b92SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx) 609cda39b92SHong Zhang { 610cda39b92SHong Zhang PetscFunctionBegin; 611cda39b92SHong Zhang PetscValidHeaderSpecific(ts,TS_COOKIE,1); 612cda39b92SHong Zhang if (Arhs) *Arhs = ts->Arhs; 613cda39b92SHong Zhang if (Alhs) *Alhs = ts->Alhs; 614cda39b92SHong Zhang if (ctx) *ctx = ts->jacP; 615cda39b92SHong Zhang PetscFunctionReturn(0); 616cda39b92SHong Zhang } 617cda39b92SHong Zhang 618cda39b92SHong Zhang #undef __FUNCT__ 6194a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian" 620d763cef2SBarry Smith /*@C 621d763cef2SBarry Smith TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 622d763cef2SBarry Smith where U_t = F(U,t), as well as the location to store the matrix. 62376f2fa84SHong Zhang Use TSSetMatrices() for linear problems. 624d763cef2SBarry Smith 625d763cef2SBarry Smith Collective on TS 626d763cef2SBarry Smith 627d763cef2SBarry Smith Input Parameters: 628d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 629d763cef2SBarry Smith . A - Jacobian matrix 630d763cef2SBarry Smith . B - preconditioner matrix (usually same as A) 631d763cef2SBarry Smith . f - the Jacobian evaluation routine 632d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 633d763cef2SBarry Smith Jacobian evaluation routine (may be PETSC_NULL) 634d763cef2SBarry Smith 635d763cef2SBarry Smith Calling sequence of func: 63687828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 637d763cef2SBarry Smith 638d763cef2SBarry Smith + t - current timestep 639d763cef2SBarry Smith . u - input vector 640d763cef2SBarry Smith . A - matrix A, where U_t = A(t)u 641d763cef2SBarry Smith . B - preconditioner matrix, usually the same as A 642d763cef2SBarry Smith . flag - flag indicating information about the preconditioner matrix 64394b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 644d763cef2SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 645d763cef2SBarry Smith 646d763cef2SBarry Smith Notes: 64794b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 648d763cef2SBarry Smith output parameter in the routine func(). Be sure to read this information! 649d763cef2SBarry Smith 650d763cef2SBarry Smith The routine func() takes Mat * as the matrix arguments rather than Mat. 651d763cef2SBarry Smith This allows the matrix evaluation routine to replace A and/or B with a 65256335db2SHong Zhang completely new matrix structure (not just different matrix elements) 653d763cef2SBarry Smith when appropriate, for instance, if the nonzero structure is changing 654d763cef2SBarry Smith throughout the global iterations. 655d763cef2SBarry Smith 656d763cef2SBarry Smith Level: beginner 657d763cef2SBarry Smith 658d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian 659d763cef2SBarry Smith 660d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(), 66176f2fa84SHong Zhang SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 662d763cef2SBarry Smith 663d763cef2SBarry Smith @*/ 66463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 665d763cef2SBarry Smith { 666d763cef2SBarry Smith PetscFunctionBegin; 6674482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 6684482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,2); 6694482741eSBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE,3); 670c9780b6fSBarry Smith PetscCheckSameComm(ts,1,A,2); 671c9780b6fSBarry Smith PetscCheckSameComm(ts,1,B,3); 672d763cef2SBarry Smith if (ts->problem_type != TS_NONLINEAR) { 67376f2fa84SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()"); 674d763cef2SBarry Smith } 675d763cef2SBarry Smith 676000e7ae3SMatthew Knepley ts->ops->rhsjacobian = f; 677d763cef2SBarry Smith ts->jacP = ctx; 6788beb423aSHong Zhang ts->Arhs = A; 679d763cef2SBarry Smith ts->B = B; 680d763cef2SBarry Smith PetscFunctionReturn(0); 681d763cef2SBarry Smith } 682d763cef2SBarry Smith 683316643e7SJed Brown 684316643e7SJed Brown #undef __FUNCT__ 685316643e7SJed Brown #define __FUNCT__ "TSSetIFunction" 686316643e7SJed Brown /*@C 687316643e7SJed Brown TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 688316643e7SJed Brown 689316643e7SJed Brown Collective on TS 690316643e7SJed Brown 691316643e7SJed Brown Input Parameters: 692316643e7SJed Brown + ts - the TS context obtained from TSCreate() 693316643e7SJed Brown . f - the function evaluation routine 694316643e7SJed Brown - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 695316643e7SJed Brown 696316643e7SJed Brown Calling sequence of f: 697316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 698316643e7SJed Brown 699316643e7SJed Brown + t - time at step/stage being solved 700316643e7SJed Brown . u - state vector 701316643e7SJed Brown . u_t - time derivative of state vector 702316643e7SJed Brown . F - function vector 703316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 704316643e7SJed Brown 705316643e7SJed Brown Important: 706316643e7SJed Brown The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 707316643e7SJed Brown 708316643e7SJed Brown Level: beginner 709316643e7SJed Brown 710316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian 711316643e7SJed Brown 712316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 713316643e7SJed Brown @*/ 714316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS ts,TSIFunction f,void *ctx) 715316643e7SJed Brown { 716316643e7SJed Brown 717316643e7SJed Brown PetscFunctionBegin; 718316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 719316643e7SJed Brown ts->ops->ifunction = f; 720316643e7SJed Brown ts->funP = ctx; 721316643e7SJed Brown PetscFunctionReturn(0); 722316643e7SJed Brown } 723316643e7SJed Brown 724316643e7SJed Brown #undef __FUNCT__ 725316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian" 726316643e7SJed Brown /*@C 727316643e7SJed Brown TSSetIJacobian - Set the function to compute the Jacobian of 728316643e7SJed Brown G(U) = F(t,U,U0+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 729316643e7SJed Brown 730316643e7SJed Brown Collective on TS 731316643e7SJed Brown 732316643e7SJed Brown Input Parameters: 733316643e7SJed Brown + ts - the TS context obtained from TSCreate() 734316643e7SJed Brown . A - Jacobian matrix 735316643e7SJed Brown . B - preconditioning matrix for A (may be same as A) 736316643e7SJed Brown . f - the Jacobian evaluation routine 737316643e7SJed Brown - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 738316643e7SJed Brown 739316643e7SJed Brown Calling sequence of f: 740316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 741316643e7SJed Brown 742316643e7SJed Brown + t - time at step/stage being solved 743316643e7SJed Brown . u - state vector 744316643e7SJed Brown . u_t - time derivative of state vector 745316643e7SJed Brown . a - shift 746316643e7SJed Brown . A - Jacobian of G(U) = F(t,U,U0+a*U), equivalent to dF/dU + a*dF/dU_t 747316643e7SJed Brown . B - preconditioning matrix for A, may be same as A 748316643e7SJed Brown . flag - flag indicating information about the preconditioner matrix 749316643e7SJed Brown structure (same as flag in KSPSetOperators()) 750316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 751316643e7SJed Brown 752316643e7SJed Brown Notes: 753316643e7SJed Brown The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 754316643e7SJed Brown 755316643e7SJed Brown Level: beginner 756316643e7SJed Brown 757316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian 758316643e7SJed Brown 759316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian() 760316643e7SJed Brown 761316643e7SJed Brown @*/ 762316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 763316643e7SJed Brown { 764316643e7SJed Brown PetscErrorCode ierr; 765316643e7SJed Brown 766316643e7SJed Brown PetscFunctionBegin; 767316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 768316643e7SJed Brown if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 769316643e7SJed Brown if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 770316643e7SJed Brown if (A) PetscCheckSameComm(ts,1,A,2); 771316643e7SJed Brown if (B) PetscCheckSameComm(ts,1,B,3); 772316643e7SJed Brown if (f) ts->ops->ijacobian = f; 773316643e7SJed Brown if (ctx) ts->jacP = ctx; 774316643e7SJed Brown if (A) { 775316643e7SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 776316643e7SJed Brown if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 777316643e7SJed Brown ts->A = A; 778316643e7SJed Brown } 779316643e7SJed Brown if (B) { 780316643e7SJed Brown ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 781316643e7SJed Brown if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);} 782316643e7SJed Brown ts->B = B; 783316643e7SJed Brown } 784316643e7SJed Brown PetscFunctionReturn(0); 785316643e7SJed Brown } 786316643e7SJed Brown 7874a2ae208SSatish Balay #undef __FUNCT__ 7884a2ae208SSatish Balay #define __FUNCT__ "TSView" 7897e2c5f70SBarry Smith /*@C 790d763cef2SBarry Smith TSView - Prints the TS data structure. 791d763cef2SBarry Smith 7924c49b128SBarry Smith Collective on TS 793d763cef2SBarry Smith 794d763cef2SBarry Smith Input Parameters: 795d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 796d763cef2SBarry Smith - viewer - visualization context 797d763cef2SBarry Smith 798d763cef2SBarry Smith Options Database Key: 799d763cef2SBarry Smith . -ts_view - calls TSView() at end of TSStep() 800d763cef2SBarry Smith 801d763cef2SBarry Smith Notes: 802d763cef2SBarry Smith The available visualization contexts include 803b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 804b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 805d763cef2SBarry Smith output where only the first processor opens 806d763cef2SBarry Smith the file. All other processors send their 807d763cef2SBarry Smith data to the first processor to print. 808d763cef2SBarry Smith 809d763cef2SBarry Smith The user can open an alternative visualization context with 810b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 811d763cef2SBarry Smith 812d763cef2SBarry Smith Level: beginner 813d763cef2SBarry Smith 814d763cef2SBarry Smith .keywords: TS, timestep, view 815d763cef2SBarry Smith 816b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 817d763cef2SBarry Smith @*/ 81863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer) 819d763cef2SBarry Smith { 820dfbe8321SBarry Smith PetscErrorCode ierr; 821a313700dSBarry Smith const TSType type; 82232077d6dSBarry Smith PetscTruth iascii,isstring; 823d763cef2SBarry Smith 824d763cef2SBarry Smith PetscFunctionBegin; 8254482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 8263050cee2SBarry Smith if (!viewer) { 8277adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 8283050cee2SBarry Smith } 8294482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 830c9780b6fSBarry Smith PetscCheckSameComm(ts,1,viewer,2); 831fd16b177SBarry Smith 83232077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 833b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 83432077d6dSBarry Smith if (iascii) { 835b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr); 836a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 837454a90a3SBarry Smith if (type) { 838b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 839184914b5SBarry Smith } else { 840b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 841184914b5SBarry Smith } 842000e7ae3SMatthew Knepley if (ts->ops->view) { 843b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 844000e7ae3SMatthew Knepley ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 845b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 846d763cef2SBarry Smith } 84777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 848a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 849d763cef2SBarry Smith if (ts->problem_type == TS_NONLINEAR) { 85077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 851d763cef2SBarry Smith } 85277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 8530f5bd95cSBarry Smith } else if (isstring) { 854a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 855b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 856d763cef2SBarry Smith } 857b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 85894b7f48cSBarry Smith if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);} 859d763cef2SBarry Smith if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 860b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 861d763cef2SBarry Smith PetscFunctionReturn(0); 862d763cef2SBarry Smith } 863d763cef2SBarry Smith 864d763cef2SBarry Smith 8654a2ae208SSatish Balay #undef __FUNCT__ 8664a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext" 867d763cef2SBarry Smith /*@C 868d763cef2SBarry Smith TSSetApplicationContext - Sets an optional user-defined context for 869d763cef2SBarry Smith the timesteppers. 870d763cef2SBarry Smith 871d763cef2SBarry Smith Collective on TS 872d763cef2SBarry Smith 873d763cef2SBarry Smith Input Parameters: 874d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 875d763cef2SBarry Smith - usrP - optional user context 876d763cef2SBarry Smith 877d763cef2SBarry Smith Level: intermediate 878d763cef2SBarry Smith 879d763cef2SBarry Smith .keywords: TS, timestep, set, application, context 880d763cef2SBarry Smith 881d763cef2SBarry Smith .seealso: TSGetApplicationContext() 882d763cef2SBarry Smith @*/ 88363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP) 884d763cef2SBarry Smith { 885d763cef2SBarry Smith PetscFunctionBegin; 8864482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 887d763cef2SBarry Smith ts->user = usrP; 888d763cef2SBarry Smith PetscFunctionReturn(0); 889d763cef2SBarry Smith } 890d763cef2SBarry Smith 8914a2ae208SSatish Balay #undef __FUNCT__ 8924a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext" 893d763cef2SBarry Smith /*@C 894d763cef2SBarry Smith TSGetApplicationContext - Gets the user-defined context for the 895d763cef2SBarry Smith timestepper. 896d763cef2SBarry Smith 897d763cef2SBarry Smith Not Collective 898d763cef2SBarry Smith 899d763cef2SBarry Smith Input Parameter: 900d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 901d763cef2SBarry Smith 902d763cef2SBarry Smith Output Parameter: 903d763cef2SBarry Smith . usrP - user context 904d763cef2SBarry Smith 905d763cef2SBarry Smith Level: intermediate 906d763cef2SBarry Smith 907d763cef2SBarry Smith .keywords: TS, timestep, get, application, context 908d763cef2SBarry Smith 909d763cef2SBarry Smith .seealso: TSSetApplicationContext() 910d763cef2SBarry Smith @*/ 91163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP) 912d763cef2SBarry Smith { 913d763cef2SBarry Smith PetscFunctionBegin; 9144482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 915d763cef2SBarry Smith *usrP = ts->user; 916d763cef2SBarry Smith PetscFunctionReturn(0); 917d763cef2SBarry Smith } 918d763cef2SBarry Smith 9194a2ae208SSatish Balay #undef __FUNCT__ 9204a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber" 921d763cef2SBarry Smith /*@ 922d763cef2SBarry Smith TSGetTimeStepNumber - Gets the current number of timesteps. 923d763cef2SBarry Smith 924d763cef2SBarry Smith Not Collective 925d763cef2SBarry Smith 926d763cef2SBarry Smith Input Parameter: 927d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 928d763cef2SBarry Smith 929d763cef2SBarry Smith Output Parameter: 930d763cef2SBarry Smith . iter - number steps so far 931d763cef2SBarry Smith 932d763cef2SBarry Smith Level: intermediate 933d763cef2SBarry Smith 934d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number 935d763cef2SBarry Smith @*/ 93663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter) 937d763cef2SBarry Smith { 938d763cef2SBarry Smith PetscFunctionBegin; 9394482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 9404482741eSBarry Smith PetscValidIntPointer(iter,2); 941d763cef2SBarry Smith *iter = ts->steps; 942d763cef2SBarry Smith PetscFunctionReturn(0); 943d763cef2SBarry Smith } 944d763cef2SBarry Smith 9454a2ae208SSatish Balay #undef __FUNCT__ 9464a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep" 947d763cef2SBarry Smith /*@ 948d763cef2SBarry Smith TSSetInitialTimeStep - Sets the initial timestep to be used, 949d763cef2SBarry Smith as well as the initial time. 950d763cef2SBarry Smith 951d763cef2SBarry Smith Collective on TS 952d763cef2SBarry Smith 953d763cef2SBarry Smith Input Parameters: 954d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 955d763cef2SBarry Smith . initial_time - the initial time 956d763cef2SBarry Smith - time_step - the size of the timestep 957d763cef2SBarry Smith 958d763cef2SBarry Smith Level: intermediate 959d763cef2SBarry Smith 960d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep() 961d763cef2SBarry Smith 962d763cef2SBarry Smith .keywords: TS, set, initial, timestep 963d763cef2SBarry Smith @*/ 96463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 965d763cef2SBarry Smith { 966d763cef2SBarry Smith PetscFunctionBegin; 9674482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 968d763cef2SBarry Smith ts->time_step = time_step; 969d763cef2SBarry Smith ts->initial_time_step = time_step; 970d763cef2SBarry Smith ts->ptime = initial_time; 971d763cef2SBarry Smith PetscFunctionReturn(0); 972d763cef2SBarry Smith } 973d763cef2SBarry Smith 9744a2ae208SSatish Balay #undef __FUNCT__ 9754a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep" 976d763cef2SBarry Smith /*@ 977d763cef2SBarry Smith TSSetTimeStep - Allows one to reset the timestep at any time, 978d763cef2SBarry Smith useful for simple pseudo-timestepping codes. 979d763cef2SBarry Smith 980d763cef2SBarry Smith Collective on TS 981d763cef2SBarry Smith 982d763cef2SBarry Smith Input Parameters: 983d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 984d763cef2SBarry Smith - time_step - the size of the timestep 985d763cef2SBarry Smith 986d763cef2SBarry Smith Level: intermediate 987d763cef2SBarry Smith 988d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 989d763cef2SBarry Smith 990d763cef2SBarry Smith .keywords: TS, set, timestep 991d763cef2SBarry Smith @*/ 99263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step) 993d763cef2SBarry Smith { 994d763cef2SBarry Smith PetscFunctionBegin; 9954482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 996d763cef2SBarry Smith ts->time_step = time_step; 997d763cef2SBarry Smith PetscFunctionReturn(0); 998d763cef2SBarry Smith } 999d763cef2SBarry Smith 10004a2ae208SSatish Balay #undef __FUNCT__ 10014a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep" 1002d763cef2SBarry Smith /*@ 1003d763cef2SBarry Smith TSGetTimeStep - Gets the current timestep size. 1004d763cef2SBarry Smith 1005d763cef2SBarry Smith Not Collective 1006d763cef2SBarry Smith 1007d763cef2SBarry Smith Input Parameter: 1008d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1009d763cef2SBarry Smith 1010d763cef2SBarry Smith Output Parameter: 1011d763cef2SBarry Smith . dt - the current timestep size 1012d763cef2SBarry Smith 1013d763cef2SBarry Smith Level: intermediate 1014d763cef2SBarry Smith 1015d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1016d763cef2SBarry Smith 1017d763cef2SBarry Smith .keywords: TS, get, timestep 1018d763cef2SBarry Smith @*/ 101963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt) 1020d763cef2SBarry Smith { 1021d763cef2SBarry Smith PetscFunctionBegin; 10224482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 10234482741eSBarry Smith PetscValidDoublePointer(dt,2); 1024d763cef2SBarry Smith *dt = ts->time_step; 1025d763cef2SBarry Smith PetscFunctionReturn(0); 1026d763cef2SBarry Smith } 1027d763cef2SBarry Smith 10284a2ae208SSatish Balay #undef __FUNCT__ 10294a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution" 1030d8e5e3e6SSatish Balay /*@ 1031d763cef2SBarry Smith TSGetSolution - Returns the solution at the present timestep. It 1032d763cef2SBarry Smith is valid to call this routine inside the function that you are evaluating 1033d763cef2SBarry Smith in order to move to the new timestep. This vector not changed until 1034d763cef2SBarry Smith the solution at the next timestep has been calculated. 1035d763cef2SBarry Smith 1036d763cef2SBarry Smith Not Collective, but Vec returned is parallel if TS is parallel 1037d763cef2SBarry Smith 1038d763cef2SBarry Smith Input Parameter: 1039d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1040d763cef2SBarry Smith 1041d763cef2SBarry Smith Output Parameter: 1042d763cef2SBarry Smith . v - the vector containing the solution 1043d763cef2SBarry Smith 1044d763cef2SBarry Smith Level: intermediate 1045d763cef2SBarry Smith 1046d763cef2SBarry Smith .seealso: TSGetTimeStep() 1047d763cef2SBarry Smith 1048d763cef2SBarry Smith .keywords: TS, timestep, get, solution 1049d763cef2SBarry Smith @*/ 105063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v) 1051d763cef2SBarry Smith { 1052d763cef2SBarry Smith PetscFunctionBegin; 10534482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 10544482741eSBarry Smith PetscValidPointer(v,2); 1055d763cef2SBarry Smith *v = ts->vec_sol_always; 1056d763cef2SBarry Smith PetscFunctionReturn(0); 1057d763cef2SBarry Smith } 1058d763cef2SBarry Smith 1059bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */ 10604a2ae208SSatish Balay #undef __FUNCT__ 1061bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType" 1062d8e5e3e6SSatish Balay /*@ 1063bdad233fSMatthew Knepley TSSetProblemType - Sets the type of problem to be solved. 1064d763cef2SBarry Smith 1065bdad233fSMatthew Knepley Not collective 1066d763cef2SBarry Smith 1067bdad233fSMatthew Knepley Input Parameters: 1068bdad233fSMatthew Knepley + ts - The TS 1069bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1070d763cef2SBarry Smith .vb 1071d763cef2SBarry Smith U_t = A U 1072d763cef2SBarry Smith U_t = A(t) U 1073d763cef2SBarry Smith U_t = F(t,U) 1074d763cef2SBarry Smith .ve 1075d763cef2SBarry Smith 1076d763cef2SBarry Smith Level: beginner 1077d763cef2SBarry Smith 1078bdad233fSMatthew Knepley .keywords: TS, problem type 1079bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1080d763cef2SBarry Smith @*/ 108163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type) 1082a7cc72afSBarry Smith { 1083d763cef2SBarry Smith PetscFunctionBegin; 10844482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1085bdad233fSMatthew Knepley ts->problem_type = type; 1086d763cef2SBarry Smith PetscFunctionReturn(0); 1087d763cef2SBarry Smith } 1088d763cef2SBarry Smith 1089bdad233fSMatthew Knepley #undef __FUNCT__ 1090bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType" 1091bdad233fSMatthew Knepley /*@C 1092bdad233fSMatthew Knepley TSGetProblemType - Gets the type of problem to be solved. 1093bdad233fSMatthew Knepley 1094bdad233fSMatthew Knepley Not collective 1095bdad233fSMatthew Knepley 1096bdad233fSMatthew Knepley Input Parameter: 1097bdad233fSMatthew Knepley . ts - The TS 1098bdad233fSMatthew Knepley 1099bdad233fSMatthew Knepley Output Parameter: 1100bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1101bdad233fSMatthew Knepley .vb 1102bdad233fSMatthew Knepley U_t = A U 1103bdad233fSMatthew Knepley U_t = A(t) U 1104bdad233fSMatthew Knepley U_t = F(t,U) 1105bdad233fSMatthew Knepley .ve 1106bdad233fSMatthew Knepley 1107bdad233fSMatthew Knepley Level: beginner 1108bdad233fSMatthew Knepley 1109bdad233fSMatthew Knepley .keywords: TS, problem type 1110bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1111bdad233fSMatthew Knepley @*/ 111263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type) 1113a7cc72afSBarry Smith { 1114bdad233fSMatthew Knepley PetscFunctionBegin; 11154482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 11164482741eSBarry Smith PetscValidIntPointer(type,2); 1117bdad233fSMatthew Knepley *type = ts->problem_type; 1118bdad233fSMatthew Knepley PetscFunctionReturn(0); 1119bdad233fSMatthew Knepley } 1120d763cef2SBarry Smith 11214a2ae208SSatish Balay #undef __FUNCT__ 11224a2ae208SSatish Balay #define __FUNCT__ "TSSetUp" 1123d763cef2SBarry Smith /*@ 1124d763cef2SBarry Smith TSSetUp - Sets up the internal data structures for the later use 1125d763cef2SBarry Smith of a timestepper. 1126d763cef2SBarry Smith 1127d763cef2SBarry Smith Collective on TS 1128d763cef2SBarry Smith 1129d763cef2SBarry Smith Input Parameter: 1130d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1131d763cef2SBarry Smith 1132d763cef2SBarry Smith Notes: 1133d763cef2SBarry Smith For basic use of the TS solvers the user need not explicitly call 1134d763cef2SBarry Smith TSSetUp(), since these actions will automatically occur during 1135d763cef2SBarry Smith the call to TSStep(). However, if one wishes to control this 1136d763cef2SBarry Smith phase separately, TSSetUp() should be called after TSCreate() 1137d763cef2SBarry Smith and optional routines of the form TSSetXXX(), but before TSStep(). 1138d763cef2SBarry Smith 1139d763cef2SBarry Smith Level: advanced 1140d763cef2SBarry Smith 1141d763cef2SBarry Smith .keywords: TS, timestep, setup 1142d763cef2SBarry Smith 1143d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy() 1144d763cef2SBarry Smith @*/ 114563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts) 1146d763cef2SBarry Smith { 1147dfbe8321SBarry Smith PetscErrorCode ierr; 1148d763cef2SBarry Smith 1149d763cef2SBarry Smith PetscFunctionBegin; 11504482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 115129bbc08cSBarry Smith if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 11527adad957SLisandro Dalcin if (!((PetscObject)ts)->type_name) { 1153d763cef2SBarry Smith ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr); 1154d763cef2SBarry Smith } 1155000e7ae3SMatthew Knepley ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1156d763cef2SBarry Smith ts->setupcalled = 1; 1157d763cef2SBarry Smith PetscFunctionReturn(0); 1158d763cef2SBarry Smith } 1159d763cef2SBarry Smith 11604a2ae208SSatish Balay #undef __FUNCT__ 11614a2ae208SSatish Balay #define __FUNCT__ "TSDestroy" 1162d8e5e3e6SSatish Balay /*@ 1163d763cef2SBarry Smith TSDestroy - Destroys the timestepper context that was created 1164d763cef2SBarry Smith with TSCreate(). 1165d763cef2SBarry Smith 1166d763cef2SBarry Smith Collective on TS 1167d763cef2SBarry Smith 1168d763cef2SBarry Smith Input Parameter: 1169d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1170d763cef2SBarry Smith 1171d763cef2SBarry Smith Level: beginner 1172d763cef2SBarry Smith 1173d763cef2SBarry Smith .keywords: TS, timestepper, destroy 1174d763cef2SBarry Smith 1175d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve() 1176d763cef2SBarry Smith @*/ 117763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts) 1178d763cef2SBarry Smith { 11796849ba73SBarry Smith PetscErrorCode ierr; 1180d763cef2SBarry Smith 1181d763cef2SBarry Smith PetscFunctionBegin; 11824482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 11837adad957SLisandro Dalcin if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0); 1184d763cef2SBarry Smith 1185be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 11860f5bd95cSBarry Smith ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 11878beb423aSHong Zhang if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr)} 118894b7f48cSBarry Smith if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);} 1189d763cef2SBarry Smith if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 11901e3347e8SBarry Smith if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);} 1191a6570f20SBarry Smith ierr = TSMonitorCancel(ts);CHKERRQ(ierr); 1192a79aaaedSSatish Balay ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1193d763cef2SBarry Smith PetscFunctionReturn(0); 1194d763cef2SBarry Smith } 1195d763cef2SBarry Smith 11964a2ae208SSatish Balay #undef __FUNCT__ 11974a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES" 1198d8e5e3e6SSatish Balay /*@ 1199d763cef2SBarry Smith TSGetSNES - Returns the SNES (nonlinear solver) associated with 1200d763cef2SBarry Smith a TS (timestepper) context. Valid only for nonlinear problems. 1201d763cef2SBarry Smith 1202d763cef2SBarry Smith Not Collective, but SNES is parallel if TS is parallel 1203d763cef2SBarry Smith 1204d763cef2SBarry Smith Input Parameter: 1205d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1206d763cef2SBarry Smith 1207d763cef2SBarry Smith Output Parameter: 1208d763cef2SBarry Smith . snes - the nonlinear solver context 1209d763cef2SBarry Smith 1210d763cef2SBarry Smith Notes: 1211d763cef2SBarry Smith The user can then directly manipulate the SNES context to set various 1212d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 121394b7f48cSBarry Smith KSP, KSP, and PC contexts as well. 1214d763cef2SBarry Smith 1215d763cef2SBarry Smith TSGetSNES() does not work for integrators that do not use SNES; in 1216d763cef2SBarry Smith this case TSGetSNES() returns PETSC_NULL in snes. 1217d763cef2SBarry Smith 1218d763cef2SBarry Smith Level: beginner 1219d763cef2SBarry Smith 1220d763cef2SBarry Smith .keywords: timestep, get, SNES 1221d763cef2SBarry Smith @*/ 122263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes) 1223d763cef2SBarry Smith { 1224d763cef2SBarry Smith PetscFunctionBegin; 12254482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 12264482741eSBarry Smith PetscValidPointer(snes,2); 1227447600ffSHong Zhang if (((PetscObject)ts)->type_name == PETSC_NULL) 1228447600ffSHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first"); 122994b7f48cSBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1230d763cef2SBarry Smith *snes = ts->snes; 1231d763cef2SBarry Smith PetscFunctionReturn(0); 1232d763cef2SBarry Smith } 1233d763cef2SBarry Smith 12344a2ae208SSatish Balay #undef __FUNCT__ 123594b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP" 1236d8e5e3e6SSatish Balay /*@ 123794b7f48cSBarry Smith TSGetKSP - Returns the KSP (linear solver) associated with 1238d763cef2SBarry Smith a TS (timestepper) context. 1239d763cef2SBarry Smith 124094b7f48cSBarry Smith Not Collective, but KSP is parallel if TS is parallel 1241d763cef2SBarry Smith 1242d763cef2SBarry Smith Input Parameter: 1243d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1244d763cef2SBarry Smith 1245d763cef2SBarry Smith Output Parameter: 124694b7f48cSBarry Smith . ksp - the nonlinear solver context 1247d763cef2SBarry Smith 1248d763cef2SBarry Smith Notes: 124994b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 1250d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 1251d763cef2SBarry Smith KSP and PC contexts as well. 1252d763cef2SBarry Smith 125394b7f48cSBarry Smith TSGetKSP() does not work for integrators that do not use KSP; 125494b7f48cSBarry Smith in this case TSGetKSP() returns PETSC_NULL in ksp. 1255d763cef2SBarry Smith 1256d763cef2SBarry Smith Level: beginner 1257d763cef2SBarry Smith 125894b7f48cSBarry Smith .keywords: timestep, get, KSP 1259d763cef2SBarry Smith @*/ 126063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp) 1261d763cef2SBarry Smith { 1262d763cef2SBarry Smith PetscFunctionBegin; 12634482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 12644482741eSBarry Smith PetscValidPointer(ksp,2); 1265988402f6SHong Zhang if (((PetscObject)ts)->type_name == PETSC_NULL) 1266988402f6SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 126729bbc08cSBarry Smith if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 126894b7f48cSBarry Smith *ksp = ts->ksp; 1269d763cef2SBarry Smith PetscFunctionReturn(0); 1270d763cef2SBarry Smith } 1271d763cef2SBarry Smith 1272d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */ 1273d763cef2SBarry Smith 12744a2ae208SSatish Balay #undef __FUNCT__ 1275adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration" 1276adb62b0dSMatthew Knepley /*@ 1277adb62b0dSMatthew Knepley TSGetDuration - Gets the maximum number of timesteps to use and 1278adb62b0dSMatthew Knepley maximum time for iteration. 1279adb62b0dSMatthew Knepley 1280adb62b0dSMatthew Knepley Collective on TS 1281adb62b0dSMatthew Knepley 1282adb62b0dSMatthew Knepley Input Parameters: 1283adb62b0dSMatthew Knepley + ts - the TS context obtained from TSCreate() 1284adb62b0dSMatthew Knepley . maxsteps - maximum number of iterations to use, or PETSC_NULL 1285adb62b0dSMatthew Knepley - maxtime - final time to iterate to, or PETSC_NULL 1286adb62b0dSMatthew Knepley 1287adb62b0dSMatthew Knepley Level: intermediate 1288adb62b0dSMatthew Knepley 1289adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time 1290adb62b0dSMatthew Knepley @*/ 129163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1292adb62b0dSMatthew Knepley { 1293adb62b0dSMatthew Knepley PetscFunctionBegin; 12944482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1295abc0a331SBarry Smith if (maxsteps) { 12964482741eSBarry Smith PetscValidIntPointer(maxsteps,2); 1297adb62b0dSMatthew Knepley *maxsteps = ts->max_steps; 1298adb62b0dSMatthew Knepley } 1299abc0a331SBarry Smith if (maxtime ) { 13004482741eSBarry Smith PetscValidScalarPointer(maxtime,3); 1301adb62b0dSMatthew Knepley *maxtime = ts->max_time; 1302adb62b0dSMatthew Knepley } 1303adb62b0dSMatthew Knepley PetscFunctionReturn(0); 1304adb62b0dSMatthew Knepley } 1305adb62b0dSMatthew Knepley 1306adb62b0dSMatthew Knepley #undef __FUNCT__ 13074a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration" 1308d763cef2SBarry Smith /*@ 1309d763cef2SBarry Smith TSSetDuration - Sets the maximum number of timesteps to use and 1310d763cef2SBarry Smith maximum time for iteration. 1311d763cef2SBarry Smith 1312d763cef2SBarry Smith Collective on TS 1313d763cef2SBarry Smith 1314d763cef2SBarry Smith Input Parameters: 1315d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1316d763cef2SBarry Smith . maxsteps - maximum number of iterations to use 1317d763cef2SBarry Smith - maxtime - final time to iterate to 1318d763cef2SBarry Smith 1319d763cef2SBarry Smith Options Database Keys: 1320d763cef2SBarry Smith . -ts_max_steps <maxsteps> - Sets maxsteps 1321d763cef2SBarry Smith . -ts_max_time <maxtime> - Sets maxtime 1322d763cef2SBarry Smith 1323d763cef2SBarry Smith Notes: 1324d763cef2SBarry Smith The default maximum number of iterations is 5000. Default time is 5.0 1325d763cef2SBarry Smith 1326d763cef2SBarry Smith Level: intermediate 1327d763cef2SBarry Smith 1328d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations 1329d763cef2SBarry Smith @*/ 133063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1331d763cef2SBarry Smith { 1332d763cef2SBarry Smith PetscFunctionBegin; 13334482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1334d763cef2SBarry Smith ts->max_steps = maxsteps; 1335d763cef2SBarry Smith ts->max_time = maxtime; 1336d763cef2SBarry Smith PetscFunctionReturn(0); 1337d763cef2SBarry Smith } 1338d763cef2SBarry Smith 13394a2ae208SSatish Balay #undef __FUNCT__ 13404a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution" 1341d763cef2SBarry Smith /*@ 1342d763cef2SBarry Smith TSSetSolution - Sets the initial solution vector 1343d763cef2SBarry Smith for use by the TS routines. 1344d763cef2SBarry Smith 1345d763cef2SBarry Smith Collective on TS and Vec 1346d763cef2SBarry Smith 1347d763cef2SBarry Smith Input Parameters: 1348d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1349d763cef2SBarry Smith - x - the solution vector 1350d763cef2SBarry Smith 1351d763cef2SBarry Smith Level: beginner 1352d763cef2SBarry Smith 1353d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions 1354d763cef2SBarry Smith @*/ 135563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x) 1356d763cef2SBarry Smith { 1357d763cef2SBarry Smith PetscFunctionBegin; 13584482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 13594482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1360d763cef2SBarry Smith ts->vec_sol = ts->vec_sol_always = x; 1361d763cef2SBarry Smith PetscFunctionReturn(0); 1362d763cef2SBarry Smith } 1363d763cef2SBarry Smith 1364e74ef692SMatthew Knepley #undef __FUNCT__ 1365e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep" 1366ac226902SBarry Smith /*@C 1367000e7ae3SMatthew Knepley TSSetPreStep - Sets the general-purpose function 1368000e7ae3SMatthew Knepley called once at the beginning of time stepping. 1369000e7ae3SMatthew Knepley 1370000e7ae3SMatthew Knepley Collective on TS 1371000e7ae3SMatthew Knepley 1372000e7ae3SMatthew Knepley Input Parameters: 1373000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1374000e7ae3SMatthew Knepley - func - The function 1375000e7ae3SMatthew Knepley 1376000e7ae3SMatthew Knepley Calling sequence of func: 1377000e7ae3SMatthew Knepley . func (TS ts); 1378000e7ae3SMatthew Knepley 1379000e7ae3SMatthew Knepley Level: intermediate 1380000e7ae3SMatthew Knepley 1381000e7ae3SMatthew Knepley .keywords: TS, timestep 1382000e7ae3SMatthew Knepley @*/ 138363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1384000e7ae3SMatthew Knepley { 1385000e7ae3SMatthew Knepley PetscFunctionBegin; 13864482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1387000e7ae3SMatthew Knepley ts->ops->prestep = func; 1388000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1389000e7ae3SMatthew Knepley } 1390000e7ae3SMatthew Knepley 1391e74ef692SMatthew Knepley #undef __FUNCT__ 1392e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep" 1393000e7ae3SMatthew Knepley /*@ 1394000e7ae3SMatthew Knepley TSDefaultPreStep - The default pre-stepping function which does nothing. 1395000e7ae3SMatthew Knepley 1396000e7ae3SMatthew Knepley Collective on TS 1397000e7ae3SMatthew Knepley 1398000e7ae3SMatthew Knepley Input Parameters: 1399000e7ae3SMatthew Knepley . ts - The TS context obtained from TSCreate() 1400000e7ae3SMatthew Knepley 1401000e7ae3SMatthew Knepley Level: developer 1402000e7ae3SMatthew Knepley 1403000e7ae3SMatthew Knepley .keywords: TS, timestep 1404000e7ae3SMatthew Knepley @*/ 140563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts) 1406000e7ae3SMatthew Knepley { 1407000e7ae3SMatthew Knepley PetscFunctionBegin; 1408000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1409000e7ae3SMatthew Knepley } 1410000e7ae3SMatthew Knepley 1411e74ef692SMatthew Knepley #undef __FUNCT__ 1412e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep" 1413ac226902SBarry Smith /*@C 1414000e7ae3SMatthew Knepley TSSetPostStep - Sets the general-purpose function 1415000e7ae3SMatthew Knepley called once at the end of time stepping. 1416000e7ae3SMatthew Knepley 1417000e7ae3SMatthew Knepley Collective on TS 1418000e7ae3SMatthew Knepley 1419000e7ae3SMatthew Knepley Input Parameters: 1420000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1421000e7ae3SMatthew Knepley - func - The function 1422000e7ae3SMatthew Knepley 1423000e7ae3SMatthew Knepley Calling sequence of func: 1424000e7ae3SMatthew Knepley . func (TS ts); 1425000e7ae3SMatthew Knepley 1426000e7ae3SMatthew Knepley Level: intermediate 1427000e7ae3SMatthew Knepley 1428000e7ae3SMatthew Knepley .keywords: TS, timestep 1429000e7ae3SMatthew Knepley @*/ 143063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1431000e7ae3SMatthew Knepley { 1432000e7ae3SMatthew Knepley PetscFunctionBegin; 14334482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1434000e7ae3SMatthew Knepley ts->ops->poststep = func; 1435000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1436000e7ae3SMatthew Knepley } 1437000e7ae3SMatthew Knepley 1438e74ef692SMatthew Knepley #undef __FUNCT__ 1439e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep" 1440000e7ae3SMatthew Knepley /*@ 1441000e7ae3SMatthew Knepley TSDefaultPostStep - The default post-stepping function which does nothing. 1442000e7ae3SMatthew Knepley 1443000e7ae3SMatthew Knepley Collective on TS 1444000e7ae3SMatthew Knepley 1445000e7ae3SMatthew Knepley Input Parameters: 1446000e7ae3SMatthew Knepley . ts - The TS context obtained from TSCreate() 1447000e7ae3SMatthew Knepley 1448000e7ae3SMatthew Knepley Level: developer 1449000e7ae3SMatthew Knepley 1450000e7ae3SMatthew Knepley .keywords: TS, timestep 1451000e7ae3SMatthew Knepley @*/ 145263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts) 1453000e7ae3SMatthew Knepley { 1454000e7ae3SMatthew Knepley PetscFunctionBegin; 1455000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1456000e7ae3SMatthew Knepley } 1457000e7ae3SMatthew Knepley 1458d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 1459d763cef2SBarry Smith 14604a2ae208SSatish Balay #undef __FUNCT__ 1461a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet" 1462d763cef2SBarry Smith /*@C 1463a6570f20SBarry Smith TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1464d763cef2SBarry Smith timestep to display the iteration's progress. 1465d763cef2SBarry Smith 1466d763cef2SBarry Smith Collective on TS 1467d763cef2SBarry Smith 1468d763cef2SBarry Smith Input Parameters: 1469d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1470d763cef2SBarry Smith . func - monitoring routine 1471329f5518SBarry Smith . mctx - [optional] user-defined context for private data for the 1472b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desired) 1473b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1474b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 1475d763cef2SBarry Smith 1476d763cef2SBarry Smith Calling sequence of func: 1477a7cc72afSBarry Smith $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1478d763cef2SBarry Smith 1479d763cef2SBarry Smith + ts - the TS context 1480d763cef2SBarry Smith . steps - iteration number 14811f06c33eSBarry Smith . time - current time 1482d763cef2SBarry Smith . x - current iterate 1483d763cef2SBarry Smith - mctx - [optional] monitoring context 1484d763cef2SBarry Smith 1485d763cef2SBarry Smith Notes: 1486d763cef2SBarry Smith This routine adds an additional monitor to the list of monitors that 1487d763cef2SBarry Smith already has been loaded. 1488d763cef2SBarry Smith 1489025f1a04SBarry Smith Fortran notes: Only a single monitor function can be set for each TS object 1490025f1a04SBarry Smith 1491d763cef2SBarry Smith Level: intermediate 1492d763cef2SBarry Smith 1493d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1494d763cef2SBarry Smith 1495a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel() 1496d763cef2SBarry Smith @*/ 1497a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1498d763cef2SBarry Smith { 1499d763cef2SBarry Smith PetscFunctionBegin; 15004482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1501d763cef2SBarry Smith if (ts->numbermonitors >= MAXTSMONITORS) { 150229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1503d763cef2SBarry Smith } 1504d763cef2SBarry Smith ts->monitor[ts->numbermonitors] = monitor; 1505329f5518SBarry Smith ts->mdestroy[ts->numbermonitors] = mdestroy; 1506d763cef2SBarry Smith ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1507d763cef2SBarry Smith PetscFunctionReturn(0); 1508d763cef2SBarry Smith } 1509d763cef2SBarry Smith 15104a2ae208SSatish Balay #undef __FUNCT__ 1511a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel" 1512d763cef2SBarry Smith /*@C 1513a6570f20SBarry Smith TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1514d763cef2SBarry Smith 1515d763cef2SBarry Smith Collective on TS 1516d763cef2SBarry Smith 1517d763cef2SBarry Smith Input Parameters: 1518d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1519d763cef2SBarry Smith 1520d763cef2SBarry Smith Notes: 1521d763cef2SBarry Smith There is no way to remove a single, specific monitor. 1522d763cef2SBarry Smith 1523d763cef2SBarry Smith Level: intermediate 1524d763cef2SBarry Smith 1525d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1526d763cef2SBarry Smith 1527a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet() 1528d763cef2SBarry Smith @*/ 1529a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts) 1530d763cef2SBarry Smith { 1531d952e501SBarry Smith PetscErrorCode ierr; 1532d952e501SBarry Smith PetscInt i; 1533d952e501SBarry Smith 1534d763cef2SBarry Smith PetscFunctionBegin; 15354482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1536d952e501SBarry Smith for (i=0; i<ts->numbermonitors; i++) { 1537d952e501SBarry Smith if (ts->mdestroy[i]) { 1538d952e501SBarry Smith ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 1539d952e501SBarry Smith } 1540d952e501SBarry Smith } 1541d763cef2SBarry Smith ts->numbermonitors = 0; 1542d763cef2SBarry Smith PetscFunctionReturn(0); 1543d763cef2SBarry Smith } 1544d763cef2SBarry Smith 15454a2ae208SSatish Balay #undef __FUNCT__ 1546a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault" 1547d8e5e3e6SSatish Balay /*@ 1548a6570f20SBarry Smith TSMonitorDefault - Sets the Default monitor 15495516499fSSatish Balay 15505516499fSSatish Balay Level: intermediate 155141251cbbSSatish Balay 15525516499fSSatish Balay .keywords: TS, set, monitor 15535516499fSSatish Balay 155441251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet() 155541251cbbSSatish Balay @*/ 1556a6570f20SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1557d763cef2SBarry Smith { 1558dfbe8321SBarry Smith PetscErrorCode ierr; 1559a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 1560d132466eSBarry Smith 1561d763cef2SBarry Smith PetscFunctionBegin; 1562a34d58ebSBarry Smith if (!ctx) { 15637adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 1564a34d58ebSBarry Smith } 1565a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1566a34d58ebSBarry Smith if (!ctx) { 1567a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 1568a34d58ebSBarry Smith } 1569d763cef2SBarry Smith PetscFunctionReturn(0); 1570d763cef2SBarry Smith } 1571d763cef2SBarry Smith 15724a2ae208SSatish Balay #undef __FUNCT__ 15734a2ae208SSatish Balay #define __FUNCT__ "TSStep" 1574d763cef2SBarry Smith /*@ 1575d763cef2SBarry Smith TSStep - Steps the requested number of timesteps. 1576d763cef2SBarry Smith 1577d763cef2SBarry Smith Collective on TS 1578d763cef2SBarry Smith 1579d763cef2SBarry Smith Input Parameter: 1580d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1581d763cef2SBarry Smith 1582d763cef2SBarry Smith Output Parameters: 1583d763cef2SBarry Smith + steps - number of iterations until termination 1584142b95e3SSatish Balay - ptime - time until termination 1585d763cef2SBarry Smith 1586d763cef2SBarry Smith Level: beginner 1587d763cef2SBarry Smith 1588d763cef2SBarry Smith .keywords: TS, timestep, solve 1589d763cef2SBarry Smith 1590d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy() 1591d763cef2SBarry Smith @*/ 159263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1593d763cef2SBarry Smith { 1594dfbe8321SBarry Smith PetscErrorCode ierr; 1595d763cef2SBarry Smith 1596d763cef2SBarry Smith PetscFunctionBegin; 15974482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1598d405a339SMatthew Knepley if (!ts->setupcalled) { 1599d405a339SMatthew Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 1600d405a339SMatthew Knepley } 1601d405a339SMatthew Knepley 1602d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1603d405a339SMatthew Knepley ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1604000e7ae3SMatthew Knepley ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1605d405a339SMatthew Knepley ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1606d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1607d405a339SMatthew Knepley 16084bb05414SBarry Smith if (!PetscPreLoadingOn) { 16097adad957SLisandro Dalcin ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1610d405a339SMatthew Knepley } 1611d763cef2SBarry Smith PetscFunctionReturn(0); 1612d763cef2SBarry Smith } 1613d763cef2SBarry Smith 16144a2ae208SSatish Balay #undef __FUNCT__ 16156a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve" 16166a4d4014SLisandro Dalcin /*@ 16176a4d4014SLisandro Dalcin TSSolve - Steps the requested number of timesteps. 16186a4d4014SLisandro Dalcin 16196a4d4014SLisandro Dalcin Collective on TS 16206a4d4014SLisandro Dalcin 16216a4d4014SLisandro Dalcin Input Parameter: 16226a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 16236a4d4014SLisandro Dalcin - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 16246a4d4014SLisandro Dalcin 16256a4d4014SLisandro Dalcin Level: beginner 16266a4d4014SLisandro Dalcin 16276a4d4014SLisandro Dalcin .keywords: TS, timestep, solve 16286a4d4014SLisandro Dalcin 16296a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep() 16306a4d4014SLisandro Dalcin @*/ 16316a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x) 16326a4d4014SLisandro Dalcin { 16336a4d4014SLisandro Dalcin PetscInt steps; 16346a4d4014SLisandro Dalcin PetscReal ptime; 16356a4d4014SLisandro Dalcin PetscErrorCode ierr; 16366a4d4014SLisandro Dalcin PetscFunctionBegin; 16376a4d4014SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_COOKIE,1); 16386a4d4014SLisandro Dalcin /* set solution vector if provided */ 16396a4d4014SLisandro Dalcin if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 16406a4d4014SLisandro Dalcin /* reset time step and iteration counters */ 16416a4d4014SLisandro Dalcin ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 16426a4d4014SLisandro Dalcin /* steps the requested number of timesteps. */ 16436a4d4014SLisandro Dalcin ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 16446a4d4014SLisandro Dalcin PetscFunctionReturn(0); 16456a4d4014SLisandro Dalcin } 16466a4d4014SLisandro Dalcin 16476a4d4014SLisandro Dalcin #undef __FUNCT__ 16484a2ae208SSatish Balay #define __FUNCT__ "TSMonitor" 1649d763cef2SBarry Smith /* 1650d763cef2SBarry Smith Runs the user provided monitor routines, if they exists. 1651d763cef2SBarry Smith */ 1652a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1653d763cef2SBarry Smith { 16546849ba73SBarry Smith PetscErrorCode ierr; 1655a7cc72afSBarry Smith PetscInt i,n = ts->numbermonitors; 1656d763cef2SBarry Smith 1657d763cef2SBarry Smith PetscFunctionBegin; 1658d763cef2SBarry Smith for (i=0; i<n; i++) { 1659142b95e3SSatish Balay ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1660d763cef2SBarry Smith } 1661d763cef2SBarry Smith PetscFunctionReturn(0); 1662d763cef2SBarry Smith } 1663d763cef2SBarry Smith 1664d763cef2SBarry Smith /* ------------------------------------------------------------------------*/ 1665d763cef2SBarry Smith 16664a2ae208SSatish Balay #undef __FUNCT__ 1667a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate" 1668d763cef2SBarry Smith /*@C 1669a6570f20SBarry Smith TSMonitorLGCreate - Creates a line graph context for use with 1670d763cef2SBarry Smith TS to monitor convergence of preconditioned residual norms. 1671d763cef2SBarry Smith 1672d763cef2SBarry Smith Collective on TS 1673d763cef2SBarry Smith 1674d763cef2SBarry Smith Input Parameters: 1675d763cef2SBarry Smith + host - the X display to open, or null for the local machine 1676d763cef2SBarry Smith . label - the title to put in the title bar 16777c922b88SBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 1678d763cef2SBarry Smith - m, n - the screen width and height in pixels 1679d763cef2SBarry Smith 1680d763cef2SBarry Smith Output Parameter: 1681d763cef2SBarry Smith . draw - the drawing context 1682d763cef2SBarry Smith 1683d763cef2SBarry Smith Options Database Key: 1684a6570f20SBarry Smith . -ts_monitor_draw - automatically sets line graph monitor 1685d763cef2SBarry Smith 1686d763cef2SBarry Smith Notes: 1687a6570f20SBarry Smith Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1688d763cef2SBarry Smith 1689d763cef2SBarry Smith Level: intermediate 1690d763cef2SBarry Smith 16917c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso 1692d763cef2SBarry Smith 1693a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet() 16947c922b88SBarry Smith 1695d763cef2SBarry Smith @*/ 1696a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1697d763cef2SBarry Smith { 1698b0a32e0cSBarry Smith PetscDraw win; 1699dfbe8321SBarry Smith PetscErrorCode ierr; 1700d763cef2SBarry Smith 1701d763cef2SBarry Smith PetscFunctionBegin; 1702b0a32e0cSBarry Smith ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1703b0a32e0cSBarry Smith ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1704b0a32e0cSBarry Smith ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1705b0a32e0cSBarry Smith ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1706d763cef2SBarry Smith 170752e6d16bSBarry Smith ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1708d763cef2SBarry Smith PetscFunctionReturn(0); 1709d763cef2SBarry Smith } 1710d763cef2SBarry Smith 17114a2ae208SSatish Balay #undef __FUNCT__ 1712a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG" 1713a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1714d763cef2SBarry Smith { 1715b0a32e0cSBarry Smith PetscDrawLG lg = (PetscDrawLG) monctx; 171687828ca2SBarry Smith PetscReal x,y = ptime; 1717dfbe8321SBarry Smith PetscErrorCode ierr; 1718d763cef2SBarry Smith 1719d763cef2SBarry Smith PetscFunctionBegin; 17207c922b88SBarry Smith if (!monctx) { 17217c922b88SBarry Smith MPI_Comm comm; 1722b0a32e0cSBarry Smith PetscViewer viewer; 17237c922b88SBarry Smith 17247c922b88SBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1725b0a32e0cSBarry Smith viewer = PETSC_VIEWER_DRAW_(comm); 1726b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 17277c922b88SBarry Smith } 17287c922b88SBarry Smith 1729b0a32e0cSBarry Smith if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 173087828ca2SBarry Smith x = (PetscReal)n; 1731b0a32e0cSBarry Smith ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1732d763cef2SBarry Smith if (n < 20 || (n % 5)) { 1733b0a32e0cSBarry Smith ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1734d763cef2SBarry Smith } 1735d763cef2SBarry Smith PetscFunctionReturn(0); 1736d763cef2SBarry Smith } 1737d763cef2SBarry Smith 17384a2ae208SSatish Balay #undef __FUNCT__ 1739a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy" 1740d763cef2SBarry Smith /*@C 1741a6570f20SBarry Smith TSMonitorLGDestroy - Destroys a line graph context that was created 1742a6570f20SBarry Smith with TSMonitorLGCreate(). 1743d763cef2SBarry Smith 1744b0a32e0cSBarry Smith Collective on PetscDrawLG 1745d763cef2SBarry Smith 1746d763cef2SBarry Smith Input Parameter: 1747d763cef2SBarry Smith . draw - the drawing context 1748d763cef2SBarry Smith 1749d763cef2SBarry Smith Level: intermediate 1750d763cef2SBarry Smith 1751d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy 1752d763cef2SBarry Smith 1753a6570f20SBarry Smith .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1754d763cef2SBarry Smith @*/ 1755a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg) 1756d763cef2SBarry Smith { 1757b0a32e0cSBarry Smith PetscDraw draw; 1758dfbe8321SBarry Smith PetscErrorCode ierr; 1759d763cef2SBarry Smith 1760d763cef2SBarry Smith PetscFunctionBegin; 1761b0a32e0cSBarry Smith ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1762b0a32e0cSBarry Smith ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1763b0a32e0cSBarry Smith ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1764d763cef2SBarry Smith PetscFunctionReturn(0); 1765d763cef2SBarry Smith } 1766d763cef2SBarry Smith 17674a2ae208SSatish Balay #undef __FUNCT__ 17684a2ae208SSatish Balay #define __FUNCT__ "TSGetTime" 1769d763cef2SBarry Smith /*@ 1770d763cef2SBarry Smith TSGetTime - Gets the current time. 1771d763cef2SBarry Smith 1772d763cef2SBarry Smith Not Collective 1773d763cef2SBarry Smith 1774d763cef2SBarry Smith Input Parameter: 1775d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1776d763cef2SBarry Smith 1777d763cef2SBarry Smith Output Parameter: 1778d763cef2SBarry Smith . t - the current time 1779d763cef2SBarry Smith 1780d763cef2SBarry Smith Level: beginner 1781d763cef2SBarry Smith 1782d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1783d763cef2SBarry Smith 1784d763cef2SBarry Smith .keywords: TS, get, time 1785d763cef2SBarry Smith @*/ 178663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t) 1787d763cef2SBarry Smith { 1788d763cef2SBarry Smith PetscFunctionBegin; 17894482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 17904482741eSBarry Smith PetscValidDoublePointer(t,2); 1791d763cef2SBarry Smith *t = ts->ptime; 1792d763cef2SBarry Smith PetscFunctionReturn(0); 1793d763cef2SBarry Smith } 1794d763cef2SBarry Smith 17954a2ae208SSatish Balay #undef __FUNCT__ 17966a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime" 17976a4d4014SLisandro Dalcin /*@ 17986a4d4014SLisandro Dalcin TSSetTime - Allows one to reset the time. 17996a4d4014SLisandro Dalcin 18006a4d4014SLisandro Dalcin Collective on TS 18016a4d4014SLisandro Dalcin 18026a4d4014SLisandro Dalcin Input Parameters: 18036a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 18046a4d4014SLisandro Dalcin - time - the time 18056a4d4014SLisandro Dalcin 18066a4d4014SLisandro Dalcin Level: intermediate 18076a4d4014SLisandro Dalcin 18086a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration() 18096a4d4014SLisandro Dalcin 18106a4d4014SLisandro Dalcin .keywords: TS, set, time 18116a4d4014SLisandro Dalcin @*/ 18126a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t) 18136a4d4014SLisandro Dalcin { 18146a4d4014SLisandro Dalcin PetscFunctionBegin; 18156a4d4014SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_COOKIE,1); 18166a4d4014SLisandro Dalcin ts->ptime = t; 18176a4d4014SLisandro Dalcin PetscFunctionReturn(0); 18186a4d4014SLisandro Dalcin } 18196a4d4014SLisandro Dalcin 18206a4d4014SLisandro Dalcin #undef __FUNCT__ 18214a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix" 1822d763cef2SBarry Smith /*@C 1823d763cef2SBarry Smith TSSetOptionsPrefix - Sets the prefix used for searching for all 1824d763cef2SBarry Smith TS options in the database. 1825d763cef2SBarry Smith 1826d763cef2SBarry Smith Collective on TS 1827d763cef2SBarry Smith 1828d763cef2SBarry Smith Input Parameter: 1829d763cef2SBarry Smith + ts - The TS context 1830d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1831d763cef2SBarry Smith 1832d763cef2SBarry Smith Notes: 1833d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1834d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1835d763cef2SBarry Smith hyphen. 1836d763cef2SBarry Smith 1837d763cef2SBarry Smith Level: advanced 1838d763cef2SBarry Smith 1839d763cef2SBarry Smith .keywords: TS, set, options, prefix, database 1840d763cef2SBarry Smith 1841d763cef2SBarry Smith .seealso: TSSetFromOptions() 1842d763cef2SBarry Smith 1843d763cef2SBarry Smith @*/ 184463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[]) 1845d763cef2SBarry Smith { 1846dfbe8321SBarry Smith PetscErrorCode ierr; 1847d763cef2SBarry Smith 1848d763cef2SBarry Smith PetscFunctionBegin; 18494482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1850d763cef2SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1851d763cef2SBarry Smith switch(ts->problem_type) { 1852d763cef2SBarry Smith case TS_NONLINEAR: 185359580b9cSBarry Smith if (ts->snes) { 1854d763cef2SBarry Smith ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 185559580b9cSBarry Smith } 1856d763cef2SBarry Smith break; 1857d763cef2SBarry Smith case TS_LINEAR: 185859580b9cSBarry Smith if (ts->ksp) { 185994b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 186059580b9cSBarry Smith } 1861d763cef2SBarry Smith break; 1862d763cef2SBarry Smith } 1863d763cef2SBarry Smith PetscFunctionReturn(0); 1864d763cef2SBarry Smith } 1865d763cef2SBarry Smith 1866d763cef2SBarry Smith 18674a2ae208SSatish Balay #undef __FUNCT__ 18684a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix" 1869d763cef2SBarry Smith /*@C 1870d763cef2SBarry Smith TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1871d763cef2SBarry Smith TS options in the database. 1872d763cef2SBarry Smith 1873d763cef2SBarry Smith Collective on TS 1874d763cef2SBarry Smith 1875d763cef2SBarry Smith Input Parameter: 1876d763cef2SBarry Smith + ts - The TS context 1877d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1878d763cef2SBarry Smith 1879d763cef2SBarry Smith Notes: 1880d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1881d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1882d763cef2SBarry Smith hyphen. 1883d763cef2SBarry Smith 1884d763cef2SBarry Smith Level: advanced 1885d763cef2SBarry Smith 1886d763cef2SBarry Smith .keywords: TS, append, options, prefix, database 1887d763cef2SBarry Smith 1888d763cef2SBarry Smith .seealso: TSGetOptionsPrefix() 1889d763cef2SBarry Smith 1890d763cef2SBarry Smith @*/ 189163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[]) 1892d763cef2SBarry Smith { 1893dfbe8321SBarry Smith PetscErrorCode ierr; 1894d763cef2SBarry Smith 1895d763cef2SBarry Smith PetscFunctionBegin; 18964482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1897d763cef2SBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1898d763cef2SBarry Smith switch(ts->problem_type) { 1899d763cef2SBarry Smith case TS_NONLINEAR: 19001ac94b3bSBarry Smith if (ts->snes) { 1901d763cef2SBarry Smith ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 19021ac94b3bSBarry Smith } 1903d763cef2SBarry Smith break; 1904d763cef2SBarry Smith case TS_LINEAR: 19051ac94b3bSBarry Smith if (ts->ksp) { 190694b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 19071ac94b3bSBarry Smith } 1908d763cef2SBarry Smith break; 1909d763cef2SBarry Smith } 1910d763cef2SBarry Smith PetscFunctionReturn(0); 1911d763cef2SBarry Smith } 1912d763cef2SBarry Smith 19134a2ae208SSatish Balay #undef __FUNCT__ 19144a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix" 1915d763cef2SBarry Smith /*@C 1916d763cef2SBarry Smith TSGetOptionsPrefix - Sets the prefix used for searching for all 1917d763cef2SBarry Smith TS options in the database. 1918d763cef2SBarry Smith 1919d763cef2SBarry Smith Not Collective 1920d763cef2SBarry Smith 1921d763cef2SBarry Smith Input Parameter: 1922d763cef2SBarry Smith . ts - The TS context 1923d763cef2SBarry Smith 1924d763cef2SBarry Smith Output Parameter: 1925d763cef2SBarry Smith . prefix - A pointer to the prefix string used 1926d763cef2SBarry Smith 1927d763cef2SBarry Smith Notes: On the fortran side, the user should pass in a string 'prifix' of 1928d763cef2SBarry Smith sufficient length to hold the prefix. 1929d763cef2SBarry Smith 1930d763cef2SBarry Smith Level: intermediate 1931d763cef2SBarry Smith 1932d763cef2SBarry Smith .keywords: TS, get, options, prefix, database 1933d763cef2SBarry Smith 1934d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix() 1935d763cef2SBarry Smith @*/ 1936e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[]) 1937d763cef2SBarry Smith { 1938dfbe8321SBarry Smith PetscErrorCode ierr; 1939d763cef2SBarry Smith 1940d763cef2SBarry Smith PetscFunctionBegin; 19414482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 19424482741eSBarry Smith PetscValidPointer(prefix,2); 1943d763cef2SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1944d763cef2SBarry Smith PetscFunctionReturn(0); 1945d763cef2SBarry Smith } 1946d763cef2SBarry Smith 19474a2ae208SSatish Balay #undef __FUNCT__ 19484a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian" 1949d763cef2SBarry Smith /*@C 1950d763cef2SBarry Smith TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1951d763cef2SBarry Smith 1952d763cef2SBarry Smith Not Collective, but parallel objects are returned if TS is parallel 1953d763cef2SBarry Smith 1954d763cef2SBarry Smith Input Parameter: 1955d763cef2SBarry Smith . ts - The TS context obtained from TSCreate() 1956d763cef2SBarry Smith 1957d763cef2SBarry Smith Output Parameters: 1958d763cef2SBarry Smith + J - The Jacobian J of F, where U_t = F(U,t) 1959d763cef2SBarry Smith . M - The preconditioner matrix, usually the same as J 1960d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine 1961d763cef2SBarry Smith 1962d763cef2SBarry Smith Notes: You can pass in PETSC_NULL for any return argument you do not need. 1963d763cef2SBarry Smith 1964d763cef2SBarry Smith Level: intermediate 1965d763cef2SBarry Smith 196626d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 1967d763cef2SBarry Smith 1968d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian 1969d763cef2SBarry Smith @*/ 197063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 1971d763cef2SBarry Smith { 1972d763cef2SBarry Smith PetscFunctionBegin; 197326d46c62SHong Zhang if (J) *J = ts->Arhs; 197426d46c62SHong Zhang if (M) *M = ts->B; 197526d46c62SHong Zhang if (ctx) *ctx = ts->jacP; 1976d763cef2SBarry Smith PetscFunctionReturn(0); 1977d763cef2SBarry Smith } 1978d763cef2SBarry Smith 19791713a123SBarry Smith #undef __FUNCT__ 1980a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution" 19811713a123SBarry Smith /*@C 1982a6570f20SBarry Smith TSMonitorSolution - Monitors progress of the TS solvers by calling 19831713a123SBarry Smith VecView() for the solution at each timestep 19841713a123SBarry Smith 19851713a123SBarry Smith Collective on TS 19861713a123SBarry Smith 19871713a123SBarry Smith Input Parameters: 19881713a123SBarry Smith + ts - the TS context 19891713a123SBarry Smith . step - current time-step 1990142b95e3SSatish Balay . ptime - current time 19911713a123SBarry Smith - dummy - either a viewer or PETSC_NULL 19921713a123SBarry Smith 19931713a123SBarry Smith Level: intermediate 19941713a123SBarry Smith 19951713a123SBarry Smith .keywords: TS, vector, monitor, view 19961713a123SBarry Smith 1997a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 19981713a123SBarry Smith @*/ 1999a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 20001713a123SBarry Smith { 2001dfbe8321SBarry Smith PetscErrorCode ierr; 20021713a123SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 20031713a123SBarry Smith 20041713a123SBarry Smith PetscFunctionBegin; 2005a34d58ebSBarry Smith if (!dummy) { 20067adad957SLisandro Dalcin viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 20071713a123SBarry Smith } 20081713a123SBarry Smith ierr = VecView(x,viewer);CHKERRQ(ierr); 20091713a123SBarry Smith PetscFunctionReturn(0); 20101713a123SBarry Smith } 20111713a123SBarry Smith 20121713a123SBarry Smith 20131713a123SBarry Smith 2014