163dd3a1aSKris Buschelman 2c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 3d763cef2SBarry Smith 4d5ba7fb7SMatthew Knepley /* Logging support */ 57087cfbeSBarry Smith PetscClassId TS_CLASSID; 6166c7f25SBarry Smith PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 7d405a339SMatthew Knepley 84a2ae208SSatish Balay #undef __FUNCT__ 9bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions" 10bdad233fSMatthew Knepley /* 11bdad233fSMatthew Knepley TSSetTypeFromOptions - Sets the type of ts from user options. 12bdad233fSMatthew Knepley 13bdad233fSMatthew Knepley Collective on TS 14bdad233fSMatthew Knepley 15bdad233fSMatthew Knepley Input Parameter: 16bdad233fSMatthew Knepley . ts - The ts 17bdad233fSMatthew Knepley 18bdad233fSMatthew Knepley Level: intermediate 19bdad233fSMatthew Knepley 20bdad233fSMatthew Knepley .keywords: TS, set, options, database, type 21bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType() 22bdad233fSMatthew Knepley */ 236849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts) 24bdad233fSMatthew Knepley { 25ace3abfcSBarry Smith PetscBool opt; 262fc52814SBarry Smith const char *defaultType; 27bdad233fSMatthew Knepley char typeName[256]; 28dfbe8321SBarry Smith PetscErrorCode ierr; 29bdad233fSMatthew Knepley 30bdad233fSMatthew Knepley PetscFunctionBegin; 317adad957SLisandro Dalcin if (((PetscObject)ts)->type_name) { 327adad957SLisandro Dalcin defaultType = ((PetscObject)ts)->type_name; 33bdad233fSMatthew Knepley } else { 349596e0b4SJed Brown defaultType = TSEULER; 35bdad233fSMatthew Knepley } 36bdad233fSMatthew Knepley 37cce0b1b2SLisandro Dalcin if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 38bdad233fSMatthew Knepley ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 39a7cc72afSBarry Smith if (opt) { 40bdad233fSMatthew Knepley ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 41bdad233fSMatthew Knepley } else { 42bdad233fSMatthew Knepley ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 43bdad233fSMatthew Knepley } 44bdad233fSMatthew Knepley PetscFunctionReturn(0); 45bdad233fSMatthew Knepley } 46bdad233fSMatthew Knepley 47bdad233fSMatthew Knepley #undef __FUNCT__ 48bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions" 49bdad233fSMatthew Knepley /*@ 50bdad233fSMatthew Knepley TSSetFromOptions - Sets various TS parameters from user options. 51bdad233fSMatthew Knepley 52bdad233fSMatthew Knepley Collective on TS 53bdad233fSMatthew Knepley 54bdad233fSMatthew Knepley Input Parameter: 55bdad233fSMatthew Knepley . ts - the TS context obtained from TSCreate() 56bdad233fSMatthew Knepley 57bdad233fSMatthew Knepley Options Database Keys: 584d91e141SJed Brown + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 59bdad233fSMatthew Knepley . -ts_max_steps maxsteps - maximum number of time-steps to take 60bdad233fSMatthew Knepley . -ts_max_time time - maximum time to compute to 61bdad233fSMatthew Knepley . -ts_dt dt - initial time step 62bdad233fSMatthew Knepley . -ts_monitor - print information at each timestep 63a6570f20SBarry Smith - -ts_monitor_draw - plot information at each timestep 64bdad233fSMatthew Knepley 65bdad233fSMatthew Knepley Level: beginner 66bdad233fSMatthew Knepley 67bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database 68bdad233fSMatthew Knepley 69a313700dSBarry Smith .seealso: TSGetType() 70bdad233fSMatthew Knepley @*/ 717087cfbeSBarry Smith PetscErrorCode TSSetFromOptions(TS ts) 72bdad233fSMatthew Knepley { 73bdad233fSMatthew Knepley PetscReal dt; 74ace3abfcSBarry Smith PetscBool opt,flg; 75dfbe8321SBarry Smith PetscErrorCode ierr; 76649052a6SBarry Smith PetscViewer monviewer; 77eabae89aSBarry Smith char monfilename[PETSC_MAX_PATH_LEN]; 78089b2837SJed Brown SNES snes; 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) { 96649052a6SBarry Smith ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 97649052a6SBarry Smith ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 98bdad233fSMatthew Knepley } 9990d69ab7SBarry Smith opt = PETSC_FALSE; 100acfcf0e5SJed Brown ierr = PetscOptionsBool("-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; 105acfcf0e5SJed Brown ierr = PetscOptionsBool("-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 } 1175d973c19SBarry Smith 1185d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1195d973c19SBarry Smith ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 120bdad233fSMatthew Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 121bdad233fSMatthew Knepley 122089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 123bdad233fSMatthew Knepley /* Handle subobject options */ 124089b2837SJed Brown if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 125089b2837SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 126bdad233fSMatthew Knepley PetscFunctionReturn(0); 127bdad233fSMatthew Knepley } 128bdad233fSMatthew Knepley 129bdad233fSMatthew Knepley #undef __FUNCT__ 130bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions" 131bdad233fSMatthew Knepley /*@ 132bdad233fSMatthew Knepley TSViewFromOptions - This function visualizes the ts based upon user options. 133bdad233fSMatthew Knepley 134bdad233fSMatthew Knepley Collective on TS 135bdad233fSMatthew Knepley 136bdad233fSMatthew Knepley Input Parameter: 137bdad233fSMatthew Knepley . ts - The ts 138bdad233fSMatthew Knepley 139bdad233fSMatthew Knepley Level: intermediate 140bdad233fSMatthew Knepley 141bdad233fSMatthew Knepley .keywords: TS, view, options, database 142bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView() 143bdad233fSMatthew Knepley @*/ 1447087cfbeSBarry Smith PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 145bdad233fSMatthew Knepley { 146bdad233fSMatthew Knepley PetscViewer viewer; 147bdad233fSMatthew Knepley PetscDraw draw; 148ace3abfcSBarry Smith PetscBool opt = PETSC_FALSE; 149e10c49a3SBarry Smith char fileName[PETSC_MAX_PATH_LEN]; 150dfbe8321SBarry Smith PetscErrorCode ierr; 151bdad233fSMatthew Knepley 152bdad233fSMatthew Knepley PetscFunctionBegin; 1537adad957SLisandro Dalcin ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 154eabae89aSBarry Smith if (opt && !PetscPreLoadingOn) { 1557adad957SLisandro Dalcin ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 156bdad233fSMatthew Knepley ierr = TSView(ts, viewer);CHKERRQ(ierr); 1576bf464f9SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 158bdad233fSMatthew Knepley } 1598e83347fSKai Germaschewski opt = PETSC_FALSE; 160acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 161a7cc72afSBarry Smith if (opt) { 1627adad957SLisandro Dalcin ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 163bdad233fSMatthew Knepley ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 164a7cc72afSBarry Smith if (title) { 1651836bdbcSSatish Balay ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 166bdad233fSMatthew Knepley } else { 167bdad233fSMatthew Knepley ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 1687adad957SLisandro Dalcin ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 169bdad233fSMatthew Knepley } 170bdad233fSMatthew Knepley ierr = TSView(ts, viewer);CHKERRQ(ierr); 171bdad233fSMatthew Knepley ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 172bdad233fSMatthew Knepley ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1736bf464f9SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 174bdad233fSMatthew Knepley } 175bdad233fSMatthew Knepley PetscFunctionReturn(0); 176bdad233fSMatthew Knepley } 177bdad233fSMatthew Knepley 178bdad233fSMatthew Knepley #undef __FUNCT__ 1794a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian" 180a7a1495cSBarry Smith /*@ 1818c385f81SBarry Smith TSComputeRHSJacobian - Computes the Jacobian matrix that has been 182a7a1495cSBarry Smith set with TSSetRHSJacobian(). 183a7a1495cSBarry Smith 184a7a1495cSBarry Smith Collective on TS and Vec 185a7a1495cSBarry Smith 186a7a1495cSBarry Smith Input Parameters: 187316643e7SJed Brown + ts - the TS context 188a7a1495cSBarry Smith . t - current timestep 189a7a1495cSBarry Smith - x - input vector 190a7a1495cSBarry Smith 191a7a1495cSBarry Smith Output Parameters: 192a7a1495cSBarry Smith + A - Jacobian matrix 193a7a1495cSBarry Smith . B - optional preconditioning matrix 194a7a1495cSBarry Smith - flag - flag indicating matrix structure 195a7a1495cSBarry Smith 196a7a1495cSBarry Smith Notes: 197a7a1495cSBarry Smith Most users should not need to explicitly call this routine, as it 198a7a1495cSBarry Smith is used internally within the nonlinear solvers. 199a7a1495cSBarry Smith 20094b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 201a7a1495cSBarry Smith flag parameter. 202a7a1495cSBarry Smith 203a7a1495cSBarry Smith Level: developer 204a7a1495cSBarry Smith 205a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix 206a7a1495cSBarry Smith 20794b7f48cSBarry Smith .seealso: TSSetRHSJacobian(), KSPSetOperators() 208a7a1495cSBarry Smith @*/ 2097087cfbeSBarry Smith PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 210a7a1495cSBarry Smith { 211dfbe8321SBarry Smith PetscErrorCode ierr; 212a7a1495cSBarry Smith 213a7a1495cSBarry Smith PetscFunctionBegin; 2140700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2150700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 216c9780b6fSBarry Smith PetscCheckSameComm(ts,1,X,3); 217000e7ae3SMatthew Knepley if (ts->ops->rhsjacobian) { 218d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 219a7a1495cSBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 220a7a1495cSBarry Smith PetscStackPush("TS user Jacobian function"); 221000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 222a7a1495cSBarry Smith PetscStackPop; 223d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 224a7a1495cSBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 2250700a824SBarry Smith PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 2260700a824SBarry Smith PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 227ef66eb69SBarry Smith } else { 228*214bc6a2SJed Brown ierr = MatZeroEntries(*A);CHKERRQ(ierr); 229*214bc6a2SJed Brown if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 230*214bc6a2SJed Brown *flg = SAME_NONZERO_PATTERN; 231ef66eb69SBarry Smith } 232a7a1495cSBarry Smith PetscFunctionReturn(0); 233a7a1495cSBarry Smith } 234a7a1495cSBarry Smith 2354a2ae208SSatish Balay #undef __FUNCT__ 2364a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction" 237316643e7SJed Brown /*@ 238d763cef2SBarry Smith TSComputeRHSFunction - Evaluates the right-hand-side function. 239d763cef2SBarry Smith 240316643e7SJed Brown Collective on TS and Vec 241316643e7SJed Brown 242316643e7SJed Brown Input Parameters: 243316643e7SJed Brown + ts - the TS context 244316643e7SJed Brown . t - current time 245316643e7SJed Brown - x - state vector 246316643e7SJed Brown 247316643e7SJed Brown Output Parameter: 248316643e7SJed Brown . y - right hand side 249316643e7SJed Brown 250316643e7SJed Brown Note: 251316643e7SJed Brown Most users should not need to explicitly call this routine, as it 252316643e7SJed Brown is used internally within the nonlinear solvers. 253316643e7SJed Brown 254316643e7SJed Brown Level: developer 255316643e7SJed Brown 256316643e7SJed Brown .keywords: TS, compute 257316643e7SJed Brown 258316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction() 259316643e7SJed Brown @*/ 260dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 261d763cef2SBarry Smith { 262dfbe8321SBarry Smith PetscErrorCode ierr; 263d763cef2SBarry Smith 264d763cef2SBarry Smith PetscFunctionBegin; 2650700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2660700a824SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2670700a824SBarry Smith PetscValidHeaderSpecific(y,VEC_CLASSID,4); 268d763cef2SBarry Smith 269d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 270000e7ae3SMatthew Knepley if (ts->ops->rhsfunction) { 271d763cef2SBarry Smith PetscStackPush("TS user right-hand-side function"); 272000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 273d763cef2SBarry Smith PetscStackPop; 274*214bc6a2SJed Brown } else { 275*214bc6a2SJed Brown ierr = VecZeroEntries(y);CHKERRQ(ierr); 276b2cd27e8SJed Brown } 277d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 278d763cef2SBarry Smith PetscFunctionReturn(0); 279d763cef2SBarry Smith } 280d763cef2SBarry Smith 2814a2ae208SSatish Balay #undef __FUNCT__ 282*214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSVec_Private" 283*214bc6a2SJed Brown static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 284*214bc6a2SJed Brown { 285*214bc6a2SJed Brown PetscErrorCode ierr; 286*214bc6a2SJed Brown 287*214bc6a2SJed Brown PetscFunctionBegin; 288*214bc6a2SJed Brown if (!ts->Frhs) { 289*214bc6a2SJed Brown ierr = VecDuplicate(ts->vec_sol,&ts->Frhs);CHKERRQ(ierr); 290*214bc6a2SJed Brown } 291*214bc6a2SJed Brown *Frhs = ts->Frhs; 292*214bc6a2SJed Brown PetscFunctionReturn(0); 293*214bc6a2SJed Brown } 294*214bc6a2SJed Brown 295*214bc6a2SJed Brown #undef __FUNCT__ 296*214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSMats_Private" 297*214bc6a2SJed Brown static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 298*214bc6a2SJed Brown { 299*214bc6a2SJed Brown PetscErrorCode ierr; 300*214bc6a2SJed Brown Mat A,B; 301*214bc6a2SJed Brown 302*214bc6a2SJed Brown PetscFunctionBegin; 303*214bc6a2SJed Brown ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 304*214bc6a2SJed Brown if (Arhs) { 305*214bc6a2SJed Brown if (!ts->Arhs) { 306*214bc6a2SJed Brown ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 307*214bc6a2SJed Brown } 308*214bc6a2SJed Brown *Arhs = ts->Arhs; 309*214bc6a2SJed Brown } 310*214bc6a2SJed Brown if (Brhs) { 311*214bc6a2SJed Brown if (!ts->Brhs) { 312*214bc6a2SJed Brown ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 313*214bc6a2SJed Brown } 314*214bc6a2SJed Brown *Brhs = ts->Brhs; 315*214bc6a2SJed Brown } 316*214bc6a2SJed Brown PetscFunctionReturn(0); 317*214bc6a2SJed Brown } 318*214bc6a2SJed Brown 319*214bc6a2SJed Brown #undef __FUNCT__ 320316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction" 321316643e7SJed Brown /*@ 322316643e7SJed Brown TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 323316643e7SJed Brown 324316643e7SJed Brown Collective on TS and Vec 325316643e7SJed Brown 326316643e7SJed Brown Input Parameters: 327316643e7SJed Brown + ts - the TS context 328316643e7SJed Brown . t - current time 329316643e7SJed Brown . X - state vector 330*214bc6a2SJed Brown . Xdot - time derivative of state vector 331*214bc6a2SJed Brown - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 332316643e7SJed Brown 333316643e7SJed Brown Output Parameter: 334316643e7SJed Brown . Y - right hand side 335316643e7SJed Brown 336316643e7SJed Brown Note: 337316643e7SJed Brown Most users should not need to explicitly call this routine, as it 338316643e7SJed Brown is used internally within the nonlinear solvers. 339316643e7SJed Brown 340316643e7SJed Brown If the user did did not write their equations in implicit form, this 341316643e7SJed Brown function recasts them in implicit form. 342316643e7SJed Brown 343316643e7SJed Brown Level: developer 344316643e7SJed Brown 345316643e7SJed Brown .keywords: TS, compute 346316643e7SJed Brown 347316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction() 348316643e7SJed Brown @*/ 349*214bc6a2SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 350316643e7SJed Brown { 351316643e7SJed Brown PetscErrorCode ierr; 352316643e7SJed Brown 353316643e7SJed Brown PetscFunctionBegin; 3540700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3550700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 3560700a824SBarry Smith PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 3570700a824SBarry Smith PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 358316643e7SJed Brown 359316643e7SJed Brown ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 360316643e7SJed Brown if (ts->ops->ifunction) { 361316643e7SJed Brown PetscStackPush("TS user implicit function"); 362316643e7SJed Brown ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 363316643e7SJed Brown PetscStackPop; 364*214bc6a2SJed Brown } 365*214bc6a2SJed Brown if (imex) { 366*214bc6a2SJed Brown if (!ts->ops->ifunction) {ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);} 367316643e7SJed Brown } else { 368*214bc6a2SJed Brown if (!ts->ops->ifunction) { 369089b2837SJed Brown ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 370ace68cafSJed Brown ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 371*214bc6a2SJed Brown } else { 372*214bc6a2SJed Brown Vec Frhs; 373*214bc6a2SJed Brown ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 374*214bc6a2SJed Brown ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 375*214bc6a2SJed Brown ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 376316643e7SJed Brown } 3774a6899ffSJed Brown } 378316643e7SJed Brown ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 379316643e7SJed Brown PetscFunctionReturn(0); 380316643e7SJed Brown } 381316643e7SJed Brown 382316643e7SJed Brown #undef __FUNCT__ 383316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian" 384316643e7SJed Brown /*@ 385316643e7SJed Brown TSComputeIJacobian - Evaluates the Jacobian of the DAE 386316643e7SJed Brown 387316643e7SJed Brown Collective on TS and Vec 388316643e7SJed Brown 389316643e7SJed Brown Input 390316643e7SJed Brown Input Parameters: 391316643e7SJed Brown + ts - the TS context 392316643e7SJed Brown . t - current timestep 393316643e7SJed Brown . X - state vector 394316643e7SJed Brown . Xdot - time derivative of state vector 395*214bc6a2SJed Brown . shift - shift to apply, see note below 396*214bc6a2SJed Brown - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 397316643e7SJed Brown 398316643e7SJed Brown Output Parameters: 399316643e7SJed Brown + A - Jacobian matrix 400316643e7SJed Brown . B - optional preconditioning matrix 401316643e7SJed Brown - flag - flag indicating matrix structure 402316643e7SJed Brown 403316643e7SJed Brown Notes: 404316643e7SJed Brown If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 405316643e7SJed Brown 406316643e7SJed Brown dF/dX + shift*dF/dXdot 407316643e7SJed Brown 408316643e7SJed Brown Most users should not need to explicitly call this routine, as it 409316643e7SJed Brown is used internally within the nonlinear solvers. 410316643e7SJed Brown 411316643e7SJed Brown Level: developer 412316643e7SJed Brown 413316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix 414316643e7SJed Brown 415316643e7SJed Brown .seealso: TSSetIJacobian() 41663495f91SJed Brown @*/ 417*214bc6a2SJed Brown PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 418316643e7SJed Brown { 419316643e7SJed Brown PetscErrorCode ierr; 420316643e7SJed Brown 421316643e7SJed Brown PetscFunctionBegin; 4220700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4230700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 4240700a824SBarry Smith PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 425316643e7SJed Brown PetscValidPointer(A,6); 4260700a824SBarry Smith PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 427316643e7SJed Brown PetscValidPointer(B,7); 4280700a824SBarry Smith PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 429316643e7SJed Brown PetscValidPointer(flg,8); 430316643e7SJed Brown 4314e684422SJed Brown *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 432316643e7SJed Brown ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 433316643e7SJed Brown if (ts->ops->ijacobian) { 434316643e7SJed Brown PetscStackPush("TS user implicit Jacobian"); 435316643e7SJed Brown ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 436316643e7SJed Brown PetscStackPop; 437*214bc6a2SJed Brown /* make sure user returned a correct Jacobian and preconditioner */ 438*214bc6a2SJed Brown PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 439*214bc6a2SJed Brown PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 4404a6899ffSJed Brown } 441*214bc6a2SJed Brown if (imex) { 442*214bc6a2SJed Brown if (!ts->ops->ijacobian) { /* system was written as Xdot = F(t,X) */ 443*214bc6a2SJed Brown ierr = MatZeroEntries(*A);CHKERRQ(ierr); 444*214bc6a2SJed Brown ierr = MatShift(*A,1.0);CHKERRQ(ierr); 445*214bc6a2SJed Brown if (*A != *B) { 446*214bc6a2SJed Brown ierr = MatZeroEntries(*B);CHKERRQ(ierr); 447*214bc6a2SJed Brown ierr = MatShift(*B,shift);CHKERRQ(ierr); 448*214bc6a2SJed Brown } 449*214bc6a2SJed Brown *flg = SAME_PRECONDITIONER; 450*214bc6a2SJed Brown } 451*214bc6a2SJed Brown } else { 452*214bc6a2SJed Brown if (!ts->ops->ijacobian) { 453*214bc6a2SJed Brown ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 454*214bc6a2SJed Brown ierr = MatScale(*A,-1);CHKERRQ(ierr); 455*214bc6a2SJed Brown ierr = MatShift(*A,shift);CHKERRQ(ierr); 456316643e7SJed Brown if (*A != *B) { 457316643e7SJed Brown ierr = MatScale(*B,-1);CHKERRQ(ierr); 458316643e7SJed Brown ierr = MatShift(*B,shift);CHKERRQ(ierr); 459316643e7SJed Brown } 460*214bc6a2SJed Brown } else if (ts->ops->rhsjacobian) { 461*214bc6a2SJed Brown Mat Arhs,Brhs; 462*214bc6a2SJed Brown MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 463*214bc6a2SJed Brown ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 464*214bc6a2SJed Brown ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 465*214bc6a2SJed Brown axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 466*214bc6a2SJed Brown ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 467*214bc6a2SJed Brown if (*A != *B) { 468*214bc6a2SJed Brown ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 469*214bc6a2SJed Brown } 470*214bc6a2SJed Brown *flg = PetscMin(*flg,flg2); 471*214bc6a2SJed Brown } 472316643e7SJed Brown } 473316643e7SJed Brown ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 474316643e7SJed Brown PetscFunctionReturn(0); 475316643e7SJed Brown } 476316643e7SJed Brown 477316643e7SJed Brown #undef __FUNCT__ 4784a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction" 479d763cef2SBarry Smith /*@C 480d763cef2SBarry Smith TSSetRHSFunction - Sets the routine for evaluating the function, 481d763cef2SBarry Smith F(t,u), where U_t = F(t,u). 482d763cef2SBarry Smith 4833f9fe445SBarry Smith Logically Collective on TS 484d763cef2SBarry Smith 485d763cef2SBarry Smith Input Parameters: 486d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 487ca94891dSJed Brown . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 488d763cef2SBarry Smith . f - routine for evaluating the right-hand-side function 489d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 490d763cef2SBarry Smith function evaluation routine (may be PETSC_NULL) 491d763cef2SBarry Smith 492d763cef2SBarry Smith Calling sequence of func: 49387828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 494d763cef2SBarry Smith 495d763cef2SBarry Smith + t - current timestep 496d763cef2SBarry Smith . u - input vector 497d763cef2SBarry Smith . F - function vector 498d763cef2SBarry Smith - ctx - [optional] user-defined function context 499d763cef2SBarry Smith 500d763cef2SBarry Smith Important: 50176f2fa84SHong Zhang The user MUST call either this routine or TSSetMatrices(). 502d763cef2SBarry Smith 503d763cef2SBarry Smith Level: beginner 504d763cef2SBarry Smith 505d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function 506d763cef2SBarry Smith 50776f2fa84SHong Zhang .seealso: TSSetMatrices() 508d763cef2SBarry Smith @*/ 509089b2837SJed Brown PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 510d763cef2SBarry Smith { 511089b2837SJed Brown PetscErrorCode ierr; 512089b2837SJed Brown SNES snes; 513d763cef2SBarry Smith 514089b2837SJed Brown PetscFunctionBegin; 5150700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 516ca94891dSJed Brown if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 5174e684422SJed Brown if (f) ts->ops->rhsfunction = f; 5184e684422SJed Brown if (ctx) ts->funP = ctx; 519089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 520089b2837SJed Brown ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 521d763cef2SBarry Smith PetscFunctionReturn(0); 522d763cef2SBarry Smith } 523d763cef2SBarry Smith 5244a2ae208SSatish Balay #undef __FUNCT__ 5254a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian" 526d763cef2SBarry Smith /*@C 527d763cef2SBarry Smith TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 528d763cef2SBarry Smith where U_t = F(U,t), as well as the location to store the matrix. 52976f2fa84SHong Zhang Use TSSetMatrices() for linear problems. 530d763cef2SBarry Smith 5313f9fe445SBarry Smith Logically Collective on TS 532d763cef2SBarry Smith 533d763cef2SBarry Smith Input Parameters: 534d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 535d763cef2SBarry Smith . A - Jacobian matrix 536d763cef2SBarry Smith . B - preconditioner matrix (usually same as A) 537d763cef2SBarry Smith . f - the Jacobian evaluation routine 538d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 539d763cef2SBarry Smith Jacobian evaluation routine (may be PETSC_NULL) 540d763cef2SBarry Smith 541d763cef2SBarry Smith Calling sequence of func: 54287828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 543d763cef2SBarry Smith 544d763cef2SBarry Smith + t - current timestep 545d763cef2SBarry Smith . u - input vector 546d763cef2SBarry Smith . A - matrix A, where U_t = A(t)u 547d763cef2SBarry Smith . B - preconditioner matrix, usually the same as A 548d763cef2SBarry Smith . flag - flag indicating information about the preconditioner matrix 54994b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 550d763cef2SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 551d763cef2SBarry Smith 552d763cef2SBarry Smith Notes: 55394b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 554d763cef2SBarry Smith output parameter in the routine func(). Be sure to read this information! 555d763cef2SBarry Smith 556d763cef2SBarry Smith The routine func() takes Mat * as the matrix arguments rather than Mat. 557d763cef2SBarry Smith This allows the matrix evaluation routine to replace A and/or B with a 55856335db2SHong Zhang completely new matrix structure (not just different matrix elements) 559d763cef2SBarry Smith when appropriate, for instance, if the nonzero structure is changing 560d763cef2SBarry Smith throughout the global iterations. 561d763cef2SBarry Smith 562d763cef2SBarry Smith Level: beginner 563d763cef2SBarry Smith 564d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian 565d763cef2SBarry Smith 566d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(), 56776f2fa84SHong Zhang SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 568d763cef2SBarry Smith 569d763cef2SBarry Smith @*/ 570089b2837SJed Brown PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 571d763cef2SBarry Smith { 572277b19d0SLisandro Dalcin PetscErrorCode ierr; 573089b2837SJed Brown SNES snes; 574277b19d0SLisandro Dalcin 575d763cef2SBarry Smith PetscFunctionBegin; 5760700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 577277b19d0SLisandro Dalcin if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 578277b19d0SLisandro Dalcin if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 579277b19d0SLisandro Dalcin if (A) PetscCheckSameComm(ts,1,A,2); 580277b19d0SLisandro Dalcin if (B) PetscCheckSameComm(ts,1,B,3); 581d763cef2SBarry Smith 582277b19d0SLisandro Dalcin if (f) ts->ops->rhsjacobian = f; 583277b19d0SLisandro Dalcin if (ctx) ts->jacP = ctx; 584089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 585089b2837SJed Brown ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 586d763cef2SBarry Smith PetscFunctionReturn(0); 587d763cef2SBarry Smith } 588d763cef2SBarry Smith 589316643e7SJed Brown 590316643e7SJed Brown #undef __FUNCT__ 591316643e7SJed Brown #define __FUNCT__ "TSSetIFunction" 592316643e7SJed Brown /*@C 593316643e7SJed Brown TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 594316643e7SJed Brown 5953f9fe445SBarry Smith Logically Collective on TS 596316643e7SJed Brown 597316643e7SJed Brown Input Parameters: 598316643e7SJed Brown + ts - the TS context obtained from TSCreate() 599ca94891dSJed Brown . r - vector to hold the residual (or PETSC_NULL to have it created internally) 600316643e7SJed Brown . f - the function evaluation routine 601316643e7SJed Brown - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 602316643e7SJed Brown 603316643e7SJed Brown Calling sequence of f: 604316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 605316643e7SJed Brown 606316643e7SJed Brown + t - time at step/stage being solved 607316643e7SJed Brown . u - state vector 608316643e7SJed Brown . u_t - time derivative of state vector 609316643e7SJed Brown . F - function vector 610316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 611316643e7SJed Brown 612316643e7SJed Brown Important: 613316643e7SJed Brown The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 614316643e7SJed Brown 615316643e7SJed Brown Level: beginner 616316643e7SJed Brown 617316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian 618316643e7SJed Brown 619316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 620316643e7SJed Brown @*/ 621089b2837SJed Brown PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 622316643e7SJed Brown { 623089b2837SJed Brown PetscErrorCode ierr; 624089b2837SJed Brown SNES snes; 625316643e7SJed Brown 626316643e7SJed Brown PetscFunctionBegin; 6270700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 628ca94891dSJed Brown if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 629089b2837SJed Brown if (f) ts->ops->ifunction = f; 630089b2837SJed Brown if (ctx) ts->funP = ctx; 631089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 632089b2837SJed Brown ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 633089b2837SJed Brown PetscFunctionReturn(0); 634089b2837SJed Brown } 635089b2837SJed Brown 636089b2837SJed Brown #undef __FUNCT__ 637089b2837SJed Brown #define __FUNCT__ "TSGetIFunction" 638089b2837SJed Brown /*@C 639089b2837SJed Brown TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 640089b2837SJed Brown 641089b2837SJed Brown Not Collective 642089b2837SJed Brown 643089b2837SJed Brown Input Parameter: 644089b2837SJed Brown . ts - the TS context 645089b2837SJed Brown 646089b2837SJed Brown Output Parameter: 647089b2837SJed Brown + r - vector to hold residual (or PETSC_NULL) 648089b2837SJed Brown . func - the function to compute residual (or PETSC_NULL) 649089b2837SJed Brown - ctx - the function context (or PETSC_NULL) 650089b2837SJed Brown 651089b2837SJed Brown Level: advanced 652089b2837SJed Brown 653089b2837SJed Brown .keywords: TS, nonlinear, get, function 654089b2837SJed Brown 655089b2837SJed Brown .seealso: TSSetIFunction(), SNESGetFunction() 656089b2837SJed Brown @*/ 657089b2837SJed Brown PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 658089b2837SJed Brown { 659089b2837SJed Brown PetscErrorCode ierr; 660089b2837SJed Brown SNES snes; 661089b2837SJed Brown 662089b2837SJed Brown PetscFunctionBegin; 663089b2837SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 664089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 665089b2837SJed Brown ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 666089b2837SJed Brown if (func) *func = ts->ops->ifunction; 667089b2837SJed Brown if (ctx) *ctx = ts->funP; 668089b2837SJed Brown PetscFunctionReturn(0); 669089b2837SJed Brown } 670089b2837SJed Brown 671089b2837SJed Brown #undef __FUNCT__ 672089b2837SJed Brown #define __FUNCT__ "TSGetRHSFunction" 673089b2837SJed Brown /*@C 674089b2837SJed Brown TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 675089b2837SJed Brown 676089b2837SJed Brown Not Collective 677089b2837SJed Brown 678089b2837SJed Brown Input Parameter: 679089b2837SJed Brown . ts - the TS context 680089b2837SJed Brown 681089b2837SJed Brown Output Parameter: 682089b2837SJed Brown + r - vector to hold computed right hand side (or PETSC_NULL) 683089b2837SJed Brown . func - the function to compute right hand side (or PETSC_NULL) 684089b2837SJed Brown - ctx - the function context (or PETSC_NULL) 685089b2837SJed Brown 686089b2837SJed Brown Level: advanced 687089b2837SJed Brown 688089b2837SJed Brown .keywords: TS, nonlinear, get, function 689089b2837SJed Brown 690089b2837SJed Brown .seealso: TSSetRhsfunction(), SNESGetFunction() 691089b2837SJed Brown @*/ 692089b2837SJed Brown PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 693089b2837SJed Brown { 694089b2837SJed Brown PetscErrorCode ierr; 695089b2837SJed Brown SNES snes; 696089b2837SJed Brown 697089b2837SJed Brown PetscFunctionBegin; 698089b2837SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 699089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 700089b2837SJed Brown ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 701089b2837SJed Brown if (func) *func = ts->ops->rhsfunction; 702089b2837SJed Brown if (ctx) *ctx = ts->funP; 703316643e7SJed Brown PetscFunctionReturn(0); 704316643e7SJed Brown } 705316643e7SJed Brown 706316643e7SJed Brown #undef __FUNCT__ 707316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian" 708316643e7SJed Brown /*@C 709a4f0a591SBarry Smith TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 710a4f0a591SBarry Smith you provided with TSSetIFunction(). 711316643e7SJed Brown 7123f9fe445SBarry Smith Logically Collective on TS 713316643e7SJed Brown 714316643e7SJed Brown Input Parameters: 715316643e7SJed Brown + ts - the TS context obtained from TSCreate() 716316643e7SJed Brown . A - Jacobian matrix 717316643e7SJed Brown . B - preconditioning matrix for A (may be same as A) 718316643e7SJed Brown . f - the Jacobian evaluation routine 719316643e7SJed Brown - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 720316643e7SJed Brown 721316643e7SJed Brown Calling sequence of f: 7221b4a444bSJed Brown $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 723316643e7SJed Brown 724316643e7SJed Brown + t - time at step/stage being solved 7251b4a444bSJed Brown . U - state vector 7261b4a444bSJed Brown . U_t - time derivative of state vector 727316643e7SJed Brown . a - shift 7281b4a444bSJed Brown . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 729316643e7SJed Brown . B - preconditioning matrix for A, may be same as A 730316643e7SJed Brown . flag - flag indicating information about the preconditioner matrix 731316643e7SJed Brown structure (same as flag in KSPSetOperators()) 732316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 733316643e7SJed Brown 734316643e7SJed Brown Notes: 735316643e7SJed Brown The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 736316643e7SJed Brown 737a4f0a591SBarry Smith The matrix dF/dU + a*dF/dU_t you provide turns out to be 738a4f0a591SBarry Smith the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 739a4f0a591SBarry Smith The time integrator internally approximates U_t by W+a*U where the positive "shift" 740a4f0a591SBarry Smith a and vector W depend on the integration method, step size, and past states. For example with 741a4f0a591SBarry Smith the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 742a4f0a591SBarry Smith W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 743a4f0a591SBarry Smith 744316643e7SJed Brown Level: beginner 745316643e7SJed Brown 746316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian 747316643e7SJed Brown 748316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian() 749316643e7SJed Brown 750316643e7SJed Brown @*/ 7517087cfbeSBarry Smith PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 752316643e7SJed Brown { 753316643e7SJed Brown PetscErrorCode ierr; 754089b2837SJed Brown SNES snes; 755316643e7SJed Brown 756316643e7SJed Brown PetscFunctionBegin; 7570700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7580700a824SBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 7590700a824SBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 760316643e7SJed Brown if (A) PetscCheckSameComm(ts,1,A,2); 761316643e7SJed Brown if (B) PetscCheckSameComm(ts,1,B,3); 762316643e7SJed Brown if (f) ts->ops->ijacobian = f; 763316643e7SJed Brown if (ctx) ts->jacP = ctx; 764089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 765089b2837SJed Brown ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 766316643e7SJed Brown PetscFunctionReturn(0); 767316643e7SJed Brown } 768316643e7SJed Brown 7694a2ae208SSatish Balay #undef __FUNCT__ 7704a2ae208SSatish Balay #define __FUNCT__ "TSView" 7717e2c5f70SBarry Smith /*@C 772d763cef2SBarry Smith TSView - Prints the TS data structure. 773d763cef2SBarry Smith 7744c49b128SBarry Smith Collective on TS 775d763cef2SBarry Smith 776d763cef2SBarry Smith Input Parameters: 777d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 778d763cef2SBarry Smith - viewer - visualization context 779d763cef2SBarry Smith 780d763cef2SBarry Smith Options Database Key: 781d763cef2SBarry Smith . -ts_view - calls TSView() at end of TSStep() 782d763cef2SBarry Smith 783d763cef2SBarry Smith Notes: 784d763cef2SBarry Smith The available visualization contexts include 785b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 786b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 787d763cef2SBarry Smith output where only the first processor opens 788d763cef2SBarry Smith the file. All other processors send their 789d763cef2SBarry Smith data to the first processor to print. 790d763cef2SBarry Smith 791d763cef2SBarry Smith The user can open an alternative visualization context with 792b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 793d763cef2SBarry Smith 794d763cef2SBarry Smith Level: beginner 795d763cef2SBarry Smith 796d763cef2SBarry Smith .keywords: TS, timestep, view 797d763cef2SBarry Smith 798b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 799d763cef2SBarry Smith @*/ 8007087cfbeSBarry Smith PetscErrorCode TSView(TS ts,PetscViewer viewer) 801d763cef2SBarry Smith { 802dfbe8321SBarry Smith PetscErrorCode ierr; 803a313700dSBarry Smith const TSType type; 804089b2837SJed Brown PetscBool iascii,isstring,isundials; 805d763cef2SBarry Smith 806d763cef2SBarry Smith PetscFunctionBegin; 8070700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8083050cee2SBarry Smith if (!viewer) { 8097adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 8103050cee2SBarry Smith } 8110700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 812c9780b6fSBarry Smith PetscCheckSameComm(ts,1,viewer,2); 813fd16b177SBarry Smith 8142692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8152692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 81632077d6dSBarry Smith if (iascii) { 817317d6ea6SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 818000e7ae3SMatthew Knepley if (ts->ops->view) { 819b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 820000e7ae3SMatthew Knepley ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 821b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 822d763cef2SBarry Smith } 82377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 824a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 825d763cef2SBarry Smith if (ts->problem_type == TS_NONLINEAR) { 82677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 827d763cef2SBarry Smith } 82877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 8290f5bd95cSBarry Smith } else if (isstring) { 830a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 831b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 832d763cef2SBarry Smith } 833b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 834089b2837SJed Brown ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 835089b2837SJed Brown if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 836b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 837d763cef2SBarry Smith PetscFunctionReturn(0); 838d763cef2SBarry Smith } 839d763cef2SBarry Smith 840d763cef2SBarry Smith 8414a2ae208SSatish Balay #undef __FUNCT__ 8424a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext" 843b07ff414SBarry Smith /*@ 844d763cef2SBarry Smith TSSetApplicationContext - Sets an optional user-defined context for 845d763cef2SBarry Smith the timesteppers. 846d763cef2SBarry Smith 8473f9fe445SBarry Smith Logically Collective on TS 848d763cef2SBarry Smith 849d763cef2SBarry Smith Input Parameters: 850d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 851d763cef2SBarry Smith - usrP - optional user context 852d763cef2SBarry Smith 853d763cef2SBarry Smith Level: intermediate 854d763cef2SBarry Smith 855d763cef2SBarry Smith .keywords: TS, timestep, set, application, context 856d763cef2SBarry Smith 857d763cef2SBarry Smith .seealso: TSGetApplicationContext() 858d763cef2SBarry Smith @*/ 8597087cfbeSBarry Smith PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 860d763cef2SBarry Smith { 861d763cef2SBarry Smith PetscFunctionBegin; 8620700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 863d763cef2SBarry Smith ts->user = usrP; 864d763cef2SBarry Smith PetscFunctionReturn(0); 865d763cef2SBarry Smith } 866d763cef2SBarry Smith 8674a2ae208SSatish Balay #undef __FUNCT__ 8684a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext" 869b07ff414SBarry Smith /*@ 870d763cef2SBarry Smith TSGetApplicationContext - Gets the user-defined context for the 871d763cef2SBarry Smith timestepper. 872d763cef2SBarry Smith 873d763cef2SBarry Smith Not Collective 874d763cef2SBarry Smith 875d763cef2SBarry Smith Input Parameter: 876d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 877d763cef2SBarry Smith 878d763cef2SBarry Smith Output Parameter: 879d763cef2SBarry Smith . usrP - user context 880d763cef2SBarry Smith 881d763cef2SBarry Smith Level: intermediate 882d763cef2SBarry Smith 883d763cef2SBarry Smith .keywords: TS, timestep, get, application, context 884d763cef2SBarry Smith 885d763cef2SBarry Smith .seealso: TSSetApplicationContext() 886d763cef2SBarry Smith @*/ 887e71120c6SJed Brown PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 888d763cef2SBarry Smith { 889d763cef2SBarry Smith PetscFunctionBegin; 8900700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 891e71120c6SJed Brown *(void**)usrP = ts->user; 892d763cef2SBarry Smith PetscFunctionReturn(0); 893d763cef2SBarry Smith } 894d763cef2SBarry Smith 8954a2ae208SSatish Balay #undef __FUNCT__ 8964a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber" 897d763cef2SBarry Smith /*@ 898d763cef2SBarry Smith TSGetTimeStepNumber - Gets the current number of timesteps. 899d763cef2SBarry Smith 900d763cef2SBarry Smith Not Collective 901d763cef2SBarry Smith 902d763cef2SBarry Smith Input Parameter: 903d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 904d763cef2SBarry Smith 905d763cef2SBarry Smith Output Parameter: 906d763cef2SBarry Smith . iter - number steps so far 907d763cef2SBarry Smith 908d763cef2SBarry Smith Level: intermediate 909d763cef2SBarry Smith 910d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number 911d763cef2SBarry Smith @*/ 9127087cfbeSBarry Smith PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 913d763cef2SBarry Smith { 914d763cef2SBarry Smith PetscFunctionBegin; 9150700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9164482741eSBarry Smith PetscValidIntPointer(iter,2); 917d763cef2SBarry Smith *iter = ts->steps; 918d763cef2SBarry Smith PetscFunctionReturn(0); 919d763cef2SBarry Smith } 920d763cef2SBarry Smith 9214a2ae208SSatish Balay #undef __FUNCT__ 9224a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep" 923d763cef2SBarry Smith /*@ 924d763cef2SBarry Smith TSSetInitialTimeStep - Sets the initial timestep to be used, 925d763cef2SBarry Smith as well as the initial time. 926d763cef2SBarry Smith 9273f9fe445SBarry Smith Logically Collective on TS 928d763cef2SBarry Smith 929d763cef2SBarry Smith Input Parameters: 930d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 931d763cef2SBarry Smith . initial_time - the initial time 932d763cef2SBarry Smith - time_step - the size of the timestep 933d763cef2SBarry Smith 934d763cef2SBarry Smith Level: intermediate 935d763cef2SBarry Smith 936d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep() 937d763cef2SBarry Smith 938d763cef2SBarry Smith .keywords: TS, set, initial, timestep 939d763cef2SBarry Smith @*/ 9407087cfbeSBarry Smith PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 941d763cef2SBarry Smith { 942d763cef2SBarry Smith PetscFunctionBegin; 9430700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 944d763cef2SBarry Smith ts->time_step = time_step; 945d763cef2SBarry Smith ts->initial_time_step = time_step; 946d763cef2SBarry Smith ts->ptime = initial_time; 947d763cef2SBarry Smith PetscFunctionReturn(0); 948d763cef2SBarry Smith } 949d763cef2SBarry Smith 9504a2ae208SSatish Balay #undef __FUNCT__ 9514a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep" 952d763cef2SBarry Smith /*@ 953d763cef2SBarry Smith TSSetTimeStep - Allows one to reset the timestep at any time, 954d763cef2SBarry Smith useful for simple pseudo-timestepping codes. 955d763cef2SBarry Smith 9563f9fe445SBarry Smith Logically Collective on TS 957d763cef2SBarry Smith 958d763cef2SBarry Smith Input Parameters: 959d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 960d763cef2SBarry Smith - time_step - the size of the timestep 961d763cef2SBarry Smith 962d763cef2SBarry Smith Level: intermediate 963d763cef2SBarry Smith 964d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 965d763cef2SBarry Smith 966d763cef2SBarry Smith .keywords: TS, set, timestep 967d763cef2SBarry Smith @*/ 9687087cfbeSBarry Smith PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 969d763cef2SBarry Smith { 970d763cef2SBarry Smith PetscFunctionBegin; 9710700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 972c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,time_step,2); 973d763cef2SBarry Smith ts->time_step = time_step; 974d763cef2SBarry Smith PetscFunctionReturn(0); 975d763cef2SBarry Smith } 976d763cef2SBarry Smith 9774a2ae208SSatish Balay #undef __FUNCT__ 9784a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep" 979d763cef2SBarry Smith /*@ 980d763cef2SBarry Smith TSGetTimeStep - Gets the current timestep size. 981d763cef2SBarry Smith 982d763cef2SBarry Smith Not Collective 983d763cef2SBarry Smith 984d763cef2SBarry Smith Input Parameter: 985d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 986d763cef2SBarry Smith 987d763cef2SBarry Smith Output Parameter: 988d763cef2SBarry Smith . dt - the current timestep size 989d763cef2SBarry Smith 990d763cef2SBarry Smith Level: intermediate 991d763cef2SBarry Smith 992d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 993d763cef2SBarry Smith 994d763cef2SBarry Smith .keywords: TS, get, timestep 995d763cef2SBarry Smith @*/ 9967087cfbeSBarry Smith PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 997d763cef2SBarry Smith { 998d763cef2SBarry Smith PetscFunctionBegin; 9990700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10004482741eSBarry Smith PetscValidDoublePointer(dt,2); 1001d763cef2SBarry Smith *dt = ts->time_step; 1002d763cef2SBarry Smith PetscFunctionReturn(0); 1003d763cef2SBarry Smith } 1004d763cef2SBarry Smith 10054a2ae208SSatish Balay #undef __FUNCT__ 10064a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution" 1007d8e5e3e6SSatish Balay /*@ 1008d763cef2SBarry Smith TSGetSolution - Returns the solution at the present timestep. It 1009d763cef2SBarry Smith is valid to call this routine inside the function that you are evaluating 1010d763cef2SBarry Smith in order to move to the new timestep. This vector not changed until 1011d763cef2SBarry Smith the solution at the next timestep has been calculated. 1012d763cef2SBarry Smith 1013d763cef2SBarry Smith Not Collective, but Vec returned is parallel if TS is parallel 1014d763cef2SBarry Smith 1015d763cef2SBarry Smith Input Parameter: 1016d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1017d763cef2SBarry Smith 1018d763cef2SBarry Smith Output Parameter: 1019d763cef2SBarry Smith . v - the vector containing the solution 1020d763cef2SBarry Smith 1021d763cef2SBarry Smith Level: intermediate 1022d763cef2SBarry Smith 1023d763cef2SBarry Smith .seealso: TSGetTimeStep() 1024d763cef2SBarry Smith 1025d763cef2SBarry Smith .keywords: TS, timestep, get, solution 1026d763cef2SBarry Smith @*/ 10277087cfbeSBarry Smith PetscErrorCode TSGetSolution(TS ts,Vec *v) 1028d763cef2SBarry Smith { 1029d763cef2SBarry Smith PetscFunctionBegin; 10300700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10314482741eSBarry Smith PetscValidPointer(v,2); 10328737fe31SLisandro Dalcin *v = ts->vec_sol; 1033d763cef2SBarry Smith PetscFunctionReturn(0); 1034d763cef2SBarry Smith } 1035d763cef2SBarry Smith 1036bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */ 10374a2ae208SSatish Balay #undef __FUNCT__ 1038bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType" 1039d8e5e3e6SSatish Balay /*@ 1040bdad233fSMatthew Knepley TSSetProblemType - Sets the type of problem to be solved. 1041d763cef2SBarry Smith 1042bdad233fSMatthew Knepley Not collective 1043d763cef2SBarry Smith 1044bdad233fSMatthew Knepley Input Parameters: 1045bdad233fSMatthew Knepley + ts - The TS 1046bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1047d763cef2SBarry Smith .vb 1048d763cef2SBarry Smith U_t = A U 1049d763cef2SBarry Smith U_t = A(t) U 1050d763cef2SBarry Smith U_t = F(t,U) 1051d763cef2SBarry Smith .ve 1052d763cef2SBarry Smith 1053d763cef2SBarry Smith Level: beginner 1054d763cef2SBarry Smith 1055bdad233fSMatthew Knepley .keywords: TS, problem type 1056bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1057d763cef2SBarry Smith @*/ 10587087cfbeSBarry Smith PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1059a7cc72afSBarry Smith { 10609e2a6581SJed Brown PetscErrorCode ierr; 10619e2a6581SJed Brown 1062d763cef2SBarry Smith PetscFunctionBegin; 10630700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1064bdad233fSMatthew Knepley ts->problem_type = type; 10659e2a6581SJed Brown if (type == TS_LINEAR) { 10669e2a6581SJed Brown SNES snes; 10679e2a6581SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 10689e2a6581SJed Brown ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 10699e2a6581SJed Brown } 1070d763cef2SBarry Smith PetscFunctionReturn(0); 1071d763cef2SBarry Smith } 1072d763cef2SBarry Smith 1073bdad233fSMatthew Knepley #undef __FUNCT__ 1074bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType" 1075bdad233fSMatthew Knepley /*@C 1076bdad233fSMatthew Knepley TSGetProblemType - Gets the type of problem to be solved. 1077bdad233fSMatthew Knepley 1078bdad233fSMatthew Knepley Not collective 1079bdad233fSMatthew Knepley 1080bdad233fSMatthew Knepley Input Parameter: 1081bdad233fSMatthew Knepley . ts - The TS 1082bdad233fSMatthew Knepley 1083bdad233fSMatthew Knepley Output Parameter: 1084bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1085bdad233fSMatthew Knepley .vb 1086089b2837SJed Brown M U_t = A U 1087089b2837SJed Brown M(t) U_t = A(t) U 1088bdad233fSMatthew Knepley U_t = F(t,U) 1089bdad233fSMatthew Knepley .ve 1090bdad233fSMatthew Knepley 1091bdad233fSMatthew Knepley Level: beginner 1092bdad233fSMatthew Knepley 1093bdad233fSMatthew Knepley .keywords: TS, problem type 1094bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1095bdad233fSMatthew Knepley @*/ 10967087cfbeSBarry Smith PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1097a7cc72afSBarry Smith { 1098bdad233fSMatthew Knepley PetscFunctionBegin; 10990700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 11004482741eSBarry Smith PetscValidIntPointer(type,2); 1101bdad233fSMatthew Knepley *type = ts->problem_type; 1102bdad233fSMatthew Knepley PetscFunctionReturn(0); 1103bdad233fSMatthew Knepley } 1104d763cef2SBarry Smith 11054a2ae208SSatish Balay #undef __FUNCT__ 11064a2ae208SSatish Balay #define __FUNCT__ "TSSetUp" 1107d763cef2SBarry Smith /*@ 1108d763cef2SBarry Smith TSSetUp - Sets up the internal data structures for the later use 1109d763cef2SBarry Smith of a timestepper. 1110d763cef2SBarry Smith 1111d763cef2SBarry Smith Collective on TS 1112d763cef2SBarry Smith 1113d763cef2SBarry Smith Input Parameter: 1114d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1115d763cef2SBarry Smith 1116d763cef2SBarry Smith Notes: 1117d763cef2SBarry Smith For basic use of the TS solvers the user need not explicitly call 1118d763cef2SBarry Smith TSSetUp(), since these actions will automatically occur during 1119d763cef2SBarry Smith the call to TSStep(). However, if one wishes to control this 1120d763cef2SBarry Smith phase separately, TSSetUp() should be called after TSCreate() 1121d763cef2SBarry Smith and optional routines of the form TSSetXXX(), but before TSStep(). 1122d763cef2SBarry Smith 1123d763cef2SBarry Smith Level: advanced 1124d763cef2SBarry Smith 1125d763cef2SBarry Smith .keywords: TS, timestep, setup 1126d763cef2SBarry Smith 1127d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy() 1128d763cef2SBarry Smith @*/ 11297087cfbeSBarry Smith PetscErrorCode TSSetUp(TS ts) 1130d763cef2SBarry Smith { 1131dfbe8321SBarry Smith PetscErrorCode ierr; 1132d763cef2SBarry Smith 1133d763cef2SBarry Smith PetscFunctionBegin; 11340700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1135277b19d0SLisandro Dalcin if (ts->setupcalled) PetscFunctionReturn(0); 1136277b19d0SLisandro Dalcin 11377adad957SLisandro Dalcin if (!((PetscObject)ts)->type_name) { 11389596e0b4SJed Brown ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1139d763cef2SBarry Smith } 1140277b19d0SLisandro Dalcin 1141277b19d0SLisandro Dalcin if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1142277b19d0SLisandro Dalcin 1143277b19d0SLisandro Dalcin if (ts->ops->setup) { 1144000e7ae3SMatthew Knepley ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1145277b19d0SLisandro Dalcin } 1146277b19d0SLisandro Dalcin 1147277b19d0SLisandro Dalcin ts->setupcalled = PETSC_TRUE; 1148277b19d0SLisandro Dalcin PetscFunctionReturn(0); 1149277b19d0SLisandro Dalcin } 1150277b19d0SLisandro Dalcin 1151277b19d0SLisandro Dalcin #undef __FUNCT__ 1152277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset" 1153277b19d0SLisandro Dalcin /*@ 1154277b19d0SLisandro Dalcin TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1155277b19d0SLisandro Dalcin 1156277b19d0SLisandro Dalcin Collective on TS 1157277b19d0SLisandro Dalcin 1158277b19d0SLisandro Dalcin Input Parameter: 1159277b19d0SLisandro Dalcin . ts - the TS context obtained from TSCreate() 1160277b19d0SLisandro Dalcin 1161277b19d0SLisandro Dalcin Level: beginner 1162277b19d0SLisandro Dalcin 1163277b19d0SLisandro Dalcin .keywords: TS, timestep, reset 1164277b19d0SLisandro Dalcin 1165277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy() 1166277b19d0SLisandro Dalcin @*/ 1167277b19d0SLisandro Dalcin PetscErrorCode TSReset(TS ts) 1168277b19d0SLisandro Dalcin { 1169277b19d0SLisandro Dalcin PetscErrorCode ierr; 1170277b19d0SLisandro Dalcin 1171277b19d0SLisandro Dalcin PetscFunctionBegin; 1172277b19d0SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1173277b19d0SLisandro Dalcin if (ts->ops->reset) { 1174277b19d0SLisandro Dalcin ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1175277b19d0SLisandro Dalcin } 1176277b19d0SLisandro Dalcin if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 11774e684422SJed Brown ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 11784e684422SJed Brown ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 1179*214bc6a2SJed Brown ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr); 11806bf464f9SBarry Smith ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1181277b19d0SLisandro Dalcin if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);} 1182277b19d0SLisandro Dalcin ts->setupcalled = PETSC_FALSE; 1183d763cef2SBarry Smith PetscFunctionReturn(0); 1184d763cef2SBarry Smith } 1185d763cef2SBarry Smith 11864a2ae208SSatish Balay #undef __FUNCT__ 11874a2ae208SSatish Balay #define __FUNCT__ "TSDestroy" 1188d8e5e3e6SSatish Balay /*@ 1189d763cef2SBarry Smith TSDestroy - Destroys the timestepper context that was created 1190d763cef2SBarry Smith with TSCreate(). 1191d763cef2SBarry Smith 1192d763cef2SBarry Smith Collective on TS 1193d763cef2SBarry Smith 1194d763cef2SBarry Smith Input Parameter: 1195d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1196d763cef2SBarry Smith 1197d763cef2SBarry Smith Level: beginner 1198d763cef2SBarry Smith 1199d763cef2SBarry Smith .keywords: TS, timestepper, destroy 1200d763cef2SBarry Smith 1201d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve() 1202d763cef2SBarry Smith @*/ 12036bf464f9SBarry Smith PetscErrorCode TSDestroy(TS *ts) 1204d763cef2SBarry Smith { 12056849ba73SBarry Smith PetscErrorCode ierr; 1206d763cef2SBarry Smith 1207d763cef2SBarry Smith PetscFunctionBegin; 12086bf464f9SBarry Smith if (!*ts) PetscFunctionReturn(0); 12096bf464f9SBarry Smith PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 12106bf464f9SBarry Smith if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1211d763cef2SBarry Smith 12126bf464f9SBarry Smith ierr = TSReset((*ts));CHKERRQ(ierr); 1213277b19d0SLisandro Dalcin 1214be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 12156bf464f9SBarry Smith ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 12166bf464f9SBarry Smith if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 12176d4c513bSLisandro Dalcin 12186bf464f9SBarry Smith ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 12196bf464f9SBarry Smith ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 12206bf464f9SBarry Smith ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 12216d4c513bSLisandro Dalcin 1222a79aaaedSSatish Balay ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1223d763cef2SBarry Smith PetscFunctionReturn(0); 1224d763cef2SBarry Smith } 1225d763cef2SBarry Smith 12264a2ae208SSatish Balay #undef __FUNCT__ 12274a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES" 1228d8e5e3e6SSatish Balay /*@ 1229d763cef2SBarry Smith TSGetSNES - Returns the SNES (nonlinear solver) associated with 1230d763cef2SBarry Smith a TS (timestepper) context. Valid only for nonlinear problems. 1231d763cef2SBarry Smith 1232d763cef2SBarry Smith Not Collective, but SNES is parallel if TS is parallel 1233d763cef2SBarry Smith 1234d763cef2SBarry Smith Input Parameter: 1235d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1236d763cef2SBarry Smith 1237d763cef2SBarry Smith Output Parameter: 1238d763cef2SBarry Smith . snes - the nonlinear solver context 1239d763cef2SBarry Smith 1240d763cef2SBarry Smith Notes: 1241d763cef2SBarry Smith The user can then directly manipulate the SNES context to set various 1242d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 124394b7f48cSBarry Smith KSP, KSP, and PC contexts as well. 1244d763cef2SBarry Smith 1245d763cef2SBarry Smith TSGetSNES() does not work for integrators that do not use SNES; in 1246d763cef2SBarry Smith this case TSGetSNES() returns PETSC_NULL in snes. 1247d763cef2SBarry Smith 1248d763cef2SBarry Smith Level: beginner 1249d763cef2SBarry Smith 1250d763cef2SBarry Smith .keywords: timestep, get, SNES 1251d763cef2SBarry Smith @*/ 12527087cfbeSBarry Smith PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1253d763cef2SBarry Smith { 1254d372ba47SLisandro Dalcin PetscErrorCode ierr; 1255d372ba47SLisandro Dalcin 1256d763cef2SBarry Smith PetscFunctionBegin; 12570700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12584482741eSBarry Smith PetscValidPointer(snes,2); 1259d372ba47SLisandro Dalcin if (!ts->snes) { 1260d372ba47SLisandro Dalcin ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1261d372ba47SLisandro Dalcin ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1262d372ba47SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 12639e2a6581SJed Brown if (ts->problem_type == TS_LINEAR) { 12649e2a6581SJed Brown ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr); 12659e2a6581SJed Brown } 1266d372ba47SLisandro Dalcin } 1267d763cef2SBarry Smith *snes = ts->snes; 1268d763cef2SBarry Smith PetscFunctionReturn(0); 1269d763cef2SBarry Smith } 1270d763cef2SBarry Smith 12714a2ae208SSatish Balay #undef __FUNCT__ 127294b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP" 1273d8e5e3e6SSatish Balay /*@ 127494b7f48cSBarry Smith TSGetKSP - Returns the KSP (linear solver) associated with 1275d763cef2SBarry Smith a TS (timestepper) context. 1276d763cef2SBarry Smith 127794b7f48cSBarry Smith Not Collective, but KSP is parallel if TS is parallel 1278d763cef2SBarry Smith 1279d763cef2SBarry Smith Input Parameter: 1280d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1281d763cef2SBarry Smith 1282d763cef2SBarry Smith Output Parameter: 128394b7f48cSBarry Smith . ksp - the nonlinear solver context 1284d763cef2SBarry Smith 1285d763cef2SBarry Smith Notes: 128694b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 1287d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 1288d763cef2SBarry Smith KSP and PC contexts as well. 1289d763cef2SBarry Smith 129094b7f48cSBarry Smith TSGetKSP() does not work for integrators that do not use KSP; 129194b7f48cSBarry Smith in this case TSGetKSP() returns PETSC_NULL in ksp. 1292d763cef2SBarry Smith 1293d763cef2SBarry Smith Level: beginner 1294d763cef2SBarry Smith 129594b7f48cSBarry Smith .keywords: timestep, get, KSP 1296d763cef2SBarry Smith @*/ 12977087cfbeSBarry Smith PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1298d763cef2SBarry Smith { 1299d372ba47SLisandro Dalcin PetscErrorCode ierr; 1300089b2837SJed Brown SNES snes; 1301d372ba47SLisandro Dalcin 1302d763cef2SBarry Smith PetscFunctionBegin; 13030700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13044482741eSBarry Smith PetscValidPointer(ksp,2); 130517186662SBarry Smith if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1306e32f2f54SBarry Smith if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1307089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1308089b2837SJed Brown ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1309d763cef2SBarry Smith PetscFunctionReturn(0); 1310d763cef2SBarry Smith } 1311d763cef2SBarry Smith 1312d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */ 1313d763cef2SBarry Smith 13144a2ae208SSatish Balay #undef __FUNCT__ 1315adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration" 1316adb62b0dSMatthew Knepley /*@ 1317adb62b0dSMatthew Knepley TSGetDuration - Gets the maximum number of timesteps to use and 1318adb62b0dSMatthew Knepley maximum time for iteration. 1319adb62b0dSMatthew Knepley 13203f9fe445SBarry Smith Not Collective 1321adb62b0dSMatthew Knepley 1322adb62b0dSMatthew Knepley Input Parameters: 1323adb62b0dSMatthew Knepley + ts - the TS context obtained from TSCreate() 1324adb62b0dSMatthew Knepley . maxsteps - maximum number of iterations to use, or PETSC_NULL 1325adb62b0dSMatthew Knepley - maxtime - final time to iterate to, or PETSC_NULL 1326adb62b0dSMatthew Knepley 1327adb62b0dSMatthew Knepley Level: intermediate 1328adb62b0dSMatthew Knepley 1329adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time 1330adb62b0dSMatthew Knepley @*/ 13317087cfbeSBarry Smith PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1332adb62b0dSMatthew Knepley { 1333adb62b0dSMatthew Knepley PetscFunctionBegin; 13340700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1335abc0a331SBarry Smith if (maxsteps) { 13364482741eSBarry Smith PetscValidIntPointer(maxsteps,2); 1337adb62b0dSMatthew Knepley *maxsteps = ts->max_steps; 1338adb62b0dSMatthew Knepley } 1339abc0a331SBarry Smith if (maxtime ) { 13404482741eSBarry Smith PetscValidScalarPointer(maxtime,3); 1341adb62b0dSMatthew Knepley *maxtime = ts->max_time; 1342adb62b0dSMatthew Knepley } 1343adb62b0dSMatthew Knepley PetscFunctionReturn(0); 1344adb62b0dSMatthew Knepley } 1345adb62b0dSMatthew Knepley 1346adb62b0dSMatthew Knepley #undef __FUNCT__ 13474a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration" 1348d763cef2SBarry Smith /*@ 1349d763cef2SBarry Smith TSSetDuration - Sets the maximum number of timesteps to use and 1350d763cef2SBarry Smith maximum time for iteration. 1351d763cef2SBarry Smith 13523f9fe445SBarry Smith Logically Collective on TS 1353d763cef2SBarry Smith 1354d763cef2SBarry Smith Input Parameters: 1355d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1356d763cef2SBarry Smith . maxsteps - maximum number of iterations to use 1357d763cef2SBarry Smith - maxtime - final time to iterate to 1358d763cef2SBarry Smith 1359d763cef2SBarry Smith Options Database Keys: 1360d763cef2SBarry Smith . -ts_max_steps <maxsteps> - Sets maxsteps 1361d763cef2SBarry Smith . -ts_max_time <maxtime> - Sets maxtime 1362d763cef2SBarry Smith 1363d763cef2SBarry Smith Notes: 1364d763cef2SBarry Smith The default maximum number of iterations is 5000. Default time is 5.0 1365d763cef2SBarry Smith 1366d763cef2SBarry Smith Level: intermediate 1367d763cef2SBarry Smith 1368d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations 1369d763cef2SBarry Smith @*/ 13707087cfbeSBarry Smith PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1371d763cef2SBarry Smith { 1372d763cef2SBarry Smith PetscFunctionBegin; 13730700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1374c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1375c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,maxtime,2); 1376d763cef2SBarry Smith ts->max_steps = maxsteps; 1377d763cef2SBarry Smith ts->max_time = maxtime; 1378d763cef2SBarry Smith PetscFunctionReturn(0); 1379d763cef2SBarry Smith } 1380d763cef2SBarry Smith 13814a2ae208SSatish Balay #undef __FUNCT__ 13824a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution" 1383d763cef2SBarry Smith /*@ 1384d763cef2SBarry Smith TSSetSolution - Sets the initial solution vector 1385d763cef2SBarry Smith for use by the TS routines. 1386d763cef2SBarry Smith 13873f9fe445SBarry Smith Logically Collective on TS and Vec 1388d763cef2SBarry Smith 1389d763cef2SBarry Smith Input Parameters: 1390d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1391d763cef2SBarry Smith - x - the solution vector 1392d763cef2SBarry Smith 1393d763cef2SBarry Smith Level: beginner 1394d763cef2SBarry Smith 1395d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions 1396d763cef2SBarry Smith @*/ 13977087cfbeSBarry Smith PetscErrorCode TSSetSolution(TS ts,Vec x) 1398d763cef2SBarry Smith { 13998737fe31SLisandro Dalcin PetscErrorCode ierr; 14008737fe31SLisandro Dalcin 1401d763cef2SBarry Smith PetscFunctionBegin; 14020700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 14030700a824SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,2); 14048737fe31SLisandro Dalcin ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 14056bf464f9SBarry Smith ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 14068737fe31SLisandro Dalcin ts->vec_sol = x; 1407d763cef2SBarry Smith PetscFunctionReturn(0); 1408d763cef2SBarry Smith } 1409d763cef2SBarry Smith 1410e74ef692SMatthew Knepley #undef __FUNCT__ 1411e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep" 1412ac226902SBarry Smith /*@C 1413000e7ae3SMatthew Knepley TSSetPreStep - Sets the general-purpose function 14143f2090d5SJed Brown called once at the beginning of each time step. 1415000e7ae3SMatthew Knepley 14163f9fe445SBarry Smith Logically Collective on TS 1417000e7ae3SMatthew Knepley 1418000e7ae3SMatthew Knepley Input Parameters: 1419000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1420000e7ae3SMatthew Knepley - func - The function 1421000e7ae3SMatthew Knepley 1422000e7ae3SMatthew Knepley Calling sequence of func: 1423000e7ae3SMatthew Knepley . func (TS ts); 1424000e7ae3SMatthew Knepley 1425000e7ae3SMatthew Knepley Level: intermediate 1426000e7ae3SMatthew Knepley 1427000e7ae3SMatthew Knepley .keywords: TS, timestep 1428000e7ae3SMatthew Knepley @*/ 14297087cfbeSBarry Smith PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1430000e7ae3SMatthew Knepley { 1431000e7ae3SMatthew Knepley PetscFunctionBegin; 14320700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1433000e7ae3SMatthew Knepley ts->ops->prestep = func; 1434000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1435000e7ae3SMatthew Knepley } 1436000e7ae3SMatthew Knepley 1437e74ef692SMatthew Knepley #undef __FUNCT__ 14383f2090d5SJed Brown #define __FUNCT__ "TSPreStep" 14393f2090d5SJed Brown /*@C 14403f2090d5SJed Brown TSPreStep - Runs the user-defined pre-step function. 14413f2090d5SJed Brown 14423f2090d5SJed Brown Collective on TS 14433f2090d5SJed Brown 14443f2090d5SJed Brown Input Parameters: 14453f2090d5SJed Brown . ts - The TS context obtained from TSCreate() 14463f2090d5SJed Brown 14473f2090d5SJed Brown Notes: 14483f2090d5SJed Brown TSPreStep() is typically used within time stepping implementations, 14493f2090d5SJed Brown so most users would not generally call this routine themselves. 14503f2090d5SJed Brown 14513f2090d5SJed Brown Level: developer 14523f2090d5SJed Brown 14533f2090d5SJed Brown .keywords: TS, timestep 14543f2090d5SJed Brown @*/ 14557087cfbeSBarry Smith PetscErrorCode TSPreStep(TS ts) 14563f2090d5SJed Brown { 14573f2090d5SJed Brown PetscErrorCode ierr; 14583f2090d5SJed Brown 14593f2090d5SJed Brown PetscFunctionBegin; 14600700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 146172ac3e02SJed Brown if (ts->ops->prestep) { 14623f2090d5SJed Brown PetscStackPush("TS PreStep function"); 14633f2090d5SJed Brown ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 14643f2090d5SJed Brown PetscStackPop; 1465312ce896SJed Brown } 14663f2090d5SJed Brown PetscFunctionReturn(0); 14673f2090d5SJed Brown } 14683f2090d5SJed Brown 14693f2090d5SJed Brown #undef __FUNCT__ 1470e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep" 1471ac226902SBarry Smith /*@C 1472000e7ae3SMatthew Knepley TSSetPostStep - Sets the general-purpose function 14733f2090d5SJed Brown called once at the end of each time step. 1474000e7ae3SMatthew Knepley 14753f9fe445SBarry Smith Logically Collective on TS 1476000e7ae3SMatthew Knepley 1477000e7ae3SMatthew Knepley Input Parameters: 1478000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1479000e7ae3SMatthew Knepley - func - The function 1480000e7ae3SMatthew Knepley 1481000e7ae3SMatthew Knepley Calling sequence of func: 1482000e7ae3SMatthew Knepley . func (TS ts); 1483000e7ae3SMatthew Knepley 1484000e7ae3SMatthew Knepley Level: intermediate 1485000e7ae3SMatthew Knepley 1486000e7ae3SMatthew Knepley .keywords: TS, timestep 1487000e7ae3SMatthew Knepley @*/ 14887087cfbeSBarry Smith PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1489000e7ae3SMatthew Knepley { 1490000e7ae3SMatthew Knepley PetscFunctionBegin; 14910700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1492000e7ae3SMatthew Knepley ts->ops->poststep = func; 1493000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1494000e7ae3SMatthew Knepley } 1495000e7ae3SMatthew Knepley 1496e74ef692SMatthew Knepley #undef __FUNCT__ 14973f2090d5SJed Brown #define __FUNCT__ "TSPostStep" 14983f2090d5SJed Brown /*@C 14993f2090d5SJed Brown TSPostStep - Runs the user-defined post-step function. 15003f2090d5SJed Brown 15013f2090d5SJed Brown Collective on TS 15023f2090d5SJed Brown 15033f2090d5SJed Brown Input Parameters: 15043f2090d5SJed Brown . ts - The TS context obtained from TSCreate() 15053f2090d5SJed Brown 15063f2090d5SJed Brown Notes: 15073f2090d5SJed Brown TSPostStep() is typically used within time stepping implementations, 15083f2090d5SJed Brown so most users would not generally call this routine themselves. 15093f2090d5SJed Brown 15103f2090d5SJed Brown Level: developer 15113f2090d5SJed Brown 15123f2090d5SJed Brown .keywords: TS, timestep 15133f2090d5SJed Brown @*/ 15147087cfbeSBarry Smith PetscErrorCode TSPostStep(TS ts) 15153f2090d5SJed Brown { 15163f2090d5SJed Brown PetscErrorCode ierr; 15173f2090d5SJed Brown 15183f2090d5SJed Brown PetscFunctionBegin; 15190700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 152072ac3e02SJed Brown if (ts->ops->poststep) { 15213f2090d5SJed Brown PetscStackPush("TS PostStep function"); 15223f2090d5SJed Brown ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 15233f2090d5SJed Brown PetscStackPop; 152472ac3e02SJed Brown } 15253f2090d5SJed Brown PetscFunctionReturn(0); 15263f2090d5SJed Brown } 15273f2090d5SJed Brown 1528d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 1529d763cef2SBarry Smith 15304a2ae208SSatish Balay #undef __FUNCT__ 1531a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet" 1532d763cef2SBarry Smith /*@C 1533a6570f20SBarry Smith TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1534d763cef2SBarry Smith timestep to display the iteration's progress. 1535d763cef2SBarry Smith 15363f9fe445SBarry Smith Logically Collective on TS 1537d763cef2SBarry Smith 1538d763cef2SBarry Smith Input Parameters: 1539d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1540d763cef2SBarry Smith . func - monitoring routine 1541329f5518SBarry Smith . mctx - [optional] user-defined context for private data for the 1542b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desired) 1543b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1544b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 1545d763cef2SBarry Smith 1546d763cef2SBarry Smith Calling sequence of func: 1547a7cc72afSBarry Smith $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1548d763cef2SBarry Smith 1549d763cef2SBarry Smith + ts - the TS context 1550d763cef2SBarry Smith . steps - iteration number 15511f06c33eSBarry Smith . time - current time 1552d763cef2SBarry Smith . x - current iterate 1553d763cef2SBarry Smith - mctx - [optional] monitoring context 1554d763cef2SBarry Smith 1555d763cef2SBarry Smith Notes: 1556d763cef2SBarry Smith This routine adds an additional monitor to the list of monitors that 1557d763cef2SBarry Smith already has been loaded. 1558d763cef2SBarry Smith 1559025f1a04SBarry Smith Fortran notes: Only a single monitor function can be set for each TS object 1560025f1a04SBarry Smith 1561d763cef2SBarry Smith Level: intermediate 1562d763cef2SBarry Smith 1563d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1564d763cef2SBarry Smith 1565a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel() 1566d763cef2SBarry Smith @*/ 1567c2efdce3SBarry Smith PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1568d763cef2SBarry Smith { 1569d763cef2SBarry Smith PetscFunctionBegin; 15700700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 157117186662SBarry Smith if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1572d763cef2SBarry Smith ts->monitor[ts->numbermonitors] = monitor; 1573329f5518SBarry Smith ts->mdestroy[ts->numbermonitors] = mdestroy; 1574d763cef2SBarry Smith ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1575d763cef2SBarry Smith PetscFunctionReturn(0); 1576d763cef2SBarry Smith } 1577d763cef2SBarry Smith 15784a2ae208SSatish Balay #undef __FUNCT__ 1579a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel" 1580d763cef2SBarry Smith /*@C 1581a6570f20SBarry Smith TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1582d763cef2SBarry Smith 15833f9fe445SBarry Smith Logically Collective on TS 1584d763cef2SBarry Smith 1585d763cef2SBarry Smith Input Parameters: 1586d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1587d763cef2SBarry Smith 1588d763cef2SBarry Smith Notes: 1589d763cef2SBarry Smith There is no way to remove a single, specific monitor. 1590d763cef2SBarry Smith 1591d763cef2SBarry Smith Level: intermediate 1592d763cef2SBarry Smith 1593d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1594d763cef2SBarry Smith 1595a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet() 1596d763cef2SBarry Smith @*/ 15977087cfbeSBarry Smith PetscErrorCode TSMonitorCancel(TS ts) 1598d763cef2SBarry Smith { 1599d952e501SBarry Smith PetscErrorCode ierr; 1600d952e501SBarry Smith PetscInt i; 1601d952e501SBarry Smith 1602d763cef2SBarry Smith PetscFunctionBegin; 16030700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1604d952e501SBarry Smith for (i=0; i<ts->numbermonitors; i++) { 1605d952e501SBarry Smith if (ts->mdestroy[i]) { 16063c4aec1bSBarry Smith ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1607d952e501SBarry Smith } 1608d952e501SBarry Smith } 1609d763cef2SBarry Smith ts->numbermonitors = 0; 1610d763cef2SBarry Smith PetscFunctionReturn(0); 1611d763cef2SBarry Smith } 1612d763cef2SBarry Smith 16134a2ae208SSatish Balay #undef __FUNCT__ 1614a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault" 1615d8e5e3e6SSatish Balay /*@ 1616a6570f20SBarry Smith TSMonitorDefault - Sets the Default monitor 16175516499fSSatish Balay 16185516499fSSatish Balay Level: intermediate 161941251cbbSSatish Balay 16205516499fSSatish Balay .keywords: TS, set, monitor 16215516499fSSatish Balay 162241251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet() 162341251cbbSSatish Balay @*/ 1624649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1625d763cef2SBarry Smith { 1626dfbe8321SBarry Smith PetscErrorCode ierr; 1627649052a6SBarry Smith PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1628d132466eSBarry Smith 1629d763cef2SBarry Smith PetscFunctionBegin; 1630649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1631649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1632649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1633d763cef2SBarry Smith PetscFunctionReturn(0); 1634d763cef2SBarry Smith } 1635d763cef2SBarry Smith 16364a2ae208SSatish Balay #undef __FUNCT__ 16374a2ae208SSatish Balay #define __FUNCT__ "TSStep" 1638d763cef2SBarry Smith /*@ 1639d763cef2SBarry Smith TSStep - Steps the requested number of timesteps. 1640d763cef2SBarry Smith 1641d763cef2SBarry Smith Collective on TS 1642d763cef2SBarry Smith 1643d763cef2SBarry Smith Input Parameter: 1644d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1645d763cef2SBarry Smith 1646d763cef2SBarry Smith Output Parameters: 1647d763cef2SBarry Smith + steps - number of iterations until termination 1648142b95e3SSatish Balay - ptime - time until termination 1649d763cef2SBarry Smith 1650d763cef2SBarry Smith Level: beginner 1651d763cef2SBarry Smith 1652d763cef2SBarry Smith .keywords: TS, timestep, solve 1653d763cef2SBarry Smith 1654d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy() 1655d763cef2SBarry Smith @*/ 16567087cfbeSBarry Smith PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1657d763cef2SBarry Smith { 1658dfbe8321SBarry Smith PetscErrorCode ierr; 1659d763cef2SBarry Smith 1660d763cef2SBarry Smith PetscFunctionBegin; 16610700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1662277b19d0SLisandro Dalcin 1663d405a339SMatthew Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 1664d405a339SMatthew Knepley 1665d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1666000e7ae3SMatthew Knepley ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1667d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1668d405a339SMatthew Knepley 16694bb05414SBarry Smith if (!PetscPreLoadingOn) { 16707adad957SLisandro Dalcin ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1671d405a339SMatthew Knepley } 1672d763cef2SBarry Smith PetscFunctionReturn(0); 1673d763cef2SBarry Smith } 1674d763cef2SBarry Smith 16754a2ae208SSatish Balay #undef __FUNCT__ 16766a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve" 16776a4d4014SLisandro Dalcin /*@ 16786a4d4014SLisandro Dalcin TSSolve - Steps the requested number of timesteps. 16796a4d4014SLisandro Dalcin 16806a4d4014SLisandro Dalcin Collective on TS 16816a4d4014SLisandro Dalcin 16826a4d4014SLisandro Dalcin Input Parameter: 16836a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 16846a4d4014SLisandro Dalcin - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 16856a4d4014SLisandro Dalcin 16866a4d4014SLisandro Dalcin Level: beginner 16876a4d4014SLisandro Dalcin 16886a4d4014SLisandro Dalcin .keywords: TS, timestep, solve 16896a4d4014SLisandro Dalcin 16906a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep() 16916a4d4014SLisandro Dalcin @*/ 16927087cfbeSBarry Smith PetscErrorCode TSSolve(TS ts, Vec x) 16936a4d4014SLisandro Dalcin { 16946a4d4014SLisandro Dalcin PetscInt steps; 16956a4d4014SLisandro Dalcin PetscReal ptime; 16966a4d4014SLisandro Dalcin PetscErrorCode ierr; 1697f22f69f0SBarry Smith 16986a4d4014SLisandro Dalcin PetscFunctionBegin; 16990700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17006a4d4014SLisandro Dalcin /* set solution vector if provided */ 17016a4d4014SLisandro Dalcin if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 17026a4d4014SLisandro Dalcin /* reset time step and iteration counters */ 17036a4d4014SLisandro Dalcin ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 17046a4d4014SLisandro Dalcin /* steps the requested number of timesteps. */ 17056a4d4014SLisandro Dalcin ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 17066a4d4014SLisandro Dalcin PetscFunctionReturn(0); 17076a4d4014SLisandro Dalcin } 17086a4d4014SLisandro Dalcin 17096a4d4014SLisandro Dalcin #undef __FUNCT__ 17104a2ae208SSatish Balay #define __FUNCT__ "TSMonitor" 1711d763cef2SBarry Smith /* 1712d763cef2SBarry Smith Runs the user provided monitor routines, if they exists. 1713d763cef2SBarry Smith */ 1714a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1715d763cef2SBarry Smith { 17166849ba73SBarry Smith PetscErrorCode ierr; 1717a7cc72afSBarry Smith PetscInt i,n = ts->numbermonitors; 1718d763cef2SBarry Smith 1719d763cef2SBarry Smith PetscFunctionBegin; 1720d763cef2SBarry Smith for (i=0; i<n; i++) { 1721142b95e3SSatish Balay ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1722d763cef2SBarry Smith } 1723d763cef2SBarry Smith PetscFunctionReturn(0); 1724d763cef2SBarry Smith } 1725d763cef2SBarry Smith 1726d763cef2SBarry Smith /* ------------------------------------------------------------------------*/ 1727d763cef2SBarry Smith 17284a2ae208SSatish Balay #undef __FUNCT__ 1729a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate" 1730d763cef2SBarry Smith /*@C 1731a6570f20SBarry Smith TSMonitorLGCreate - Creates a line graph context for use with 1732d763cef2SBarry Smith TS to monitor convergence of preconditioned residual norms. 1733d763cef2SBarry Smith 1734d763cef2SBarry Smith Collective on TS 1735d763cef2SBarry Smith 1736d763cef2SBarry Smith Input Parameters: 1737d763cef2SBarry Smith + host - the X display to open, or null for the local machine 1738d763cef2SBarry Smith . label - the title to put in the title bar 17397c922b88SBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 1740d763cef2SBarry Smith - m, n - the screen width and height in pixels 1741d763cef2SBarry Smith 1742d763cef2SBarry Smith Output Parameter: 1743d763cef2SBarry Smith . draw - the drawing context 1744d763cef2SBarry Smith 1745d763cef2SBarry Smith Options Database Key: 1746a6570f20SBarry Smith . -ts_monitor_draw - automatically sets line graph monitor 1747d763cef2SBarry Smith 1748d763cef2SBarry Smith Notes: 1749a6570f20SBarry Smith Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1750d763cef2SBarry Smith 1751d763cef2SBarry Smith Level: intermediate 1752d763cef2SBarry Smith 17537c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso 1754d763cef2SBarry Smith 1755a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet() 17567c922b88SBarry Smith 1757d763cef2SBarry Smith @*/ 17587087cfbeSBarry Smith PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1759d763cef2SBarry Smith { 1760b0a32e0cSBarry Smith PetscDraw win; 1761dfbe8321SBarry Smith PetscErrorCode ierr; 1762d763cef2SBarry Smith 1763d763cef2SBarry Smith PetscFunctionBegin; 1764b0a32e0cSBarry Smith ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1765b0a32e0cSBarry Smith ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1766b0a32e0cSBarry Smith ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1767b0a32e0cSBarry Smith ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1768d763cef2SBarry Smith 176952e6d16bSBarry Smith ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1770d763cef2SBarry Smith PetscFunctionReturn(0); 1771d763cef2SBarry Smith } 1772d763cef2SBarry Smith 17734a2ae208SSatish Balay #undef __FUNCT__ 1774a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG" 1775a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1776d763cef2SBarry Smith { 1777b0a32e0cSBarry Smith PetscDrawLG lg = (PetscDrawLG) monctx; 177887828ca2SBarry Smith PetscReal x,y = ptime; 1779dfbe8321SBarry Smith PetscErrorCode ierr; 1780d763cef2SBarry Smith 1781d763cef2SBarry Smith PetscFunctionBegin; 17827c922b88SBarry Smith if (!monctx) { 17837c922b88SBarry Smith MPI_Comm comm; 1784b0a32e0cSBarry Smith PetscViewer viewer; 17857c922b88SBarry Smith 17867c922b88SBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1787b0a32e0cSBarry Smith viewer = PETSC_VIEWER_DRAW_(comm); 1788b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 17897c922b88SBarry Smith } 17907c922b88SBarry Smith 1791b0a32e0cSBarry Smith if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 179287828ca2SBarry Smith x = (PetscReal)n; 1793b0a32e0cSBarry Smith ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1794d763cef2SBarry Smith if (n < 20 || (n % 5)) { 1795b0a32e0cSBarry Smith ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1796d763cef2SBarry Smith } 1797d763cef2SBarry Smith PetscFunctionReturn(0); 1798d763cef2SBarry Smith } 1799d763cef2SBarry Smith 18004a2ae208SSatish Balay #undef __FUNCT__ 1801a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy" 1802d763cef2SBarry Smith /*@C 1803a6570f20SBarry Smith TSMonitorLGDestroy - Destroys a line graph context that was created 1804a6570f20SBarry Smith with TSMonitorLGCreate(). 1805d763cef2SBarry Smith 1806b0a32e0cSBarry Smith Collective on PetscDrawLG 1807d763cef2SBarry Smith 1808d763cef2SBarry Smith Input Parameter: 1809d763cef2SBarry Smith . draw - the drawing context 1810d763cef2SBarry Smith 1811d763cef2SBarry Smith Level: intermediate 1812d763cef2SBarry Smith 1813d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy 1814d763cef2SBarry Smith 1815a6570f20SBarry Smith .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1816d763cef2SBarry Smith @*/ 18176bf464f9SBarry Smith PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 1818d763cef2SBarry Smith { 1819b0a32e0cSBarry Smith PetscDraw draw; 1820dfbe8321SBarry Smith PetscErrorCode ierr; 1821d763cef2SBarry Smith 1822d763cef2SBarry Smith PetscFunctionBegin; 18236bf464f9SBarry Smith ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 18246bf464f9SBarry Smith ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 1825b0a32e0cSBarry Smith ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1826d763cef2SBarry Smith PetscFunctionReturn(0); 1827d763cef2SBarry Smith } 1828d763cef2SBarry Smith 18294a2ae208SSatish Balay #undef __FUNCT__ 18304a2ae208SSatish Balay #define __FUNCT__ "TSGetTime" 1831d763cef2SBarry Smith /*@ 1832d763cef2SBarry Smith TSGetTime - Gets the current time. 1833d763cef2SBarry Smith 1834d763cef2SBarry Smith Not Collective 1835d763cef2SBarry Smith 1836d763cef2SBarry Smith Input Parameter: 1837d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1838d763cef2SBarry Smith 1839d763cef2SBarry Smith Output Parameter: 1840d763cef2SBarry Smith . t - the current time 1841d763cef2SBarry Smith 1842d763cef2SBarry Smith Level: beginner 1843d763cef2SBarry Smith 1844d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1845d763cef2SBarry Smith 1846d763cef2SBarry Smith .keywords: TS, get, time 1847d763cef2SBarry Smith @*/ 18487087cfbeSBarry Smith PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1849d763cef2SBarry Smith { 1850d763cef2SBarry Smith PetscFunctionBegin; 18510700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 18524482741eSBarry Smith PetscValidDoublePointer(t,2); 1853d763cef2SBarry Smith *t = ts->ptime; 1854d763cef2SBarry Smith PetscFunctionReturn(0); 1855d763cef2SBarry Smith } 1856d763cef2SBarry Smith 18574a2ae208SSatish Balay #undef __FUNCT__ 18586a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime" 18596a4d4014SLisandro Dalcin /*@ 18606a4d4014SLisandro Dalcin TSSetTime - Allows one to reset the time. 18616a4d4014SLisandro Dalcin 18623f9fe445SBarry Smith Logically Collective on TS 18636a4d4014SLisandro Dalcin 18646a4d4014SLisandro Dalcin Input Parameters: 18656a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 18666a4d4014SLisandro Dalcin - time - the time 18676a4d4014SLisandro Dalcin 18686a4d4014SLisandro Dalcin Level: intermediate 18696a4d4014SLisandro Dalcin 18706a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration() 18716a4d4014SLisandro Dalcin 18726a4d4014SLisandro Dalcin .keywords: TS, set, time 18736a4d4014SLisandro Dalcin @*/ 18747087cfbeSBarry Smith PetscErrorCode TSSetTime(TS ts, PetscReal t) 18756a4d4014SLisandro Dalcin { 18766a4d4014SLisandro Dalcin PetscFunctionBegin; 18770700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1878c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,t,2); 18796a4d4014SLisandro Dalcin ts->ptime = t; 18806a4d4014SLisandro Dalcin PetscFunctionReturn(0); 18816a4d4014SLisandro Dalcin } 18826a4d4014SLisandro Dalcin 18836a4d4014SLisandro Dalcin #undef __FUNCT__ 18844a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix" 1885d763cef2SBarry Smith /*@C 1886d763cef2SBarry Smith TSSetOptionsPrefix - Sets the prefix used for searching for all 1887d763cef2SBarry Smith TS options in the database. 1888d763cef2SBarry Smith 18893f9fe445SBarry Smith Logically Collective on TS 1890d763cef2SBarry Smith 1891d763cef2SBarry Smith Input Parameter: 1892d763cef2SBarry Smith + ts - The TS context 1893d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1894d763cef2SBarry Smith 1895d763cef2SBarry Smith Notes: 1896d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1897d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1898d763cef2SBarry Smith hyphen. 1899d763cef2SBarry Smith 1900d763cef2SBarry Smith Level: advanced 1901d763cef2SBarry Smith 1902d763cef2SBarry Smith .keywords: TS, set, options, prefix, database 1903d763cef2SBarry Smith 1904d763cef2SBarry Smith .seealso: TSSetFromOptions() 1905d763cef2SBarry Smith 1906d763cef2SBarry Smith @*/ 19077087cfbeSBarry Smith PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1908d763cef2SBarry Smith { 1909dfbe8321SBarry Smith PetscErrorCode ierr; 1910089b2837SJed Brown SNES snes; 1911d763cef2SBarry Smith 1912d763cef2SBarry Smith PetscFunctionBegin; 19130700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1914d763cef2SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1915089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1916089b2837SJed Brown ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1917d763cef2SBarry Smith PetscFunctionReturn(0); 1918d763cef2SBarry Smith } 1919d763cef2SBarry Smith 1920d763cef2SBarry Smith 19214a2ae208SSatish Balay #undef __FUNCT__ 19224a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix" 1923d763cef2SBarry Smith /*@C 1924d763cef2SBarry Smith TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1925d763cef2SBarry Smith TS options in the database. 1926d763cef2SBarry Smith 19273f9fe445SBarry Smith Logically Collective on TS 1928d763cef2SBarry Smith 1929d763cef2SBarry Smith Input Parameter: 1930d763cef2SBarry Smith + ts - The TS context 1931d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1932d763cef2SBarry Smith 1933d763cef2SBarry Smith Notes: 1934d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1935d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1936d763cef2SBarry Smith hyphen. 1937d763cef2SBarry Smith 1938d763cef2SBarry Smith Level: advanced 1939d763cef2SBarry Smith 1940d763cef2SBarry Smith .keywords: TS, append, options, prefix, database 1941d763cef2SBarry Smith 1942d763cef2SBarry Smith .seealso: TSGetOptionsPrefix() 1943d763cef2SBarry Smith 1944d763cef2SBarry Smith @*/ 19457087cfbeSBarry Smith PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1946d763cef2SBarry Smith { 1947dfbe8321SBarry Smith PetscErrorCode ierr; 1948089b2837SJed Brown SNES snes; 1949d763cef2SBarry Smith 1950d763cef2SBarry Smith PetscFunctionBegin; 19510700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1952d763cef2SBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1953089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1954089b2837SJed Brown ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1955d763cef2SBarry Smith PetscFunctionReturn(0); 1956d763cef2SBarry Smith } 1957d763cef2SBarry Smith 19584a2ae208SSatish Balay #undef __FUNCT__ 19594a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix" 1960d763cef2SBarry Smith /*@C 1961d763cef2SBarry Smith TSGetOptionsPrefix - Sets the prefix used for searching for all 1962d763cef2SBarry Smith TS options in the database. 1963d763cef2SBarry Smith 1964d763cef2SBarry Smith Not Collective 1965d763cef2SBarry Smith 1966d763cef2SBarry Smith Input Parameter: 1967d763cef2SBarry Smith . ts - The TS context 1968d763cef2SBarry Smith 1969d763cef2SBarry Smith Output Parameter: 1970d763cef2SBarry Smith . prefix - A pointer to the prefix string used 1971d763cef2SBarry Smith 1972d763cef2SBarry Smith Notes: On the fortran side, the user should pass in a string 'prifix' of 1973d763cef2SBarry Smith sufficient length to hold the prefix. 1974d763cef2SBarry Smith 1975d763cef2SBarry Smith Level: intermediate 1976d763cef2SBarry Smith 1977d763cef2SBarry Smith .keywords: TS, get, options, prefix, database 1978d763cef2SBarry Smith 1979d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix() 1980d763cef2SBarry Smith @*/ 19817087cfbeSBarry Smith PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 1982d763cef2SBarry Smith { 1983dfbe8321SBarry Smith PetscErrorCode ierr; 1984d763cef2SBarry Smith 1985d763cef2SBarry Smith PetscFunctionBegin; 19860700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 19874482741eSBarry Smith PetscValidPointer(prefix,2); 1988d763cef2SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1989d763cef2SBarry Smith PetscFunctionReturn(0); 1990d763cef2SBarry Smith } 1991d763cef2SBarry Smith 19924a2ae208SSatish Balay #undef __FUNCT__ 19934a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian" 1994d763cef2SBarry Smith /*@C 1995d763cef2SBarry Smith TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1996d763cef2SBarry Smith 1997d763cef2SBarry Smith Not Collective, but parallel objects are returned if TS is parallel 1998d763cef2SBarry Smith 1999d763cef2SBarry Smith Input Parameter: 2000d763cef2SBarry Smith . ts - The TS context obtained from TSCreate() 2001d763cef2SBarry Smith 2002d763cef2SBarry Smith Output Parameters: 2003d763cef2SBarry Smith + J - The Jacobian J of F, where U_t = F(U,t) 2004d763cef2SBarry Smith . M - The preconditioner matrix, usually the same as J 2005089b2837SJed Brown . func - Function to compute the Jacobian of the RHS 2006d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine 2007d763cef2SBarry Smith 2008d763cef2SBarry Smith Notes: You can pass in PETSC_NULL for any return argument you do not need. 2009d763cef2SBarry Smith 2010d763cef2SBarry Smith Level: intermediate 2011d763cef2SBarry Smith 201226d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2013d763cef2SBarry Smith 2014d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian 2015d763cef2SBarry Smith @*/ 2016089b2837SJed Brown PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2017d763cef2SBarry Smith { 2018089b2837SJed Brown PetscErrorCode ierr; 2019089b2837SJed Brown SNES snes; 2020089b2837SJed Brown 2021d763cef2SBarry Smith PetscFunctionBegin; 2022089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2023089b2837SJed Brown ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2024089b2837SJed Brown if (func) *func = ts->ops->rhsjacobian; 202526d46c62SHong Zhang if (ctx) *ctx = ts->jacP; 2026d763cef2SBarry Smith PetscFunctionReturn(0); 2027d763cef2SBarry Smith } 2028d763cef2SBarry Smith 20291713a123SBarry Smith #undef __FUNCT__ 20302eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian" 20312eca1d9cSJed Brown /*@C 20322eca1d9cSJed Brown TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 20332eca1d9cSJed Brown 20342eca1d9cSJed Brown Not Collective, but parallel objects are returned if TS is parallel 20352eca1d9cSJed Brown 20362eca1d9cSJed Brown Input Parameter: 20372eca1d9cSJed Brown . ts - The TS context obtained from TSCreate() 20382eca1d9cSJed Brown 20392eca1d9cSJed Brown Output Parameters: 20402eca1d9cSJed Brown + A - The Jacobian of F(t,U,U_t) 20412eca1d9cSJed Brown . B - The preconditioner matrix, often the same as A 20422eca1d9cSJed Brown . f - The function to compute the matrices 20432eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine 20442eca1d9cSJed Brown 20452eca1d9cSJed Brown Notes: You can pass in PETSC_NULL for any return argument you do not need. 20462eca1d9cSJed Brown 20472eca1d9cSJed Brown Level: advanced 20482eca1d9cSJed Brown 20492eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 20502eca1d9cSJed Brown 20512eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian 20522eca1d9cSJed Brown @*/ 20537087cfbeSBarry Smith PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 20542eca1d9cSJed Brown { 2055089b2837SJed Brown PetscErrorCode ierr; 2056089b2837SJed Brown SNES snes; 2057089b2837SJed Brown 20582eca1d9cSJed Brown PetscFunctionBegin; 2059089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2060089b2837SJed Brown ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 20612eca1d9cSJed Brown if (f) *f = ts->ops->ijacobian; 20622eca1d9cSJed Brown if (ctx) *ctx = ts->jacP; 20632eca1d9cSJed Brown PetscFunctionReturn(0); 20642eca1d9cSJed Brown } 20652eca1d9cSJed Brown 20662eca1d9cSJed Brown #undef __FUNCT__ 2067a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution" 20681713a123SBarry Smith /*@C 2069a6570f20SBarry Smith TSMonitorSolution - Monitors progress of the TS solvers by calling 20701713a123SBarry Smith VecView() for the solution at each timestep 20711713a123SBarry Smith 20721713a123SBarry Smith Collective on TS 20731713a123SBarry Smith 20741713a123SBarry Smith Input Parameters: 20751713a123SBarry Smith + ts - the TS context 20761713a123SBarry Smith . step - current time-step 2077142b95e3SSatish Balay . ptime - current time 20781713a123SBarry Smith - dummy - either a viewer or PETSC_NULL 20791713a123SBarry Smith 20801713a123SBarry Smith Level: intermediate 20811713a123SBarry Smith 20821713a123SBarry Smith .keywords: TS, vector, monitor, view 20831713a123SBarry Smith 2084a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 20851713a123SBarry Smith @*/ 20867087cfbeSBarry Smith PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 20871713a123SBarry Smith { 2088dfbe8321SBarry Smith PetscErrorCode ierr; 20891713a123SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 20901713a123SBarry Smith 20911713a123SBarry Smith PetscFunctionBegin; 2092a34d58ebSBarry Smith if (!dummy) { 20937adad957SLisandro Dalcin viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 20941713a123SBarry Smith } 20951713a123SBarry Smith ierr = VecView(x,viewer);CHKERRQ(ierr); 20961713a123SBarry Smith PetscFunctionReturn(0); 20971713a123SBarry Smith } 20981713a123SBarry Smith 20991713a123SBarry Smith 21006c699258SBarry Smith #undef __FUNCT__ 21016c699258SBarry Smith #define __FUNCT__ "TSSetDM" 21026c699258SBarry Smith /*@ 21036c699258SBarry Smith TSSetDM - Sets the DM that may be used by some preconditioners 21046c699258SBarry Smith 21053f9fe445SBarry Smith Logically Collective on TS and DM 21066c699258SBarry Smith 21076c699258SBarry Smith Input Parameters: 21086c699258SBarry Smith + ts - the preconditioner context 21096c699258SBarry Smith - dm - the dm 21106c699258SBarry Smith 21116c699258SBarry Smith Level: intermediate 21126c699258SBarry Smith 21136c699258SBarry Smith 21146c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 21156c699258SBarry Smith @*/ 21167087cfbeSBarry Smith PetscErrorCode TSSetDM(TS ts,DM dm) 21176c699258SBarry Smith { 21186c699258SBarry Smith PetscErrorCode ierr; 2119089b2837SJed Brown SNES snes; 21206c699258SBarry Smith 21216c699258SBarry Smith PetscFunctionBegin; 21220700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 212370663e4aSLisandro Dalcin ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 21246bf464f9SBarry Smith ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 21256c699258SBarry Smith ts->dm = dm; 2126089b2837SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2127089b2837SJed Brown ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 21286c699258SBarry Smith PetscFunctionReturn(0); 21296c699258SBarry Smith } 21306c699258SBarry Smith 21316c699258SBarry Smith #undef __FUNCT__ 21326c699258SBarry Smith #define __FUNCT__ "TSGetDM" 21336c699258SBarry Smith /*@ 21346c699258SBarry Smith TSGetDM - Gets the DM that may be used by some preconditioners 21356c699258SBarry Smith 21363f9fe445SBarry Smith Not Collective 21376c699258SBarry Smith 21386c699258SBarry Smith Input Parameter: 21396c699258SBarry Smith . ts - the preconditioner context 21406c699258SBarry Smith 21416c699258SBarry Smith Output Parameter: 21426c699258SBarry Smith . dm - the dm 21436c699258SBarry Smith 21446c699258SBarry Smith Level: intermediate 21456c699258SBarry Smith 21466c699258SBarry Smith 21476c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 21486c699258SBarry Smith @*/ 21497087cfbeSBarry Smith PetscErrorCode TSGetDM(TS ts,DM *dm) 21506c699258SBarry Smith { 21516c699258SBarry Smith PetscFunctionBegin; 21520700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 21536c699258SBarry Smith *dm = ts->dm; 21546c699258SBarry Smith PetscFunctionReturn(0); 21556c699258SBarry Smith } 21561713a123SBarry Smith 21570f5c6efeSJed Brown #undef __FUNCT__ 21580f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction" 21590f5c6efeSJed Brown /*@ 21600f5c6efeSJed Brown SNESTSFormFunction - Function to evaluate nonlinear residual 21610f5c6efeSJed Brown 21623f9fe445SBarry Smith Logically Collective on SNES 21630f5c6efeSJed Brown 21640f5c6efeSJed Brown Input Parameter: 2165d42a1c89SJed Brown + snes - nonlinear solver 21660f5c6efeSJed Brown . X - the current state at which to evaluate the residual 2167d42a1c89SJed Brown - ctx - user context, must be a TS 21680f5c6efeSJed Brown 21690f5c6efeSJed Brown Output Parameter: 21700f5c6efeSJed Brown . F - the nonlinear residual 21710f5c6efeSJed Brown 21720f5c6efeSJed Brown Notes: 21730f5c6efeSJed Brown This function is not normally called by users and is automatically registered with the SNES used by TS. 21740f5c6efeSJed Brown It is most frequently passed to MatFDColoringSetFunction(). 21750f5c6efeSJed Brown 21760f5c6efeSJed Brown Level: advanced 21770f5c6efeSJed Brown 21780f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction() 21790f5c6efeSJed Brown @*/ 21807087cfbeSBarry Smith PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 21810f5c6efeSJed Brown { 21820f5c6efeSJed Brown TS ts = (TS)ctx; 21830f5c6efeSJed Brown PetscErrorCode ierr; 21840f5c6efeSJed Brown 21850f5c6efeSJed Brown PetscFunctionBegin; 21860f5c6efeSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 21870f5c6efeSJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 21880f5c6efeSJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 21890f5c6efeSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,4); 21900f5c6efeSJed Brown ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 21910f5c6efeSJed Brown PetscFunctionReturn(0); 21920f5c6efeSJed Brown } 21930f5c6efeSJed Brown 21940f5c6efeSJed Brown #undef __FUNCT__ 21950f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian" 21960f5c6efeSJed Brown /*@ 21970f5c6efeSJed Brown SNESTSFormJacobian - Function to evaluate the Jacobian 21980f5c6efeSJed Brown 21990f5c6efeSJed Brown Collective on SNES 22000f5c6efeSJed Brown 22010f5c6efeSJed Brown Input Parameter: 22020f5c6efeSJed Brown + snes - nonlinear solver 22030f5c6efeSJed Brown . X - the current state at which to evaluate the residual 22040f5c6efeSJed Brown - ctx - user context, must be a TS 22050f5c6efeSJed Brown 22060f5c6efeSJed Brown Output Parameter: 22070f5c6efeSJed Brown + A - the Jacobian 22080f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A) 22090f5c6efeSJed Brown - flag - indicates any structure change in the matrix 22100f5c6efeSJed Brown 22110f5c6efeSJed Brown Notes: 22120f5c6efeSJed Brown This function is not normally called by users and is automatically registered with the SNES used by TS. 22130f5c6efeSJed Brown 22140f5c6efeSJed Brown Level: developer 22150f5c6efeSJed Brown 22160f5c6efeSJed Brown .seealso: SNESSetJacobian() 22170f5c6efeSJed Brown @*/ 22187087cfbeSBarry Smith PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 22190f5c6efeSJed Brown { 22200f5c6efeSJed Brown TS ts = (TS)ctx; 22210f5c6efeSJed Brown PetscErrorCode ierr; 22220f5c6efeSJed Brown 22230f5c6efeSJed Brown PetscFunctionBegin; 22240f5c6efeSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 22250f5c6efeSJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 22260f5c6efeSJed Brown PetscValidPointer(A,3); 22270f5c6efeSJed Brown PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 22280f5c6efeSJed Brown PetscValidPointer(B,4); 22290f5c6efeSJed Brown PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 22300f5c6efeSJed Brown PetscValidPointer(flag,5); 22310f5c6efeSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,6); 22320f5c6efeSJed Brown ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 22330f5c6efeSJed Brown PetscFunctionReturn(0); 22340f5c6efeSJed Brown } 2235325fc9f4SBarry Smith 2236325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 2237c6db04a5SJed Brown #include <mex.h> 2238325fc9f4SBarry Smith 2239325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2240325fc9f4SBarry Smith 2241325fc9f4SBarry Smith #undef __FUNCT__ 2242325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab" 2243325fc9f4SBarry Smith /* 2244325fc9f4SBarry Smith TSComputeFunction_Matlab - Calls the function that has been set with 2245325fc9f4SBarry Smith TSSetFunctionMatlab(). 2246325fc9f4SBarry Smith 2247325fc9f4SBarry Smith Collective on TS 2248325fc9f4SBarry Smith 2249325fc9f4SBarry Smith Input Parameters: 2250325fc9f4SBarry Smith + snes - the TS context 2251325fc9f4SBarry Smith - x - input vector 2252325fc9f4SBarry Smith 2253325fc9f4SBarry Smith Output Parameter: 2254325fc9f4SBarry Smith . y - function vector, as set by TSSetFunction() 2255325fc9f4SBarry Smith 2256325fc9f4SBarry Smith Notes: 2257325fc9f4SBarry Smith TSComputeFunction() is typically used within nonlinear solvers 2258325fc9f4SBarry Smith implementations, so most users would not generally call this routine 2259325fc9f4SBarry Smith themselves. 2260325fc9f4SBarry Smith 2261325fc9f4SBarry Smith Level: developer 2262325fc9f4SBarry Smith 2263325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function 2264325fc9f4SBarry Smith 2265325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction() 2266325fc9f4SBarry Smith */ 22677087cfbeSBarry Smith PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2268325fc9f4SBarry Smith { 2269325fc9f4SBarry Smith PetscErrorCode ierr; 2270325fc9f4SBarry Smith TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2271325fc9f4SBarry Smith int nlhs = 1,nrhs = 7; 2272325fc9f4SBarry Smith mxArray *plhs[1],*prhs[7]; 2273325fc9f4SBarry Smith long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2274325fc9f4SBarry Smith 2275325fc9f4SBarry Smith PetscFunctionBegin; 2276325fc9f4SBarry Smith PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2277325fc9f4SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2278325fc9f4SBarry Smith PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2279325fc9f4SBarry Smith PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2280325fc9f4SBarry Smith PetscCheckSameComm(snes,1,x,3); 2281325fc9f4SBarry Smith PetscCheckSameComm(snes,1,y,5); 2282325fc9f4SBarry Smith 2283325fc9f4SBarry Smith ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2284325fc9f4SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2285880f3077SBarry Smith ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2286325fc9f4SBarry Smith ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2287325fc9f4SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 2288325fc9f4SBarry Smith prhs[1] = mxCreateDoubleScalar(time); 2289325fc9f4SBarry Smith prhs[2] = mxCreateDoubleScalar((double)lx); 2290325fc9f4SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lxdot); 2291325fc9f4SBarry Smith prhs[4] = mxCreateDoubleScalar((double)ly); 2292325fc9f4SBarry Smith prhs[5] = mxCreateString(sctx->funcname); 2293325fc9f4SBarry Smith prhs[6] = sctx->ctx; 2294325fc9f4SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2295325fc9f4SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2296325fc9f4SBarry Smith mxDestroyArray(prhs[0]); 2297325fc9f4SBarry Smith mxDestroyArray(prhs[1]); 2298325fc9f4SBarry Smith mxDestroyArray(prhs[2]); 2299325fc9f4SBarry Smith mxDestroyArray(prhs[3]); 2300325fc9f4SBarry Smith mxDestroyArray(prhs[4]); 2301325fc9f4SBarry Smith mxDestroyArray(prhs[5]); 2302325fc9f4SBarry Smith mxDestroyArray(plhs[0]); 2303325fc9f4SBarry Smith PetscFunctionReturn(0); 2304325fc9f4SBarry Smith } 2305325fc9f4SBarry Smith 2306325fc9f4SBarry Smith 2307325fc9f4SBarry Smith #undef __FUNCT__ 2308325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab" 2309325fc9f4SBarry Smith /* 2310325fc9f4SBarry Smith TSSetFunctionMatlab - Sets the function evaluation routine and function 2311325fc9f4SBarry Smith vector for use by the TS routines in solving ODEs 2312e3c5b3baSBarry Smith equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2313325fc9f4SBarry Smith 2314325fc9f4SBarry Smith Logically Collective on TS 2315325fc9f4SBarry Smith 2316325fc9f4SBarry Smith Input Parameters: 2317325fc9f4SBarry Smith + ts - the TS context 2318325fc9f4SBarry Smith - func - function evaluation routine 2319325fc9f4SBarry Smith 2320325fc9f4SBarry Smith Calling sequence of func: 2321325fc9f4SBarry Smith $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2322325fc9f4SBarry Smith 2323325fc9f4SBarry Smith Level: beginner 2324325fc9f4SBarry Smith 2325325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function 2326325fc9f4SBarry Smith 2327325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2328325fc9f4SBarry Smith */ 23297087cfbeSBarry Smith PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2330325fc9f4SBarry Smith { 2331325fc9f4SBarry Smith PetscErrorCode ierr; 2332325fc9f4SBarry Smith TSMatlabContext *sctx; 2333325fc9f4SBarry Smith 2334325fc9f4SBarry Smith PetscFunctionBegin; 2335325fc9f4SBarry Smith /* currently sctx is memory bleed */ 2336325fc9f4SBarry Smith ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2337325fc9f4SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2338325fc9f4SBarry Smith /* 2339325fc9f4SBarry Smith This should work, but it doesn't 2340325fc9f4SBarry Smith sctx->ctx = ctx; 2341325fc9f4SBarry Smith mexMakeArrayPersistent(sctx->ctx); 2342325fc9f4SBarry Smith */ 2343325fc9f4SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 2344325fc9f4SBarry Smith ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2345325fc9f4SBarry Smith PetscFunctionReturn(0); 2346325fc9f4SBarry Smith } 2347325fc9f4SBarry Smith 2348325fc9f4SBarry Smith #undef __FUNCT__ 2349325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab" 2350325fc9f4SBarry Smith /* 2351325fc9f4SBarry Smith TSComputeJacobian_Matlab - Calls the function that has been set with 2352325fc9f4SBarry Smith TSSetJacobianMatlab(). 2353325fc9f4SBarry Smith 2354325fc9f4SBarry Smith Collective on TS 2355325fc9f4SBarry Smith 2356325fc9f4SBarry Smith Input Parameters: 2357325fc9f4SBarry Smith + snes - the TS context 2358325fc9f4SBarry Smith . x - input vector 2359325fc9f4SBarry Smith . A, B - the matrices 2360325fc9f4SBarry Smith - ctx - user context 2361325fc9f4SBarry Smith 2362325fc9f4SBarry Smith Output Parameter: 2363325fc9f4SBarry Smith . flag - structure of the matrix 2364325fc9f4SBarry Smith 2365325fc9f4SBarry Smith Level: developer 2366325fc9f4SBarry Smith 2367325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function 2368325fc9f4SBarry Smith 2369325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction() 2370325fc9f4SBarry Smith @*/ 23717087cfbeSBarry Smith PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2372325fc9f4SBarry Smith { 2373325fc9f4SBarry Smith PetscErrorCode ierr; 2374325fc9f4SBarry Smith TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2375325fc9f4SBarry Smith int nlhs = 2,nrhs = 9; 2376325fc9f4SBarry Smith mxArray *plhs[2],*prhs[9]; 2377325fc9f4SBarry Smith long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2378325fc9f4SBarry Smith 2379325fc9f4SBarry Smith PetscFunctionBegin; 2380325fc9f4SBarry Smith PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2381325fc9f4SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2382325fc9f4SBarry Smith 2383325fc9f4SBarry Smith /* call Matlab function in ctx with arguments x and y */ 2384325fc9f4SBarry Smith 2385325fc9f4SBarry Smith ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2386325fc9f4SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2387325fc9f4SBarry Smith ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2388325fc9f4SBarry Smith ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2389325fc9f4SBarry Smith ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2390325fc9f4SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 2391325fc9f4SBarry Smith prhs[1] = mxCreateDoubleScalar((double)time); 2392325fc9f4SBarry Smith prhs[2] = mxCreateDoubleScalar((double)lx); 2393325fc9f4SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lxdot); 2394325fc9f4SBarry Smith prhs[4] = mxCreateDoubleScalar((double)shift); 2395325fc9f4SBarry Smith prhs[5] = mxCreateDoubleScalar((double)lA); 2396325fc9f4SBarry Smith prhs[6] = mxCreateDoubleScalar((double)lB); 2397325fc9f4SBarry Smith prhs[7] = mxCreateString(sctx->funcname); 2398325fc9f4SBarry Smith prhs[8] = sctx->ctx; 2399325fc9f4SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2400325fc9f4SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2401325fc9f4SBarry Smith *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2402325fc9f4SBarry Smith mxDestroyArray(prhs[0]); 2403325fc9f4SBarry Smith mxDestroyArray(prhs[1]); 2404325fc9f4SBarry Smith mxDestroyArray(prhs[2]); 2405325fc9f4SBarry Smith mxDestroyArray(prhs[3]); 2406325fc9f4SBarry Smith mxDestroyArray(prhs[4]); 2407325fc9f4SBarry Smith mxDestroyArray(prhs[5]); 2408325fc9f4SBarry Smith mxDestroyArray(prhs[6]); 2409325fc9f4SBarry Smith mxDestroyArray(prhs[7]); 2410325fc9f4SBarry Smith mxDestroyArray(plhs[0]); 2411325fc9f4SBarry Smith mxDestroyArray(plhs[1]); 2412325fc9f4SBarry Smith PetscFunctionReturn(0); 2413325fc9f4SBarry Smith } 2414325fc9f4SBarry Smith 2415325fc9f4SBarry Smith 2416325fc9f4SBarry Smith #undef __FUNCT__ 2417325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab" 2418325fc9f4SBarry Smith /* 2419325fc9f4SBarry Smith TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2420e3c5b3baSBarry Smith vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function 2421325fc9f4SBarry Smith 2422325fc9f4SBarry Smith Logically Collective on TS 2423325fc9f4SBarry Smith 2424325fc9f4SBarry Smith Input Parameters: 2425325fc9f4SBarry Smith + snes - the TS context 2426325fc9f4SBarry Smith . A,B - Jacobian matrices 2427325fc9f4SBarry Smith . func - function evaluation routine 2428325fc9f4SBarry Smith - ctx - user context 2429325fc9f4SBarry Smith 2430325fc9f4SBarry Smith Calling sequence of func: 2431325fc9f4SBarry Smith $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2432325fc9f4SBarry Smith 2433325fc9f4SBarry Smith 2434325fc9f4SBarry Smith Level: developer 2435325fc9f4SBarry Smith 2436325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function 2437325fc9f4SBarry Smith 2438325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2439325fc9f4SBarry Smith */ 24407087cfbeSBarry Smith PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2441325fc9f4SBarry Smith { 2442325fc9f4SBarry Smith PetscErrorCode ierr; 2443325fc9f4SBarry Smith TSMatlabContext *sctx; 2444325fc9f4SBarry Smith 2445325fc9f4SBarry Smith PetscFunctionBegin; 2446325fc9f4SBarry Smith /* currently sctx is memory bleed */ 2447325fc9f4SBarry Smith ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2448325fc9f4SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2449325fc9f4SBarry Smith /* 2450325fc9f4SBarry Smith This should work, but it doesn't 2451325fc9f4SBarry Smith sctx->ctx = ctx; 2452325fc9f4SBarry Smith mexMakeArrayPersistent(sctx->ctx); 2453325fc9f4SBarry Smith */ 2454325fc9f4SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 2455325fc9f4SBarry Smith ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2456325fc9f4SBarry Smith PetscFunctionReturn(0); 2457325fc9f4SBarry Smith } 2458325fc9f4SBarry Smith 2459b5b1a830SBarry Smith #undef __FUNCT__ 2460b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab" 2461b5b1a830SBarry Smith /* 2462b5b1a830SBarry Smith TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2463b5b1a830SBarry Smith 2464b5b1a830SBarry Smith Collective on TS 2465b5b1a830SBarry Smith 2466b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction() 2467b5b1a830SBarry Smith @*/ 24687087cfbeSBarry Smith PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2469b5b1a830SBarry Smith { 2470b5b1a830SBarry Smith PetscErrorCode ierr; 2471b5b1a830SBarry Smith TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2472a530c242SBarry Smith int nlhs = 1,nrhs = 6; 2473b5b1a830SBarry Smith mxArray *plhs[1],*prhs[6]; 2474b5b1a830SBarry Smith long long int lx = 0,ls = 0; 2475b5b1a830SBarry Smith 2476b5b1a830SBarry Smith PetscFunctionBegin; 2477b5b1a830SBarry Smith PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2478b5b1a830SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2479b5b1a830SBarry Smith 2480b5b1a830SBarry Smith ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2481b5b1a830SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2482b5b1a830SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 2483b5b1a830SBarry Smith prhs[1] = mxCreateDoubleScalar((double)it); 2484b5b1a830SBarry Smith prhs[2] = mxCreateDoubleScalar((double)time); 2485b5b1a830SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lx); 2486b5b1a830SBarry Smith prhs[4] = mxCreateString(sctx->funcname); 2487b5b1a830SBarry Smith prhs[5] = sctx->ctx; 2488b5b1a830SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2489b5b1a830SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2490b5b1a830SBarry Smith mxDestroyArray(prhs[0]); 2491b5b1a830SBarry Smith mxDestroyArray(prhs[1]); 2492b5b1a830SBarry Smith mxDestroyArray(prhs[2]); 2493b5b1a830SBarry Smith mxDestroyArray(prhs[3]); 2494b5b1a830SBarry Smith mxDestroyArray(prhs[4]); 2495b5b1a830SBarry Smith mxDestroyArray(plhs[0]); 2496b5b1a830SBarry Smith PetscFunctionReturn(0); 2497b5b1a830SBarry Smith } 2498b5b1a830SBarry Smith 2499b5b1a830SBarry Smith 2500b5b1a830SBarry Smith #undef __FUNCT__ 2501b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab" 2502b5b1a830SBarry Smith /* 2503b5b1a830SBarry Smith TSMonitorSetMatlab - Sets the monitor function from Matlab 2504b5b1a830SBarry Smith 2505b5b1a830SBarry Smith Level: developer 2506b5b1a830SBarry Smith 2507b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function 2508b5b1a830SBarry Smith 2509b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2510b5b1a830SBarry Smith */ 25117087cfbeSBarry Smith PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2512b5b1a830SBarry Smith { 2513b5b1a830SBarry Smith PetscErrorCode ierr; 2514b5b1a830SBarry Smith TSMatlabContext *sctx; 2515b5b1a830SBarry Smith 2516b5b1a830SBarry Smith PetscFunctionBegin; 2517b5b1a830SBarry Smith /* currently sctx is memory bleed */ 2518b5b1a830SBarry Smith ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2519b5b1a830SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2520b5b1a830SBarry Smith /* 2521b5b1a830SBarry Smith This should work, but it doesn't 2522b5b1a830SBarry Smith sctx->ctx = ctx; 2523b5b1a830SBarry Smith mexMakeArrayPersistent(sctx->ctx); 2524b5b1a830SBarry Smith */ 2525b5b1a830SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 2526b5b1a830SBarry Smith ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2527b5b1a830SBarry Smith PetscFunctionReturn(0); 2528b5b1a830SBarry Smith } 2529b5b1a830SBarry Smith 2530325fc9f4SBarry Smith #endif 2531