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 if (ts->ops->rhsfunction) { 363316643e7SJed Brown PetscStackPush("TS user right-hand-side function"); 364316643e7SJed Brown ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr); 365316643e7SJed Brown PetscStackPop; 366316643e7SJed Brown } else { 367316643e7SJed Brown if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 368316643e7SJed Brown MatStructure flg; 369316643e7SJed Brown PetscStackPush("TS user right-hand-side matrix function"); 370316643e7SJed Brown ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 371316643e7SJed Brown PetscStackPop; 372316643e7SJed Brown } 373316643e7SJed Brown ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr); 374316643e7SJed Brown } 375316643e7SJed Brown 376316643e7SJed Brown /* Convert to implicit form */ 377ace68cafSJed Brown ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 378316643e7SJed Brown } 379316643e7SJed Brown ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 380316643e7SJed Brown PetscFunctionReturn(0); 381316643e7SJed Brown } 382316643e7SJed Brown 383316643e7SJed Brown #undef __FUNCT__ 384316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian" 385316643e7SJed Brown /*@ 386316643e7SJed Brown TSComputeIJacobian - Evaluates the Jacobian of the DAE 387316643e7SJed Brown 388316643e7SJed Brown Collective on TS and Vec 389316643e7SJed Brown 390316643e7SJed Brown Input 391316643e7SJed Brown Input Parameters: 392316643e7SJed Brown + ts - the TS context 393316643e7SJed Brown . t - current timestep 394316643e7SJed Brown . X - state vector 395316643e7SJed Brown . Xdot - time derivative of state vector 396316643e7SJed Brown - shift - shift to apply, see note below 397316643e7SJed Brown 398316643e7SJed Brown Output Parameters: 399316643e7SJed Brown + A - Jacobian matrix 400316643e7SJed Brown . B - optional preconditioning matrix 401316643e7SJed Brown - flag - flag indicating matrix structure 402316643e7SJed Brown 403316643e7SJed Brown Notes: 404316643e7SJed Brown If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 405316643e7SJed Brown 406316643e7SJed Brown dF/dX + shift*dF/dXdot 407316643e7SJed Brown 408316643e7SJed Brown Most users should not need to explicitly call this routine, as it 409316643e7SJed Brown is used internally within the nonlinear solvers. 410316643e7SJed Brown 411316643e7SJed Brown TSComputeIJacobian() is valid only for TS_NONLINEAR 412316643e7SJed Brown 413316643e7SJed Brown Level: developer 414316643e7SJed Brown 415316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix 416316643e7SJed Brown 417316643e7SJed Brown .seealso: TSSetIJacobian() 41863495f91SJed Brown @*/ 419316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg) 420316643e7SJed Brown { 421316643e7SJed Brown PetscErrorCode ierr; 422316643e7SJed Brown 423316643e7SJed Brown PetscFunctionBegin; 424316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 425316643e7SJed Brown PetscValidHeaderSpecific(X,VEC_COOKIE,3); 426316643e7SJed Brown PetscValidHeaderSpecific(Xdot,VEC_COOKIE,4); 427316643e7SJed Brown PetscValidPointer(A,6); 428316643e7SJed Brown PetscValidHeaderSpecific(*A,MAT_COOKIE,6); 429316643e7SJed Brown PetscValidPointer(B,7); 430316643e7SJed Brown PetscValidHeaderSpecific(*B,MAT_COOKIE,7); 431316643e7SJed Brown PetscValidPointer(flg,8); 432316643e7SJed Brown 433316643e7SJed Brown ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 434316643e7SJed Brown if (ts->ops->ijacobian) { 435316643e7SJed Brown PetscStackPush("TS user implicit Jacobian"); 436316643e7SJed Brown ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 437316643e7SJed Brown PetscStackPop; 438316643e7SJed Brown } else { 439316643e7SJed Brown if (ts->ops->rhsjacobian) { 440316643e7SJed Brown PetscStackPush("TS user right-hand-side Jacobian"); 441316643e7SJed Brown ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 442316643e7SJed Brown PetscStackPop; 443316643e7SJed Brown } else { 444316643e7SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 445316643e7SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 446316643e7SJed Brown if (*A != *B) { 447316643e7SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 448316643e7SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 449316643e7SJed Brown } 450316643e7SJed Brown } 451316643e7SJed Brown 452316643e7SJed Brown /* Convert to implicit form */ 453316643e7SJed Brown /* inefficient because these operations will normally traverse all matrix elements twice */ 454316643e7SJed Brown ierr = MatScale(*A,-1);CHKERRQ(ierr); 455316643e7SJed Brown ierr = MatShift(*A,shift);CHKERRQ(ierr); 456316643e7SJed Brown if (*A != *B) { 457316643e7SJed Brown ierr = MatScale(*B,-1);CHKERRQ(ierr); 458316643e7SJed Brown ierr = MatShift(*B,shift);CHKERRQ(ierr); 459316643e7SJed Brown } 460316643e7SJed Brown } 461316643e7SJed Brown ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 462316643e7SJed Brown PetscFunctionReturn(0); 463316643e7SJed Brown } 464316643e7SJed Brown 465316643e7SJed Brown #undef __FUNCT__ 4664a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction" 467d763cef2SBarry Smith /*@C 468d763cef2SBarry Smith TSSetRHSFunction - Sets the routine for evaluating the function, 469d763cef2SBarry Smith F(t,u), where U_t = F(t,u). 470d763cef2SBarry Smith 471d763cef2SBarry Smith Collective on TS 472d763cef2SBarry Smith 473d763cef2SBarry Smith Input Parameters: 474d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 475d763cef2SBarry Smith . f - routine for evaluating the right-hand-side function 476d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 477d763cef2SBarry Smith function evaluation routine (may be PETSC_NULL) 478d763cef2SBarry Smith 479d763cef2SBarry Smith Calling sequence of func: 48087828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 481d763cef2SBarry Smith 482d763cef2SBarry Smith + t - current timestep 483d763cef2SBarry Smith . u - input vector 484d763cef2SBarry Smith . F - function vector 485d763cef2SBarry Smith - ctx - [optional] user-defined function context 486d763cef2SBarry Smith 487d763cef2SBarry Smith Important: 48876f2fa84SHong Zhang The user MUST call either this routine or TSSetMatrices(). 489d763cef2SBarry Smith 490d763cef2SBarry Smith Level: beginner 491d763cef2SBarry Smith 492d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function 493d763cef2SBarry Smith 49476f2fa84SHong Zhang .seealso: TSSetMatrices() 495d763cef2SBarry Smith @*/ 49663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 497d763cef2SBarry Smith { 498d763cef2SBarry Smith PetscFunctionBegin; 499d763cef2SBarry Smith 5004482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 501d763cef2SBarry Smith if (ts->problem_type == TS_LINEAR) { 50229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem"); 503d763cef2SBarry Smith } 504000e7ae3SMatthew Knepley ts->ops->rhsfunction = f; 505d763cef2SBarry Smith ts->funP = ctx; 506d763cef2SBarry Smith PetscFunctionReturn(0); 507d763cef2SBarry Smith } 508d763cef2SBarry Smith 5094a2ae208SSatish Balay #undef __FUNCT__ 51095f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices" 51195f0b562SHong Zhang /*@C 51295f0b562SHong Zhang TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs, 51395f0b562SHong Zhang where Alhs(t) U_t = Arhs(t) U. 51495f0b562SHong Zhang 51595f0b562SHong Zhang Collective on TS 51695f0b562SHong Zhang 51795f0b562SHong Zhang Input Parameters: 51895f0b562SHong Zhang + ts - the TS context obtained from TSCreate() 51995f0b562SHong Zhang . Arhs - matrix 52095f0b562SHong Zhang . frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 52195f0b562SHong Zhang if Arhs is not a function of t. 52295f0b562SHong Zhang . Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix. 52395f0b562SHong Zhang . flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 52495f0b562SHong Zhang if Alhs is not a function of t. 52595f0b562SHong Zhang . flag - flag indicating information about the matrix structure of Arhs and Alhs. 52695f0b562SHong Zhang The available options are 52795f0b562SHong Zhang SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs 52895f0b562SHong Zhang DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs 52995f0b562SHong Zhang - ctx - [optional] user-defined context for private data for the 53095f0b562SHong Zhang matrix evaluation routine (may be PETSC_NULL) 53195f0b562SHong Zhang 53295f0b562SHong Zhang Calling sequence of func: 53395f0b562SHong Zhang $ func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx); 53495f0b562SHong Zhang 53595f0b562SHong Zhang + t - current timestep 53695f0b562SHong Zhang . A - matrix A, where U_t = A(t) U 53795f0b562SHong Zhang . B - preconditioner matrix, usually the same as A 53895f0b562SHong Zhang . flag - flag indicating information about the preconditioner matrix 53995f0b562SHong Zhang structure (same as flag in KSPSetOperators()) 54095f0b562SHong Zhang - ctx - [optional] user-defined context for matrix evaluation routine 54195f0b562SHong Zhang 54295f0b562SHong Zhang Notes: 54395f0b562SHong Zhang The routine func() takes Mat* as the matrix arguments rather than Mat. 54495f0b562SHong Zhang This allows the matrix evaluation routine to replace Arhs or Alhs with a 54595f0b562SHong Zhang completely new new matrix structure (not just different matrix elements) 54695f0b562SHong Zhang when appropriate, for instance, if the nonzero structure is changing 54795f0b562SHong Zhang throughout the global iterations. 54895f0b562SHong Zhang 54995f0b562SHong Zhang Important: 55095f0b562SHong Zhang The user MUST call either this routine or TSSetRHSFunction(). 55195f0b562SHong Zhang 55295f0b562SHong Zhang Level: beginner 55395f0b562SHong Zhang 55495f0b562SHong Zhang .keywords: TS, timestep, set, matrix 55595f0b562SHong Zhang 55695f0b562SHong Zhang .seealso: TSSetRHSFunction() 55795f0b562SHong Zhang @*/ 55895f0b562SHong 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) 55995f0b562SHong Zhang { 56095f0b562SHong Zhang PetscFunctionBegin; 56195f0b562SHong Zhang PetscValidHeaderSpecific(ts,TS_COOKIE,1); 56292af4f6aSHong Zhang if (Arhs){ 56395f0b562SHong Zhang PetscValidHeaderSpecific(Arhs,MAT_COOKIE,2); 56495f0b562SHong Zhang PetscCheckSameComm(ts,1,Arhs,2); 56595f0b562SHong Zhang ts->Arhs = Arhs; 56692af4f6aSHong Zhang ts->ops->rhsmatrix = frhs; 56792af4f6aSHong Zhang } 56892af4f6aSHong Zhang if (Alhs){ 56992af4f6aSHong Zhang PetscValidHeaderSpecific(Alhs,MAT_COOKIE,4); 57092af4f6aSHong Zhang PetscCheckSameComm(ts,1,Alhs,4); 57195f0b562SHong Zhang ts->Alhs = Alhs; 57292af4f6aSHong Zhang ts->ops->lhsmatrix = flhs; 57392af4f6aSHong Zhang } 57492af4f6aSHong Zhang 57592af4f6aSHong Zhang ts->jacP = ctx; 57695f0b562SHong Zhang ts->matflg = flag; 57795f0b562SHong Zhang PetscFunctionReturn(0); 57895f0b562SHong Zhang } 579d763cef2SBarry Smith 580aa644b49SHong Zhang #undef __FUNCT__ 581cda39b92SHong Zhang #define __FUNCT__ "TSGetMatrices" 582cda39b92SHong Zhang /*@C 583cda39b92SHong Zhang TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep, 584cda39b92SHong Zhang where Alhs(t) U_t = Arhs(t) U. 585cda39b92SHong Zhang 586cda39b92SHong Zhang Not Collective, but parallel objects are returned if TS is parallel 587cda39b92SHong Zhang 588cda39b92SHong Zhang Input Parameter: 589cda39b92SHong Zhang . ts - The TS context obtained from TSCreate() 590cda39b92SHong Zhang 591cda39b92SHong Zhang Output Parameters: 592cda39b92SHong Zhang + Arhs - The right-hand side matrix 593cda39b92SHong Zhang . Alhs - The left-hand side matrix 594cda39b92SHong Zhang - ctx - User-defined context for matrix evaluation routine 595cda39b92SHong Zhang 596cda39b92SHong Zhang Notes: You can pass in PETSC_NULL for any return argument you do not need. 597cda39b92SHong Zhang 598cda39b92SHong Zhang Level: intermediate 599cda39b92SHong Zhang 600cda39b92SHong Zhang .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 601cda39b92SHong Zhang 602cda39b92SHong Zhang .keywords: TS, timestep, get, matrix 603cda39b92SHong Zhang 604cda39b92SHong Zhang @*/ 605cda39b92SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx) 606cda39b92SHong Zhang { 607cda39b92SHong Zhang PetscFunctionBegin; 608cda39b92SHong Zhang PetscValidHeaderSpecific(ts,TS_COOKIE,1); 609cda39b92SHong Zhang if (Arhs) *Arhs = ts->Arhs; 610cda39b92SHong Zhang if (Alhs) *Alhs = ts->Alhs; 611cda39b92SHong Zhang if (ctx) *ctx = ts->jacP; 612cda39b92SHong Zhang PetscFunctionReturn(0); 613cda39b92SHong Zhang } 614cda39b92SHong Zhang 615cda39b92SHong Zhang #undef __FUNCT__ 6164a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian" 617d763cef2SBarry Smith /*@C 618d763cef2SBarry Smith TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 619d763cef2SBarry Smith where U_t = F(U,t), as well as the location to store the matrix. 62076f2fa84SHong Zhang Use TSSetMatrices() for linear problems. 621d763cef2SBarry Smith 622d763cef2SBarry Smith Collective on TS 623d763cef2SBarry Smith 624d763cef2SBarry Smith Input Parameters: 625d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 626d763cef2SBarry Smith . A - Jacobian matrix 627d763cef2SBarry Smith . B - preconditioner matrix (usually same as A) 628d763cef2SBarry Smith . f - the Jacobian evaluation routine 629d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 630d763cef2SBarry Smith Jacobian evaluation routine (may be PETSC_NULL) 631d763cef2SBarry Smith 632d763cef2SBarry Smith Calling sequence of func: 63387828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 634d763cef2SBarry Smith 635d763cef2SBarry Smith + t - current timestep 636d763cef2SBarry Smith . u - input vector 637d763cef2SBarry Smith . A - matrix A, where U_t = A(t)u 638d763cef2SBarry Smith . B - preconditioner matrix, usually the same as A 639d763cef2SBarry Smith . flag - flag indicating information about the preconditioner matrix 64094b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 641d763cef2SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 642d763cef2SBarry Smith 643d763cef2SBarry Smith Notes: 64494b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 645d763cef2SBarry Smith output parameter in the routine func(). Be sure to read this information! 646d763cef2SBarry Smith 647d763cef2SBarry Smith The routine func() takes Mat * as the matrix arguments rather than Mat. 648d763cef2SBarry Smith This allows the matrix evaluation routine to replace A and/or B with a 64956335db2SHong Zhang completely new matrix structure (not just different matrix elements) 650d763cef2SBarry Smith when appropriate, for instance, if the nonzero structure is changing 651d763cef2SBarry Smith throughout the global iterations. 652d763cef2SBarry Smith 653d763cef2SBarry Smith Level: beginner 654d763cef2SBarry Smith 655d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian 656d763cef2SBarry Smith 657d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(), 65876f2fa84SHong Zhang SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 659d763cef2SBarry Smith 660d763cef2SBarry Smith @*/ 66163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 662d763cef2SBarry Smith { 663d763cef2SBarry Smith PetscFunctionBegin; 6644482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 6654482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,2); 6664482741eSBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE,3); 667c9780b6fSBarry Smith PetscCheckSameComm(ts,1,A,2); 668c9780b6fSBarry Smith PetscCheckSameComm(ts,1,B,3); 669d763cef2SBarry Smith if (ts->problem_type != TS_NONLINEAR) { 67076f2fa84SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()"); 671d763cef2SBarry Smith } 672d763cef2SBarry Smith 673000e7ae3SMatthew Knepley ts->ops->rhsjacobian = f; 674d763cef2SBarry Smith ts->jacP = ctx; 6758beb423aSHong Zhang ts->Arhs = A; 676d763cef2SBarry Smith ts->B = B; 677d763cef2SBarry Smith PetscFunctionReturn(0); 678d763cef2SBarry Smith } 679d763cef2SBarry Smith 680316643e7SJed Brown 681316643e7SJed Brown #undef __FUNCT__ 682316643e7SJed Brown #define __FUNCT__ "TSSetIFunction" 683316643e7SJed Brown /*@C 684316643e7SJed Brown TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 685316643e7SJed Brown 686316643e7SJed Brown Collective on TS 687316643e7SJed Brown 688316643e7SJed Brown Input Parameters: 689316643e7SJed Brown + ts - the TS context obtained from TSCreate() 690316643e7SJed Brown . f - the function evaluation routine 691316643e7SJed Brown - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 692316643e7SJed Brown 693316643e7SJed Brown Calling sequence of f: 694316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 695316643e7SJed Brown 696316643e7SJed Brown + t - time at step/stage being solved 697316643e7SJed Brown . u - state vector 698316643e7SJed Brown . u_t - time derivative of state vector 699316643e7SJed Brown . F - function vector 700316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 701316643e7SJed Brown 702316643e7SJed Brown Important: 703316643e7SJed Brown The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 704316643e7SJed Brown 705316643e7SJed Brown Level: beginner 706316643e7SJed Brown 707316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian 708316643e7SJed Brown 709316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 710316643e7SJed Brown @*/ 711316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS ts,TSIFunction f,void *ctx) 712316643e7SJed Brown { 713316643e7SJed Brown 714316643e7SJed Brown PetscFunctionBegin; 715316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 716316643e7SJed Brown ts->ops->ifunction = f; 717316643e7SJed Brown ts->funP = ctx; 718316643e7SJed Brown PetscFunctionReturn(0); 719316643e7SJed Brown } 720316643e7SJed Brown 721316643e7SJed Brown #undef __FUNCT__ 722316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian" 723316643e7SJed Brown /*@C 724316643e7SJed Brown TSSetIJacobian - Set the function to compute the Jacobian of 725316643e7SJed Brown G(U) = F(t,U,U0+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 726316643e7SJed Brown 727316643e7SJed Brown Collective on TS 728316643e7SJed Brown 729316643e7SJed Brown Input Parameters: 730316643e7SJed Brown + ts - the TS context obtained from TSCreate() 731316643e7SJed Brown . A - Jacobian matrix 732316643e7SJed Brown . B - preconditioning matrix for A (may be same as A) 733316643e7SJed Brown . f - the Jacobian evaluation routine 734316643e7SJed Brown - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 735316643e7SJed Brown 736316643e7SJed Brown Calling sequence of f: 737316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 738316643e7SJed Brown 739316643e7SJed Brown + t - time at step/stage being solved 740316643e7SJed Brown . u - state vector 741316643e7SJed Brown . u_t - time derivative of state vector 742316643e7SJed Brown . a - shift 743316643e7SJed Brown . A - Jacobian of G(U) = F(t,U,U0+a*U), equivalent to dF/dU + a*dF/dU_t 744316643e7SJed Brown . B - preconditioning matrix for A, may be same as A 745316643e7SJed Brown . flag - flag indicating information about the preconditioner matrix 746316643e7SJed Brown structure (same as flag in KSPSetOperators()) 747316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 748316643e7SJed Brown 749316643e7SJed Brown Notes: 750316643e7SJed Brown The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 751316643e7SJed Brown 752316643e7SJed Brown Level: beginner 753316643e7SJed Brown 754316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian 755316643e7SJed Brown 756316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian() 757316643e7SJed Brown 758316643e7SJed Brown @*/ 759316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 760316643e7SJed Brown { 761316643e7SJed Brown PetscErrorCode ierr; 762316643e7SJed Brown 763316643e7SJed Brown PetscFunctionBegin; 764316643e7SJed Brown PetscValidHeaderSpecific(ts,TS_COOKIE,1); 765316643e7SJed Brown if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 766316643e7SJed Brown if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 767316643e7SJed Brown if (A) PetscCheckSameComm(ts,1,A,2); 768316643e7SJed Brown if (B) PetscCheckSameComm(ts,1,B,3); 769316643e7SJed Brown if (f) ts->ops->ijacobian = f; 770316643e7SJed Brown if (ctx) ts->jacP = ctx; 771316643e7SJed Brown if (A) { 772316643e7SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 773316643e7SJed Brown if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 774316643e7SJed Brown ts->A = A; 775316643e7SJed Brown } 776*dc3f620dSJed Brown #if 0 777*dc3f620dSJed Brown /* The sane and consistent alternative */ 778316643e7SJed Brown if (B) { 779316643e7SJed Brown ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 780316643e7SJed Brown if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);} 781316643e7SJed Brown ts->B = B; 782316643e7SJed Brown } 783*dc3f620dSJed Brown #else 784*dc3f620dSJed Brown /* Don't reference B because TSDestroy() doesn't destroy it. These ownership semantics are awkward and inconsistent. */ 785*dc3f620dSJed Brown if (B) ts->B = B; 786*dc3f620dSJed Brown #endif 787316643e7SJed Brown PetscFunctionReturn(0); 788316643e7SJed Brown } 789316643e7SJed Brown 7904a2ae208SSatish Balay #undef __FUNCT__ 7914a2ae208SSatish Balay #define __FUNCT__ "TSView" 7927e2c5f70SBarry Smith /*@C 793d763cef2SBarry Smith TSView - Prints the TS data structure. 794d763cef2SBarry Smith 7954c49b128SBarry Smith Collective on TS 796d763cef2SBarry Smith 797d763cef2SBarry Smith Input Parameters: 798d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 799d763cef2SBarry Smith - viewer - visualization context 800d763cef2SBarry Smith 801d763cef2SBarry Smith Options Database Key: 802d763cef2SBarry Smith . -ts_view - calls TSView() at end of TSStep() 803d763cef2SBarry Smith 804d763cef2SBarry Smith Notes: 805d763cef2SBarry Smith The available visualization contexts include 806b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 807b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 808d763cef2SBarry Smith output where only the first processor opens 809d763cef2SBarry Smith the file. All other processors send their 810d763cef2SBarry Smith data to the first processor to print. 811d763cef2SBarry Smith 812d763cef2SBarry Smith The user can open an alternative visualization context with 813b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 814d763cef2SBarry Smith 815d763cef2SBarry Smith Level: beginner 816d763cef2SBarry Smith 817d763cef2SBarry Smith .keywords: TS, timestep, view 818d763cef2SBarry Smith 819b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 820d763cef2SBarry Smith @*/ 82163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer) 822d763cef2SBarry Smith { 823dfbe8321SBarry Smith PetscErrorCode ierr; 824a313700dSBarry Smith const TSType type; 82532077d6dSBarry Smith PetscTruth iascii,isstring; 826d763cef2SBarry Smith 827d763cef2SBarry Smith PetscFunctionBegin; 8284482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 8293050cee2SBarry Smith if (!viewer) { 8307adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 8313050cee2SBarry Smith } 8324482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 833c9780b6fSBarry Smith PetscCheckSameComm(ts,1,viewer,2); 834fd16b177SBarry Smith 83532077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 836b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 83732077d6dSBarry Smith if (iascii) { 838b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr); 839a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 840454a90a3SBarry Smith if (type) { 841b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 842184914b5SBarry Smith } else { 843b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 844184914b5SBarry Smith } 845000e7ae3SMatthew Knepley if (ts->ops->view) { 846b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 847000e7ae3SMatthew Knepley ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 848b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 849d763cef2SBarry Smith } 85077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 851a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 852d763cef2SBarry Smith if (ts->problem_type == TS_NONLINEAR) { 85377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 854d763cef2SBarry Smith } 85577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 8560f5bd95cSBarry Smith } else if (isstring) { 857a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 858b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 859d763cef2SBarry Smith } 860b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 86194b7f48cSBarry Smith if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);} 862d763cef2SBarry Smith if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 863b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 864d763cef2SBarry Smith PetscFunctionReturn(0); 865d763cef2SBarry Smith } 866d763cef2SBarry Smith 867d763cef2SBarry Smith 8684a2ae208SSatish Balay #undef __FUNCT__ 8694a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext" 870d763cef2SBarry Smith /*@C 871d763cef2SBarry Smith TSSetApplicationContext - Sets an optional user-defined context for 872d763cef2SBarry Smith the timesteppers. 873d763cef2SBarry Smith 874d763cef2SBarry Smith Collective on TS 875d763cef2SBarry Smith 876d763cef2SBarry Smith Input Parameters: 877d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 878d763cef2SBarry Smith - usrP - optional user context 879d763cef2SBarry Smith 880d763cef2SBarry Smith Level: intermediate 881d763cef2SBarry Smith 882d763cef2SBarry Smith .keywords: TS, timestep, set, application, context 883d763cef2SBarry Smith 884d763cef2SBarry Smith .seealso: TSGetApplicationContext() 885d763cef2SBarry Smith @*/ 88663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP) 887d763cef2SBarry Smith { 888d763cef2SBarry Smith PetscFunctionBegin; 8894482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 890d763cef2SBarry Smith ts->user = usrP; 891d763cef2SBarry Smith PetscFunctionReturn(0); 892d763cef2SBarry Smith } 893d763cef2SBarry Smith 8944a2ae208SSatish Balay #undef __FUNCT__ 8954a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext" 896d763cef2SBarry Smith /*@C 897d763cef2SBarry Smith TSGetApplicationContext - Gets the user-defined context for the 898d763cef2SBarry Smith timestepper. 899d763cef2SBarry Smith 900d763cef2SBarry Smith Not Collective 901d763cef2SBarry Smith 902d763cef2SBarry Smith Input Parameter: 903d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 904d763cef2SBarry Smith 905d763cef2SBarry Smith Output Parameter: 906d763cef2SBarry Smith . usrP - user context 907d763cef2SBarry Smith 908d763cef2SBarry Smith Level: intermediate 909d763cef2SBarry Smith 910d763cef2SBarry Smith .keywords: TS, timestep, get, application, context 911d763cef2SBarry Smith 912d763cef2SBarry Smith .seealso: TSSetApplicationContext() 913d763cef2SBarry Smith @*/ 91463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP) 915d763cef2SBarry Smith { 916d763cef2SBarry Smith PetscFunctionBegin; 9174482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 918d763cef2SBarry Smith *usrP = ts->user; 919d763cef2SBarry Smith PetscFunctionReturn(0); 920d763cef2SBarry Smith } 921d763cef2SBarry Smith 9224a2ae208SSatish Balay #undef __FUNCT__ 9234a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber" 924d763cef2SBarry Smith /*@ 925d763cef2SBarry Smith TSGetTimeStepNumber - Gets the current number of timesteps. 926d763cef2SBarry Smith 927d763cef2SBarry Smith Not Collective 928d763cef2SBarry Smith 929d763cef2SBarry Smith Input Parameter: 930d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 931d763cef2SBarry Smith 932d763cef2SBarry Smith Output Parameter: 933d763cef2SBarry Smith . iter - number steps so far 934d763cef2SBarry Smith 935d763cef2SBarry Smith Level: intermediate 936d763cef2SBarry Smith 937d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number 938d763cef2SBarry Smith @*/ 93963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter) 940d763cef2SBarry Smith { 941d763cef2SBarry Smith PetscFunctionBegin; 9424482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 9434482741eSBarry Smith PetscValidIntPointer(iter,2); 944d763cef2SBarry Smith *iter = ts->steps; 945d763cef2SBarry Smith PetscFunctionReturn(0); 946d763cef2SBarry Smith } 947d763cef2SBarry Smith 9484a2ae208SSatish Balay #undef __FUNCT__ 9494a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep" 950d763cef2SBarry Smith /*@ 951d763cef2SBarry Smith TSSetInitialTimeStep - Sets the initial timestep to be used, 952d763cef2SBarry Smith as well as the initial time. 953d763cef2SBarry Smith 954d763cef2SBarry Smith Collective on TS 955d763cef2SBarry Smith 956d763cef2SBarry Smith Input Parameters: 957d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 958d763cef2SBarry Smith . initial_time - the initial time 959d763cef2SBarry Smith - time_step - the size of the timestep 960d763cef2SBarry Smith 961d763cef2SBarry Smith Level: intermediate 962d763cef2SBarry Smith 963d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep() 964d763cef2SBarry Smith 965d763cef2SBarry Smith .keywords: TS, set, initial, timestep 966d763cef2SBarry Smith @*/ 96763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 968d763cef2SBarry Smith { 969d763cef2SBarry Smith PetscFunctionBegin; 9704482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 971d763cef2SBarry Smith ts->time_step = time_step; 972d763cef2SBarry Smith ts->initial_time_step = time_step; 973d763cef2SBarry Smith ts->ptime = initial_time; 974d763cef2SBarry Smith PetscFunctionReturn(0); 975d763cef2SBarry Smith } 976d763cef2SBarry Smith 9774a2ae208SSatish Balay #undef __FUNCT__ 9784a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep" 979d763cef2SBarry Smith /*@ 980d763cef2SBarry Smith TSSetTimeStep - Allows one to reset the timestep at any time, 981d763cef2SBarry Smith useful for simple pseudo-timestepping codes. 982d763cef2SBarry Smith 983d763cef2SBarry Smith Collective on TS 984d763cef2SBarry Smith 985d763cef2SBarry Smith Input Parameters: 986d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 987d763cef2SBarry Smith - time_step - the size of the timestep 988d763cef2SBarry Smith 989d763cef2SBarry Smith Level: intermediate 990d763cef2SBarry Smith 991d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 992d763cef2SBarry Smith 993d763cef2SBarry Smith .keywords: TS, set, timestep 994d763cef2SBarry Smith @*/ 99563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step) 996d763cef2SBarry Smith { 997d763cef2SBarry Smith PetscFunctionBegin; 9984482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 999d763cef2SBarry Smith ts->time_step = time_step; 1000d763cef2SBarry Smith PetscFunctionReturn(0); 1001d763cef2SBarry Smith } 1002d763cef2SBarry Smith 10034a2ae208SSatish Balay #undef __FUNCT__ 10044a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep" 1005d763cef2SBarry Smith /*@ 1006d763cef2SBarry Smith TSGetTimeStep - Gets the current timestep size. 1007d763cef2SBarry Smith 1008d763cef2SBarry Smith Not Collective 1009d763cef2SBarry Smith 1010d763cef2SBarry Smith Input Parameter: 1011d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1012d763cef2SBarry Smith 1013d763cef2SBarry Smith Output Parameter: 1014d763cef2SBarry Smith . dt - the current timestep size 1015d763cef2SBarry Smith 1016d763cef2SBarry Smith Level: intermediate 1017d763cef2SBarry Smith 1018d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1019d763cef2SBarry Smith 1020d763cef2SBarry Smith .keywords: TS, get, timestep 1021d763cef2SBarry Smith @*/ 102263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt) 1023d763cef2SBarry Smith { 1024d763cef2SBarry Smith PetscFunctionBegin; 10254482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 10264482741eSBarry Smith PetscValidDoublePointer(dt,2); 1027d763cef2SBarry Smith *dt = ts->time_step; 1028d763cef2SBarry Smith PetscFunctionReturn(0); 1029d763cef2SBarry Smith } 1030d763cef2SBarry Smith 10314a2ae208SSatish Balay #undef __FUNCT__ 10324a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution" 1033d8e5e3e6SSatish Balay /*@ 1034d763cef2SBarry Smith TSGetSolution - Returns the solution at the present timestep. It 1035d763cef2SBarry Smith is valid to call this routine inside the function that you are evaluating 1036d763cef2SBarry Smith in order to move to the new timestep. This vector not changed until 1037d763cef2SBarry Smith the solution at the next timestep has been calculated. 1038d763cef2SBarry Smith 1039d763cef2SBarry Smith Not Collective, but Vec returned is parallel if TS is parallel 1040d763cef2SBarry Smith 1041d763cef2SBarry Smith Input Parameter: 1042d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1043d763cef2SBarry Smith 1044d763cef2SBarry Smith Output Parameter: 1045d763cef2SBarry Smith . v - the vector containing the solution 1046d763cef2SBarry Smith 1047d763cef2SBarry Smith Level: intermediate 1048d763cef2SBarry Smith 1049d763cef2SBarry Smith .seealso: TSGetTimeStep() 1050d763cef2SBarry Smith 1051d763cef2SBarry Smith .keywords: TS, timestep, get, solution 1052d763cef2SBarry Smith @*/ 105363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v) 1054d763cef2SBarry Smith { 1055d763cef2SBarry Smith PetscFunctionBegin; 10564482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 10574482741eSBarry Smith PetscValidPointer(v,2); 1058d763cef2SBarry Smith *v = ts->vec_sol_always; 1059d763cef2SBarry Smith PetscFunctionReturn(0); 1060d763cef2SBarry Smith } 1061d763cef2SBarry Smith 1062bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */ 10634a2ae208SSatish Balay #undef __FUNCT__ 1064bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType" 1065d8e5e3e6SSatish Balay /*@ 1066bdad233fSMatthew Knepley TSSetProblemType - Sets the type of problem to be solved. 1067d763cef2SBarry Smith 1068bdad233fSMatthew Knepley Not collective 1069d763cef2SBarry Smith 1070bdad233fSMatthew Knepley Input Parameters: 1071bdad233fSMatthew Knepley + ts - The TS 1072bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1073d763cef2SBarry Smith .vb 1074d763cef2SBarry Smith U_t = A U 1075d763cef2SBarry Smith U_t = A(t) U 1076d763cef2SBarry Smith U_t = F(t,U) 1077d763cef2SBarry Smith .ve 1078d763cef2SBarry Smith 1079d763cef2SBarry Smith Level: beginner 1080d763cef2SBarry Smith 1081bdad233fSMatthew Knepley .keywords: TS, problem type 1082bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1083d763cef2SBarry Smith @*/ 108463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type) 1085a7cc72afSBarry Smith { 1086d763cef2SBarry Smith PetscFunctionBegin; 10874482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1088bdad233fSMatthew Knepley ts->problem_type = type; 1089d763cef2SBarry Smith PetscFunctionReturn(0); 1090d763cef2SBarry Smith } 1091d763cef2SBarry Smith 1092bdad233fSMatthew Knepley #undef __FUNCT__ 1093bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType" 1094bdad233fSMatthew Knepley /*@C 1095bdad233fSMatthew Knepley TSGetProblemType - Gets the type of problem to be solved. 1096bdad233fSMatthew Knepley 1097bdad233fSMatthew Knepley Not collective 1098bdad233fSMatthew Knepley 1099bdad233fSMatthew Knepley Input Parameter: 1100bdad233fSMatthew Knepley . ts - The TS 1101bdad233fSMatthew Knepley 1102bdad233fSMatthew Knepley Output Parameter: 1103bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1104bdad233fSMatthew Knepley .vb 1105bdad233fSMatthew Knepley U_t = A U 1106bdad233fSMatthew Knepley U_t = A(t) U 1107bdad233fSMatthew Knepley U_t = F(t,U) 1108bdad233fSMatthew Knepley .ve 1109bdad233fSMatthew Knepley 1110bdad233fSMatthew Knepley Level: beginner 1111bdad233fSMatthew Knepley 1112bdad233fSMatthew Knepley .keywords: TS, problem type 1113bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1114bdad233fSMatthew Knepley @*/ 111563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type) 1116a7cc72afSBarry Smith { 1117bdad233fSMatthew Knepley PetscFunctionBegin; 11184482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 11194482741eSBarry Smith PetscValidIntPointer(type,2); 1120bdad233fSMatthew Knepley *type = ts->problem_type; 1121bdad233fSMatthew Knepley PetscFunctionReturn(0); 1122bdad233fSMatthew Knepley } 1123d763cef2SBarry Smith 11244a2ae208SSatish Balay #undef __FUNCT__ 11254a2ae208SSatish Balay #define __FUNCT__ "TSSetUp" 1126d763cef2SBarry Smith /*@ 1127d763cef2SBarry Smith TSSetUp - Sets up the internal data structures for the later use 1128d763cef2SBarry Smith of a timestepper. 1129d763cef2SBarry Smith 1130d763cef2SBarry Smith Collective on TS 1131d763cef2SBarry Smith 1132d763cef2SBarry Smith Input Parameter: 1133d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1134d763cef2SBarry Smith 1135d763cef2SBarry Smith Notes: 1136d763cef2SBarry Smith For basic use of the TS solvers the user need not explicitly call 1137d763cef2SBarry Smith TSSetUp(), since these actions will automatically occur during 1138d763cef2SBarry Smith the call to TSStep(). However, if one wishes to control this 1139d763cef2SBarry Smith phase separately, TSSetUp() should be called after TSCreate() 1140d763cef2SBarry Smith and optional routines of the form TSSetXXX(), but before TSStep(). 1141d763cef2SBarry Smith 1142d763cef2SBarry Smith Level: advanced 1143d763cef2SBarry Smith 1144d763cef2SBarry Smith .keywords: TS, timestep, setup 1145d763cef2SBarry Smith 1146d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy() 1147d763cef2SBarry Smith @*/ 114863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts) 1149d763cef2SBarry Smith { 1150dfbe8321SBarry Smith PetscErrorCode ierr; 1151d763cef2SBarry Smith 1152d763cef2SBarry Smith PetscFunctionBegin; 11534482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 115429bbc08cSBarry Smith if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 11557adad957SLisandro Dalcin if (!((PetscObject)ts)->type_name) { 1156d763cef2SBarry Smith ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr); 1157d763cef2SBarry Smith } 1158000e7ae3SMatthew Knepley ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1159d763cef2SBarry Smith ts->setupcalled = 1; 1160d763cef2SBarry Smith PetscFunctionReturn(0); 1161d763cef2SBarry Smith } 1162d763cef2SBarry Smith 11634a2ae208SSatish Balay #undef __FUNCT__ 11644a2ae208SSatish Balay #define __FUNCT__ "TSDestroy" 1165d8e5e3e6SSatish Balay /*@ 1166d763cef2SBarry Smith TSDestroy - Destroys the timestepper context that was created 1167d763cef2SBarry Smith with TSCreate(). 1168d763cef2SBarry Smith 1169d763cef2SBarry Smith Collective on TS 1170d763cef2SBarry Smith 1171d763cef2SBarry Smith Input Parameter: 1172d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1173d763cef2SBarry Smith 1174d763cef2SBarry Smith Level: beginner 1175d763cef2SBarry Smith 1176d763cef2SBarry Smith .keywords: TS, timestepper, destroy 1177d763cef2SBarry Smith 1178d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve() 1179d763cef2SBarry Smith @*/ 118063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts) 1181d763cef2SBarry Smith { 11826849ba73SBarry Smith PetscErrorCode ierr; 1183d763cef2SBarry Smith 1184d763cef2SBarry Smith PetscFunctionBegin; 11854482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 11867adad957SLisandro Dalcin if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0); 1187d763cef2SBarry Smith 1188be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 11890f5bd95cSBarry Smith ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 11908beb423aSHong Zhang if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr)} 119194b7f48cSBarry Smith if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);} 1192d763cef2SBarry Smith if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 11931e3347e8SBarry Smith if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);} 1194a6570f20SBarry Smith ierr = TSMonitorCancel(ts);CHKERRQ(ierr); 1195a79aaaedSSatish Balay ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1196d763cef2SBarry Smith PetscFunctionReturn(0); 1197d763cef2SBarry Smith } 1198d763cef2SBarry Smith 11994a2ae208SSatish Balay #undef __FUNCT__ 12004a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES" 1201d8e5e3e6SSatish Balay /*@ 1202d763cef2SBarry Smith TSGetSNES - Returns the SNES (nonlinear solver) associated with 1203d763cef2SBarry Smith a TS (timestepper) context. Valid only for nonlinear problems. 1204d763cef2SBarry Smith 1205d763cef2SBarry Smith Not Collective, but SNES is parallel if TS is parallel 1206d763cef2SBarry Smith 1207d763cef2SBarry Smith Input Parameter: 1208d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1209d763cef2SBarry Smith 1210d763cef2SBarry Smith Output Parameter: 1211d763cef2SBarry Smith . snes - the nonlinear solver context 1212d763cef2SBarry Smith 1213d763cef2SBarry Smith Notes: 1214d763cef2SBarry Smith The user can then directly manipulate the SNES context to set various 1215d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 121694b7f48cSBarry Smith KSP, KSP, and PC contexts as well. 1217d763cef2SBarry Smith 1218d763cef2SBarry Smith TSGetSNES() does not work for integrators that do not use SNES; in 1219d763cef2SBarry Smith this case TSGetSNES() returns PETSC_NULL in snes. 1220d763cef2SBarry Smith 1221d763cef2SBarry Smith Level: beginner 1222d763cef2SBarry Smith 1223d763cef2SBarry Smith .keywords: timestep, get, SNES 1224d763cef2SBarry Smith @*/ 122563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes) 1226d763cef2SBarry Smith { 1227d763cef2SBarry Smith PetscFunctionBegin; 12284482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 12294482741eSBarry Smith PetscValidPointer(snes,2); 1230447600ffSHong Zhang if (((PetscObject)ts)->type_name == PETSC_NULL) 1231447600ffSHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first"); 123294b7f48cSBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1233d763cef2SBarry Smith *snes = ts->snes; 1234d763cef2SBarry Smith PetscFunctionReturn(0); 1235d763cef2SBarry Smith } 1236d763cef2SBarry Smith 12374a2ae208SSatish Balay #undef __FUNCT__ 123894b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP" 1239d8e5e3e6SSatish Balay /*@ 124094b7f48cSBarry Smith TSGetKSP - Returns the KSP (linear solver) associated with 1241d763cef2SBarry Smith a TS (timestepper) context. 1242d763cef2SBarry Smith 124394b7f48cSBarry Smith Not Collective, but KSP is parallel if TS is parallel 1244d763cef2SBarry Smith 1245d763cef2SBarry Smith Input Parameter: 1246d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1247d763cef2SBarry Smith 1248d763cef2SBarry Smith Output Parameter: 124994b7f48cSBarry Smith . ksp - the nonlinear solver context 1250d763cef2SBarry Smith 1251d763cef2SBarry Smith Notes: 125294b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 1253d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 1254d763cef2SBarry Smith KSP and PC contexts as well. 1255d763cef2SBarry Smith 125694b7f48cSBarry Smith TSGetKSP() does not work for integrators that do not use KSP; 125794b7f48cSBarry Smith in this case TSGetKSP() returns PETSC_NULL in ksp. 1258d763cef2SBarry Smith 1259d763cef2SBarry Smith Level: beginner 1260d763cef2SBarry Smith 126194b7f48cSBarry Smith .keywords: timestep, get, KSP 1262d763cef2SBarry Smith @*/ 126363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp) 1264d763cef2SBarry Smith { 1265d763cef2SBarry Smith PetscFunctionBegin; 12664482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 12674482741eSBarry Smith PetscValidPointer(ksp,2); 1268988402f6SHong Zhang if (((PetscObject)ts)->type_name == PETSC_NULL) 1269988402f6SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 127029bbc08cSBarry Smith if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 127194b7f48cSBarry Smith *ksp = ts->ksp; 1272d763cef2SBarry Smith PetscFunctionReturn(0); 1273d763cef2SBarry Smith } 1274d763cef2SBarry Smith 1275d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */ 1276d763cef2SBarry Smith 12774a2ae208SSatish Balay #undef __FUNCT__ 1278adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration" 1279adb62b0dSMatthew Knepley /*@ 1280adb62b0dSMatthew Knepley TSGetDuration - Gets the maximum number of timesteps to use and 1281adb62b0dSMatthew Knepley maximum time for iteration. 1282adb62b0dSMatthew Knepley 1283adb62b0dSMatthew Knepley Collective on TS 1284adb62b0dSMatthew Knepley 1285adb62b0dSMatthew Knepley Input Parameters: 1286adb62b0dSMatthew Knepley + ts - the TS context obtained from TSCreate() 1287adb62b0dSMatthew Knepley . maxsteps - maximum number of iterations to use, or PETSC_NULL 1288adb62b0dSMatthew Knepley - maxtime - final time to iterate to, or PETSC_NULL 1289adb62b0dSMatthew Knepley 1290adb62b0dSMatthew Knepley Level: intermediate 1291adb62b0dSMatthew Knepley 1292adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time 1293adb62b0dSMatthew Knepley @*/ 129463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1295adb62b0dSMatthew Knepley { 1296adb62b0dSMatthew Knepley PetscFunctionBegin; 12974482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1298abc0a331SBarry Smith if (maxsteps) { 12994482741eSBarry Smith PetscValidIntPointer(maxsteps,2); 1300adb62b0dSMatthew Knepley *maxsteps = ts->max_steps; 1301adb62b0dSMatthew Knepley } 1302abc0a331SBarry Smith if (maxtime ) { 13034482741eSBarry Smith PetscValidScalarPointer(maxtime,3); 1304adb62b0dSMatthew Knepley *maxtime = ts->max_time; 1305adb62b0dSMatthew Knepley } 1306adb62b0dSMatthew Knepley PetscFunctionReturn(0); 1307adb62b0dSMatthew Knepley } 1308adb62b0dSMatthew Knepley 1309adb62b0dSMatthew Knepley #undef __FUNCT__ 13104a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration" 1311d763cef2SBarry Smith /*@ 1312d763cef2SBarry Smith TSSetDuration - Sets the maximum number of timesteps to use and 1313d763cef2SBarry Smith maximum time for iteration. 1314d763cef2SBarry Smith 1315d763cef2SBarry Smith Collective on TS 1316d763cef2SBarry Smith 1317d763cef2SBarry Smith Input Parameters: 1318d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1319d763cef2SBarry Smith . maxsteps - maximum number of iterations to use 1320d763cef2SBarry Smith - maxtime - final time to iterate to 1321d763cef2SBarry Smith 1322d763cef2SBarry Smith Options Database Keys: 1323d763cef2SBarry Smith . -ts_max_steps <maxsteps> - Sets maxsteps 1324d763cef2SBarry Smith . -ts_max_time <maxtime> - Sets maxtime 1325d763cef2SBarry Smith 1326d763cef2SBarry Smith Notes: 1327d763cef2SBarry Smith The default maximum number of iterations is 5000. Default time is 5.0 1328d763cef2SBarry Smith 1329d763cef2SBarry Smith Level: intermediate 1330d763cef2SBarry Smith 1331d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations 1332d763cef2SBarry Smith @*/ 133363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1334d763cef2SBarry Smith { 1335d763cef2SBarry Smith PetscFunctionBegin; 13364482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1337d763cef2SBarry Smith ts->max_steps = maxsteps; 1338d763cef2SBarry Smith ts->max_time = maxtime; 1339d763cef2SBarry Smith PetscFunctionReturn(0); 1340d763cef2SBarry Smith } 1341d763cef2SBarry Smith 13424a2ae208SSatish Balay #undef __FUNCT__ 13434a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution" 1344d763cef2SBarry Smith /*@ 1345d763cef2SBarry Smith TSSetSolution - Sets the initial solution vector 1346d763cef2SBarry Smith for use by the TS routines. 1347d763cef2SBarry Smith 1348d763cef2SBarry Smith Collective on TS and Vec 1349d763cef2SBarry Smith 1350d763cef2SBarry Smith Input Parameters: 1351d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1352d763cef2SBarry Smith - x - the solution vector 1353d763cef2SBarry Smith 1354d763cef2SBarry Smith Level: beginner 1355d763cef2SBarry Smith 1356d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions 1357d763cef2SBarry Smith @*/ 135863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x) 1359d763cef2SBarry Smith { 1360d763cef2SBarry Smith PetscFunctionBegin; 13614482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 13624482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1363d763cef2SBarry Smith ts->vec_sol = ts->vec_sol_always = x; 1364d763cef2SBarry Smith PetscFunctionReturn(0); 1365d763cef2SBarry Smith } 1366d763cef2SBarry Smith 1367e74ef692SMatthew Knepley #undef __FUNCT__ 1368e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep" 1369ac226902SBarry Smith /*@C 1370000e7ae3SMatthew Knepley TSSetPreStep - Sets the general-purpose function 1371000e7ae3SMatthew Knepley called once at the beginning of time stepping. 1372000e7ae3SMatthew Knepley 1373000e7ae3SMatthew Knepley Collective on TS 1374000e7ae3SMatthew Knepley 1375000e7ae3SMatthew Knepley Input Parameters: 1376000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1377000e7ae3SMatthew Knepley - func - The function 1378000e7ae3SMatthew Knepley 1379000e7ae3SMatthew Knepley Calling sequence of func: 1380000e7ae3SMatthew Knepley . func (TS ts); 1381000e7ae3SMatthew Knepley 1382000e7ae3SMatthew Knepley Level: intermediate 1383000e7ae3SMatthew Knepley 1384000e7ae3SMatthew Knepley .keywords: TS, timestep 1385000e7ae3SMatthew Knepley @*/ 138663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1387000e7ae3SMatthew Knepley { 1388000e7ae3SMatthew Knepley PetscFunctionBegin; 13894482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1390000e7ae3SMatthew Knepley ts->ops->prestep = func; 1391000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1392000e7ae3SMatthew Knepley } 1393000e7ae3SMatthew Knepley 1394e74ef692SMatthew Knepley #undef __FUNCT__ 1395e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep" 1396000e7ae3SMatthew Knepley /*@ 1397000e7ae3SMatthew Knepley TSDefaultPreStep - The default pre-stepping function which does nothing. 1398000e7ae3SMatthew Knepley 1399000e7ae3SMatthew Knepley Collective on TS 1400000e7ae3SMatthew Knepley 1401000e7ae3SMatthew Knepley Input Parameters: 1402000e7ae3SMatthew Knepley . ts - The TS context obtained from TSCreate() 1403000e7ae3SMatthew Knepley 1404000e7ae3SMatthew Knepley Level: developer 1405000e7ae3SMatthew Knepley 1406000e7ae3SMatthew Knepley .keywords: TS, timestep 1407000e7ae3SMatthew Knepley @*/ 140863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts) 1409000e7ae3SMatthew Knepley { 1410000e7ae3SMatthew Knepley PetscFunctionBegin; 1411000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1412000e7ae3SMatthew Knepley } 1413000e7ae3SMatthew Knepley 1414e74ef692SMatthew Knepley #undef __FUNCT__ 1415e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep" 1416ac226902SBarry Smith /*@C 1417000e7ae3SMatthew Knepley TSSetPostStep - Sets the general-purpose function 1418000e7ae3SMatthew Knepley called once at the end of time stepping. 1419000e7ae3SMatthew Knepley 1420000e7ae3SMatthew Knepley Collective on TS 1421000e7ae3SMatthew Knepley 1422000e7ae3SMatthew Knepley Input Parameters: 1423000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1424000e7ae3SMatthew Knepley - func - The function 1425000e7ae3SMatthew Knepley 1426000e7ae3SMatthew Knepley Calling sequence of func: 1427000e7ae3SMatthew Knepley . func (TS ts); 1428000e7ae3SMatthew Knepley 1429000e7ae3SMatthew Knepley Level: intermediate 1430000e7ae3SMatthew Knepley 1431000e7ae3SMatthew Knepley .keywords: TS, timestep 1432000e7ae3SMatthew Knepley @*/ 143363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1434000e7ae3SMatthew Knepley { 1435000e7ae3SMatthew Knepley PetscFunctionBegin; 14364482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1437000e7ae3SMatthew Knepley ts->ops->poststep = func; 1438000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1439000e7ae3SMatthew Knepley } 1440000e7ae3SMatthew Knepley 1441e74ef692SMatthew Knepley #undef __FUNCT__ 1442e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep" 1443000e7ae3SMatthew Knepley /*@ 1444000e7ae3SMatthew Knepley TSDefaultPostStep - The default post-stepping function which does nothing. 1445000e7ae3SMatthew Knepley 1446000e7ae3SMatthew Knepley Collective on TS 1447000e7ae3SMatthew Knepley 1448000e7ae3SMatthew Knepley Input Parameters: 1449000e7ae3SMatthew Knepley . ts - The TS context obtained from TSCreate() 1450000e7ae3SMatthew Knepley 1451000e7ae3SMatthew Knepley Level: developer 1452000e7ae3SMatthew Knepley 1453000e7ae3SMatthew Knepley .keywords: TS, timestep 1454000e7ae3SMatthew Knepley @*/ 145563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts) 1456000e7ae3SMatthew Knepley { 1457000e7ae3SMatthew Knepley PetscFunctionBegin; 1458000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1459000e7ae3SMatthew Knepley } 1460000e7ae3SMatthew Knepley 1461d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 1462d763cef2SBarry Smith 14634a2ae208SSatish Balay #undef __FUNCT__ 1464a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet" 1465d763cef2SBarry Smith /*@C 1466a6570f20SBarry Smith TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1467d763cef2SBarry Smith timestep to display the iteration's progress. 1468d763cef2SBarry Smith 1469d763cef2SBarry Smith Collective on TS 1470d763cef2SBarry Smith 1471d763cef2SBarry Smith Input Parameters: 1472d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1473d763cef2SBarry Smith . func - monitoring routine 1474329f5518SBarry Smith . mctx - [optional] user-defined context for private data for the 1475b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desired) 1476b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1477b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 1478d763cef2SBarry Smith 1479d763cef2SBarry Smith Calling sequence of func: 1480a7cc72afSBarry Smith $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1481d763cef2SBarry Smith 1482d763cef2SBarry Smith + ts - the TS context 1483d763cef2SBarry Smith . steps - iteration number 14841f06c33eSBarry Smith . time - current time 1485d763cef2SBarry Smith . x - current iterate 1486d763cef2SBarry Smith - mctx - [optional] monitoring context 1487d763cef2SBarry Smith 1488d763cef2SBarry Smith Notes: 1489d763cef2SBarry Smith This routine adds an additional monitor to the list of monitors that 1490d763cef2SBarry Smith already has been loaded. 1491d763cef2SBarry Smith 1492025f1a04SBarry Smith Fortran notes: Only a single monitor function can be set for each TS object 1493025f1a04SBarry Smith 1494d763cef2SBarry Smith Level: intermediate 1495d763cef2SBarry Smith 1496d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1497d763cef2SBarry Smith 1498a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel() 1499d763cef2SBarry Smith @*/ 1500a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1501d763cef2SBarry Smith { 1502d763cef2SBarry Smith PetscFunctionBegin; 15034482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1504d763cef2SBarry Smith if (ts->numbermonitors >= MAXTSMONITORS) { 150529bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1506d763cef2SBarry Smith } 1507d763cef2SBarry Smith ts->monitor[ts->numbermonitors] = monitor; 1508329f5518SBarry Smith ts->mdestroy[ts->numbermonitors] = mdestroy; 1509d763cef2SBarry Smith ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1510d763cef2SBarry Smith PetscFunctionReturn(0); 1511d763cef2SBarry Smith } 1512d763cef2SBarry Smith 15134a2ae208SSatish Balay #undef __FUNCT__ 1514a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel" 1515d763cef2SBarry Smith /*@C 1516a6570f20SBarry Smith TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1517d763cef2SBarry Smith 1518d763cef2SBarry Smith Collective on TS 1519d763cef2SBarry Smith 1520d763cef2SBarry Smith Input Parameters: 1521d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1522d763cef2SBarry Smith 1523d763cef2SBarry Smith Notes: 1524d763cef2SBarry Smith There is no way to remove a single, specific monitor. 1525d763cef2SBarry Smith 1526d763cef2SBarry Smith Level: intermediate 1527d763cef2SBarry Smith 1528d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1529d763cef2SBarry Smith 1530a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet() 1531d763cef2SBarry Smith @*/ 1532a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts) 1533d763cef2SBarry Smith { 1534d952e501SBarry Smith PetscErrorCode ierr; 1535d952e501SBarry Smith PetscInt i; 1536d952e501SBarry Smith 1537d763cef2SBarry Smith PetscFunctionBegin; 15384482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1539d952e501SBarry Smith for (i=0; i<ts->numbermonitors; i++) { 1540d952e501SBarry Smith if (ts->mdestroy[i]) { 1541d952e501SBarry Smith ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 1542d952e501SBarry Smith } 1543d952e501SBarry Smith } 1544d763cef2SBarry Smith ts->numbermonitors = 0; 1545d763cef2SBarry Smith PetscFunctionReturn(0); 1546d763cef2SBarry Smith } 1547d763cef2SBarry Smith 15484a2ae208SSatish Balay #undef __FUNCT__ 1549a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault" 1550d8e5e3e6SSatish Balay /*@ 1551a6570f20SBarry Smith TSMonitorDefault - Sets the Default monitor 15525516499fSSatish Balay 15535516499fSSatish Balay Level: intermediate 155441251cbbSSatish Balay 15555516499fSSatish Balay .keywords: TS, set, monitor 15565516499fSSatish Balay 155741251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet() 155841251cbbSSatish Balay @*/ 1559a6570f20SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1560d763cef2SBarry Smith { 1561dfbe8321SBarry Smith PetscErrorCode ierr; 1562a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 1563d132466eSBarry Smith 1564d763cef2SBarry Smith PetscFunctionBegin; 1565a34d58ebSBarry Smith if (!ctx) { 15667adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 1567a34d58ebSBarry Smith } 1568a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1569a34d58ebSBarry Smith if (!ctx) { 1570a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 1571a34d58ebSBarry Smith } 1572d763cef2SBarry Smith PetscFunctionReturn(0); 1573d763cef2SBarry Smith } 1574d763cef2SBarry Smith 15754a2ae208SSatish Balay #undef __FUNCT__ 15764a2ae208SSatish Balay #define __FUNCT__ "TSStep" 1577d763cef2SBarry Smith /*@ 1578d763cef2SBarry Smith TSStep - Steps the requested number of timesteps. 1579d763cef2SBarry Smith 1580d763cef2SBarry Smith Collective on TS 1581d763cef2SBarry Smith 1582d763cef2SBarry Smith Input Parameter: 1583d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1584d763cef2SBarry Smith 1585d763cef2SBarry Smith Output Parameters: 1586d763cef2SBarry Smith + steps - number of iterations until termination 1587142b95e3SSatish Balay - ptime - time until termination 1588d763cef2SBarry Smith 1589d763cef2SBarry Smith Level: beginner 1590d763cef2SBarry Smith 1591d763cef2SBarry Smith .keywords: TS, timestep, solve 1592d763cef2SBarry Smith 1593d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy() 1594d763cef2SBarry Smith @*/ 159563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1596d763cef2SBarry Smith { 1597dfbe8321SBarry Smith PetscErrorCode ierr; 1598d763cef2SBarry Smith 1599d763cef2SBarry Smith PetscFunctionBegin; 16004482741eSBarry Smith PetscValidHeaderSpecific(ts, TS_COOKIE,1); 1601d405a339SMatthew Knepley if (!ts->setupcalled) { 1602d405a339SMatthew Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 1603d405a339SMatthew Knepley } 1604d405a339SMatthew Knepley 1605d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1606d405a339SMatthew Knepley ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1607000e7ae3SMatthew Knepley ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1608d405a339SMatthew Knepley ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1609d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1610d405a339SMatthew Knepley 16114bb05414SBarry Smith if (!PetscPreLoadingOn) { 16127adad957SLisandro Dalcin ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1613d405a339SMatthew Knepley } 1614d763cef2SBarry Smith PetscFunctionReturn(0); 1615d763cef2SBarry Smith } 1616d763cef2SBarry Smith 16174a2ae208SSatish Balay #undef __FUNCT__ 16186a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve" 16196a4d4014SLisandro Dalcin /*@ 16206a4d4014SLisandro Dalcin TSSolve - Steps the requested number of timesteps. 16216a4d4014SLisandro Dalcin 16226a4d4014SLisandro Dalcin Collective on TS 16236a4d4014SLisandro Dalcin 16246a4d4014SLisandro Dalcin Input Parameter: 16256a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 16266a4d4014SLisandro Dalcin - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 16276a4d4014SLisandro Dalcin 16286a4d4014SLisandro Dalcin Level: beginner 16296a4d4014SLisandro Dalcin 16306a4d4014SLisandro Dalcin .keywords: TS, timestep, solve 16316a4d4014SLisandro Dalcin 16326a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep() 16336a4d4014SLisandro Dalcin @*/ 16346a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x) 16356a4d4014SLisandro Dalcin { 16366a4d4014SLisandro Dalcin PetscInt steps; 16376a4d4014SLisandro Dalcin PetscReal ptime; 16386a4d4014SLisandro Dalcin PetscErrorCode ierr; 16396a4d4014SLisandro Dalcin PetscFunctionBegin; 16406a4d4014SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_COOKIE,1); 16416a4d4014SLisandro Dalcin /* set solution vector if provided */ 16426a4d4014SLisandro Dalcin if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 16436a4d4014SLisandro Dalcin /* reset time step and iteration counters */ 16446a4d4014SLisandro Dalcin ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 16456a4d4014SLisandro Dalcin /* steps the requested number of timesteps. */ 16466a4d4014SLisandro Dalcin ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 16476a4d4014SLisandro Dalcin PetscFunctionReturn(0); 16486a4d4014SLisandro Dalcin } 16496a4d4014SLisandro Dalcin 16506a4d4014SLisandro Dalcin #undef __FUNCT__ 16514a2ae208SSatish Balay #define __FUNCT__ "TSMonitor" 1652d763cef2SBarry Smith /* 1653d763cef2SBarry Smith Runs the user provided monitor routines, if they exists. 1654d763cef2SBarry Smith */ 1655a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1656d763cef2SBarry Smith { 16576849ba73SBarry Smith PetscErrorCode ierr; 1658a7cc72afSBarry Smith PetscInt i,n = ts->numbermonitors; 1659d763cef2SBarry Smith 1660d763cef2SBarry Smith PetscFunctionBegin; 1661d763cef2SBarry Smith for (i=0; i<n; i++) { 1662142b95e3SSatish Balay ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1663d763cef2SBarry Smith } 1664d763cef2SBarry Smith PetscFunctionReturn(0); 1665d763cef2SBarry Smith } 1666d763cef2SBarry Smith 1667d763cef2SBarry Smith /* ------------------------------------------------------------------------*/ 1668d763cef2SBarry Smith 16694a2ae208SSatish Balay #undef __FUNCT__ 1670a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate" 1671d763cef2SBarry Smith /*@C 1672a6570f20SBarry Smith TSMonitorLGCreate - Creates a line graph context for use with 1673d763cef2SBarry Smith TS to monitor convergence of preconditioned residual norms. 1674d763cef2SBarry Smith 1675d763cef2SBarry Smith Collective on TS 1676d763cef2SBarry Smith 1677d763cef2SBarry Smith Input Parameters: 1678d763cef2SBarry Smith + host - the X display to open, or null for the local machine 1679d763cef2SBarry Smith . label - the title to put in the title bar 16807c922b88SBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 1681d763cef2SBarry Smith - m, n - the screen width and height in pixels 1682d763cef2SBarry Smith 1683d763cef2SBarry Smith Output Parameter: 1684d763cef2SBarry Smith . draw - the drawing context 1685d763cef2SBarry Smith 1686d763cef2SBarry Smith Options Database Key: 1687a6570f20SBarry Smith . -ts_monitor_draw - automatically sets line graph monitor 1688d763cef2SBarry Smith 1689d763cef2SBarry Smith Notes: 1690a6570f20SBarry Smith Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1691d763cef2SBarry Smith 1692d763cef2SBarry Smith Level: intermediate 1693d763cef2SBarry Smith 16947c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso 1695d763cef2SBarry Smith 1696a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet() 16977c922b88SBarry Smith 1698d763cef2SBarry Smith @*/ 1699a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1700d763cef2SBarry Smith { 1701b0a32e0cSBarry Smith PetscDraw win; 1702dfbe8321SBarry Smith PetscErrorCode ierr; 1703d763cef2SBarry Smith 1704d763cef2SBarry Smith PetscFunctionBegin; 1705b0a32e0cSBarry Smith ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1706b0a32e0cSBarry Smith ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1707b0a32e0cSBarry Smith ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1708b0a32e0cSBarry Smith ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1709d763cef2SBarry Smith 171052e6d16bSBarry Smith ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1711d763cef2SBarry Smith PetscFunctionReturn(0); 1712d763cef2SBarry Smith } 1713d763cef2SBarry Smith 17144a2ae208SSatish Balay #undef __FUNCT__ 1715a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG" 1716a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1717d763cef2SBarry Smith { 1718b0a32e0cSBarry Smith PetscDrawLG lg = (PetscDrawLG) monctx; 171987828ca2SBarry Smith PetscReal x,y = ptime; 1720dfbe8321SBarry Smith PetscErrorCode ierr; 1721d763cef2SBarry Smith 1722d763cef2SBarry Smith PetscFunctionBegin; 17237c922b88SBarry Smith if (!monctx) { 17247c922b88SBarry Smith MPI_Comm comm; 1725b0a32e0cSBarry Smith PetscViewer viewer; 17267c922b88SBarry Smith 17277c922b88SBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1728b0a32e0cSBarry Smith viewer = PETSC_VIEWER_DRAW_(comm); 1729b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 17307c922b88SBarry Smith } 17317c922b88SBarry Smith 1732b0a32e0cSBarry Smith if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 173387828ca2SBarry Smith x = (PetscReal)n; 1734b0a32e0cSBarry Smith ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1735d763cef2SBarry Smith if (n < 20 || (n % 5)) { 1736b0a32e0cSBarry Smith ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1737d763cef2SBarry Smith } 1738d763cef2SBarry Smith PetscFunctionReturn(0); 1739d763cef2SBarry Smith } 1740d763cef2SBarry Smith 17414a2ae208SSatish Balay #undef __FUNCT__ 1742a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy" 1743d763cef2SBarry Smith /*@C 1744a6570f20SBarry Smith TSMonitorLGDestroy - Destroys a line graph context that was created 1745a6570f20SBarry Smith with TSMonitorLGCreate(). 1746d763cef2SBarry Smith 1747b0a32e0cSBarry Smith Collective on PetscDrawLG 1748d763cef2SBarry Smith 1749d763cef2SBarry Smith Input Parameter: 1750d763cef2SBarry Smith . draw - the drawing context 1751d763cef2SBarry Smith 1752d763cef2SBarry Smith Level: intermediate 1753d763cef2SBarry Smith 1754d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy 1755d763cef2SBarry Smith 1756a6570f20SBarry Smith .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1757d763cef2SBarry Smith @*/ 1758a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg) 1759d763cef2SBarry Smith { 1760b0a32e0cSBarry Smith PetscDraw draw; 1761dfbe8321SBarry Smith PetscErrorCode ierr; 1762d763cef2SBarry Smith 1763d763cef2SBarry Smith PetscFunctionBegin; 1764b0a32e0cSBarry Smith ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1765b0a32e0cSBarry Smith ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1766b0a32e0cSBarry Smith ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1767d763cef2SBarry Smith PetscFunctionReturn(0); 1768d763cef2SBarry Smith } 1769d763cef2SBarry Smith 17704a2ae208SSatish Balay #undef __FUNCT__ 17714a2ae208SSatish Balay #define __FUNCT__ "TSGetTime" 1772d763cef2SBarry Smith /*@ 1773d763cef2SBarry Smith TSGetTime - Gets the current time. 1774d763cef2SBarry Smith 1775d763cef2SBarry Smith Not Collective 1776d763cef2SBarry Smith 1777d763cef2SBarry Smith Input Parameter: 1778d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1779d763cef2SBarry Smith 1780d763cef2SBarry Smith Output Parameter: 1781d763cef2SBarry Smith . t - the current time 1782d763cef2SBarry Smith 1783d763cef2SBarry Smith Level: beginner 1784d763cef2SBarry Smith 1785d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1786d763cef2SBarry Smith 1787d763cef2SBarry Smith .keywords: TS, get, time 1788d763cef2SBarry Smith @*/ 178963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t) 1790d763cef2SBarry Smith { 1791d763cef2SBarry Smith PetscFunctionBegin; 17924482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 17934482741eSBarry Smith PetscValidDoublePointer(t,2); 1794d763cef2SBarry Smith *t = ts->ptime; 1795d763cef2SBarry Smith PetscFunctionReturn(0); 1796d763cef2SBarry Smith } 1797d763cef2SBarry Smith 17984a2ae208SSatish Balay #undef __FUNCT__ 17996a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime" 18006a4d4014SLisandro Dalcin /*@ 18016a4d4014SLisandro Dalcin TSSetTime - Allows one to reset the time. 18026a4d4014SLisandro Dalcin 18036a4d4014SLisandro Dalcin Collective on TS 18046a4d4014SLisandro Dalcin 18056a4d4014SLisandro Dalcin Input Parameters: 18066a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 18076a4d4014SLisandro Dalcin - time - the time 18086a4d4014SLisandro Dalcin 18096a4d4014SLisandro Dalcin Level: intermediate 18106a4d4014SLisandro Dalcin 18116a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration() 18126a4d4014SLisandro Dalcin 18136a4d4014SLisandro Dalcin .keywords: TS, set, time 18146a4d4014SLisandro Dalcin @*/ 18156a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t) 18166a4d4014SLisandro Dalcin { 18176a4d4014SLisandro Dalcin PetscFunctionBegin; 18186a4d4014SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_COOKIE,1); 18196a4d4014SLisandro Dalcin ts->ptime = t; 18206a4d4014SLisandro Dalcin PetscFunctionReturn(0); 18216a4d4014SLisandro Dalcin } 18226a4d4014SLisandro Dalcin 18236a4d4014SLisandro Dalcin #undef __FUNCT__ 18244a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix" 1825d763cef2SBarry Smith /*@C 1826d763cef2SBarry Smith TSSetOptionsPrefix - Sets the prefix used for searching for all 1827d763cef2SBarry Smith TS options in the database. 1828d763cef2SBarry Smith 1829d763cef2SBarry Smith Collective on TS 1830d763cef2SBarry Smith 1831d763cef2SBarry Smith Input Parameter: 1832d763cef2SBarry Smith + ts - The TS context 1833d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1834d763cef2SBarry Smith 1835d763cef2SBarry Smith Notes: 1836d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1837d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1838d763cef2SBarry Smith hyphen. 1839d763cef2SBarry Smith 1840d763cef2SBarry Smith Level: advanced 1841d763cef2SBarry Smith 1842d763cef2SBarry Smith .keywords: TS, set, options, prefix, database 1843d763cef2SBarry Smith 1844d763cef2SBarry Smith .seealso: TSSetFromOptions() 1845d763cef2SBarry Smith 1846d763cef2SBarry Smith @*/ 184763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[]) 1848d763cef2SBarry Smith { 1849dfbe8321SBarry Smith PetscErrorCode ierr; 1850d763cef2SBarry Smith 1851d763cef2SBarry Smith PetscFunctionBegin; 18524482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1853d763cef2SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1854d763cef2SBarry Smith switch(ts->problem_type) { 1855d763cef2SBarry Smith case TS_NONLINEAR: 185659580b9cSBarry Smith if (ts->snes) { 1857d763cef2SBarry Smith ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 185859580b9cSBarry Smith } 1859d763cef2SBarry Smith break; 1860d763cef2SBarry Smith case TS_LINEAR: 186159580b9cSBarry Smith if (ts->ksp) { 186294b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 186359580b9cSBarry Smith } 1864d763cef2SBarry Smith break; 1865d763cef2SBarry Smith } 1866d763cef2SBarry Smith PetscFunctionReturn(0); 1867d763cef2SBarry Smith } 1868d763cef2SBarry Smith 1869d763cef2SBarry Smith 18704a2ae208SSatish Balay #undef __FUNCT__ 18714a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix" 1872d763cef2SBarry Smith /*@C 1873d763cef2SBarry Smith TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1874d763cef2SBarry Smith TS options in the database. 1875d763cef2SBarry Smith 1876d763cef2SBarry Smith Collective on TS 1877d763cef2SBarry Smith 1878d763cef2SBarry Smith Input Parameter: 1879d763cef2SBarry Smith + ts - The TS context 1880d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1881d763cef2SBarry Smith 1882d763cef2SBarry Smith Notes: 1883d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1884d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1885d763cef2SBarry Smith hyphen. 1886d763cef2SBarry Smith 1887d763cef2SBarry Smith Level: advanced 1888d763cef2SBarry Smith 1889d763cef2SBarry Smith .keywords: TS, append, options, prefix, database 1890d763cef2SBarry Smith 1891d763cef2SBarry Smith .seealso: TSGetOptionsPrefix() 1892d763cef2SBarry Smith 1893d763cef2SBarry Smith @*/ 189463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[]) 1895d763cef2SBarry Smith { 1896dfbe8321SBarry Smith PetscErrorCode ierr; 1897d763cef2SBarry Smith 1898d763cef2SBarry Smith PetscFunctionBegin; 18994482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 1900d763cef2SBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1901d763cef2SBarry Smith switch(ts->problem_type) { 1902d763cef2SBarry Smith case TS_NONLINEAR: 19031ac94b3bSBarry Smith if (ts->snes) { 1904d763cef2SBarry Smith ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 19051ac94b3bSBarry Smith } 1906d763cef2SBarry Smith break; 1907d763cef2SBarry Smith case TS_LINEAR: 19081ac94b3bSBarry Smith if (ts->ksp) { 190994b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 19101ac94b3bSBarry Smith } 1911d763cef2SBarry Smith break; 1912d763cef2SBarry Smith } 1913d763cef2SBarry Smith PetscFunctionReturn(0); 1914d763cef2SBarry Smith } 1915d763cef2SBarry Smith 19164a2ae208SSatish Balay #undef __FUNCT__ 19174a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix" 1918d763cef2SBarry Smith /*@C 1919d763cef2SBarry Smith TSGetOptionsPrefix - Sets the prefix used for searching for all 1920d763cef2SBarry Smith TS options in the database. 1921d763cef2SBarry Smith 1922d763cef2SBarry Smith Not Collective 1923d763cef2SBarry Smith 1924d763cef2SBarry Smith Input Parameter: 1925d763cef2SBarry Smith . ts - The TS context 1926d763cef2SBarry Smith 1927d763cef2SBarry Smith Output Parameter: 1928d763cef2SBarry Smith . prefix - A pointer to the prefix string used 1929d763cef2SBarry Smith 1930d763cef2SBarry Smith Notes: On the fortran side, the user should pass in a string 'prifix' of 1931d763cef2SBarry Smith sufficient length to hold the prefix. 1932d763cef2SBarry Smith 1933d763cef2SBarry Smith Level: intermediate 1934d763cef2SBarry Smith 1935d763cef2SBarry Smith .keywords: TS, get, options, prefix, database 1936d763cef2SBarry Smith 1937d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix() 1938d763cef2SBarry Smith @*/ 1939e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[]) 1940d763cef2SBarry Smith { 1941dfbe8321SBarry Smith PetscErrorCode ierr; 1942d763cef2SBarry Smith 1943d763cef2SBarry Smith PetscFunctionBegin; 19444482741eSBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE,1); 19454482741eSBarry Smith PetscValidPointer(prefix,2); 1946d763cef2SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1947d763cef2SBarry Smith PetscFunctionReturn(0); 1948d763cef2SBarry Smith } 1949d763cef2SBarry Smith 19504a2ae208SSatish Balay #undef __FUNCT__ 19514a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian" 1952d763cef2SBarry Smith /*@C 1953d763cef2SBarry Smith TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1954d763cef2SBarry Smith 1955d763cef2SBarry Smith Not Collective, but parallel objects are returned if TS is parallel 1956d763cef2SBarry Smith 1957d763cef2SBarry Smith Input Parameter: 1958d763cef2SBarry Smith . ts - The TS context obtained from TSCreate() 1959d763cef2SBarry Smith 1960d763cef2SBarry Smith Output Parameters: 1961d763cef2SBarry Smith + J - The Jacobian J of F, where U_t = F(U,t) 1962d763cef2SBarry Smith . M - The preconditioner matrix, usually the same as J 1963d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine 1964d763cef2SBarry Smith 1965d763cef2SBarry Smith Notes: You can pass in PETSC_NULL for any return argument you do not need. 1966d763cef2SBarry Smith 1967d763cef2SBarry Smith Level: intermediate 1968d763cef2SBarry Smith 196926d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 1970d763cef2SBarry Smith 1971d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian 1972d763cef2SBarry Smith @*/ 197363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 1974d763cef2SBarry Smith { 1975d763cef2SBarry Smith PetscFunctionBegin; 197626d46c62SHong Zhang if (J) *J = ts->Arhs; 197726d46c62SHong Zhang if (M) *M = ts->B; 197826d46c62SHong Zhang if (ctx) *ctx = ts->jacP; 1979d763cef2SBarry Smith PetscFunctionReturn(0); 1980d763cef2SBarry Smith } 1981d763cef2SBarry Smith 19821713a123SBarry Smith #undef __FUNCT__ 1983a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution" 19841713a123SBarry Smith /*@C 1985a6570f20SBarry Smith TSMonitorSolution - Monitors progress of the TS solvers by calling 19861713a123SBarry Smith VecView() for the solution at each timestep 19871713a123SBarry Smith 19881713a123SBarry Smith Collective on TS 19891713a123SBarry Smith 19901713a123SBarry Smith Input Parameters: 19911713a123SBarry Smith + ts - the TS context 19921713a123SBarry Smith . step - current time-step 1993142b95e3SSatish Balay . ptime - current time 19941713a123SBarry Smith - dummy - either a viewer or PETSC_NULL 19951713a123SBarry Smith 19961713a123SBarry Smith Level: intermediate 19971713a123SBarry Smith 19981713a123SBarry Smith .keywords: TS, vector, monitor, view 19991713a123SBarry Smith 2000a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 20011713a123SBarry Smith @*/ 2002a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 20031713a123SBarry Smith { 2004dfbe8321SBarry Smith PetscErrorCode ierr; 20051713a123SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 20061713a123SBarry Smith 20071713a123SBarry Smith PetscFunctionBegin; 2008a34d58ebSBarry Smith if (!dummy) { 20097adad957SLisandro Dalcin viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 20101713a123SBarry Smith } 20111713a123SBarry Smith ierr = VecView(x,viewer);CHKERRQ(ierr); 20121713a123SBarry Smith PetscFunctionReturn(0); 20131713a123SBarry Smith } 20141713a123SBarry Smith 20151713a123SBarry Smith 20161713a123SBarry Smith 2017