163dd3a1aSKris Buschelman #define PETSCTS_DLL 263dd3a1aSKris Buschelman 37c4f633dSBarry Smith #include "private/tsimpl.h" /*I "petscts.h" I*/ 4d763cef2SBarry Smith 5d5ba7fb7SMatthew Knepley /* Logging support */ 60700a824SBarry Smith PetscClassId PETSCTS_DLLEXPORT TS_CLASSID; 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 { 359596e0b4SJed Brown defaultType = TSEULER; 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: 594d91e141SJed Brown + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 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; 810700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,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: 140e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,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; 2330700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2340700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 235c9780b6fSBarry Smith PetscCheckSameComm(ts,1,X,3); 23617186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only"); 237000e7ae3SMatthew Knepley if (ts->ops->rhsjacobian) { 238d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 239a7a1495cSBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 240a7a1495cSBarry Smith PetscStackPush("TS user Jacobian function"); 241000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 242a7a1495cSBarry Smith PetscStackPop; 243d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 244a7a1495cSBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 2450700a824SBarry Smith PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 2460700a824SBarry Smith PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 247ef66eb69SBarry Smith } else { 248ef66eb69SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 249ef66eb69SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 250ef66eb69SBarry Smith if (*A != *B) { 251ef66eb69SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252ef66eb69SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 253ef66eb69SBarry Smith } 254ef66eb69SBarry Smith } 255a7a1495cSBarry Smith PetscFunctionReturn(0); 256a7a1495cSBarry Smith } 257a7a1495cSBarry Smith 2584a2ae208SSatish Balay #undef __FUNCT__ 2594a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction" 260316643e7SJed Brown /*@ 261d763cef2SBarry Smith TSComputeRHSFunction - Evaluates the right-hand-side function. 262d763cef2SBarry Smith 263316643e7SJed Brown Collective on TS and Vec 264316643e7SJed Brown 265316643e7SJed Brown Input Parameters: 266316643e7SJed Brown + ts - the TS context 267316643e7SJed Brown . t - current time 268316643e7SJed Brown - x - state vector 269316643e7SJed Brown 270316643e7SJed Brown Output Parameter: 271316643e7SJed Brown . y - right hand side 272316643e7SJed Brown 273316643e7SJed Brown Note: 274316643e7SJed Brown Most users should not need to explicitly call this routine, as it 275316643e7SJed Brown is used internally within the nonlinear solvers. 276316643e7SJed Brown 277316643e7SJed Brown If the user did not provide a function but merely a matrix, 278d763cef2SBarry Smith this routine applies the matrix. 279316643e7SJed Brown 280316643e7SJed Brown Level: developer 281316643e7SJed Brown 282316643e7SJed Brown .keywords: TS, compute 283316643e7SJed Brown 284316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction() 285316643e7SJed Brown @*/ 286dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 287d763cef2SBarry Smith { 288dfbe8321SBarry Smith PetscErrorCode ierr; 289d763cef2SBarry Smith 290d763cef2SBarry Smith PetscFunctionBegin; 2910700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2920700a824SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2930700a824SBarry Smith PetscValidHeaderSpecific(y,VEC_CLASSID,4); 294d763cef2SBarry Smith 295d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 296000e7ae3SMatthew Knepley if (ts->ops->rhsfunction) { 297d763cef2SBarry Smith PetscStackPush("TS user right-hand-side function"); 298000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 299d763cef2SBarry Smith PetscStackPop; 3007c922b88SBarry Smith } else { 301000e7ae3SMatthew Knepley if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 302d763cef2SBarry Smith MatStructure flg; 303d763cef2SBarry Smith PetscStackPush("TS user right-hand-side matrix function"); 3048beb423aSHong Zhang ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 305d763cef2SBarry Smith PetscStackPop; 306d763cef2SBarry Smith } 3078beb423aSHong Zhang ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr); 3087c922b88SBarry Smith } 309d763cef2SBarry Smith 310d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 311d763cef2SBarry Smith 312d763cef2SBarry Smith PetscFunctionReturn(0); 313d763cef2SBarry Smith } 314d763cef2SBarry Smith 3154a2ae208SSatish Balay #undef __FUNCT__ 316316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction" 317316643e7SJed Brown /*@ 318316643e7SJed Brown TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 319316643e7SJed Brown 320316643e7SJed Brown Collective on TS and Vec 321316643e7SJed Brown 322316643e7SJed Brown Input Parameters: 323316643e7SJed Brown + ts - the TS context 324316643e7SJed Brown . t - current time 325316643e7SJed Brown . X - state vector 326316643e7SJed Brown - Xdot - time derivative of state vector 327316643e7SJed Brown 328316643e7SJed Brown Output Parameter: 329316643e7SJed Brown . Y - right hand side 330316643e7SJed Brown 331316643e7SJed Brown Note: 332316643e7SJed Brown Most users should not need to explicitly call this routine, as it 333316643e7SJed Brown is used internally within the nonlinear solvers. 334316643e7SJed Brown 335316643e7SJed Brown If the user did did not write their equations in implicit form, this 336316643e7SJed Brown function recasts them in implicit form. 337316643e7SJed Brown 338316643e7SJed Brown Level: developer 339316643e7SJed Brown 340316643e7SJed Brown .keywords: TS, compute 341316643e7SJed Brown 342316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction() 343316643e7SJed Brown @*/ 344316643e7SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y) 345316643e7SJed Brown { 346316643e7SJed Brown PetscErrorCode ierr; 347316643e7SJed Brown 348316643e7SJed Brown PetscFunctionBegin; 3490700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3500700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 3510700a824SBarry Smith PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 3520700a824SBarry Smith PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 353316643e7SJed Brown 354316643e7SJed Brown ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 355316643e7SJed Brown if (ts->ops->ifunction) { 356316643e7SJed Brown PetscStackPush("TS user implicit function"); 357316643e7SJed Brown ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 358316643e7SJed Brown PetscStackPop; 359316643e7SJed Brown } else { 360316643e7SJed Brown if (ts->ops->rhsfunction) { 361316643e7SJed Brown PetscStackPush("TS user right-hand-side function"); 362316643e7SJed Brown ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr); 363316643e7SJed Brown PetscStackPop; 364316643e7SJed Brown } else { 365316643e7SJed Brown if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 366316643e7SJed Brown MatStructure flg; 3674a6899ffSJed Brown /* Note: flg is not being used. 3684a6899ffSJed Brown For it to be useful, we'd have to cache it and then apply it in TSComputeIJacobian. 3694a6899ffSJed Brown */ 370316643e7SJed Brown PetscStackPush("TS user right-hand-side matrix function"); 371316643e7SJed Brown ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 372316643e7SJed Brown PetscStackPop; 373316643e7SJed Brown } 374316643e7SJed Brown ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr); 375316643e7SJed Brown } 376316643e7SJed Brown 3774a6899ffSJed Brown /* Convert to implicit form: F(X,Xdot) = Alhs * Xdot - Frhs(X) */ 3784a6899ffSJed Brown if (ts->Alhs) { 3794a6899ffSJed Brown if (ts->ops->lhsmatrix) { 3804a6899ffSJed Brown MatStructure flg; 3814a6899ffSJed Brown PetscStackPush("TS user left-hand-side matrix function"); 3824a6899ffSJed Brown ierr = (*ts->ops->lhsmatrix)(ts,t,&ts->Alhs,PETSC_NULL,&flg,ts->jacP);CHKERRQ(ierr); 3834a6899ffSJed Brown PetscStackPop; 3844a6899ffSJed Brown } 3854a6899ffSJed Brown ierr = VecScale(Y,-1.);CHKERRQ(ierr); 3864a6899ffSJed Brown ierr = MatMultAdd(ts->Alhs,Xdot,Y,Y);CHKERRQ(ierr); 3874a6899ffSJed Brown } else { 388ace68cafSJed Brown ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 389316643e7SJed Brown } 3904a6899ffSJed Brown } 391316643e7SJed Brown ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 392316643e7SJed Brown PetscFunctionReturn(0); 393316643e7SJed Brown } 394316643e7SJed Brown 395316643e7SJed Brown #undef __FUNCT__ 396316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian" 397316643e7SJed Brown /*@ 398316643e7SJed Brown TSComputeIJacobian - Evaluates the Jacobian of the DAE 399316643e7SJed Brown 400316643e7SJed Brown Collective on TS and Vec 401316643e7SJed Brown 402316643e7SJed Brown Input 403316643e7SJed Brown Input Parameters: 404316643e7SJed Brown + ts - the TS context 405316643e7SJed Brown . t - current timestep 406316643e7SJed Brown . X - state vector 407316643e7SJed Brown . Xdot - time derivative of state vector 408316643e7SJed Brown - shift - shift to apply, see note below 409316643e7SJed Brown 410316643e7SJed Brown Output Parameters: 411316643e7SJed Brown + A - Jacobian matrix 412316643e7SJed Brown . B - optional preconditioning matrix 413316643e7SJed Brown - flag - flag indicating matrix structure 414316643e7SJed Brown 415316643e7SJed Brown Notes: 416316643e7SJed Brown If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 417316643e7SJed Brown 418316643e7SJed Brown dF/dX + shift*dF/dXdot 419316643e7SJed Brown 420316643e7SJed Brown Most users should not need to explicitly call this routine, as it 421316643e7SJed Brown is used internally within the nonlinear solvers. 422316643e7SJed Brown 423316643e7SJed Brown TSComputeIJacobian() is valid only for TS_NONLINEAR 424316643e7SJed Brown 425316643e7SJed Brown Level: developer 426316643e7SJed Brown 427316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix 428316643e7SJed Brown 429316643e7SJed Brown .seealso: TSSetIJacobian() 43063495f91SJed Brown @*/ 431316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg) 432316643e7SJed Brown { 433316643e7SJed Brown PetscErrorCode ierr; 434316643e7SJed Brown 435316643e7SJed Brown PetscFunctionBegin; 4360700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4370700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 4380700a824SBarry Smith PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 439316643e7SJed Brown PetscValidPointer(A,6); 4400700a824SBarry Smith PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 441316643e7SJed Brown PetscValidPointer(B,7); 4420700a824SBarry Smith PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 443316643e7SJed Brown PetscValidPointer(flg,8); 444316643e7SJed Brown 4454a6899ffSJed Brown *flg = SAME_NONZERO_PATTERN; /* In case it we're solving a linear problem in which case it wouldn't get initialized below. */ 446316643e7SJed Brown ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 447316643e7SJed Brown if (ts->ops->ijacobian) { 448316643e7SJed Brown PetscStackPush("TS user implicit Jacobian"); 449316643e7SJed Brown ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 450316643e7SJed Brown PetscStackPop; 451316643e7SJed Brown } else { 452316643e7SJed Brown if (ts->ops->rhsjacobian) { 453316643e7SJed Brown PetscStackPush("TS user right-hand-side Jacobian"); 454316643e7SJed Brown ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 455316643e7SJed Brown PetscStackPop; 456316643e7SJed Brown } else { 457316643e7SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 458316643e7SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 459316643e7SJed Brown if (*A != *B) { 460316643e7SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 461316643e7SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 462316643e7SJed Brown } 463316643e7SJed Brown } 464316643e7SJed Brown 465316643e7SJed Brown /* Convert to implicit form */ 466316643e7SJed Brown /* inefficient because these operations will normally traverse all matrix elements twice */ 467316643e7SJed Brown ierr = MatScale(*A,-1);CHKERRQ(ierr); 4684a6899ffSJed Brown if (ts->Alhs) { 4694a6899ffSJed Brown ierr = MatAXPY(*A,shift,ts->Alhs,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 4704a6899ffSJed Brown } else { 471316643e7SJed Brown ierr = MatShift(*A,shift);CHKERRQ(ierr); 4724a6899ffSJed Brown } 473316643e7SJed Brown if (*A != *B) { 474316643e7SJed Brown ierr = MatScale(*B,-1);CHKERRQ(ierr); 475316643e7SJed Brown ierr = MatShift(*B,shift);CHKERRQ(ierr); 476316643e7SJed Brown } 477316643e7SJed Brown } 478316643e7SJed Brown ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 479316643e7SJed Brown PetscFunctionReturn(0); 480316643e7SJed Brown } 481316643e7SJed Brown 482316643e7SJed Brown #undef __FUNCT__ 4834a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction" 484d763cef2SBarry Smith /*@C 485d763cef2SBarry Smith TSSetRHSFunction - Sets the routine for evaluating the function, 486d763cef2SBarry Smith F(t,u), where U_t = F(t,u). 487d763cef2SBarry Smith 488d763cef2SBarry Smith Collective on TS 489d763cef2SBarry Smith 490d763cef2SBarry Smith Input Parameters: 491d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 492d763cef2SBarry Smith . f - routine for evaluating the right-hand-side function 493d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 494d763cef2SBarry Smith function evaluation routine (may be PETSC_NULL) 495d763cef2SBarry Smith 496d763cef2SBarry Smith Calling sequence of func: 49787828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 498d763cef2SBarry Smith 499d763cef2SBarry Smith + t - current timestep 500d763cef2SBarry Smith . u - input vector 501d763cef2SBarry Smith . F - function vector 502d763cef2SBarry Smith - ctx - [optional] user-defined function context 503d763cef2SBarry Smith 504d763cef2SBarry Smith Important: 50576f2fa84SHong Zhang The user MUST call either this routine or TSSetMatrices(). 506d763cef2SBarry Smith 507d763cef2SBarry Smith Level: beginner 508d763cef2SBarry Smith 509d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function 510d763cef2SBarry Smith 51176f2fa84SHong Zhang .seealso: TSSetMatrices() 512d763cef2SBarry Smith @*/ 51363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 514d763cef2SBarry Smith { 515d763cef2SBarry Smith PetscFunctionBegin; 516d763cef2SBarry Smith 5170700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 51817186662SBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem"); 519000e7ae3SMatthew Knepley ts->ops->rhsfunction = f; 520d763cef2SBarry Smith ts->funP = ctx; 521d763cef2SBarry Smith PetscFunctionReturn(0); 522d763cef2SBarry Smith } 523d763cef2SBarry Smith 5244a2ae208SSatish Balay #undef __FUNCT__ 52595f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices" 52695f0b562SHong Zhang /*@C 52795f0b562SHong Zhang TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs, 52895f0b562SHong Zhang where Alhs(t) U_t = Arhs(t) U. 52995f0b562SHong Zhang 53095f0b562SHong Zhang Collective on TS 53195f0b562SHong Zhang 53295f0b562SHong Zhang Input Parameters: 53395f0b562SHong Zhang + ts - the TS context obtained from TSCreate() 53495f0b562SHong Zhang . Arhs - matrix 53595f0b562SHong Zhang . frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 53695f0b562SHong Zhang if Arhs is not a function of t. 53795f0b562SHong Zhang . Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix. 53895f0b562SHong Zhang . flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 53995f0b562SHong Zhang if Alhs is not a function of t. 54095f0b562SHong Zhang . flag - flag indicating information about the matrix structure of Arhs and Alhs. 54195f0b562SHong Zhang The available options are 54295f0b562SHong Zhang SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs 54395f0b562SHong Zhang DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs 54495f0b562SHong Zhang - ctx - [optional] user-defined context for private data for the 54595f0b562SHong Zhang matrix evaluation routine (may be PETSC_NULL) 54695f0b562SHong Zhang 54795f0b562SHong Zhang Calling sequence of func: 54895f0b562SHong Zhang $ func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx); 54995f0b562SHong Zhang 55095f0b562SHong Zhang + t - current timestep 55195f0b562SHong Zhang . A - matrix A, where U_t = A(t) U 55295f0b562SHong Zhang . B - preconditioner matrix, usually the same as A 55395f0b562SHong Zhang . flag - flag indicating information about the preconditioner matrix 55495f0b562SHong Zhang structure (same as flag in KSPSetOperators()) 55595f0b562SHong Zhang - ctx - [optional] user-defined context for matrix evaluation routine 55695f0b562SHong Zhang 55795f0b562SHong Zhang Notes: 55895f0b562SHong Zhang The routine func() takes Mat* as the matrix arguments rather than Mat. 55995f0b562SHong Zhang This allows the matrix evaluation routine to replace Arhs or Alhs with a 56095f0b562SHong Zhang completely new new matrix structure (not just different matrix elements) 56195f0b562SHong Zhang when appropriate, for instance, if the nonzero structure is changing 56295f0b562SHong Zhang throughout the global iterations. 56395f0b562SHong Zhang 56495f0b562SHong Zhang Important: 56595f0b562SHong Zhang The user MUST call either this routine or TSSetRHSFunction(). 56695f0b562SHong Zhang 56795f0b562SHong Zhang Level: beginner 56895f0b562SHong Zhang 56995f0b562SHong Zhang .keywords: TS, timestep, set, matrix 57095f0b562SHong Zhang 57195f0b562SHong Zhang .seealso: TSSetRHSFunction() 57295f0b562SHong Zhang @*/ 57395f0b562SHong 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) 57495f0b562SHong Zhang { 57595f0b562SHong Zhang PetscFunctionBegin; 5760700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 57792af4f6aSHong Zhang if (Arhs){ 5780700a824SBarry Smith PetscValidHeaderSpecific(Arhs,MAT_CLASSID,2); 57995f0b562SHong Zhang PetscCheckSameComm(ts,1,Arhs,2); 58095f0b562SHong Zhang ts->Arhs = Arhs; 58192af4f6aSHong Zhang ts->ops->rhsmatrix = frhs; 58292af4f6aSHong Zhang } 58392af4f6aSHong Zhang if (Alhs){ 5840700a824SBarry Smith PetscValidHeaderSpecific(Alhs,MAT_CLASSID,4); 58592af4f6aSHong Zhang PetscCheckSameComm(ts,1,Alhs,4); 58695f0b562SHong Zhang ts->Alhs = Alhs; 58792af4f6aSHong Zhang ts->ops->lhsmatrix = flhs; 58892af4f6aSHong Zhang } 58992af4f6aSHong Zhang 59092af4f6aSHong Zhang ts->jacP = ctx; 59195f0b562SHong Zhang ts->matflg = flag; 59295f0b562SHong Zhang PetscFunctionReturn(0); 59395f0b562SHong Zhang } 594d763cef2SBarry Smith 595aa644b49SHong Zhang #undef __FUNCT__ 596cda39b92SHong Zhang #define __FUNCT__ "TSGetMatrices" 597cda39b92SHong Zhang /*@C 598cda39b92SHong Zhang TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep, 599cda39b92SHong Zhang where Alhs(t) U_t = Arhs(t) U. 600cda39b92SHong Zhang 601cda39b92SHong Zhang Not Collective, but parallel objects are returned if TS is parallel 602cda39b92SHong Zhang 603cda39b92SHong Zhang Input Parameter: 604cda39b92SHong Zhang . ts - The TS context obtained from TSCreate() 605cda39b92SHong Zhang 606cda39b92SHong Zhang Output Parameters: 607cda39b92SHong Zhang + Arhs - The right-hand side matrix 608cda39b92SHong Zhang . Alhs - The left-hand side matrix 609cda39b92SHong Zhang - ctx - User-defined context for matrix evaluation routine 610cda39b92SHong Zhang 611cda39b92SHong Zhang Notes: You can pass in PETSC_NULL for any return argument you do not need. 612cda39b92SHong Zhang 613cda39b92SHong Zhang Level: intermediate 614cda39b92SHong Zhang 615cda39b92SHong Zhang .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 616cda39b92SHong Zhang 617cda39b92SHong Zhang .keywords: TS, timestep, get, matrix 618cda39b92SHong Zhang 619cda39b92SHong Zhang @*/ 620cda39b92SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx) 621cda39b92SHong Zhang { 622cda39b92SHong Zhang PetscFunctionBegin; 6230700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 624cda39b92SHong Zhang if (Arhs) *Arhs = ts->Arhs; 625cda39b92SHong Zhang if (Alhs) *Alhs = ts->Alhs; 626cda39b92SHong Zhang if (ctx) *ctx = ts->jacP; 627cda39b92SHong Zhang PetscFunctionReturn(0); 628cda39b92SHong Zhang } 629cda39b92SHong Zhang 630cda39b92SHong Zhang #undef __FUNCT__ 6314a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian" 632d763cef2SBarry Smith /*@C 633d763cef2SBarry Smith TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 634d763cef2SBarry Smith where U_t = F(U,t), as well as the location to store the matrix. 63576f2fa84SHong Zhang Use TSSetMatrices() for linear problems. 636d763cef2SBarry Smith 637d763cef2SBarry Smith Collective on TS 638d763cef2SBarry Smith 639d763cef2SBarry Smith Input Parameters: 640d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 641d763cef2SBarry Smith . A - Jacobian matrix 642d763cef2SBarry Smith . B - preconditioner matrix (usually same as A) 643d763cef2SBarry Smith . f - the Jacobian evaluation routine 644d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 645d763cef2SBarry Smith Jacobian evaluation routine (may be PETSC_NULL) 646d763cef2SBarry Smith 647d763cef2SBarry Smith Calling sequence of func: 64887828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 649d763cef2SBarry Smith 650d763cef2SBarry Smith + t - current timestep 651d763cef2SBarry Smith . u - input vector 652d763cef2SBarry Smith . A - matrix A, where U_t = A(t)u 653d763cef2SBarry Smith . B - preconditioner matrix, usually the same as A 654d763cef2SBarry Smith . flag - flag indicating information about the preconditioner matrix 65594b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 656d763cef2SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 657d763cef2SBarry Smith 658d763cef2SBarry Smith Notes: 65994b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 660d763cef2SBarry Smith output parameter in the routine func(). Be sure to read this information! 661d763cef2SBarry Smith 662d763cef2SBarry Smith The routine func() takes Mat * as the matrix arguments rather than Mat. 663d763cef2SBarry Smith This allows the matrix evaluation routine to replace A and/or B with a 66456335db2SHong Zhang completely new matrix structure (not just different matrix elements) 665d763cef2SBarry Smith when appropriate, for instance, if the nonzero structure is changing 666d763cef2SBarry Smith throughout the global iterations. 667d763cef2SBarry Smith 668d763cef2SBarry Smith Level: beginner 669d763cef2SBarry Smith 670d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian 671d763cef2SBarry Smith 672d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(), 67376f2fa84SHong Zhang SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 674d763cef2SBarry Smith 675d763cef2SBarry Smith @*/ 67663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 677d763cef2SBarry Smith { 678d763cef2SBarry Smith PetscFunctionBegin; 6790700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6800700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,2); 6810700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,3); 682c9780b6fSBarry Smith PetscCheckSameComm(ts,1,A,2); 683c9780b6fSBarry Smith PetscCheckSameComm(ts,1,B,3); 68417186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()"); 685d763cef2SBarry Smith 686000e7ae3SMatthew Knepley ts->ops->rhsjacobian = f; 687d763cef2SBarry Smith ts->jacP = ctx; 6888beb423aSHong Zhang ts->Arhs = A; 689d763cef2SBarry Smith ts->B = B; 690d763cef2SBarry Smith PetscFunctionReturn(0); 691d763cef2SBarry Smith } 692d763cef2SBarry Smith 693316643e7SJed Brown 694316643e7SJed Brown #undef __FUNCT__ 695316643e7SJed Brown #define __FUNCT__ "TSSetIFunction" 696316643e7SJed Brown /*@C 697316643e7SJed Brown TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 698316643e7SJed Brown 699316643e7SJed Brown Collective on TS 700316643e7SJed Brown 701316643e7SJed Brown Input Parameters: 702316643e7SJed Brown + ts - the TS context obtained from TSCreate() 703316643e7SJed Brown . f - the function evaluation routine 704316643e7SJed Brown - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 705316643e7SJed Brown 706316643e7SJed Brown Calling sequence of f: 707316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 708316643e7SJed Brown 709316643e7SJed Brown + t - time at step/stage being solved 710316643e7SJed Brown . u - state vector 711316643e7SJed Brown . u_t - time derivative of state vector 712316643e7SJed Brown . F - function vector 713316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 714316643e7SJed Brown 715316643e7SJed Brown Important: 716316643e7SJed Brown The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 717316643e7SJed Brown 718316643e7SJed Brown Level: beginner 719316643e7SJed Brown 720316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian 721316643e7SJed Brown 722316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 723316643e7SJed Brown @*/ 724316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS ts,TSIFunction f,void *ctx) 725316643e7SJed Brown { 726316643e7SJed Brown 727316643e7SJed Brown PetscFunctionBegin; 7280700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 729316643e7SJed Brown ts->ops->ifunction = f; 730316643e7SJed Brown ts->funP = ctx; 731316643e7SJed Brown PetscFunctionReturn(0); 732316643e7SJed Brown } 733316643e7SJed Brown 734316643e7SJed Brown #undef __FUNCT__ 735316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian" 736316643e7SJed Brown /*@C 737316643e7SJed Brown TSSetIJacobian - Set the function to compute the Jacobian of 738316643e7SJed Brown G(U) = F(t,U,U0+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 739316643e7SJed Brown 740316643e7SJed Brown Collective on TS 741316643e7SJed Brown 742316643e7SJed Brown Input Parameters: 743316643e7SJed Brown + ts - the TS context obtained from TSCreate() 744316643e7SJed Brown . A - Jacobian matrix 745316643e7SJed Brown . B - preconditioning matrix for A (may be same as A) 746316643e7SJed Brown . f - the Jacobian evaluation routine 747316643e7SJed Brown - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 748316643e7SJed Brown 749316643e7SJed Brown Calling sequence of f: 750316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 751316643e7SJed Brown 752316643e7SJed Brown + t - time at step/stage being solved 753316643e7SJed Brown . u - state vector 754316643e7SJed Brown . u_t - time derivative of state vector 755316643e7SJed Brown . a - shift 756316643e7SJed Brown . A - Jacobian of G(U) = F(t,U,U0+a*U), equivalent to dF/dU + a*dF/dU_t 757316643e7SJed Brown . B - preconditioning matrix for A, may be same as A 758316643e7SJed Brown . flag - flag indicating information about the preconditioner matrix 759316643e7SJed Brown structure (same as flag in KSPSetOperators()) 760316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 761316643e7SJed Brown 762316643e7SJed Brown Notes: 763316643e7SJed Brown The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 764316643e7SJed Brown 765316643e7SJed Brown Level: beginner 766316643e7SJed Brown 767316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian 768316643e7SJed Brown 769316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian() 770316643e7SJed Brown 771316643e7SJed Brown @*/ 772316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 773316643e7SJed Brown { 774316643e7SJed Brown PetscErrorCode ierr; 775316643e7SJed Brown 776316643e7SJed Brown PetscFunctionBegin; 7770700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7780700a824SBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 7790700a824SBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 780316643e7SJed Brown if (A) PetscCheckSameComm(ts,1,A,2); 781316643e7SJed Brown if (B) PetscCheckSameComm(ts,1,B,3); 782316643e7SJed Brown if (f) ts->ops->ijacobian = f; 783316643e7SJed Brown if (ctx) ts->jacP = ctx; 784316643e7SJed Brown if (A) { 785316643e7SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 786316643e7SJed Brown if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 787316643e7SJed Brown ts->A = A; 788316643e7SJed Brown } 789dc3f620dSJed Brown #if 0 790dc3f620dSJed Brown /* The sane and consistent alternative */ 791316643e7SJed Brown if (B) { 792316643e7SJed Brown ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 793316643e7SJed Brown if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);} 794316643e7SJed Brown ts->B = B; 795316643e7SJed Brown } 796dc3f620dSJed Brown #else 797dc3f620dSJed Brown /* Don't reference B because TSDestroy() doesn't destroy it. These ownership semantics are awkward and inconsistent. */ 798dc3f620dSJed Brown if (B) ts->B = B; 799dc3f620dSJed Brown #endif 800316643e7SJed Brown PetscFunctionReturn(0); 801316643e7SJed Brown } 802316643e7SJed Brown 8034a2ae208SSatish Balay #undef __FUNCT__ 8044a2ae208SSatish Balay #define __FUNCT__ "TSView" 8057e2c5f70SBarry Smith /*@C 806d763cef2SBarry Smith TSView - Prints the TS data structure. 807d763cef2SBarry Smith 8084c49b128SBarry Smith Collective on TS 809d763cef2SBarry Smith 810d763cef2SBarry Smith Input Parameters: 811d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 812d763cef2SBarry Smith - viewer - visualization context 813d763cef2SBarry Smith 814d763cef2SBarry Smith Options Database Key: 815d763cef2SBarry Smith . -ts_view - calls TSView() at end of TSStep() 816d763cef2SBarry Smith 817d763cef2SBarry Smith Notes: 818d763cef2SBarry Smith The available visualization contexts include 819b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 820b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 821d763cef2SBarry Smith output where only the first processor opens 822d763cef2SBarry Smith the file. All other processors send their 823d763cef2SBarry Smith data to the first processor to print. 824d763cef2SBarry Smith 825d763cef2SBarry Smith The user can open an alternative visualization context with 826b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 827d763cef2SBarry Smith 828d763cef2SBarry Smith Level: beginner 829d763cef2SBarry Smith 830d763cef2SBarry Smith .keywords: TS, timestep, view 831d763cef2SBarry Smith 832b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 833d763cef2SBarry Smith @*/ 83463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer) 835d763cef2SBarry Smith { 836dfbe8321SBarry Smith PetscErrorCode ierr; 837a313700dSBarry Smith const TSType type; 83832077d6dSBarry Smith PetscTruth iascii,isstring; 839d763cef2SBarry Smith 840d763cef2SBarry Smith PetscFunctionBegin; 8410700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8423050cee2SBarry Smith if (!viewer) { 8437adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 8443050cee2SBarry Smith } 8450700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 846c9780b6fSBarry Smith PetscCheckSameComm(ts,1,viewer,2); 847fd16b177SBarry Smith 84832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 849b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 85032077d6dSBarry Smith if (iascii) { 851b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr); 852a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 853454a90a3SBarry Smith if (type) { 854b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 855184914b5SBarry Smith } else { 856b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 857184914b5SBarry Smith } 858000e7ae3SMatthew Knepley if (ts->ops->view) { 859b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 860000e7ae3SMatthew Knepley ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 861b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 862d763cef2SBarry Smith } 86377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 864a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 865d763cef2SBarry Smith if (ts->problem_type == TS_NONLINEAR) { 86677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 867d763cef2SBarry Smith } 86877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 8690f5bd95cSBarry Smith } else if (isstring) { 870a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 871b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 872d763cef2SBarry Smith } 873b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 87494b7f48cSBarry Smith if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);} 875d763cef2SBarry Smith if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 876b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 877d763cef2SBarry Smith PetscFunctionReturn(0); 878d763cef2SBarry Smith } 879d763cef2SBarry Smith 880d763cef2SBarry Smith 8814a2ae208SSatish Balay #undef __FUNCT__ 8824a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext" 883d763cef2SBarry Smith /*@C 884d763cef2SBarry Smith TSSetApplicationContext - Sets an optional user-defined context for 885d763cef2SBarry Smith the timesteppers. 886d763cef2SBarry Smith 887d763cef2SBarry Smith Collective on TS 888d763cef2SBarry Smith 889d763cef2SBarry Smith Input Parameters: 890d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 891d763cef2SBarry Smith - usrP - optional user context 892d763cef2SBarry Smith 893d763cef2SBarry Smith Level: intermediate 894d763cef2SBarry Smith 895d763cef2SBarry Smith .keywords: TS, timestep, set, application, context 896d763cef2SBarry Smith 897d763cef2SBarry Smith .seealso: TSGetApplicationContext() 898d763cef2SBarry Smith @*/ 89963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP) 900d763cef2SBarry Smith { 901d763cef2SBarry Smith PetscFunctionBegin; 9020700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 903d763cef2SBarry Smith ts->user = usrP; 904d763cef2SBarry Smith PetscFunctionReturn(0); 905d763cef2SBarry Smith } 906d763cef2SBarry Smith 9074a2ae208SSatish Balay #undef __FUNCT__ 9084a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext" 909d763cef2SBarry Smith /*@C 910d763cef2SBarry Smith TSGetApplicationContext - Gets the user-defined context for the 911d763cef2SBarry Smith timestepper. 912d763cef2SBarry Smith 913d763cef2SBarry Smith Not Collective 914d763cef2SBarry Smith 915d763cef2SBarry Smith Input Parameter: 916d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 917d763cef2SBarry Smith 918d763cef2SBarry Smith Output Parameter: 919d763cef2SBarry Smith . usrP - user context 920d763cef2SBarry Smith 921d763cef2SBarry Smith Level: intermediate 922d763cef2SBarry Smith 923d763cef2SBarry Smith .keywords: TS, timestep, get, application, context 924d763cef2SBarry Smith 925d763cef2SBarry Smith .seealso: TSSetApplicationContext() 926d763cef2SBarry Smith @*/ 92763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP) 928d763cef2SBarry Smith { 929d763cef2SBarry Smith PetscFunctionBegin; 9300700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 931d763cef2SBarry Smith *usrP = ts->user; 932d763cef2SBarry Smith PetscFunctionReturn(0); 933d763cef2SBarry Smith } 934d763cef2SBarry Smith 9354a2ae208SSatish Balay #undef __FUNCT__ 9364a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber" 937d763cef2SBarry Smith /*@ 938d763cef2SBarry Smith TSGetTimeStepNumber - Gets the current number of timesteps. 939d763cef2SBarry Smith 940d763cef2SBarry Smith Not Collective 941d763cef2SBarry Smith 942d763cef2SBarry Smith Input Parameter: 943d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 944d763cef2SBarry Smith 945d763cef2SBarry Smith Output Parameter: 946d763cef2SBarry Smith . iter - number steps so far 947d763cef2SBarry Smith 948d763cef2SBarry Smith Level: intermediate 949d763cef2SBarry Smith 950d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number 951d763cef2SBarry Smith @*/ 95263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter) 953d763cef2SBarry Smith { 954d763cef2SBarry Smith PetscFunctionBegin; 9550700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9564482741eSBarry Smith PetscValidIntPointer(iter,2); 957d763cef2SBarry Smith *iter = ts->steps; 958d763cef2SBarry Smith PetscFunctionReturn(0); 959d763cef2SBarry Smith } 960d763cef2SBarry Smith 9614a2ae208SSatish Balay #undef __FUNCT__ 9624a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep" 963d763cef2SBarry Smith /*@ 964d763cef2SBarry Smith TSSetInitialTimeStep - Sets the initial timestep to be used, 965d763cef2SBarry Smith as well as the initial time. 966d763cef2SBarry Smith 967d763cef2SBarry Smith Collective on TS 968d763cef2SBarry Smith 969d763cef2SBarry Smith Input Parameters: 970d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 971d763cef2SBarry Smith . initial_time - the initial time 972d763cef2SBarry Smith - time_step - the size of the timestep 973d763cef2SBarry Smith 974d763cef2SBarry Smith Level: intermediate 975d763cef2SBarry Smith 976d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep() 977d763cef2SBarry Smith 978d763cef2SBarry Smith .keywords: TS, set, initial, timestep 979d763cef2SBarry Smith @*/ 98063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 981d763cef2SBarry Smith { 982d763cef2SBarry Smith PetscFunctionBegin; 9830700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 984d763cef2SBarry Smith ts->time_step = time_step; 985d763cef2SBarry Smith ts->initial_time_step = time_step; 986d763cef2SBarry Smith ts->ptime = initial_time; 987d763cef2SBarry Smith PetscFunctionReturn(0); 988d763cef2SBarry Smith } 989d763cef2SBarry Smith 9904a2ae208SSatish Balay #undef __FUNCT__ 9914a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep" 992d763cef2SBarry Smith /*@ 993d763cef2SBarry Smith TSSetTimeStep - Allows one to reset the timestep at any time, 994d763cef2SBarry Smith useful for simple pseudo-timestepping codes. 995d763cef2SBarry Smith 996d763cef2SBarry Smith Collective on TS 997d763cef2SBarry Smith 998d763cef2SBarry Smith Input Parameters: 999d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1000d763cef2SBarry Smith - time_step - the size of the timestep 1001d763cef2SBarry Smith 1002d763cef2SBarry Smith Level: intermediate 1003d763cef2SBarry Smith 1004d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1005d763cef2SBarry Smith 1006d763cef2SBarry Smith .keywords: TS, set, timestep 1007d763cef2SBarry Smith @*/ 100863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step) 1009d763cef2SBarry Smith { 1010d763cef2SBarry Smith PetscFunctionBegin; 10110700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1012d763cef2SBarry Smith ts->time_step = time_step; 1013d763cef2SBarry Smith PetscFunctionReturn(0); 1014d763cef2SBarry Smith } 1015d763cef2SBarry Smith 10164a2ae208SSatish Balay #undef __FUNCT__ 10174a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep" 1018d763cef2SBarry Smith /*@ 1019d763cef2SBarry Smith TSGetTimeStep - Gets the current timestep size. 1020d763cef2SBarry Smith 1021d763cef2SBarry Smith Not Collective 1022d763cef2SBarry Smith 1023d763cef2SBarry Smith Input Parameter: 1024d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1025d763cef2SBarry Smith 1026d763cef2SBarry Smith Output Parameter: 1027d763cef2SBarry Smith . dt - the current timestep size 1028d763cef2SBarry Smith 1029d763cef2SBarry Smith Level: intermediate 1030d763cef2SBarry Smith 1031d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1032d763cef2SBarry Smith 1033d763cef2SBarry Smith .keywords: TS, get, timestep 1034d763cef2SBarry Smith @*/ 103563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt) 1036d763cef2SBarry Smith { 1037d763cef2SBarry Smith PetscFunctionBegin; 10380700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10394482741eSBarry Smith PetscValidDoublePointer(dt,2); 1040d763cef2SBarry Smith *dt = ts->time_step; 1041d763cef2SBarry Smith PetscFunctionReturn(0); 1042d763cef2SBarry Smith } 1043d763cef2SBarry Smith 10444a2ae208SSatish Balay #undef __FUNCT__ 10454a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution" 1046d8e5e3e6SSatish Balay /*@ 1047d763cef2SBarry Smith TSGetSolution - Returns the solution at the present timestep. It 1048d763cef2SBarry Smith is valid to call this routine inside the function that you are evaluating 1049d763cef2SBarry Smith in order to move to the new timestep. This vector not changed until 1050d763cef2SBarry Smith the solution at the next timestep has been calculated. 1051d763cef2SBarry Smith 1052d763cef2SBarry Smith Not Collective, but Vec returned is parallel if TS is parallel 1053d763cef2SBarry Smith 1054d763cef2SBarry Smith Input Parameter: 1055d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1056d763cef2SBarry Smith 1057d763cef2SBarry Smith Output Parameter: 1058d763cef2SBarry Smith . v - the vector containing the solution 1059d763cef2SBarry Smith 1060d763cef2SBarry Smith Level: intermediate 1061d763cef2SBarry Smith 1062d763cef2SBarry Smith .seealso: TSGetTimeStep() 1063d763cef2SBarry Smith 1064d763cef2SBarry Smith .keywords: TS, timestep, get, solution 1065d763cef2SBarry Smith @*/ 106663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v) 1067d763cef2SBarry Smith { 1068d763cef2SBarry Smith PetscFunctionBegin; 10690700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10704482741eSBarry Smith PetscValidPointer(v,2); 1071d763cef2SBarry Smith *v = ts->vec_sol_always; 1072d763cef2SBarry Smith PetscFunctionReturn(0); 1073d763cef2SBarry Smith } 1074d763cef2SBarry Smith 1075bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */ 10764a2ae208SSatish Balay #undef __FUNCT__ 1077bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType" 1078d8e5e3e6SSatish Balay /*@ 1079bdad233fSMatthew Knepley TSSetProblemType - Sets the type of problem to be solved. 1080d763cef2SBarry Smith 1081bdad233fSMatthew Knepley Not collective 1082d763cef2SBarry Smith 1083bdad233fSMatthew Knepley Input Parameters: 1084bdad233fSMatthew Knepley + ts - The TS 1085bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1086d763cef2SBarry Smith .vb 1087d763cef2SBarry Smith U_t = A U 1088d763cef2SBarry Smith U_t = A(t) U 1089d763cef2SBarry Smith U_t = F(t,U) 1090d763cef2SBarry Smith .ve 1091d763cef2SBarry Smith 1092d763cef2SBarry Smith Level: beginner 1093d763cef2SBarry Smith 1094bdad233fSMatthew Knepley .keywords: TS, problem type 1095bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1096d763cef2SBarry Smith @*/ 109763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type) 1098a7cc72afSBarry Smith { 1099d763cef2SBarry Smith PetscFunctionBegin; 11000700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1101bdad233fSMatthew Knepley ts->problem_type = type; 1102d763cef2SBarry Smith PetscFunctionReturn(0); 1103d763cef2SBarry Smith } 1104d763cef2SBarry Smith 1105bdad233fSMatthew Knepley #undef __FUNCT__ 1106bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType" 1107bdad233fSMatthew Knepley /*@C 1108bdad233fSMatthew Knepley TSGetProblemType - Gets the type of problem to be solved. 1109bdad233fSMatthew Knepley 1110bdad233fSMatthew Knepley Not collective 1111bdad233fSMatthew Knepley 1112bdad233fSMatthew Knepley Input Parameter: 1113bdad233fSMatthew Knepley . ts - The TS 1114bdad233fSMatthew Knepley 1115bdad233fSMatthew Knepley Output Parameter: 1116bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1117bdad233fSMatthew Knepley .vb 1118bdad233fSMatthew Knepley U_t = A U 1119bdad233fSMatthew Knepley U_t = A(t) U 1120bdad233fSMatthew Knepley U_t = F(t,U) 1121bdad233fSMatthew Knepley .ve 1122bdad233fSMatthew Knepley 1123bdad233fSMatthew Knepley Level: beginner 1124bdad233fSMatthew Knepley 1125bdad233fSMatthew Knepley .keywords: TS, problem type 1126bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1127bdad233fSMatthew Knepley @*/ 112863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type) 1129a7cc72afSBarry Smith { 1130bdad233fSMatthew Knepley PetscFunctionBegin; 11310700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 11324482741eSBarry Smith PetscValidIntPointer(type,2); 1133bdad233fSMatthew Knepley *type = ts->problem_type; 1134bdad233fSMatthew Knepley PetscFunctionReturn(0); 1135bdad233fSMatthew Knepley } 1136d763cef2SBarry Smith 11374a2ae208SSatish Balay #undef __FUNCT__ 11384a2ae208SSatish Balay #define __FUNCT__ "TSSetUp" 1139d763cef2SBarry Smith /*@ 1140d763cef2SBarry Smith TSSetUp - Sets up the internal data structures for the later use 1141d763cef2SBarry Smith of a timestepper. 1142d763cef2SBarry Smith 1143d763cef2SBarry Smith Collective on TS 1144d763cef2SBarry Smith 1145d763cef2SBarry Smith Input Parameter: 1146d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1147d763cef2SBarry Smith 1148d763cef2SBarry Smith Notes: 1149d763cef2SBarry Smith For basic use of the TS solvers the user need not explicitly call 1150d763cef2SBarry Smith TSSetUp(), since these actions will automatically occur during 1151d763cef2SBarry Smith the call to TSStep(). However, if one wishes to control this 1152d763cef2SBarry Smith phase separately, TSSetUp() should be called after TSCreate() 1153d763cef2SBarry Smith and optional routines of the form TSSetXXX(), but before TSStep(). 1154d763cef2SBarry Smith 1155d763cef2SBarry Smith Level: advanced 1156d763cef2SBarry Smith 1157d763cef2SBarry Smith .keywords: TS, timestep, setup 1158d763cef2SBarry Smith 1159d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy() 1160d763cef2SBarry Smith @*/ 116163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts) 1162d763cef2SBarry Smith { 1163dfbe8321SBarry Smith PetscErrorCode ierr; 1164d763cef2SBarry Smith 1165d763cef2SBarry Smith PetscFunctionBegin; 11660700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1167e32f2f54SBarry Smith if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 11687adad957SLisandro Dalcin if (!((PetscObject)ts)->type_name) { 11699596e0b4SJed Brown ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1170d763cef2SBarry Smith } 1171000e7ae3SMatthew Knepley ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1172d763cef2SBarry Smith ts->setupcalled = 1; 1173d763cef2SBarry Smith PetscFunctionReturn(0); 1174d763cef2SBarry Smith } 1175d763cef2SBarry Smith 11764a2ae208SSatish Balay #undef __FUNCT__ 11774a2ae208SSatish Balay #define __FUNCT__ "TSDestroy" 1178d8e5e3e6SSatish Balay /*@ 1179d763cef2SBarry Smith TSDestroy - Destroys the timestepper context that was created 1180d763cef2SBarry Smith with TSCreate(). 1181d763cef2SBarry Smith 1182d763cef2SBarry Smith Collective on TS 1183d763cef2SBarry Smith 1184d763cef2SBarry Smith Input Parameter: 1185d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1186d763cef2SBarry Smith 1187d763cef2SBarry Smith Level: beginner 1188d763cef2SBarry Smith 1189d763cef2SBarry Smith .keywords: TS, timestepper, destroy 1190d763cef2SBarry Smith 1191d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve() 1192d763cef2SBarry Smith @*/ 119363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts) 1194d763cef2SBarry Smith { 11956849ba73SBarry Smith PetscErrorCode ierr; 1196d763cef2SBarry Smith 1197d763cef2SBarry Smith PetscFunctionBegin; 11980700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 11997adad957SLisandro Dalcin if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0); 1200d763cef2SBarry Smith 1201be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 12020f5bd95cSBarry Smith ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 12036c699258SBarry Smith 12046c699258SBarry Smith if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);} 1205cb9801acSJed Brown if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 120694b7f48cSBarry Smith if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);} 1207d763cef2SBarry Smith if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 12081e3347e8SBarry Smith if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);} 1209a6570f20SBarry Smith ierr = TSMonitorCancel(ts);CHKERRQ(ierr); 1210a79aaaedSSatish Balay ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1211d763cef2SBarry Smith PetscFunctionReturn(0); 1212d763cef2SBarry Smith } 1213d763cef2SBarry Smith 12144a2ae208SSatish Balay #undef __FUNCT__ 12154a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES" 1216d8e5e3e6SSatish Balay /*@ 1217d763cef2SBarry Smith TSGetSNES - Returns the SNES (nonlinear solver) associated with 1218d763cef2SBarry Smith a TS (timestepper) context. Valid only for nonlinear problems. 1219d763cef2SBarry Smith 1220d763cef2SBarry Smith Not Collective, but SNES is parallel if TS is parallel 1221d763cef2SBarry Smith 1222d763cef2SBarry Smith Input Parameter: 1223d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1224d763cef2SBarry Smith 1225d763cef2SBarry Smith Output Parameter: 1226d763cef2SBarry Smith . snes - the nonlinear solver context 1227d763cef2SBarry Smith 1228d763cef2SBarry Smith Notes: 1229d763cef2SBarry Smith The user can then directly manipulate the SNES context to set various 1230d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 123194b7f48cSBarry Smith KSP, KSP, and PC contexts as well. 1232d763cef2SBarry Smith 1233d763cef2SBarry Smith TSGetSNES() does not work for integrators that do not use SNES; in 1234d763cef2SBarry Smith this case TSGetSNES() returns PETSC_NULL in snes. 1235d763cef2SBarry Smith 1236d763cef2SBarry Smith Level: beginner 1237d763cef2SBarry Smith 1238d763cef2SBarry Smith .keywords: timestep, get, SNES 1239d763cef2SBarry Smith @*/ 124063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes) 1241d763cef2SBarry Smith { 1242d763cef2SBarry Smith PetscFunctionBegin; 12430700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12444482741eSBarry Smith PetscValidPointer(snes,2); 124517186662SBarry Smith if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first"); 1246e32f2f54SBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1247d763cef2SBarry Smith *snes = ts->snes; 1248d763cef2SBarry Smith PetscFunctionReturn(0); 1249d763cef2SBarry Smith } 1250d763cef2SBarry Smith 12514a2ae208SSatish Balay #undef __FUNCT__ 125294b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP" 1253d8e5e3e6SSatish Balay /*@ 125494b7f48cSBarry Smith TSGetKSP - Returns the KSP (linear solver) associated with 1255d763cef2SBarry Smith a TS (timestepper) context. 1256d763cef2SBarry Smith 125794b7f48cSBarry Smith Not Collective, but KSP is parallel if TS is parallel 1258d763cef2SBarry Smith 1259d763cef2SBarry Smith Input Parameter: 1260d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1261d763cef2SBarry Smith 1262d763cef2SBarry Smith Output Parameter: 126394b7f48cSBarry Smith . ksp - the nonlinear solver context 1264d763cef2SBarry Smith 1265d763cef2SBarry Smith Notes: 126694b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 1267d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 1268d763cef2SBarry Smith KSP and PC contexts as well. 1269d763cef2SBarry Smith 127094b7f48cSBarry Smith TSGetKSP() does not work for integrators that do not use KSP; 127194b7f48cSBarry Smith in this case TSGetKSP() returns PETSC_NULL in ksp. 1272d763cef2SBarry Smith 1273d763cef2SBarry Smith Level: beginner 1274d763cef2SBarry Smith 127594b7f48cSBarry Smith .keywords: timestep, get, KSP 1276d763cef2SBarry Smith @*/ 127763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp) 1278d763cef2SBarry Smith { 1279d763cef2SBarry Smith PetscFunctionBegin; 12800700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12814482741eSBarry Smith PetscValidPointer(ksp,2); 128217186662SBarry Smith if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1283e32f2f54SBarry Smith if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 128494b7f48cSBarry Smith *ksp = ts->ksp; 1285d763cef2SBarry Smith PetscFunctionReturn(0); 1286d763cef2SBarry Smith } 1287d763cef2SBarry Smith 1288d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */ 1289d763cef2SBarry Smith 12904a2ae208SSatish Balay #undef __FUNCT__ 1291adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration" 1292adb62b0dSMatthew Knepley /*@ 1293adb62b0dSMatthew Knepley TSGetDuration - Gets the maximum number of timesteps to use and 1294adb62b0dSMatthew Knepley maximum time for iteration. 1295adb62b0dSMatthew Knepley 1296adb62b0dSMatthew Knepley Collective on TS 1297adb62b0dSMatthew Knepley 1298adb62b0dSMatthew Knepley Input Parameters: 1299adb62b0dSMatthew Knepley + ts - the TS context obtained from TSCreate() 1300adb62b0dSMatthew Knepley . maxsteps - maximum number of iterations to use, or PETSC_NULL 1301adb62b0dSMatthew Knepley - maxtime - final time to iterate to, or PETSC_NULL 1302adb62b0dSMatthew Knepley 1303adb62b0dSMatthew Knepley Level: intermediate 1304adb62b0dSMatthew Knepley 1305adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time 1306adb62b0dSMatthew Knepley @*/ 130763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1308adb62b0dSMatthew Knepley { 1309adb62b0dSMatthew Knepley PetscFunctionBegin; 13100700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1311abc0a331SBarry Smith if (maxsteps) { 13124482741eSBarry Smith PetscValidIntPointer(maxsteps,2); 1313adb62b0dSMatthew Knepley *maxsteps = ts->max_steps; 1314adb62b0dSMatthew Knepley } 1315abc0a331SBarry Smith if (maxtime ) { 13164482741eSBarry Smith PetscValidScalarPointer(maxtime,3); 1317adb62b0dSMatthew Knepley *maxtime = ts->max_time; 1318adb62b0dSMatthew Knepley } 1319adb62b0dSMatthew Knepley PetscFunctionReturn(0); 1320adb62b0dSMatthew Knepley } 1321adb62b0dSMatthew Knepley 1322adb62b0dSMatthew Knepley #undef __FUNCT__ 13234a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration" 1324d763cef2SBarry Smith /*@ 1325d763cef2SBarry Smith TSSetDuration - Sets the maximum number of timesteps to use and 1326d763cef2SBarry Smith maximum time for iteration. 1327d763cef2SBarry Smith 1328d763cef2SBarry Smith Collective on TS 1329d763cef2SBarry Smith 1330d763cef2SBarry Smith Input Parameters: 1331d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1332d763cef2SBarry Smith . maxsteps - maximum number of iterations to use 1333d763cef2SBarry Smith - maxtime - final time to iterate to 1334d763cef2SBarry Smith 1335d763cef2SBarry Smith Options Database Keys: 1336d763cef2SBarry Smith . -ts_max_steps <maxsteps> - Sets maxsteps 1337d763cef2SBarry Smith . -ts_max_time <maxtime> - Sets maxtime 1338d763cef2SBarry Smith 1339d763cef2SBarry Smith Notes: 1340d763cef2SBarry Smith The default maximum number of iterations is 5000. Default time is 5.0 1341d763cef2SBarry Smith 1342d763cef2SBarry Smith Level: intermediate 1343d763cef2SBarry Smith 1344d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations 1345d763cef2SBarry Smith @*/ 134663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1347d763cef2SBarry Smith { 1348d763cef2SBarry Smith PetscFunctionBegin; 13490700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1350d763cef2SBarry Smith ts->max_steps = maxsteps; 1351d763cef2SBarry Smith ts->max_time = maxtime; 1352d763cef2SBarry Smith PetscFunctionReturn(0); 1353d763cef2SBarry Smith } 1354d763cef2SBarry Smith 13554a2ae208SSatish Balay #undef __FUNCT__ 13564a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution" 1357d763cef2SBarry Smith /*@ 1358d763cef2SBarry Smith TSSetSolution - Sets the initial solution vector 1359d763cef2SBarry Smith for use by the TS routines. 1360d763cef2SBarry Smith 1361d763cef2SBarry Smith Collective on TS and Vec 1362d763cef2SBarry Smith 1363d763cef2SBarry Smith Input Parameters: 1364d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1365d763cef2SBarry Smith - x - the solution vector 1366d763cef2SBarry Smith 1367d763cef2SBarry Smith Level: beginner 1368d763cef2SBarry Smith 1369d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions 1370d763cef2SBarry Smith @*/ 137163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x) 1372d763cef2SBarry Smith { 1373d763cef2SBarry Smith PetscFunctionBegin; 13740700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13750700a824SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1376d763cef2SBarry Smith ts->vec_sol = ts->vec_sol_always = x; 1377d763cef2SBarry Smith PetscFunctionReturn(0); 1378d763cef2SBarry Smith } 1379d763cef2SBarry Smith 1380e74ef692SMatthew Knepley #undef __FUNCT__ 1381e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep" 1382ac226902SBarry Smith /*@C 1383000e7ae3SMatthew Knepley TSSetPreStep - Sets the general-purpose function 13843f2090d5SJed Brown called once at the beginning of each time step. 1385000e7ae3SMatthew Knepley 1386000e7ae3SMatthew Knepley Collective on TS 1387000e7ae3SMatthew Knepley 1388000e7ae3SMatthew Knepley Input Parameters: 1389000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1390000e7ae3SMatthew Knepley - func - The function 1391000e7ae3SMatthew Knepley 1392000e7ae3SMatthew Knepley Calling sequence of func: 1393000e7ae3SMatthew Knepley . func (TS ts); 1394000e7ae3SMatthew Knepley 1395000e7ae3SMatthew Knepley Level: intermediate 1396000e7ae3SMatthew Knepley 1397000e7ae3SMatthew Knepley .keywords: TS, timestep 1398000e7ae3SMatthew Knepley @*/ 139963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1400000e7ae3SMatthew Knepley { 1401000e7ae3SMatthew Knepley PetscFunctionBegin; 14020700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1403000e7ae3SMatthew Knepley ts->ops->prestep = func; 1404000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1405000e7ae3SMatthew Knepley } 1406000e7ae3SMatthew Knepley 1407e74ef692SMatthew Knepley #undef __FUNCT__ 14083f2090d5SJed Brown #define __FUNCT__ "TSPreStep" 14093f2090d5SJed Brown /*@C 14103f2090d5SJed Brown TSPreStep - Runs the user-defined pre-step function. 14113f2090d5SJed Brown 14123f2090d5SJed Brown Collective on TS 14133f2090d5SJed Brown 14143f2090d5SJed Brown Input Parameters: 14153f2090d5SJed Brown . ts - The TS context obtained from TSCreate() 14163f2090d5SJed Brown 14173f2090d5SJed Brown Notes: 14183f2090d5SJed Brown TSPreStep() is typically used within time stepping implementations, 14193f2090d5SJed Brown so most users would not generally call this routine themselves. 14203f2090d5SJed Brown 14213f2090d5SJed Brown Level: developer 14223f2090d5SJed Brown 14233f2090d5SJed Brown .keywords: TS, timestep 14243f2090d5SJed Brown @*/ 14253f2090d5SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSPreStep(TS ts) 14263f2090d5SJed Brown { 14273f2090d5SJed Brown PetscErrorCode ierr; 14283f2090d5SJed Brown 14293f2090d5SJed Brown PetscFunctionBegin; 14300700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 143172ac3e02SJed Brown if (ts->ops->prestep) { 14323f2090d5SJed Brown PetscStackPush("TS PreStep function"); 14333f2090d5SJed Brown CHKMEMQ; 14343f2090d5SJed Brown ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 14353f2090d5SJed Brown CHKMEMQ; 14363f2090d5SJed Brown PetscStackPop; 1437312ce896SJed Brown } 14383f2090d5SJed Brown PetscFunctionReturn(0); 14393f2090d5SJed Brown } 14403f2090d5SJed Brown 14413f2090d5SJed Brown #undef __FUNCT__ 1442e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep" 1443000e7ae3SMatthew Knepley /*@ 1444000e7ae3SMatthew Knepley TSDefaultPreStep - The default pre-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 TSDefaultPreStep(TS ts) 1456000e7ae3SMatthew Knepley { 1457000e7ae3SMatthew Knepley PetscFunctionBegin; 1458000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1459000e7ae3SMatthew Knepley } 1460000e7ae3SMatthew Knepley 1461e74ef692SMatthew Knepley #undef __FUNCT__ 1462e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep" 1463ac226902SBarry Smith /*@C 1464000e7ae3SMatthew Knepley TSSetPostStep - Sets the general-purpose function 14653f2090d5SJed Brown called once at the end of each time step. 1466000e7ae3SMatthew Knepley 1467000e7ae3SMatthew Knepley Collective on TS 1468000e7ae3SMatthew Knepley 1469000e7ae3SMatthew Knepley Input Parameters: 1470000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1471000e7ae3SMatthew Knepley - func - The function 1472000e7ae3SMatthew Knepley 1473000e7ae3SMatthew Knepley Calling sequence of func: 1474000e7ae3SMatthew Knepley . func (TS ts); 1475000e7ae3SMatthew Knepley 1476000e7ae3SMatthew Knepley Level: intermediate 1477000e7ae3SMatthew Knepley 1478000e7ae3SMatthew Knepley .keywords: TS, timestep 1479000e7ae3SMatthew Knepley @*/ 148063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1481000e7ae3SMatthew Knepley { 1482000e7ae3SMatthew Knepley PetscFunctionBegin; 14830700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1484000e7ae3SMatthew Knepley ts->ops->poststep = func; 1485000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1486000e7ae3SMatthew Knepley } 1487000e7ae3SMatthew Knepley 1488e74ef692SMatthew Knepley #undef __FUNCT__ 14893f2090d5SJed Brown #define __FUNCT__ "TSPostStep" 14903f2090d5SJed Brown /*@C 14913f2090d5SJed Brown TSPostStep - Runs the user-defined post-step function. 14923f2090d5SJed Brown 14933f2090d5SJed Brown Collective on TS 14943f2090d5SJed Brown 14953f2090d5SJed Brown Input Parameters: 14963f2090d5SJed Brown . ts - The TS context obtained from TSCreate() 14973f2090d5SJed Brown 14983f2090d5SJed Brown Notes: 14993f2090d5SJed Brown TSPostStep() is typically used within time stepping implementations, 15003f2090d5SJed Brown so most users would not generally call this routine themselves. 15013f2090d5SJed Brown 15023f2090d5SJed Brown Level: developer 15033f2090d5SJed Brown 15043f2090d5SJed Brown .keywords: TS, timestep 15053f2090d5SJed Brown @*/ 15063f2090d5SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSPostStep(TS ts) 15073f2090d5SJed Brown { 15083f2090d5SJed Brown PetscErrorCode ierr; 15093f2090d5SJed Brown 15103f2090d5SJed Brown PetscFunctionBegin; 15110700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 151272ac3e02SJed Brown if (ts->ops->poststep) { 15133f2090d5SJed Brown PetscStackPush("TS PostStep function"); 15143f2090d5SJed Brown CHKMEMQ; 15153f2090d5SJed Brown ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 15163f2090d5SJed Brown CHKMEMQ; 15173f2090d5SJed Brown PetscStackPop; 151872ac3e02SJed Brown } 15193f2090d5SJed Brown PetscFunctionReturn(0); 15203f2090d5SJed Brown } 15213f2090d5SJed Brown 15223f2090d5SJed Brown #undef __FUNCT__ 1523e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep" 1524000e7ae3SMatthew Knepley /*@ 1525000e7ae3SMatthew Knepley TSDefaultPostStep - The default post-stepping function which does nothing. 1526000e7ae3SMatthew Knepley 1527000e7ae3SMatthew Knepley Collective on TS 1528000e7ae3SMatthew Knepley 1529000e7ae3SMatthew Knepley Input Parameters: 1530000e7ae3SMatthew Knepley . ts - The TS context obtained from TSCreate() 1531000e7ae3SMatthew Knepley 1532000e7ae3SMatthew Knepley Level: developer 1533000e7ae3SMatthew Knepley 1534000e7ae3SMatthew Knepley .keywords: TS, timestep 1535000e7ae3SMatthew Knepley @*/ 153663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts) 1537000e7ae3SMatthew Knepley { 1538000e7ae3SMatthew Knepley PetscFunctionBegin; 1539000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1540000e7ae3SMatthew Knepley } 1541000e7ae3SMatthew Knepley 1542d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 1543d763cef2SBarry Smith 15444a2ae208SSatish Balay #undef __FUNCT__ 1545a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet" 1546d763cef2SBarry Smith /*@C 1547a6570f20SBarry Smith TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1548d763cef2SBarry Smith timestep to display the iteration's progress. 1549d763cef2SBarry Smith 1550d763cef2SBarry Smith Collective on TS 1551d763cef2SBarry Smith 1552d763cef2SBarry Smith Input Parameters: 1553d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1554d763cef2SBarry Smith . func - monitoring routine 1555329f5518SBarry Smith . mctx - [optional] user-defined context for private data for the 1556b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desired) 1557b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1558b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 1559d763cef2SBarry Smith 1560d763cef2SBarry Smith Calling sequence of func: 1561a7cc72afSBarry Smith $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1562d763cef2SBarry Smith 1563d763cef2SBarry Smith + ts - the TS context 1564d763cef2SBarry Smith . steps - iteration number 15651f06c33eSBarry Smith . time - current time 1566d763cef2SBarry Smith . x - current iterate 1567d763cef2SBarry Smith - mctx - [optional] monitoring context 1568d763cef2SBarry Smith 1569d763cef2SBarry Smith Notes: 1570d763cef2SBarry Smith This routine adds an additional monitor to the list of monitors that 1571d763cef2SBarry Smith already has been loaded. 1572d763cef2SBarry Smith 1573025f1a04SBarry Smith Fortran notes: Only a single monitor function can be set for each TS object 1574025f1a04SBarry Smith 1575d763cef2SBarry Smith Level: intermediate 1576d763cef2SBarry Smith 1577d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1578d763cef2SBarry Smith 1579a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel() 1580d763cef2SBarry Smith @*/ 1581a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1582d763cef2SBarry Smith { 1583d763cef2SBarry Smith PetscFunctionBegin; 15840700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 158517186662SBarry Smith if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1586d763cef2SBarry Smith ts->monitor[ts->numbermonitors] = monitor; 1587329f5518SBarry Smith ts->mdestroy[ts->numbermonitors] = mdestroy; 1588d763cef2SBarry Smith ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1589d763cef2SBarry Smith PetscFunctionReturn(0); 1590d763cef2SBarry Smith } 1591d763cef2SBarry Smith 15924a2ae208SSatish Balay #undef __FUNCT__ 1593a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel" 1594d763cef2SBarry Smith /*@C 1595a6570f20SBarry Smith TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1596d763cef2SBarry Smith 1597d763cef2SBarry Smith Collective on TS 1598d763cef2SBarry Smith 1599d763cef2SBarry Smith Input Parameters: 1600d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1601d763cef2SBarry Smith 1602d763cef2SBarry Smith Notes: 1603d763cef2SBarry Smith There is no way to remove a single, specific monitor. 1604d763cef2SBarry Smith 1605d763cef2SBarry Smith Level: intermediate 1606d763cef2SBarry Smith 1607d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1608d763cef2SBarry Smith 1609a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet() 1610d763cef2SBarry Smith @*/ 1611a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts) 1612d763cef2SBarry Smith { 1613d952e501SBarry Smith PetscErrorCode ierr; 1614d952e501SBarry Smith PetscInt i; 1615d952e501SBarry Smith 1616d763cef2SBarry Smith PetscFunctionBegin; 16170700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1618d952e501SBarry Smith for (i=0; i<ts->numbermonitors; i++) { 1619d952e501SBarry Smith if (ts->mdestroy[i]) { 1620d952e501SBarry Smith ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 1621d952e501SBarry Smith } 1622d952e501SBarry Smith } 1623d763cef2SBarry Smith ts->numbermonitors = 0; 1624d763cef2SBarry Smith PetscFunctionReturn(0); 1625d763cef2SBarry Smith } 1626d763cef2SBarry Smith 16274a2ae208SSatish Balay #undef __FUNCT__ 1628a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault" 1629d8e5e3e6SSatish Balay /*@ 1630a6570f20SBarry Smith TSMonitorDefault - Sets the Default monitor 16315516499fSSatish Balay 16325516499fSSatish Balay Level: intermediate 163341251cbbSSatish Balay 16345516499fSSatish Balay .keywords: TS, set, monitor 16355516499fSSatish Balay 163641251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet() 163741251cbbSSatish Balay @*/ 1638a6570f20SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1639d763cef2SBarry Smith { 1640dfbe8321SBarry Smith PetscErrorCode ierr; 1641a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 1642d132466eSBarry Smith 1643d763cef2SBarry Smith PetscFunctionBegin; 1644a34d58ebSBarry Smith if (!ctx) { 16457adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 1646a34d58ebSBarry Smith } 1647f22f69f0SBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1648a34d58ebSBarry Smith if (!ctx) { 1649a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 1650a34d58ebSBarry Smith } 1651d763cef2SBarry Smith PetscFunctionReturn(0); 1652d763cef2SBarry Smith } 1653d763cef2SBarry Smith 16544a2ae208SSatish Balay #undef __FUNCT__ 16554a2ae208SSatish Balay #define __FUNCT__ "TSStep" 1656d763cef2SBarry Smith /*@ 1657d763cef2SBarry Smith TSStep - Steps the requested number of timesteps. 1658d763cef2SBarry Smith 1659d763cef2SBarry Smith Collective on TS 1660d763cef2SBarry Smith 1661d763cef2SBarry Smith Input Parameter: 1662d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1663d763cef2SBarry Smith 1664d763cef2SBarry Smith Output Parameters: 1665d763cef2SBarry Smith + steps - number of iterations until termination 1666142b95e3SSatish Balay - ptime - time until termination 1667d763cef2SBarry Smith 1668d763cef2SBarry Smith Level: beginner 1669d763cef2SBarry Smith 1670d763cef2SBarry Smith .keywords: TS, timestep, solve 1671d763cef2SBarry Smith 1672d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy() 1673d763cef2SBarry Smith @*/ 167463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1675d763cef2SBarry Smith { 1676dfbe8321SBarry Smith PetscErrorCode ierr; 1677d763cef2SBarry Smith 1678d763cef2SBarry Smith PetscFunctionBegin; 16790700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1680d405a339SMatthew Knepley if (!ts->setupcalled) { 1681d405a339SMatthew Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 1682d405a339SMatthew Knepley } 1683d405a339SMatthew Knepley 1684d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1685000e7ae3SMatthew Knepley ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1686d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1687d405a339SMatthew Knepley 16884bb05414SBarry Smith if (!PetscPreLoadingOn) { 16897adad957SLisandro Dalcin ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1690d405a339SMatthew Knepley } 1691d763cef2SBarry Smith PetscFunctionReturn(0); 1692d763cef2SBarry Smith } 1693d763cef2SBarry Smith 16944a2ae208SSatish Balay #undef __FUNCT__ 16956a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve" 16966a4d4014SLisandro Dalcin /*@ 16976a4d4014SLisandro Dalcin TSSolve - Steps the requested number of timesteps. 16986a4d4014SLisandro Dalcin 16996a4d4014SLisandro Dalcin Collective on TS 17006a4d4014SLisandro Dalcin 17016a4d4014SLisandro Dalcin Input Parameter: 17026a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 17036a4d4014SLisandro Dalcin - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 17046a4d4014SLisandro Dalcin 17056a4d4014SLisandro Dalcin Level: beginner 17066a4d4014SLisandro Dalcin 17076a4d4014SLisandro Dalcin .keywords: TS, timestep, solve 17086a4d4014SLisandro Dalcin 17096a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep() 17106a4d4014SLisandro Dalcin @*/ 17116a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x) 17126a4d4014SLisandro Dalcin { 17136a4d4014SLisandro Dalcin PetscInt steps; 17146a4d4014SLisandro Dalcin PetscReal ptime; 17156a4d4014SLisandro Dalcin PetscErrorCode ierr; 1716f22f69f0SBarry Smith 17176a4d4014SLisandro Dalcin PetscFunctionBegin; 17180700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17196a4d4014SLisandro Dalcin /* set solution vector if provided */ 17206a4d4014SLisandro Dalcin if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 17216a4d4014SLisandro Dalcin /* reset time step and iteration counters */ 17226a4d4014SLisandro Dalcin ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 17236a4d4014SLisandro Dalcin /* steps the requested number of timesteps. */ 17246a4d4014SLisandro Dalcin ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 17256a4d4014SLisandro Dalcin PetscFunctionReturn(0); 17266a4d4014SLisandro Dalcin } 17276a4d4014SLisandro Dalcin 17286a4d4014SLisandro Dalcin #undef __FUNCT__ 17294a2ae208SSatish Balay #define __FUNCT__ "TSMonitor" 1730d763cef2SBarry Smith /* 1731d763cef2SBarry Smith Runs the user provided monitor routines, if they exists. 1732d763cef2SBarry Smith */ 1733a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1734d763cef2SBarry Smith { 17356849ba73SBarry Smith PetscErrorCode ierr; 1736a7cc72afSBarry Smith PetscInt i,n = ts->numbermonitors; 1737d763cef2SBarry Smith 1738d763cef2SBarry Smith PetscFunctionBegin; 1739d763cef2SBarry Smith for (i=0; i<n; i++) { 1740142b95e3SSatish Balay ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1741d763cef2SBarry Smith } 1742d763cef2SBarry Smith PetscFunctionReturn(0); 1743d763cef2SBarry Smith } 1744d763cef2SBarry Smith 1745d763cef2SBarry Smith /* ------------------------------------------------------------------------*/ 1746d763cef2SBarry Smith 17474a2ae208SSatish Balay #undef __FUNCT__ 1748a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate" 1749d763cef2SBarry Smith /*@C 1750a6570f20SBarry Smith TSMonitorLGCreate - Creates a line graph context for use with 1751d763cef2SBarry Smith TS to monitor convergence of preconditioned residual norms. 1752d763cef2SBarry Smith 1753d763cef2SBarry Smith Collective on TS 1754d763cef2SBarry Smith 1755d763cef2SBarry Smith Input Parameters: 1756d763cef2SBarry Smith + host - the X display to open, or null for the local machine 1757d763cef2SBarry Smith . label - the title to put in the title bar 17587c922b88SBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 1759d763cef2SBarry Smith - m, n - the screen width and height in pixels 1760d763cef2SBarry Smith 1761d763cef2SBarry Smith Output Parameter: 1762d763cef2SBarry Smith . draw - the drawing context 1763d763cef2SBarry Smith 1764d763cef2SBarry Smith Options Database Key: 1765a6570f20SBarry Smith . -ts_monitor_draw - automatically sets line graph monitor 1766d763cef2SBarry Smith 1767d763cef2SBarry Smith Notes: 1768a6570f20SBarry Smith Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1769d763cef2SBarry Smith 1770d763cef2SBarry Smith Level: intermediate 1771d763cef2SBarry Smith 17727c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso 1773d763cef2SBarry Smith 1774a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet() 17757c922b88SBarry Smith 1776d763cef2SBarry Smith @*/ 1777a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1778d763cef2SBarry Smith { 1779b0a32e0cSBarry Smith PetscDraw win; 1780dfbe8321SBarry Smith PetscErrorCode ierr; 1781d763cef2SBarry Smith 1782d763cef2SBarry Smith PetscFunctionBegin; 1783b0a32e0cSBarry Smith ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1784b0a32e0cSBarry Smith ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1785b0a32e0cSBarry Smith ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1786b0a32e0cSBarry Smith ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1787d763cef2SBarry Smith 178852e6d16bSBarry Smith ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1789d763cef2SBarry Smith PetscFunctionReturn(0); 1790d763cef2SBarry Smith } 1791d763cef2SBarry Smith 17924a2ae208SSatish Balay #undef __FUNCT__ 1793a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG" 1794a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1795d763cef2SBarry Smith { 1796b0a32e0cSBarry Smith PetscDrawLG lg = (PetscDrawLG) monctx; 179787828ca2SBarry Smith PetscReal x,y = ptime; 1798dfbe8321SBarry Smith PetscErrorCode ierr; 1799d763cef2SBarry Smith 1800d763cef2SBarry Smith PetscFunctionBegin; 18017c922b88SBarry Smith if (!monctx) { 18027c922b88SBarry Smith MPI_Comm comm; 1803b0a32e0cSBarry Smith PetscViewer viewer; 18047c922b88SBarry Smith 18057c922b88SBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1806b0a32e0cSBarry Smith viewer = PETSC_VIEWER_DRAW_(comm); 1807b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 18087c922b88SBarry Smith } 18097c922b88SBarry Smith 1810b0a32e0cSBarry Smith if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 181187828ca2SBarry Smith x = (PetscReal)n; 1812b0a32e0cSBarry Smith ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1813d763cef2SBarry Smith if (n < 20 || (n % 5)) { 1814b0a32e0cSBarry Smith ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1815d763cef2SBarry Smith } 1816d763cef2SBarry Smith PetscFunctionReturn(0); 1817d763cef2SBarry Smith } 1818d763cef2SBarry Smith 18194a2ae208SSatish Balay #undef __FUNCT__ 1820a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy" 1821d763cef2SBarry Smith /*@C 1822a6570f20SBarry Smith TSMonitorLGDestroy - Destroys a line graph context that was created 1823a6570f20SBarry Smith with TSMonitorLGCreate(). 1824d763cef2SBarry Smith 1825b0a32e0cSBarry Smith Collective on PetscDrawLG 1826d763cef2SBarry Smith 1827d763cef2SBarry Smith Input Parameter: 1828d763cef2SBarry Smith . draw - the drawing context 1829d763cef2SBarry Smith 1830d763cef2SBarry Smith Level: intermediate 1831d763cef2SBarry Smith 1832d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy 1833d763cef2SBarry Smith 1834a6570f20SBarry Smith .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1835d763cef2SBarry Smith @*/ 1836a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg) 1837d763cef2SBarry Smith { 1838b0a32e0cSBarry Smith PetscDraw draw; 1839dfbe8321SBarry Smith PetscErrorCode ierr; 1840d763cef2SBarry Smith 1841d763cef2SBarry Smith PetscFunctionBegin; 1842b0a32e0cSBarry Smith ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1843b0a32e0cSBarry Smith ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1844b0a32e0cSBarry Smith ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1845d763cef2SBarry Smith PetscFunctionReturn(0); 1846d763cef2SBarry Smith } 1847d763cef2SBarry Smith 18484a2ae208SSatish Balay #undef __FUNCT__ 18494a2ae208SSatish Balay #define __FUNCT__ "TSGetTime" 1850d763cef2SBarry Smith /*@ 1851d763cef2SBarry Smith TSGetTime - Gets the current time. 1852d763cef2SBarry Smith 1853d763cef2SBarry Smith Not Collective 1854d763cef2SBarry Smith 1855d763cef2SBarry Smith Input Parameter: 1856d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1857d763cef2SBarry Smith 1858d763cef2SBarry Smith Output Parameter: 1859d763cef2SBarry Smith . t - the current time 1860d763cef2SBarry Smith 1861d763cef2SBarry Smith Level: beginner 1862d763cef2SBarry Smith 1863d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1864d763cef2SBarry Smith 1865d763cef2SBarry Smith .keywords: TS, get, time 1866d763cef2SBarry Smith @*/ 186763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t) 1868d763cef2SBarry Smith { 1869d763cef2SBarry Smith PetscFunctionBegin; 18700700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 18714482741eSBarry Smith PetscValidDoublePointer(t,2); 1872d763cef2SBarry Smith *t = ts->ptime; 1873d763cef2SBarry Smith PetscFunctionReturn(0); 1874d763cef2SBarry Smith } 1875d763cef2SBarry Smith 18764a2ae208SSatish Balay #undef __FUNCT__ 18776a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime" 18786a4d4014SLisandro Dalcin /*@ 18796a4d4014SLisandro Dalcin TSSetTime - Allows one to reset the time. 18806a4d4014SLisandro Dalcin 18816a4d4014SLisandro Dalcin Collective on TS 18826a4d4014SLisandro Dalcin 18836a4d4014SLisandro Dalcin Input Parameters: 18846a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 18856a4d4014SLisandro Dalcin - time - the time 18866a4d4014SLisandro Dalcin 18876a4d4014SLisandro Dalcin Level: intermediate 18886a4d4014SLisandro Dalcin 18896a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration() 18906a4d4014SLisandro Dalcin 18916a4d4014SLisandro Dalcin .keywords: TS, set, time 18926a4d4014SLisandro Dalcin @*/ 18936a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t) 18946a4d4014SLisandro Dalcin { 18956a4d4014SLisandro Dalcin PetscFunctionBegin; 18960700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 18976a4d4014SLisandro Dalcin ts->ptime = t; 18986a4d4014SLisandro Dalcin PetscFunctionReturn(0); 18996a4d4014SLisandro Dalcin } 19006a4d4014SLisandro Dalcin 19016a4d4014SLisandro Dalcin #undef __FUNCT__ 19024a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix" 1903d763cef2SBarry Smith /*@C 1904d763cef2SBarry Smith TSSetOptionsPrefix - Sets the prefix used for searching for all 1905d763cef2SBarry Smith TS options in the database. 1906d763cef2SBarry Smith 1907d763cef2SBarry Smith Collective on TS 1908d763cef2SBarry Smith 1909d763cef2SBarry Smith Input Parameter: 1910d763cef2SBarry Smith + ts - The TS context 1911d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1912d763cef2SBarry Smith 1913d763cef2SBarry Smith Notes: 1914d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1915d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1916d763cef2SBarry Smith hyphen. 1917d763cef2SBarry Smith 1918d763cef2SBarry Smith Level: advanced 1919d763cef2SBarry Smith 1920d763cef2SBarry Smith .keywords: TS, set, options, prefix, database 1921d763cef2SBarry Smith 1922d763cef2SBarry Smith .seealso: TSSetFromOptions() 1923d763cef2SBarry Smith 1924d763cef2SBarry Smith @*/ 192563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[]) 1926d763cef2SBarry Smith { 1927dfbe8321SBarry Smith PetscErrorCode ierr; 1928d763cef2SBarry Smith 1929d763cef2SBarry Smith PetscFunctionBegin; 19300700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1931d763cef2SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1932d763cef2SBarry Smith switch(ts->problem_type) { 1933d763cef2SBarry Smith case TS_NONLINEAR: 193459580b9cSBarry Smith if (ts->snes) { 1935d763cef2SBarry Smith ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 193659580b9cSBarry Smith } 1937d763cef2SBarry Smith break; 1938d763cef2SBarry Smith case TS_LINEAR: 193959580b9cSBarry Smith if (ts->ksp) { 194094b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 194159580b9cSBarry Smith } 1942d763cef2SBarry Smith break; 1943d763cef2SBarry Smith } 1944d763cef2SBarry Smith PetscFunctionReturn(0); 1945d763cef2SBarry Smith } 1946d763cef2SBarry Smith 1947d763cef2SBarry Smith 19484a2ae208SSatish Balay #undef __FUNCT__ 19494a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix" 1950d763cef2SBarry Smith /*@C 1951d763cef2SBarry Smith TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1952d763cef2SBarry Smith TS options in the database. 1953d763cef2SBarry Smith 1954d763cef2SBarry Smith Collective on TS 1955d763cef2SBarry Smith 1956d763cef2SBarry Smith Input Parameter: 1957d763cef2SBarry Smith + ts - The TS context 1958d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1959d763cef2SBarry Smith 1960d763cef2SBarry Smith Notes: 1961d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1962d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1963d763cef2SBarry Smith hyphen. 1964d763cef2SBarry Smith 1965d763cef2SBarry Smith Level: advanced 1966d763cef2SBarry Smith 1967d763cef2SBarry Smith .keywords: TS, append, options, prefix, database 1968d763cef2SBarry Smith 1969d763cef2SBarry Smith .seealso: TSGetOptionsPrefix() 1970d763cef2SBarry Smith 1971d763cef2SBarry Smith @*/ 197263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[]) 1973d763cef2SBarry Smith { 1974dfbe8321SBarry Smith PetscErrorCode ierr; 1975d763cef2SBarry Smith 1976d763cef2SBarry Smith PetscFunctionBegin; 19770700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1978d763cef2SBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1979d763cef2SBarry Smith switch(ts->problem_type) { 1980d763cef2SBarry Smith case TS_NONLINEAR: 19811ac94b3bSBarry Smith if (ts->snes) { 1982d763cef2SBarry Smith ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 19831ac94b3bSBarry Smith } 1984d763cef2SBarry Smith break; 1985d763cef2SBarry Smith case TS_LINEAR: 19861ac94b3bSBarry Smith if (ts->ksp) { 198794b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 19881ac94b3bSBarry Smith } 1989d763cef2SBarry Smith break; 1990d763cef2SBarry Smith } 1991d763cef2SBarry Smith PetscFunctionReturn(0); 1992d763cef2SBarry Smith } 1993d763cef2SBarry Smith 19944a2ae208SSatish Balay #undef __FUNCT__ 19954a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix" 1996d763cef2SBarry Smith /*@C 1997d763cef2SBarry Smith TSGetOptionsPrefix - Sets the prefix used for searching for all 1998d763cef2SBarry Smith TS options in the database. 1999d763cef2SBarry Smith 2000d763cef2SBarry Smith Not Collective 2001d763cef2SBarry Smith 2002d763cef2SBarry Smith Input Parameter: 2003d763cef2SBarry Smith . ts - The TS context 2004d763cef2SBarry Smith 2005d763cef2SBarry Smith Output Parameter: 2006d763cef2SBarry Smith . prefix - A pointer to the prefix string used 2007d763cef2SBarry Smith 2008d763cef2SBarry Smith Notes: On the fortran side, the user should pass in a string 'prifix' of 2009d763cef2SBarry Smith sufficient length to hold the prefix. 2010d763cef2SBarry Smith 2011d763cef2SBarry Smith Level: intermediate 2012d763cef2SBarry Smith 2013d763cef2SBarry Smith .keywords: TS, get, options, prefix, database 2014d763cef2SBarry Smith 2015d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix() 2016d763cef2SBarry Smith @*/ 2017e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[]) 2018d763cef2SBarry Smith { 2019dfbe8321SBarry Smith PetscErrorCode ierr; 2020d763cef2SBarry Smith 2021d763cef2SBarry Smith PetscFunctionBegin; 20220700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 20234482741eSBarry Smith PetscValidPointer(prefix,2); 2024d763cef2SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2025d763cef2SBarry Smith PetscFunctionReturn(0); 2026d763cef2SBarry Smith } 2027d763cef2SBarry Smith 20284a2ae208SSatish Balay #undef __FUNCT__ 20294a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian" 2030d763cef2SBarry Smith /*@C 2031d763cef2SBarry Smith TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2032d763cef2SBarry Smith 2033d763cef2SBarry Smith Not Collective, but parallel objects are returned if TS is parallel 2034d763cef2SBarry Smith 2035d763cef2SBarry Smith Input Parameter: 2036d763cef2SBarry Smith . ts - The TS context obtained from TSCreate() 2037d763cef2SBarry Smith 2038d763cef2SBarry Smith Output Parameters: 2039d763cef2SBarry Smith + J - The Jacobian J of F, where U_t = F(U,t) 2040d763cef2SBarry Smith . M - The preconditioner matrix, usually the same as J 2041d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine 2042d763cef2SBarry Smith 2043d763cef2SBarry Smith Notes: You can pass in PETSC_NULL for any return argument you do not need. 2044d763cef2SBarry Smith 2045d763cef2SBarry Smith Level: intermediate 2046d763cef2SBarry Smith 204726d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2048d763cef2SBarry Smith 2049d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian 2050d763cef2SBarry Smith @*/ 205163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 2052d763cef2SBarry Smith { 2053d763cef2SBarry Smith PetscFunctionBegin; 205426d46c62SHong Zhang if (J) *J = ts->Arhs; 205526d46c62SHong Zhang if (M) *M = ts->B; 205626d46c62SHong Zhang if (ctx) *ctx = ts->jacP; 2057d763cef2SBarry Smith PetscFunctionReturn(0); 2058d763cef2SBarry Smith } 2059d763cef2SBarry Smith 20601713a123SBarry Smith #undef __FUNCT__ 20612eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian" 20622eca1d9cSJed Brown /*@C 20632eca1d9cSJed Brown TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 20642eca1d9cSJed Brown 20652eca1d9cSJed Brown Not Collective, but parallel objects are returned if TS is parallel 20662eca1d9cSJed Brown 20672eca1d9cSJed Brown Input Parameter: 20682eca1d9cSJed Brown . ts - The TS context obtained from TSCreate() 20692eca1d9cSJed Brown 20702eca1d9cSJed Brown Output Parameters: 20712eca1d9cSJed Brown + A - The Jacobian of F(t,U,U_t) 20722eca1d9cSJed Brown . B - The preconditioner matrix, often the same as A 20732eca1d9cSJed Brown . f - The function to compute the matrices 20742eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine 20752eca1d9cSJed Brown 20762eca1d9cSJed Brown Notes: You can pass in PETSC_NULL for any return argument you do not need. 20772eca1d9cSJed Brown 20782eca1d9cSJed Brown Level: advanced 20792eca1d9cSJed Brown 20802eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 20812eca1d9cSJed Brown 20822eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian 20832eca1d9cSJed Brown @*/ 20842eca1d9cSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 20852eca1d9cSJed Brown { 20862eca1d9cSJed Brown PetscFunctionBegin; 20872eca1d9cSJed Brown if (A) *A = ts->A; 20882eca1d9cSJed Brown if (B) *B = ts->B; 20892eca1d9cSJed Brown if (f) *f = ts->ops->ijacobian; 20902eca1d9cSJed Brown if (ctx) *ctx = ts->jacP; 20912eca1d9cSJed Brown PetscFunctionReturn(0); 20922eca1d9cSJed Brown } 20932eca1d9cSJed Brown 20942eca1d9cSJed Brown #undef __FUNCT__ 2095a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution" 20961713a123SBarry Smith /*@C 2097a6570f20SBarry Smith TSMonitorSolution - Monitors progress of the TS solvers by calling 20981713a123SBarry Smith VecView() for the solution at each timestep 20991713a123SBarry Smith 21001713a123SBarry Smith Collective on TS 21011713a123SBarry Smith 21021713a123SBarry Smith Input Parameters: 21031713a123SBarry Smith + ts - the TS context 21041713a123SBarry Smith . step - current time-step 2105142b95e3SSatish Balay . ptime - current time 21061713a123SBarry Smith - dummy - either a viewer or PETSC_NULL 21071713a123SBarry Smith 21081713a123SBarry Smith Level: intermediate 21091713a123SBarry Smith 21101713a123SBarry Smith .keywords: TS, vector, monitor, view 21111713a123SBarry Smith 2112a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 21131713a123SBarry Smith @*/ 2114a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 21151713a123SBarry Smith { 2116dfbe8321SBarry Smith PetscErrorCode ierr; 21171713a123SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 21181713a123SBarry Smith 21191713a123SBarry Smith PetscFunctionBegin; 2120a34d58ebSBarry Smith if (!dummy) { 21217adad957SLisandro Dalcin viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 21221713a123SBarry Smith } 21231713a123SBarry Smith ierr = VecView(x,viewer);CHKERRQ(ierr); 21241713a123SBarry Smith PetscFunctionReturn(0); 21251713a123SBarry Smith } 21261713a123SBarry Smith 21271713a123SBarry Smith 21286c699258SBarry Smith #undef __FUNCT__ 21296c699258SBarry Smith #define __FUNCT__ "TSSetDM" 21306c699258SBarry Smith /*@ 21316c699258SBarry Smith TSSetDM - Sets the DM that may be used by some preconditioners 21326c699258SBarry Smith 21336c699258SBarry Smith Collective on TS 21346c699258SBarry Smith 21356c699258SBarry Smith Input Parameters: 21366c699258SBarry Smith + ts - the preconditioner context 21376c699258SBarry Smith - dm - the dm 21386c699258SBarry Smith 21396c699258SBarry Smith Level: intermediate 21406c699258SBarry Smith 21416c699258SBarry Smith 21426c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 21436c699258SBarry Smith @*/ 21446c699258SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSSetDM(TS ts,DM dm) 21456c699258SBarry Smith { 21466c699258SBarry Smith PetscErrorCode ierr; 21476c699258SBarry Smith 21486c699258SBarry Smith PetscFunctionBegin; 21490700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2150*70663e4aSLisandro Dalcin ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 21516c699258SBarry Smith if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);} 21526c699258SBarry Smith ts->dm = dm; 2153*70663e4aSLisandro Dalcin if (ts->snes) { 2154*70663e4aSLisandro Dalcin ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr); 2155*70663e4aSLisandro Dalcin } 2156*70663e4aSLisandro Dalcin if (ts->ksp) { 2157*70663e4aSLisandro Dalcin ierr = KSPSetDM(ts->ksp,dm);CHKERRQ(ierr); 2158*70663e4aSLisandro Dalcin ierr = KSPSetDMActive(ts->ksp,PETSC_FALSE);CHKERRQ(ierr); 2159*70663e4aSLisandro Dalcin } 21606c699258SBarry Smith PetscFunctionReturn(0); 21616c699258SBarry Smith } 21626c699258SBarry Smith 21636c699258SBarry Smith #undef __FUNCT__ 21646c699258SBarry Smith #define __FUNCT__ "TSGetDM" 21656c699258SBarry Smith /*@ 21666c699258SBarry Smith TSGetDM - Gets the DM that may be used by some preconditioners 21676c699258SBarry Smith 21686c699258SBarry Smith Collective on TS 21696c699258SBarry Smith 21706c699258SBarry Smith Input Parameter: 21716c699258SBarry Smith . ts - the preconditioner context 21726c699258SBarry Smith 21736c699258SBarry Smith Output Parameter: 21746c699258SBarry Smith . dm - the dm 21756c699258SBarry Smith 21766c699258SBarry Smith Level: intermediate 21776c699258SBarry Smith 21786c699258SBarry Smith 21796c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 21806c699258SBarry Smith @*/ 21816c699258SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetDM(TS ts,DM *dm) 21826c699258SBarry Smith { 21836c699258SBarry Smith PetscFunctionBegin; 21840700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 21856c699258SBarry Smith *dm = ts->dm; 21866c699258SBarry Smith PetscFunctionReturn(0); 21876c699258SBarry Smith } 21881713a123SBarry Smith 21890f5c6efeSJed Brown #undef __FUNCT__ 21900f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction" 21910f5c6efeSJed Brown /*@ 21920f5c6efeSJed Brown SNESTSFormFunction - Function to evaluate nonlinear residual 21930f5c6efeSJed Brown 21940f5c6efeSJed Brown Collective on SNES 21950f5c6efeSJed Brown 21960f5c6efeSJed Brown Input Parameter: 2197d42a1c89SJed Brown + snes - nonlinear solver 21980f5c6efeSJed Brown . X - the current state at which to evaluate the residual 2199d42a1c89SJed Brown - ctx - user context, must be a TS 22000f5c6efeSJed Brown 22010f5c6efeSJed Brown Output Parameter: 22020f5c6efeSJed Brown . F - the nonlinear residual 22030f5c6efeSJed Brown 22040f5c6efeSJed Brown Notes: 22050f5c6efeSJed Brown This function is not normally called by users and is automatically registered with the SNES used by TS. 22060f5c6efeSJed Brown It is most frequently passed to MatFDColoringSetFunction(). 22070f5c6efeSJed Brown 22080f5c6efeSJed Brown Level: advanced 22090f5c6efeSJed Brown 22100f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction() 22110f5c6efeSJed Brown @*/ 22120f5c6efeSJed Brown PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 22130f5c6efeSJed Brown { 22140f5c6efeSJed Brown TS ts = (TS)ctx; 22150f5c6efeSJed Brown PetscErrorCode ierr; 22160f5c6efeSJed Brown 22170f5c6efeSJed Brown PetscFunctionBegin; 22180f5c6efeSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 22190f5c6efeSJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 22200f5c6efeSJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 22210f5c6efeSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,4); 22220f5c6efeSJed Brown ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 22230f5c6efeSJed Brown PetscFunctionReturn(0); 22240f5c6efeSJed Brown } 22250f5c6efeSJed Brown 22260f5c6efeSJed Brown #undef __FUNCT__ 22270f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian" 22280f5c6efeSJed Brown /*@ 22290f5c6efeSJed Brown SNESTSFormJacobian - Function to evaluate the Jacobian 22300f5c6efeSJed Brown 22310f5c6efeSJed Brown Collective on SNES 22320f5c6efeSJed Brown 22330f5c6efeSJed Brown Input Parameter: 22340f5c6efeSJed Brown + snes - nonlinear solver 22350f5c6efeSJed Brown . X - the current state at which to evaluate the residual 22360f5c6efeSJed Brown - ctx - user context, must be a TS 22370f5c6efeSJed Brown 22380f5c6efeSJed Brown Output Parameter: 22390f5c6efeSJed Brown + A - the Jacobian 22400f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A) 22410f5c6efeSJed Brown - flag - indicates any structure change in the matrix 22420f5c6efeSJed Brown 22430f5c6efeSJed Brown Notes: 22440f5c6efeSJed Brown This function is not normally called by users and is automatically registered with the SNES used by TS. 22450f5c6efeSJed Brown 22460f5c6efeSJed Brown Level: developer 22470f5c6efeSJed Brown 22480f5c6efeSJed Brown .seealso: SNESSetJacobian() 22490f5c6efeSJed Brown @*/ 22500f5c6efeSJed Brown PetscErrorCode PETSCTS_DLLEXPORT SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 22510f5c6efeSJed Brown { 22520f5c6efeSJed Brown TS ts = (TS)ctx; 22530f5c6efeSJed Brown PetscErrorCode ierr; 22540f5c6efeSJed Brown 22550f5c6efeSJed Brown PetscFunctionBegin; 22560f5c6efeSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 22570f5c6efeSJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 22580f5c6efeSJed Brown PetscValidPointer(A,3); 22590f5c6efeSJed Brown PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 22600f5c6efeSJed Brown PetscValidPointer(B,4); 22610f5c6efeSJed Brown PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 22620f5c6efeSJed Brown PetscValidPointer(flag,5); 22630f5c6efeSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,6); 22640f5c6efeSJed Brown ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 22650f5c6efeSJed Brown PetscFunctionReturn(0); 22660f5c6efeSJed Brown } 2267