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]; 78bdad233fSMatthew Knepley 79bdad233fSMatthew Knepley PetscFunctionBegin; 800700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 817adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr); 82bdad233fSMatthew Knepley 83bdad233fSMatthew Knepley /* Handle generic TS options */ 84bdad233fSMatthew Knepley ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 85bdad233fSMatthew Knepley ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 86bdad233fSMatthew Knepley ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr); 87bdad233fSMatthew Knepley ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr); 88a7cc72afSBarry Smith if (opt) { 89bdad233fSMatthew Knepley ts->initial_time_step = ts->time_step = dt; 90bdad233fSMatthew Knepley } 91bdad233fSMatthew Knepley 92bdad233fSMatthew Knepley /* Monitor options */ 93a6570f20SBarry Smith ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 94eabae89aSBarry Smith if (flg) { 95649052a6SBarry Smith ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 96649052a6SBarry Smith ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 97bdad233fSMatthew Knepley } 98*5180491cSLisandro Dalcin ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 99*5180491cSLisandro Dalcin if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);} 100*5180491cSLisandro Dalcin 10190d69ab7SBarry Smith opt = PETSC_FALSE; 102acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 103a7cc72afSBarry Smith if (opt) { 104a6570f20SBarry Smith ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 105bdad233fSMatthew Knepley } 10690d69ab7SBarry Smith opt = PETSC_FALSE; 107acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 108a7cc72afSBarry Smith if (opt) { 109a6570f20SBarry Smith ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 110bdad233fSMatthew Knepley } 111bdad233fSMatthew Knepley 112bdad233fSMatthew Knepley /* Handle TS type options */ 113bdad233fSMatthew Knepley ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 114bdad233fSMatthew Knepley 115bdad233fSMatthew Knepley /* Handle specific TS options */ 116abc0a331SBarry Smith if (ts->ops->setfromoptions) { 117bdad233fSMatthew Knepley ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 118bdad233fSMatthew Knepley } 1195d973c19SBarry Smith 1205d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1215d973c19SBarry Smith ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 122bdad233fSMatthew Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 123bdad233fSMatthew Knepley 124bdad233fSMatthew Knepley /* Handle subobject options */ 125bdad233fSMatthew Knepley switch(ts->problem_type) { 126156fc9a6SMatthew Knepley /* Should check for implicit/explicit */ 127bdad233fSMatthew Knepley case TS_LINEAR: 128abc0a331SBarry Smith if (ts->ksp) { 1298beb423aSHong Zhang ierr = KSPSetOperators(ts->ksp,ts->Arhs,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 13094b7f48cSBarry Smith ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr); 131156fc9a6SMatthew Knepley } 132bdad233fSMatthew Knepley break; 133bdad233fSMatthew Knepley case TS_NONLINEAR: 134abc0a331SBarry Smith if (ts->snes) { 1357c236d22SBarry Smith /* this is a bit of a hack, but it gets the matrix information into SNES earlier 1367c236d22SBarry Smith so that SNES and KSP have more information to pick reasonable defaults 1377c0b301bSJed Brown before they allow users to set options 1387c0b301bSJed Brown * If ts->A has been set at this point, we are probably using the implicit form 1397c0b301bSJed Brown and Arhs will never be used. */ 1407c0b301bSJed Brown ierr = SNESSetJacobian(ts->snes,ts->A?ts->A:ts->Arhs,ts->B,0,ts);CHKERRQ(ierr); 141bdad233fSMatthew Knepley ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 142156fc9a6SMatthew Knepley } 143bdad233fSMatthew Knepley break; 144bdad233fSMatthew Knepley default: 145e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type); 146bdad233fSMatthew Knepley } 147bdad233fSMatthew Knepley 148bdad233fSMatthew Knepley PetscFunctionReturn(0); 149bdad233fSMatthew Knepley } 150bdad233fSMatthew Knepley 151bdad233fSMatthew Knepley #undef __FUNCT__ 152bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions" 153bdad233fSMatthew Knepley /*@ 154bdad233fSMatthew Knepley TSViewFromOptions - This function visualizes the ts based upon user options. 155bdad233fSMatthew Knepley 156bdad233fSMatthew Knepley Collective on TS 157bdad233fSMatthew Knepley 158bdad233fSMatthew Knepley Input Parameter: 159bdad233fSMatthew Knepley . ts - The ts 160bdad233fSMatthew Knepley 161bdad233fSMatthew Knepley Level: intermediate 162bdad233fSMatthew Knepley 163bdad233fSMatthew Knepley .keywords: TS, view, options, database 164bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView() 165bdad233fSMatthew Knepley @*/ 1667087cfbeSBarry Smith PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 167bdad233fSMatthew Knepley { 168bdad233fSMatthew Knepley PetscViewer viewer; 169bdad233fSMatthew Knepley PetscDraw draw; 170ace3abfcSBarry Smith PetscBool opt = PETSC_FALSE; 171e10c49a3SBarry Smith char fileName[PETSC_MAX_PATH_LEN]; 172dfbe8321SBarry Smith PetscErrorCode ierr; 173bdad233fSMatthew Knepley 174bdad233fSMatthew Knepley PetscFunctionBegin; 1757adad957SLisandro Dalcin ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 176eabae89aSBarry Smith if (opt && !PetscPreLoadingOn) { 1777adad957SLisandro Dalcin ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 178bdad233fSMatthew Knepley ierr = TSView(ts, viewer);CHKERRQ(ierr); 1796bf464f9SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 180bdad233fSMatthew Knepley } 1818e83347fSKai Germaschewski opt = PETSC_FALSE; 182acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 183a7cc72afSBarry Smith if (opt) { 1847adad957SLisandro Dalcin ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 185bdad233fSMatthew Knepley ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 186a7cc72afSBarry Smith if (title) { 1871836bdbcSSatish Balay ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 188bdad233fSMatthew Knepley } else { 189bdad233fSMatthew Knepley ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 1907adad957SLisandro Dalcin ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 191bdad233fSMatthew Knepley } 192bdad233fSMatthew Knepley ierr = TSView(ts, viewer);CHKERRQ(ierr); 193bdad233fSMatthew Knepley ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 194bdad233fSMatthew Knepley ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1956bf464f9SBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 196bdad233fSMatthew Knepley } 197bdad233fSMatthew Knepley PetscFunctionReturn(0); 198bdad233fSMatthew Knepley } 199bdad233fSMatthew Knepley 200bdad233fSMatthew Knepley #undef __FUNCT__ 2014a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian" 202a7a1495cSBarry Smith /*@ 2038c385f81SBarry Smith TSComputeRHSJacobian - Computes the Jacobian matrix that has been 204a7a1495cSBarry Smith set with TSSetRHSJacobian(). 205a7a1495cSBarry Smith 206a7a1495cSBarry Smith Collective on TS and Vec 207a7a1495cSBarry Smith 208a7a1495cSBarry Smith Input Parameters: 209316643e7SJed Brown + ts - the TS context 210a7a1495cSBarry Smith . t - current timestep 211a7a1495cSBarry Smith - x - input vector 212a7a1495cSBarry Smith 213a7a1495cSBarry Smith Output Parameters: 214a7a1495cSBarry Smith + A - Jacobian matrix 215a7a1495cSBarry Smith . B - optional preconditioning matrix 216a7a1495cSBarry Smith - flag - flag indicating matrix structure 217a7a1495cSBarry Smith 218a7a1495cSBarry Smith Notes: 219a7a1495cSBarry Smith Most users should not need to explicitly call this routine, as it 220a7a1495cSBarry Smith is used internally within the nonlinear solvers. 221a7a1495cSBarry Smith 22294b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 223a7a1495cSBarry Smith flag parameter. 224a7a1495cSBarry Smith 225a7a1495cSBarry Smith TSComputeJacobian() is valid only for TS_NONLINEAR 226a7a1495cSBarry Smith 227a7a1495cSBarry Smith Level: developer 228a7a1495cSBarry Smith 229a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix 230a7a1495cSBarry Smith 23194b7f48cSBarry Smith .seealso: TSSetRHSJacobian(), KSPSetOperators() 232a7a1495cSBarry Smith @*/ 2337087cfbeSBarry Smith PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 234a7a1495cSBarry Smith { 235dfbe8321SBarry Smith PetscErrorCode ierr; 236a7a1495cSBarry Smith 237a7a1495cSBarry Smith PetscFunctionBegin; 2380700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2390700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 240c9780b6fSBarry Smith PetscCheckSameComm(ts,1,X,3); 24117186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only"); 242000e7ae3SMatthew Knepley if (ts->ops->rhsjacobian) { 243d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 244a7a1495cSBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 245a7a1495cSBarry Smith PetscStackPush("TS user Jacobian function"); 246000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 247a7a1495cSBarry Smith PetscStackPop; 248d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 249a7a1495cSBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 2500700a824SBarry Smith PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 2510700a824SBarry Smith PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 252ef66eb69SBarry Smith } else { 253ef66eb69SBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 254ef66eb69SBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255ef66eb69SBarry Smith if (*A != *B) { 256ef66eb69SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 257ef66eb69SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 258ef66eb69SBarry Smith } 259ef66eb69SBarry Smith } 260a7a1495cSBarry Smith PetscFunctionReturn(0); 261a7a1495cSBarry Smith } 262a7a1495cSBarry Smith 2634a2ae208SSatish Balay #undef __FUNCT__ 2644a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction" 265316643e7SJed Brown /*@ 266d763cef2SBarry Smith TSComputeRHSFunction - Evaluates the right-hand-side function. 267d763cef2SBarry Smith 268316643e7SJed Brown Collective on TS and Vec 269316643e7SJed Brown 270316643e7SJed Brown Input Parameters: 271316643e7SJed Brown + ts - the TS context 272316643e7SJed Brown . t - current time 273316643e7SJed Brown - x - state vector 274316643e7SJed Brown 275316643e7SJed Brown Output Parameter: 276316643e7SJed Brown . y - right hand side 277316643e7SJed Brown 278316643e7SJed Brown Note: 279316643e7SJed Brown Most users should not need to explicitly call this routine, as it 280316643e7SJed Brown is used internally within the nonlinear solvers. 281316643e7SJed Brown 282316643e7SJed Brown If the user did not provide a function but merely a matrix, 283d763cef2SBarry Smith this routine applies the matrix. 284316643e7SJed Brown 285316643e7SJed Brown Level: developer 286316643e7SJed Brown 287316643e7SJed Brown .keywords: TS, compute 288316643e7SJed Brown 289316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction() 290316643e7SJed Brown @*/ 291dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 292d763cef2SBarry Smith { 293dfbe8321SBarry Smith PetscErrorCode ierr; 294d763cef2SBarry Smith 295d763cef2SBarry Smith PetscFunctionBegin; 2960700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2970700a824SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2980700a824SBarry Smith PetscValidHeaderSpecific(y,VEC_CLASSID,4); 299d763cef2SBarry Smith 300d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 301000e7ae3SMatthew Knepley if (ts->ops->rhsfunction) { 302d763cef2SBarry Smith PetscStackPush("TS user right-hand-side function"); 303000e7ae3SMatthew Knepley ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 304d763cef2SBarry Smith PetscStackPop; 305d16dab79SJed Brown } else if (ts->Arhs) { 306000e7ae3SMatthew Knepley if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 307d763cef2SBarry Smith MatStructure flg; 308d763cef2SBarry Smith PetscStackPush("TS user right-hand-side matrix function"); 3098beb423aSHong Zhang ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 310d763cef2SBarry Smith PetscStackPop; 311d763cef2SBarry Smith } 3128beb423aSHong Zhang ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr); 313d16dab79SJed Brown } else SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"No RHS provided, must call TSSetRHSFunction() or TSSetMatrices()"); 314d763cef2SBarry Smith 315d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 316d763cef2SBarry Smith 317d763cef2SBarry Smith PetscFunctionReturn(0); 318d763cef2SBarry Smith } 319d763cef2SBarry Smith 3204a2ae208SSatish Balay #undef __FUNCT__ 321316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction" 322316643e7SJed Brown /*@ 323316643e7SJed Brown TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 324316643e7SJed Brown 325316643e7SJed Brown Collective on TS and Vec 326316643e7SJed Brown 327316643e7SJed Brown Input Parameters: 328316643e7SJed Brown + ts - the TS context 329316643e7SJed Brown . t - current time 330316643e7SJed Brown . X - state vector 331316643e7SJed Brown - Xdot - time derivative of state vector 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 @*/ 349316643e7SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y) 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; 364316643e7SJed Brown } else { 365316643e7SJed Brown if (ts->ops->rhsfunction) { 366316643e7SJed Brown PetscStackPush("TS user right-hand-side function"); 367316643e7SJed Brown ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr); 368316643e7SJed Brown PetscStackPop; 369316643e7SJed Brown } else { 370316643e7SJed Brown if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */ 371316643e7SJed Brown MatStructure flg; 3724a6899ffSJed Brown /* Note: flg is not being used. 3734a6899ffSJed Brown For it to be useful, we'd have to cache it and then apply it in TSComputeIJacobian. 3744a6899ffSJed Brown */ 375316643e7SJed Brown PetscStackPush("TS user right-hand-side matrix function"); 376316643e7SJed Brown ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 377316643e7SJed Brown PetscStackPop; 378316643e7SJed Brown } 379316643e7SJed Brown ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr); 380316643e7SJed Brown } 381316643e7SJed Brown 3824a6899ffSJed Brown /* Convert to implicit form: F(X,Xdot) = Alhs * Xdot - Frhs(X) */ 3834a6899ffSJed Brown if (ts->Alhs) { 3844a6899ffSJed Brown if (ts->ops->lhsmatrix) { 3854a6899ffSJed Brown MatStructure flg; 3864a6899ffSJed Brown PetscStackPush("TS user left-hand-side matrix function"); 3874a6899ffSJed Brown ierr = (*ts->ops->lhsmatrix)(ts,t,&ts->Alhs,PETSC_NULL,&flg,ts->jacP);CHKERRQ(ierr); 3884a6899ffSJed Brown PetscStackPop; 3894a6899ffSJed Brown } 3904a6899ffSJed Brown ierr = VecScale(Y,-1.);CHKERRQ(ierr); 3914a6899ffSJed Brown ierr = MatMultAdd(ts->Alhs,Xdot,Y,Y);CHKERRQ(ierr); 3924a6899ffSJed Brown } else { 393ace68cafSJed Brown ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 394316643e7SJed Brown } 3954a6899ffSJed Brown } 396316643e7SJed Brown ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 397316643e7SJed Brown PetscFunctionReturn(0); 398316643e7SJed Brown } 399316643e7SJed Brown 400316643e7SJed Brown #undef __FUNCT__ 401316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian" 402316643e7SJed Brown /*@ 403316643e7SJed Brown TSComputeIJacobian - Evaluates the Jacobian of the DAE 404316643e7SJed Brown 405316643e7SJed Brown Collective on TS and Vec 406316643e7SJed Brown 407316643e7SJed Brown Input 408316643e7SJed Brown Input Parameters: 409316643e7SJed Brown + ts - the TS context 410316643e7SJed Brown . t - current timestep 411316643e7SJed Brown . X - state vector 412316643e7SJed Brown . Xdot - time derivative of state vector 413316643e7SJed Brown - shift - shift to apply, see note below 414316643e7SJed Brown 415316643e7SJed Brown Output Parameters: 416316643e7SJed Brown + A - Jacobian matrix 417316643e7SJed Brown . B - optional preconditioning matrix 418316643e7SJed Brown - flag - flag indicating matrix structure 419316643e7SJed Brown 420316643e7SJed Brown Notes: 421316643e7SJed Brown If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 422316643e7SJed Brown 423316643e7SJed Brown dF/dX + shift*dF/dXdot 424316643e7SJed Brown 425316643e7SJed Brown Most users should not need to explicitly call this routine, as it 426316643e7SJed Brown is used internally within the nonlinear solvers. 427316643e7SJed Brown 428316643e7SJed Brown TSComputeIJacobian() is valid only for TS_NONLINEAR 429316643e7SJed Brown 430316643e7SJed Brown Level: developer 431316643e7SJed Brown 432316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix 433316643e7SJed Brown 434316643e7SJed Brown .seealso: TSSetIJacobian() 43563495f91SJed Brown @*/ 4367087cfbeSBarry Smith PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg) 437316643e7SJed Brown { 438316643e7SJed Brown PetscErrorCode ierr; 439316643e7SJed Brown 440316643e7SJed Brown PetscFunctionBegin; 4410700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4420700a824SBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,3); 4430700a824SBarry Smith PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 444316643e7SJed Brown PetscValidPointer(A,6); 4450700a824SBarry Smith PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 446316643e7SJed Brown PetscValidPointer(B,7); 4470700a824SBarry Smith PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 448316643e7SJed Brown PetscValidPointer(flg,8); 449316643e7SJed Brown 4504a6899ffSJed Brown *flg = SAME_NONZERO_PATTERN; /* In case it we're solving a linear problem in which case it wouldn't get initialized below. */ 451316643e7SJed Brown ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 452316643e7SJed Brown if (ts->ops->ijacobian) { 453316643e7SJed Brown PetscStackPush("TS user implicit Jacobian"); 454316643e7SJed Brown ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 455316643e7SJed Brown PetscStackPop; 456316643e7SJed Brown } else { 457316643e7SJed Brown if (ts->ops->rhsjacobian) { 458316643e7SJed Brown PetscStackPush("TS user right-hand-side Jacobian"); 459316643e7SJed Brown ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 460316643e7SJed Brown PetscStackPop; 461316643e7SJed Brown } else { 462316643e7SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 463316643e7SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 464316643e7SJed Brown if (*A != *B) { 465316643e7SJed Brown ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 466316643e7SJed Brown ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 467316643e7SJed Brown } 468316643e7SJed Brown } 469316643e7SJed Brown 470316643e7SJed Brown /* Convert to implicit form */ 471316643e7SJed Brown /* inefficient because these operations will normally traverse all matrix elements twice */ 472316643e7SJed Brown ierr = MatScale(*A,-1);CHKERRQ(ierr); 4734a6899ffSJed Brown if (ts->Alhs) { 4744a6899ffSJed Brown ierr = MatAXPY(*A,shift,ts->Alhs,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 4754a6899ffSJed Brown } else { 476316643e7SJed Brown ierr = MatShift(*A,shift);CHKERRQ(ierr); 4774a6899ffSJed Brown } 478316643e7SJed Brown if (*A != *B) { 479316643e7SJed Brown ierr = MatScale(*B,-1);CHKERRQ(ierr); 480316643e7SJed Brown ierr = MatShift(*B,shift);CHKERRQ(ierr); 481316643e7SJed Brown } 482316643e7SJed Brown } 483316643e7SJed Brown ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 484316643e7SJed Brown PetscFunctionReturn(0); 485316643e7SJed Brown } 486316643e7SJed Brown 487316643e7SJed Brown #undef __FUNCT__ 4884a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction" 489d763cef2SBarry Smith /*@C 490d763cef2SBarry Smith TSSetRHSFunction - Sets the routine for evaluating the function, 491d763cef2SBarry Smith F(t,u), where U_t = F(t,u). 492d763cef2SBarry Smith 4933f9fe445SBarry Smith Logically Collective on TS 494d763cef2SBarry Smith 495d763cef2SBarry Smith Input Parameters: 496d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 497d763cef2SBarry Smith . f - routine for evaluating the right-hand-side function 498d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 499d763cef2SBarry Smith function evaluation routine (may be PETSC_NULL) 500d763cef2SBarry Smith 501d763cef2SBarry Smith Calling sequence of func: 50287828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 503d763cef2SBarry Smith 504d763cef2SBarry Smith + t - current timestep 505d763cef2SBarry Smith . u - input vector 506d763cef2SBarry Smith . F - function vector 507d763cef2SBarry Smith - ctx - [optional] user-defined function context 508d763cef2SBarry Smith 509d763cef2SBarry Smith Important: 51076f2fa84SHong Zhang The user MUST call either this routine or TSSetMatrices(). 511d763cef2SBarry Smith 512d763cef2SBarry Smith Level: beginner 513d763cef2SBarry Smith 514d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function 515d763cef2SBarry Smith 51676f2fa84SHong Zhang .seealso: TSSetMatrices() 517d763cef2SBarry Smith @*/ 5187087cfbeSBarry Smith PetscErrorCode TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 519d763cef2SBarry Smith { 520d763cef2SBarry Smith PetscFunctionBegin; 521d763cef2SBarry Smith 5220700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 52317186662SBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem"); 524000e7ae3SMatthew Knepley ts->ops->rhsfunction = f; 525d763cef2SBarry Smith ts->funP = ctx; 526d763cef2SBarry Smith PetscFunctionReturn(0); 527d763cef2SBarry Smith } 528d763cef2SBarry Smith 5294a2ae208SSatish Balay #undef __FUNCT__ 53095f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices" 53195f0b562SHong Zhang /*@C 53295f0b562SHong Zhang TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs, 53395f0b562SHong Zhang where Alhs(t) U_t = Arhs(t) U. 53495f0b562SHong Zhang 5353f9fe445SBarry Smith Logically Collective on TS 53695f0b562SHong Zhang 53795f0b562SHong Zhang Input Parameters: 53895f0b562SHong Zhang + ts - the TS context obtained from TSCreate() 53995f0b562SHong Zhang . Arhs - matrix 54095f0b562SHong Zhang . frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 54195f0b562SHong Zhang if Arhs is not a function of t. 54295f0b562SHong Zhang . Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix. 54395f0b562SHong Zhang . flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 54495f0b562SHong Zhang if Alhs is not a function of t. 54595f0b562SHong Zhang . flag - flag indicating information about the matrix structure of Arhs and Alhs. 54695f0b562SHong Zhang The available options are 54795f0b562SHong Zhang SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs 54895f0b562SHong Zhang DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs 54995f0b562SHong Zhang - ctx - [optional] user-defined context for private data for the 55095f0b562SHong Zhang matrix evaluation routine (may be PETSC_NULL) 55195f0b562SHong Zhang 55295f0b562SHong Zhang Calling sequence of func: 55395f0b562SHong Zhang $ func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx); 55495f0b562SHong Zhang 55595f0b562SHong Zhang + t - current timestep 55695f0b562SHong Zhang . A - matrix A, where U_t = A(t) U 55795f0b562SHong Zhang . B - preconditioner matrix, usually the same as A 55895f0b562SHong Zhang . flag - flag indicating information about the preconditioner matrix 55995f0b562SHong Zhang structure (same as flag in KSPSetOperators()) 56095f0b562SHong Zhang - ctx - [optional] user-defined context for matrix evaluation routine 56195f0b562SHong Zhang 56295f0b562SHong Zhang Notes: 56395f0b562SHong Zhang The routine func() takes Mat* as the matrix arguments rather than Mat. 56495f0b562SHong Zhang This allows the matrix evaluation routine to replace Arhs or Alhs with a 56595f0b562SHong Zhang completely new new matrix structure (not just different matrix elements) 56695f0b562SHong Zhang when appropriate, for instance, if the nonzero structure is changing 56795f0b562SHong Zhang throughout the global iterations. 56895f0b562SHong Zhang 56995f0b562SHong Zhang Important: 57095f0b562SHong Zhang The user MUST call either this routine or TSSetRHSFunction(). 57195f0b562SHong Zhang 57295f0b562SHong Zhang Level: beginner 57395f0b562SHong Zhang 57495f0b562SHong Zhang .keywords: TS, timestep, set, matrix 57595f0b562SHong Zhang 57695f0b562SHong Zhang .seealso: TSSetRHSFunction() 57795f0b562SHong Zhang @*/ 5787087cfbeSBarry Smith PetscErrorCode TSSetMatrices(TS ts,Mat Arhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx) 57995f0b562SHong Zhang { 580277b19d0SLisandro Dalcin PetscErrorCode ierr; 581277b19d0SLisandro Dalcin 58295f0b562SHong Zhang PetscFunctionBegin; 5830700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 584277b19d0SLisandro Dalcin if (frhs) ts->ops->rhsmatrix = frhs; 585277b19d0SLisandro Dalcin if (flhs) ts->ops->lhsmatrix = flhs; 586277b19d0SLisandro Dalcin if (ctx) ts->jacP = ctx; 58792af4f6aSHong Zhang if (Arhs){ 5880700a824SBarry Smith PetscValidHeaderSpecific(Arhs,MAT_CLASSID,2); 58995f0b562SHong Zhang PetscCheckSameComm(ts,1,Arhs,2); 590277b19d0SLisandro Dalcin ierr = PetscObjectReference((PetscObject)Arhs);CHKERRQ(ierr); 5916bf464f9SBarry Smith ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 59295f0b562SHong Zhang ts->Arhs = Arhs; 59392af4f6aSHong Zhang } 59492af4f6aSHong Zhang if (Alhs){ 5950700a824SBarry Smith PetscValidHeaderSpecific(Alhs,MAT_CLASSID,4); 59692af4f6aSHong Zhang PetscCheckSameComm(ts,1,Alhs,4); 597277b19d0SLisandro Dalcin ierr = PetscObjectReference((PetscObject)Alhs);CHKERRQ(ierr); 5986bf464f9SBarry Smith ierr = MatDestroy(&ts->Alhs);CHKERRQ(ierr); 59995f0b562SHong Zhang ts->Alhs = Alhs; 60092af4f6aSHong Zhang } 60195f0b562SHong Zhang ts->matflg = flag; 60295f0b562SHong Zhang PetscFunctionReturn(0); 60395f0b562SHong Zhang } 604d763cef2SBarry Smith 605aa644b49SHong Zhang #undef __FUNCT__ 606cda39b92SHong Zhang #define __FUNCT__ "TSGetMatrices" 607cda39b92SHong Zhang /*@C 608cda39b92SHong Zhang TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep, 609cda39b92SHong Zhang where Alhs(t) U_t = Arhs(t) U. 610cda39b92SHong Zhang 611cda39b92SHong Zhang Not Collective, but parallel objects are returned if TS is parallel 612cda39b92SHong Zhang 613cda39b92SHong Zhang Input Parameter: 614cda39b92SHong Zhang . ts - The TS context obtained from TSCreate() 615cda39b92SHong Zhang 616cda39b92SHong Zhang Output Parameters: 617cda39b92SHong Zhang + Arhs - The right-hand side matrix 618cda39b92SHong Zhang . Alhs - The left-hand side matrix 619cda39b92SHong Zhang - ctx - User-defined context for matrix evaluation routine 620cda39b92SHong Zhang 621cda39b92SHong Zhang Notes: You can pass in PETSC_NULL for any return argument you do not need. 622cda39b92SHong Zhang 623cda39b92SHong Zhang Level: intermediate 624cda39b92SHong Zhang 625cda39b92SHong Zhang .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 626cda39b92SHong Zhang 627cda39b92SHong Zhang .keywords: TS, timestep, get, matrix 628cda39b92SHong Zhang 629cda39b92SHong Zhang @*/ 6307087cfbeSBarry Smith PetscErrorCode TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx) 631cda39b92SHong Zhang { 632cda39b92SHong Zhang PetscFunctionBegin; 6330700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 634cda39b92SHong Zhang if (Arhs) *Arhs = ts->Arhs; 635cda39b92SHong Zhang if (Alhs) *Alhs = ts->Alhs; 636cda39b92SHong Zhang if (ctx) *ctx = ts->jacP; 637cda39b92SHong Zhang PetscFunctionReturn(0); 638cda39b92SHong Zhang } 639cda39b92SHong Zhang 640cda39b92SHong Zhang #undef __FUNCT__ 6414a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian" 642d763cef2SBarry Smith /*@C 643d763cef2SBarry Smith TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 644d763cef2SBarry Smith where U_t = F(U,t), as well as the location to store the matrix. 64576f2fa84SHong Zhang Use TSSetMatrices() for linear problems. 646d763cef2SBarry Smith 6473f9fe445SBarry Smith Logically Collective on TS 648d763cef2SBarry Smith 649d763cef2SBarry Smith Input Parameters: 650d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 651d763cef2SBarry Smith . A - Jacobian matrix 652d763cef2SBarry Smith . B - preconditioner matrix (usually same as A) 653d763cef2SBarry Smith . f - the Jacobian evaluation routine 654d763cef2SBarry Smith - ctx - [optional] user-defined context for private data for the 655d763cef2SBarry Smith Jacobian evaluation routine (may be PETSC_NULL) 656d763cef2SBarry Smith 657d763cef2SBarry Smith Calling sequence of func: 65887828ca2SBarry Smith $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 659d763cef2SBarry Smith 660d763cef2SBarry Smith + t - current timestep 661d763cef2SBarry Smith . u - input vector 662d763cef2SBarry Smith . A - matrix A, where U_t = A(t)u 663d763cef2SBarry Smith . B - preconditioner matrix, usually the same as A 664d763cef2SBarry Smith . flag - flag indicating information about the preconditioner matrix 66594b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 666d763cef2SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 667d763cef2SBarry Smith 668d763cef2SBarry Smith Notes: 66994b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 670d763cef2SBarry Smith output parameter in the routine func(). Be sure to read this information! 671d763cef2SBarry Smith 672d763cef2SBarry Smith The routine func() takes Mat * as the matrix arguments rather than Mat. 673d763cef2SBarry Smith This allows the matrix evaluation routine to replace A and/or B with a 67456335db2SHong Zhang completely new matrix structure (not just different matrix elements) 675d763cef2SBarry Smith when appropriate, for instance, if the nonzero structure is changing 676d763cef2SBarry Smith throughout the global iterations. 677d763cef2SBarry Smith 678d763cef2SBarry Smith Level: beginner 679d763cef2SBarry Smith 680d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian 681d763cef2SBarry Smith 682d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(), 68376f2fa84SHong Zhang SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 684d763cef2SBarry Smith 685d763cef2SBarry Smith @*/ 6867087cfbeSBarry Smith PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 687d763cef2SBarry Smith { 688277b19d0SLisandro Dalcin PetscErrorCode ierr; 689277b19d0SLisandro Dalcin 690d763cef2SBarry Smith PetscFunctionBegin; 6910700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 692277b19d0SLisandro Dalcin if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 693277b19d0SLisandro Dalcin if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 694277b19d0SLisandro Dalcin if (A) PetscCheckSameComm(ts,1,A,2); 695277b19d0SLisandro Dalcin if (B) PetscCheckSameComm(ts,1,B,3); 69617186662SBarry Smith if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()"); 697d763cef2SBarry Smith 698277b19d0SLisandro Dalcin if (f) ts->ops->rhsjacobian = f; 699277b19d0SLisandro Dalcin if (ctx) ts->jacP = ctx; 700277b19d0SLisandro Dalcin if (A) { 701277b19d0SLisandro Dalcin ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 7026bf464f9SBarry Smith ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 7038beb423aSHong Zhang ts->Arhs = A; 704277b19d0SLisandro Dalcin } 705277b19d0SLisandro Dalcin if (B) { 706277b19d0SLisandro Dalcin ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 7076bf464f9SBarry Smith ierr = MatDestroy(&ts->B);CHKERRQ(ierr); 708d763cef2SBarry Smith ts->B = B; 709277b19d0SLisandro Dalcin } 710d763cef2SBarry Smith PetscFunctionReturn(0); 711d763cef2SBarry Smith } 712d763cef2SBarry Smith 713316643e7SJed Brown 714316643e7SJed Brown #undef __FUNCT__ 715316643e7SJed Brown #define __FUNCT__ "TSSetIFunction" 716316643e7SJed Brown /*@C 717316643e7SJed Brown TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 718316643e7SJed Brown 7193f9fe445SBarry Smith Logically Collective on TS 720316643e7SJed Brown 721316643e7SJed Brown Input Parameters: 722316643e7SJed Brown + ts - the TS context obtained from TSCreate() 723316643e7SJed Brown . f - the function evaluation routine 724316643e7SJed Brown - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 725316643e7SJed Brown 726316643e7SJed Brown Calling sequence of f: 727316643e7SJed Brown $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 728316643e7SJed Brown 729316643e7SJed Brown + t - time at step/stage being solved 730316643e7SJed Brown . u - state vector 731316643e7SJed Brown . u_t - time derivative of state vector 732316643e7SJed Brown . F - function vector 733316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 734316643e7SJed Brown 735316643e7SJed Brown Important: 736316643e7SJed Brown The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 737316643e7SJed Brown 738316643e7SJed Brown Level: beginner 739316643e7SJed Brown 740316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian 741316643e7SJed Brown 742316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 743316643e7SJed Brown @*/ 7447087cfbeSBarry Smith PetscErrorCode TSSetIFunction(TS ts,TSIFunction f,void *ctx) 745316643e7SJed Brown { 746316643e7SJed Brown 747316643e7SJed Brown PetscFunctionBegin; 7480700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 749316643e7SJed Brown ts->ops->ifunction = f; 750316643e7SJed Brown ts->funP = ctx; 751316643e7SJed Brown PetscFunctionReturn(0); 752316643e7SJed Brown } 753316643e7SJed Brown 754316643e7SJed Brown #undef __FUNCT__ 755316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian" 756316643e7SJed Brown /*@C 757a4f0a591SBarry Smith TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 758a4f0a591SBarry Smith you provided with TSSetIFunction(). 759316643e7SJed Brown 7603f9fe445SBarry Smith Logically Collective on TS 761316643e7SJed Brown 762316643e7SJed Brown Input Parameters: 763316643e7SJed Brown + ts - the TS context obtained from TSCreate() 764316643e7SJed Brown . A - Jacobian matrix 765316643e7SJed Brown . B - preconditioning matrix for A (may be same as A) 766316643e7SJed Brown . f - the Jacobian evaluation routine 767316643e7SJed Brown - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 768316643e7SJed Brown 769316643e7SJed Brown Calling sequence of f: 7701b4a444bSJed Brown $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 771316643e7SJed Brown 772316643e7SJed Brown + t - time at step/stage being solved 7731b4a444bSJed Brown . U - state vector 7741b4a444bSJed Brown . U_t - time derivative of state vector 775316643e7SJed Brown . a - shift 7761b4a444bSJed Brown . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 777316643e7SJed Brown . B - preconditioning matrix for A, may be same as A 778316643e7SJed Brown . flag - flag indicating information about the preconditioner matrix 779316643e7SJed Brown structure (same as flag in KSPSetOperators()) 780316643e7SJed Brown - ctx - [optional] user-defined context for matrix evaluation routine 781316643e7SJed Brown 782316643e7SJed Brown Notes: 783316643e7SJed Brown The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 784316643e7SJed Brown 785a4f0a591SBarry Smith The matrix dF/dU + a*dF/dU_t you provide turns out to be 786a4f0a591SBarry 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. 787a4f0a591SBarry Smith The time integrator internally approximates U_t by W+a*U where the positive "shift" 788a4f0a591SBarry Smith a and vector W depend on the integration method, step size, and past states. For example with 789a4f0a591SBarry Smith the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 790a4f0a591SBarry Smith W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 791a4f0a591SBarry Smith 792316643e7SJed Brown Level: beginner 793316643e7SJed Brown 794316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian 795316643e7SJed Brown 796316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian() 797316643e7SJed Brown 798316643e7SJed Brown @*/ 7997087cfbeSBarry Smith PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 800316643e7SJed Brown { 801316643e7SJed Brown PetscErrorCode ierr; 802316643e7SJed Brown 803316643e7SJed Brown PetscFunctionBegin; 8040700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8050700a824SBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 8060700a824SBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 807316643e7SJed Brown if (A) PetscCheckSameComm(ts,1,A,2); 808316643e7SJed Brown if (B) PetscCheckSameComm(ts,1,B,3); 809316643e7SJed Brown if (f) ts->ops->ijacobian = f; 810316643e7SJed Brown if (ctx) ts->jacP = ctx; 811316643e7SJed Brown if (A) { 812316643e7SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8136bf464f9SBarry Smith ierr = MatDestroy(&ts->A);CHKERRQ(ierr); 814316643e7SJed Brown ts->A = A; 815316643e7SJed Brown } 816316643e7SJed Brown if (B) { 817316643e7SJed Brown ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 8186bf464f9SBarry Smith ierr = MatDestroy(&ts->B);CHKERRQ(ierr); 819316643e7SJed Brown ts->B = B; 820316643e7SJed Brown } 821316643e7SJed Brown PetscFunctionReturn(0); 822316643e7SJed Brown } 823316643e7SJed Brown 8244a2ae208SSatish Balay #undef __FUNCT__ 8254a2ae208SSatish Balay #define __FUNCT__ "TSView" 8267e2c5f70SBarry Smith /*@C 827d763cef2SBarry Smith TSView - Prints the TS data structure. 828d763cef2SBarry Smith 8294c49b128SBarry Smith Collective on TS 830d763cef2SBarry Smith 831d763cef2SBarry Smith Input Parameters: 832d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 833d763cef2SBarry Smith - viewer - visualization context 834d763cef2SBarry Smith 835d763cef2SBarry Smith Options Database Key: 836d763cef2SBarry Smith . -ts_view - calls TSView() at end of TSStep() 837d763cef2SBarry Smith 838d763cef2SBarry Smith Notes: 839d763cef2SBarry Smith The available visualization contexts include 840b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 841b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 842d763cef2SBarry Smith output where only the first processor opens 843d763cef2SBarry Smith the file. All other processors send their 844d763cef2SBarry Smith data to the first processor to print. 845d763cef2SBarry Smith 846d763cef2SBarry Smith The user can open an alternative visualization context with 847b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 848d763cef2SBarry Smith 849d763cef2SBarry Smith Level: beginner 850d763cef2SBarry Smith 851d763cef2SBarry Smith .keywords: TS, timestep, view 852d763cef2SBarry Smith 853b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 854d763cef2SBarry Smith @*/ 8557087cfbeSBarry Smith PetscErrorCode TSView(TS ts,PetscViewer viewer) 856d763cef2SBarry Smith { 857dfbe8321SBarry Smith PetscErrorCode ierr; 858a313700dSBarry Smith const TSType type; 859ace3abfcSBarry Smith PetscBool iascii,isstring; 860d763cef2SBarry Smith 861d763cef2SBarry Smith PetscFunctionBegin; 8620700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8633050cee2SBarry Smith if (!viewer) { 8647adad957SLisandro Dalcin ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 8653050cee2SBarry Smith } 8660700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 867c9780b6fSBarry Smith PetscCheckSameComm(ts,1,viewer,2); 868fd16b177SBarry Smith 8692692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8702692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 87132077d6dSBarry Smith if (iascii) { 872317d6ea6SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 873000e7ae3SMatthew Knepley if (ts->ops->view) { 874b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 875000e7ae3SMatthew Knepley ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 876b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 877d763cef2SBarry Smith } 87877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 879a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 880d763cef2SBarry Smith if (ts->problem_type == TS_NONLINEAR) { 88177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 882d763cef2SBarry Smith } 88377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 8840f5bd95cSBarry Smith } else if (isstring) { 885a313700dSBarry Smith ierr = TSGetType(ts,&type);CHKERRQ(ierr); 886b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 887d763cef2SBarry Smith } 888b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 88994b7f48cSBarry Smith if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);} 890d763cef2SBarry Smith if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 891b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 892d763cef2SBarry Smith PetscFunctionReturn(0); 893d763cef2SBarry Smith } 894d763cef2SBarry Smith 895d763cef2SBarry Smith 8964a2ae208SSatish Balay #undef __FUNCT__ 8974a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext" 898b07ff414SBarry Smith /*@ 899d763cef2SBarry Smith TSSetApplicationContext - Sets an optional user-defined context for 900d763cef2SBarry Smith the timesteppers. 901d763cef2SBarry Smith 9023f9fe445SBarry Smith Logically Collective on TS 903d763cef2SBarry Smith 904d763cef2SBarry Smith Input Parameters: 905d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 906d763cef2SBarry Smith - usrP - optional user context 907d763cef2SBarry Smith 908d763cef2SBarry Smith Level: intermediate 909d763cef2SBarry Smith 910d763cef2SBarry Smith .keywords: TS, timestep, set, application, context 911d763cef2SBarry Smith 912d763cef2SBarry Smith .seealso: TSGetApplicationContext() 913d763cef2SBarry Smith @*/ 9147087cfbeSBarry Smith PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 915d763cef2SBarry Smith { 916d763cef2SBarry Smith PetscFunctionBegin; 9170700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 918d763cef2SBarry Smith ts->user = usrP; 919d763cef2SBarry Smith PetscFunctionReturn(0); 920d763cef2SBarry Smith } 921d763cef2SBarry Smith 9224a2ae208SSatish Balay #undef __FUNCT__ 9234a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext" 924b07ff414SBarry Smith /*@ 925d763cef2SBarry Smith TSGetApplicationContext - Gets the user-defined context for the 926d763cef2SBarry Smith timestepper. 927d763cef2SBarry Smith 928d763cef2SBarry Smith Not Collective 929d763cef2SBarry Smith 930d763cef2SBarry Smith Input Parameter: 931d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 932d763cef2SBarry Smith 933d763cef2SBarry Smith Output Parameter: 934d763cef2SBarry Smith . usrP - user context 935d763cef2SBarry Smith 936d763cef2SBarry Smith Level: intermediate 937d763cef2SBarry Smith 938d763cef2SBarry Smith .keywords: TS, timestep, get, application, context 939d763cef2SBarry Smith 940d763cef2SBarry Smith .seealso: TSSetApplicationContext() 941d763cef2SBarry Smith @*/ 942e71120c6SJed Brown PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 943d763cef2SBarry Smith { 944d763cef2SBarry Smith PetscFunctionBegin; 9450700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 946e71120c6SJed Brown *(void**)usrP = ts->user; 947d763cef2SBarry Smith PetscFunctionReturn(0); 948d763cef2SBarry Smith } 949d763cef2SBarry Smith 9504a2ae208SSatish Balay #undef __FUNCT__ 9514a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber" 952d763cef2SBarry Smith /*@ 953d763cef2SBarry Smith TSGetTimeStepNumber - Gets the current number of timesteps. 954d763cef2SBarry Smith 955d763cef2SBarry Smith Not Collective 956d763cef2SBarry Smith 957d763cef2SBarry Smith Input Parameter: 958d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 959d763cef2SBarry Smith 960d763cef2SBarry Smith Output Parameter: 961d763cef2SBarry Smith . iter - number steps so far 962d763cef2SBarry Smith 963d763cef2SBarry Smith Level: intermediate 964d763cef2SBarry Smith 965d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number 966d763cef2SBarry Smith @*/ 9677087cfbeSBarry Smith PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 968d763cef2SBarry Smith { 969d763cef2SBarry Smith PetscFunctionBegin; 9700700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9714482741eSBarry Smith PetscValidIntPointer(iter,2); 972d763cef2SBarry Smith *iter = ts->steps; 973d763cef2SBarry Smith PetscFunctionReturn(0); 974d763cef2SBarry Smith } 975d763cef2SBarry Smith 9764a2ae208SSatish Balay #undef __FUNCT__ 9774a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep" 978d763cef2SBarry Smith /*@ 979d763cef2SBarry Smith TSSetInitialTimeStep - Sets the initial timestep to be used, 980d763cef2SBarry Smith as well as the initial time. 981d763cef2SBarry Smith 9823f9fe445SBarry Smith Logically Collective on TS 983d763cef2SBarry Smith 984d763cef2SBarry Smith Input Parameters: 985d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 986d763cef2SBarry Smith . initial_time - the initial time 987d763cef2SBarry Smith - time_step - the size of the timestep 988d763cef2SBarry Smith 989d763cef2SBarry Smith Level: intermediate 990d763cef2SBarry Smith 991d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep() 992d763cef2SBarry Smith 993d763cef2SBarry Smith .keywords: TS, set, initial, timestep 994d763cef2SBarry Smith @*/ 9957087cfbeSBarry Smith PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 996d763cef2SBarry Smith { 997d763cef2SBarry Smith PetscFunctionBegin; 9980700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 999d763cef2SBarry Smith ts->time_step = time_step; 1000d763cef2SBarry Smith ts->initial_time_step = time_step; 1001d763cef2SBarry Smith ts->ptime = initial_time; 1002d763cef2SBarry Smith PetscFunctionReturn(0); 1003d763cef2SBarry Smith } 1004d763cef2SBarry Smith 10054a2ae208SSatish Balay #undef __FUNCT__ 10064a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep" 1007d763cef2SBarry Smith /*@ 1008d763cef2SBarry Smith TSSetTimeStep - Allows one to reset the timestep at any time, 1009d763cef2SBarry Smith useful for simple pseudo-timestepping codes. 1010d763cef2SBarry Smith 10113f9fe445SBarry Smith Logically Collective on TS 1012d763cef2SBarry Smith 1013d763cef2SBarry Smith Input Parameters: 1014d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1015d763cef2SBarry Smith - time_step - the size of the timestep 1016d763cef2SBarry Smith 1017d763cef2SBarry Smith Level: intermediate 1018d763cef2SBarry Smith 1019d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1020d763cef2SBarry Smith 1021d763cef2SBarry Smith .keywords: TS, set, timestep 1022d763cef2SBarry Smith @*/ 10237087cfbeSBarry Smith PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 1024d763cef2SBarry Smith { 1025d763cef2SBarry Smith PetscFunctionBegin; 10260700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1027c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,time_step,2); 1028d763cef2SBarry Smith ts->time_step = time_step; 1029d763cef2SBarry Smith PetscFunctionReturn(0); 1030d763cef2SBarry Smith } 1031d763cef2SBarry Smith 10324a2ae208SSatish Balay #undef __FUNCT__ 10334a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep" 1034d763cef2SBarry Smith /*@ 1035d763cef2SBarry Smith TSGetTimeStep - Gets the current timestep size. 1036d763cef2SBarry Smith 1037d763cef2SBarry Smith Not Collective 1038d763cef2SBarry Smith 1039d763cef2SBarry Smith Input Parameter: 1040d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1041d763cef2SBarry Smith 1042d763cef2SBarry Smith Output Parameter: 1043d763cef2SBarry Smith . dt - the current timestep size 1044d763cef2SBarry Smith 1045d763cef2SBarry Smith Level: intermediate 1046d763cef2SBarry Smith 1047d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1048d763cef2SBarry Smith 1049d763cef2SBarry Smith .keywords: TS, get, timestep 1050d763cef2SBarry Smith @*/ 10517087cfbeSBarry Smith PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1052d763cef2SBarry Smith { 1053d763cef2SBarry Smith PetscFunctionBegin; 10540700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10554482741eSBarry Smith PetscValidDoublePointer(dt,2); 1056d763cef2SBarry Smith *dt = ts->time_step; 1057d763cef2SBarry Smith PetscFunctionReturn(0); 1058d763cef2SBarry Smith } 1059d763cef2SBarry Smith 10604a2ae208SSatish Balay #undef __FUNCT__ 10614a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution" 1062d8e5e3e6SSatish Balay /*@ 1063d763cef2SBarry Smith TSGetSolution - Returns the solution at the present timestep. It 1064d763cef2SBarry Smith is valid to call this routine inside the function that you are evaluating 1065d763cef2SBarry Smith in order to move to the new timestep. This vector not changed until 1066d763cef2SBarry Smith the solution at the next timestep has been calculated. 1067d763cef2SBarry Smith 1068d763cef2SBarry Smith Not Collective, but Vec returned is parallel if TS is parallel 1069d763cef2SBarry Smith 1070d763cef2SBarry Smith Input Parameter: 1071d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1072d763cef2SBarry Smith 1073d763cef2SBarry Smith Output Parameter: 1074d763cef2SBarry Smith . v - the vector containing the solution 1075d763cef2SBarry Smith 1076d763cef2SBarry Smith Level: intermediate 1077d763cef2SBarry Smith 1078d763cef2SBarry Smith .seealso: TSGetTimeStep() 1079d763cef2SBarry Smith 1080d763cef2SBarry Smith .keywords: TS, timestep, get, solution 1081d763cef2SBarry Smith @*/ 10827087cfbeSBarry Smith PetscErrorCode TSGetSolution(TS ts,Vec *v) 1083d763cef2SBarry Smith { 1084d763cef2SBarry Smith PetscFunctionBegin; 10850700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10864482741eSBarry Smith PetscValidPointer(v,2); 10878737fe31SLisandro Dalcin *v = ts->vec_sol; 1088d763cef2SBarry Smith PetscFunctionReturn(0); 1089d763cef2SBarry Smith } 1090d763cef2SBarry Smith 1091bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */ 10924a2ae208SSatish Balay #undef __FUNCT__ 1093bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType" 1094d8e5e3e6SSatish Balay /*@ 1095bdad233fSMatthew Knepley TSSetProblemType - Sets the type of problem to be solved. 1096d763cef2SBarry Smith 1097bdad233fSMatthew Knepley Not collective 1098d763cef2SBarry Smith 1099bdad233fSMatthew Knepley Input Parameters: 1100bdad233fSMatthew Knepley + ts - The TS 1101bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1102d763cef2SBarry Smith .vb 1103d763cef2SBarry Smith U_t = A U 1104d763cef2SBarry Smith U_t = A(t) U 1105d763cef2SBarry Smith U_t = F(t,U) 1106d763cef2SBarry Smith .ve 1107d763cef2SBarry Smith 1108d763cef2SBarry Smith Level: beginner 1109d763cef2SBarry Smith 1110bdad233fSMatthew Knepley .keywords: TS, problem type 1111bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1112d763cef2SBarry Smith @*/ 11137087cfbeSBarry Smith PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1114a7cc72afSBarry Smith { 1115d763cef2SBarry Smith PetscFunctionBegin; 11160700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1117bdad233fSMatthew Knepley ts->problem_type = type; 1118d763cef2SBarry Smith PetscFunctionReturn(0); 1119d763cef2SBarry Smith } 1120d763cef2SBarry Smith 1121bdad233fSMatthew Knepley #undef __FUNCT__ 1122bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType" 1123bdad233fSMatthew Knepley /*@C 1124bdad233fSMatthew Knepley TSGetProblemType - Gets the type of problem to be solved. 1125bdad233fSMatthew Knepley 1126bdad233fSMatthew Knepley Not collective 1127bdad233fSMatthew Knepley 1128bdad233fSMatthew Knepley Input Parameter: 1129bdad233fSMatthew Knepley . ts - The TS 1130bdad233fSMatthew Knepley 1131bdad233fSMatthew Knepley Output Parameter: 1132bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1133bdad233fSMatthew Knepley .vb 1134bdad233fSMatthew Knepley U_t = A U 1135bdad233fSMatthew Knepley U_t = A(t) U 1136bdad233fSMatthew Knepley U_t = F(t,U) 1137bdad233fSMatthew Knepley .ve 1138bdad233fSMatthew Knepley 1139bdad233fSMatthew Knepley Level: beginner 1140bdad233fSMatthew Knepley 1141bdad233fSMatthew Knepley .keywords: TS, problem type 1142bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS 1143bdad233fSMatthew Knepley @*/ 11447087cfbeSBarry Smith PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1145a7cc72afSBarry Smith { 1146bdad233fSMatthew Knepley PetscFunctionBegin; 11470700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 11484482741eSBarry Smith PetscValidIntPointer(type,2); 1149bdad233fSMatthew Knepley *type = ts->problem_type; 1150bdad233fSMatthew Knepley PetscFunctionReturn(0); 1151bdad233fSMatthew Knepley } 1152d763cef2SBarry Smith 11534a2ae208SSatish Balay #undef __FUNCT__ 11544a2ae208SSatish Balay #define __FUNCT__ "TSSetUp" 1155d763cef2SBarry Smith /*@ 1156d763cef2SBarry Smith TSSetUp - Sets up the internal data structures for the later use 1157d763cef2SBarry Smith of a timestepper. 1158d763cef2SBarry Smith 1159d763cef2SBarry Smith Collective on TS 1160d763cef2SBarry Smith 1161d763cef2SBarry Smith Input Parameter: 1162d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1163d763cef2SBarry Smith 1164d763cef2SBarry Smith Notes: 1165d763cef2SBarry Smith For basic use of the TS solvers the user need not explicitly call 1166d763cef2SBarry Smith TSSetUp(), since these actions will automatically occur during 1167d763cef2SBarry Smith the call to TSStep(). However, if one wishes to control this 1168d763cef2SBarry Smith phase separately, TSSetUp() should be called after TSCreate() 1169d763cef2SBarry Smith and optional routines of the form TSSetXXX(), but before TSStep(). 1170d763cef2SBarry Smith 1171d763cef2SBarry Smith Level: advanced 1172d763cef2SBarry Smith 1173d763cef2SBarry Smith .keywords: TS, timestep, setup 1174d763cef2SBarry Smith 1175d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy() 1176d763cef2SBarry Smith @*/ 11777087cfbeSBarry Smith PetscErrorCode TSSetUp(TS ts) 1178d763cef2SBarry Smith { 1179dfbe8321SBarry Smith PetscErrorCode ierr; 1180d763cef2SBarry Smith 1181d763cef2SBarry Smith PetscFunctionBegin; 11820700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1183277b19d0SLisandro Dalcin if (ts->setupcalled) PetscFunctionReturn(0); 1184277b19d0SLisandro Dalcin 11857adad957SLisandro Dalcin if (!((PetscObject)ts)->type_name) { 11869596e0b4SJed Brown ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1187d763cef2SBarry Smith } 1188277b19d0SLisandro Dalcin 1189277b19d0SLisandro Dalcin if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1190277b19d0SLisandro Dalcin 1191277b19d0SLisandro Dalcin if (ts->ops->setup) { 1192000e7ae3SMatthew Knepley ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1193277b19d0SLisandro Dalcin } 1194277b19d0SLisandro Dalcin 1195277b19d0SLisandro Dalcin ts->setupcalled = PETSC_TRUE; 1196277b19d0SLisandro Dalcin PetscFunctionReturn(0); 1197277b19d0SLisandro Dalcin } 1198277b19d0SLisandro Dalcin 1199277b19d0SLisandro Dalcin #undef __FUNCT__ 1200277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset" 1201277b19d0SLisandro Dalcin /*@ 1202277b19d0SLisandro Dalcin TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1203277b19d0SLisandro Dalcin 1204277b19d0SLisandro Dalcin Collective on TS 1205277b19d0SLisandro Dalcin 1206277b19d0SLisandro Dalcin Input Parameter: 1207277b19d0SLisandro Dalcin . ts - the TS context obtained from TSCreate() 1208277b19d0SLisandro Dalcin 1209277b19d0SLisandro Dalcin Level: beginner 1210277b19d0SLisandro Dalcin 1211277b19d0SLisandro Dalcin .keywords: TS, timestep, reset 1212277b19d0SLisandro Dalcin 1213277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy() 1214277b19d0SLisandro Dalcin @*/ 1215277b19d0SLisandro Dalcin PetscErrorCode TSReset(TS ts) 1216277b19d0SLisandro Dalcin { 1217277b19d0SLisandro Dalcin PetscErrorCode ierr; 1218277b19d0SLisandro Dalcin 1219277b19d0SLisandro Dalcin PetscFunctionBegin; 1220277b19d0SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1221277b19d0SLisandro Dalcin if (ts->ops->reset) { 1222277b19d0SLisandro Dalcin ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1223277b19d0SLisandro Dalcin } 1224277b19d0SLisandro Dalcin if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1225277b19d0SLisandro Dalcin if (ts->ksp) {ierr = KSPReset(ts->ksp);CHKERRQ(ierr);} 12266bf464f9SBarry Smith ierr = MatDestroy(&ts->A);CHKERRQ(ierr); 12276bf464f9SBarry Smith ierr = MatDestroy(&ts->B);CHKERRQ(ierr); 12286bf464f9SBarry Smith ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 12296bf464f9SBarry Smith ierr = MatDestroy(&ts->Alhs);CHKERRQ(ierr); 12306bf464f9SBarry Smith ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1231277b19d0SLisandro Dalcin if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);} 1232277b19d0SLisandro Dalcin ts->setupcalled = PETSC_FALSE; 1233d763cef2SBarry Smith PetscFunctionReturn(0); 1234d763cef2SBarry Smith } 1235d763cef2SBarry Smith 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "TSDestroy" 1238d8e5e3e6SSatish Balay /*@ 1239d763cef2SBarry Smith TSDestroy - Destroys the timestepper context that was created 1240d763cef2SBarry Smith with TSCreate(). 1241d763cef2SBarry Smith 1242d763cef2SBarry Smith Collective on TS 1243d763cef2SBarry Smith 1244d763cef2SBarry Smith Input Parameter: 1245d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1246d763cef2SBarry Smith 1247d763cef2SBarry Smith Level: beginner 1248d763cef2SBarry Smith 1249d763cef2SBarry Smith .keywords: TS, timestepper, destroy 1250d763cef2SBarry Smith 1251d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve() 1252d763cef2SBarry Smith @*/ 12536bf464f9SBarry Smith PetscErrorCode TSDestroy(TS *ts) 1254d763cef2SBarry Smith { 12556849ba73SBarry Smith PetscErrorCode ierr; 1256d763cef2SBarry Smith 1257d763cef2SBarry Smith PetscFunctionBegin; 12586bf464f9SBarry Smith if (!*ts) PetscFunctionReturn(0); 12596bf464f9SBarry Smith PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 12606bf464f9SBarry Smith if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1261d763cef2SBarry Smith 12626bf464f9SBarry Smith ierr = TSReset((*ts));CHKERRQ(ierr); 1263277b19d0SLisandro Dalcin 1264be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 12656bf464f9SBarry Smith ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 12666bf464f9SBarry Smith if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 12676d4c513bSLisandro Dalcin 12686bf464f9SBarry Smith ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 12696bf464f9SBarry Smith ierr = KSPDestroy(&(*ts)->ksp);CHKERRQ(ierr); 12706bf464f9SBarry Smith ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 12716bf464f9SBarry Smith ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 12726d4c513bSLisandro Dalcin 1273a79aaaedSSatish Balay ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1274d763cef2SBarry Smith PetscFunctionReturn(0); 1275d763cef2SBarry Smith } 1276d763cef2SBarry Smith 12774a2ae208SSatish Balay #undef __FUNCT__ 12784a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES" 1279d8e5e3e6SSatish Balay /*@ 1280d763cef2SBarry Smith TSGetSNES - Returns the SNES (nonlinear solver) associated with 1281d763cef2SBarry Smith a TS (timestepper) context. Valid only for nonlinear problems. 1282d763cef2SBarry Smith 1283d763cef2SBarry Smith Not Collective, but SNES is parallel if TS is parallel 1284d763cef2SBarry Smith 1285d763cef2SBarry Smith Input Parameter: 1286d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1287d763cef2SBarry Smith 1288d763cef2SBarry Smith Output Parameter: 1289d763cef2SBarry Smith . snes - the nonlinear solver context 1290d763cef2SBarry Smith 1291d763cef2SBarry Smith Notes: 1292d763cef2SBarry Smith The user can then directly manipulate the SNES context to set various 1293d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 129494b7f48cSBarry Smith KSP, KSP, and PC contexts as well. 1295d763cef2SBarry Smith 1296d763cef2SBarry Smith TSGetSNES() does not work for integrators that do not use SNES; in 1297d763cef2SBarry Smith this case TSGetSNES() returns PETSC_NULL in snes. 1298d763cef2SBarry Smith 1299d763cef2SBarry Smith Level: beginner 1300d763cef2SBarry Smith 1301d763cef2SBarry Smith .keywords: timestep, get, SNES 1302d763cef2SBarry Smith @*/ 13037087cfbeSBarry Smith PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1304d763cef2SBarry Smith { 1305d372ba47SLisandro Dalcin PetscErrorCode ierr; 1306d372ba47SLisandro Dalcin 1307d763cef2SBarry Smith PetscFunctionBegin; 13080700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13094482741eSBarry Smith PetscValidPointer(snes,2); 131017186662SBarry Smith if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first"); 1311d372ba47SLisandro Dalcin if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1312d372ba47SLisandro Dalcin if (!ts->snes) { 1313d372ba47SLisandro Dalcin ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1314d372ba47SLisandro Dalcin ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1315d372ba47SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1316d372ba47SLisandro Dalcin } 1317d763cef2SBarry Smith *snes = ts->snes; 1318d763cef2SBarry Smith PetscFunctionReturn(0); 1319d763cef2SBarry Smith } 1320d763cef2SBarry Smith 13214a2ae208SSatish Balay #undef __FUNCT__ 132294b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP" 1323d8e5e3e6SSatish Balay /*@ 132494b7f48cSBarry Smith TSGetKSP - Returns the KSP (linear solver) associated with 1325d763cef2SBarry Smith a TS (timestepper) context. 1326d763cef2SBarry Smith 132794b7f48cSBarry Smith Not Collective, but KSP is parallel if TS is parallel 1328d763cef2SBarry Smith 1329d763cef2SBarry Smith Input Parameter: 1330d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1331d763cef2SBarry Smith 1332d763cef2SBarry Smith Output Parameter: 133394b7f48cSBarry Smith . ksp - the nonlinear solver context 1334d763cef2SBarry Smith 1335d763cef2SBarry Smith Notes: 133694b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 1337d763cef2SBarry Smith options, etc. Likewise, the user can then extract and manipulate the 1338d763cef2SBarry Smith KSP and PC contexts as well. 1339d763cef2SBarry Smith 134094b7f48cSBarry Smith TSGetKSP() does not work for integrators that do not use KSP; 134194b7f48cSBarry Smith in this case TSGetKSP() returns PETSC_NULL in ksp. 1342d763cef2SBarry Smith 1343d763cef2SBarry Smith Level: beginner 1344d763cef2SBarry Smith 134594b7f48cSBarry Smith .keywords: timestep, get, KSP 1346d763cef2SBarry Smith @*/ 13477087cfbeSBarry Smith PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1348d763cef2SBarry Smith { 1349d372ba47SLisandro Dalcin PetscErrorCode ierr; 1350d372ba47SLisandro Dalcin 1351d763cef2SBarry Smith PetscFunctionBegin; 13520700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13534482741eSBarry Smith PetscValidPointer(ksp,2); 135417186662SBarry Smith if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1355e32f2f54SBarry Smith if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1356d372ba47SLisandro Dalcin if (!ts->ksp) { 1357d372ba47SLisandro Dalcin ierr = KSPCreate(((PetscObject)ts)->comm,&ts->ksp);CHKERRQ(ierr); 1358d372ba47SLisandro Dalcin ierr = PetscLogObjectParent(ts,ts->ksp);CHKERRQ(ierr); 1359d372ba47SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ts->ksp,(PetscObject)ts,1);CHKERRQ(ierr); 1360d372ba47SLisandro Dalcin } 136194b7f48cSBarry Smith *ksp = ts->ksp; 1362d763cef2SBarry Smith PetscFunctionReturn(0); 1363d763cef2SBarry Smith } 1364d763cef2SBarry Smith 1365d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */ 1366d763cef2SBarry Smith 13674a2ae208SSatish Balay #undef __FUNCT__ 1368adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration" 1369adb62b0dSMatthew Knepley /*@ 1370adb62b0dSMatthew Knepley TSGetDuration - Gets the maximum number of timesteps to use and 1371adb62b0dSMatthew Knepley maximum time for iteration. 1372adb62b0dSMatthew Knepley 13733f9fe445SBarry Smith Not Collective 1374adb62b0dSMatthew Knepley 1375adb62b0dSMatthew Knepley Input Parameters: 1376adb62b0dSMatthew Knepley + ts - the TS context obtained from TSCreate() 1377adb62b0dSMatthew Knepley . maxsteps - maximum number of iterations to use, or PETSC_NULL 1378adb62b0dSMatthew Knepley - maxtime - final time to iterate to, or PETSC_NULL 1379adb62b0dSMatthew Knepley 1380adb62b0dSMatthew Knepley Level: intermediate 1381adb62b0dSMatthew Knepley 1382adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time 1383adb62b0dSMatthew Knepley @*/ 13847087cfbeSBarry Smith PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1385adb62b0dSMatthew Knepley { 1386adb62b0dSMatthew Knepley PetscFunctionBegin; 13870700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1388abc0a331SBarry Smith if (maxsteps) { 13894482741eSBarry Smith PetscValidIntPointer(maxsteps,2); 1390adb62b0dSMatthew Knepley *maxsteps = ts->max_steps; 1391adb62b0dSMatthew Knepley } 1392abc0a331SBarry Smith if (maxtime ) { 13934482741eSBarry Smith PetscValidScalarPointer(maxtime,3); 1394adb62b0dSMatthew Knepley *maxtime = ts->max_time; 1395adb62b0dSMatthew Knepley } 1396adb62b0dSMatthew Knepley PetscFunctionReturn(0); 1397adb62b0dSMatthew Knepley } 1398adb62b0dSMatthew Knepley 1399adb62b0dSMatthew Knepley #undef __FUNCT__ 14004a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration" 1401d763cef2SBarry Smith /*@ 1402d763cef2SBarry Smith TSSetDuration - Sets the maximum number of timesteps to use and 1403d763cef2SBarry Smith maximum time for iteration. 1404d763cef2SBarry Smith 14053f9fe445SBarry Smith Logically Collective on TS 1406d763cef2SBarry Smith 1407d763cef2SBarry Smith Input Parameters: 1408d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1409d763cef2SBarry Smith . maxsteps - maximum number of iterations to use 1410d763cef2SBarry Smith - maxtime - final time to iterate to 1411d763cef2SBarry Smith 1412d763cef2SBarry Smith Options Database Keys: 1413d763cef2SBarry Smith . -ts_max_steps <maxsteps> - Sets maxsteps 1414d763cef2SBarry Smith . -ts_max_time <maxtime> - Sets maxtime 1415d763cef2SBarry Smith 1416d763cef2SBarry Smith Notes: 1417d763cef2SBarry Smith The default maximum number of iterations is 5000. Default time is 5.0 1418d763cef2SBarry Smith 1419d763cef2SBarry Smith Level: intermediate 1420d763cef2SBarry Smith 1421d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations 1422d763cef2SBarry Smith @*/ 14237087cfbeSBarry Smith PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1424d763cef2SBarry Smith { 1425d763cef2SBarry Smith PetscFunctionBegin; 14260700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1427c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1428c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,maxtime,2); 1429d763cef2SBarry Smith ts->max_steps = maxsteps; 1430d763cef2SBarry Smith ts->max_time = maxtime; 1431d763cef2SBarry Smith PetscFunctionReturn(0); 1432d763cef2SBarry Smith } 1433d763cef2SBarry Smith 14344a2ae208SSatish Balay #undef __FUNCT__ 14354a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution" 1436d763cef2SBarry Smith /*@ 1437d763cef2SBarry Smith TSSetSolution - Sets the initial solution vector 1438d763cef2SBarry Smith for use by the TS routines. 1439d763cef2SBarry Smith 14403f9fe445SBarry Smith Logically Collective on TS and Vec 1441d763cef2SBarry Smith 1442d763cef2SBarry Smith Input Parameters: 1443d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1444d763cef2SBarry Smith - x - the solution vector 1445d763cef2SBarry Smith 1446d763cef2SBarry Smith Level: beginner 1447d763cef2SBarry Smith 1448d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions 1449d763cef2SBarry Smith @*/ 14507087cfbeSBarry Smith PetscErrorCode TSSetSolution(TS ts,Vec x) 1451d763cef2SBarry Smith { 14528737fe31SLisandro Dalcin PetscErrorCode ierr; 14538737fe31SLisandro Dalcin 1454d763cef2SBarry Smith PetscFunctionBegin; 14550700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 14560700a824SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,2); 14578737fe31SLisandro Dalcin ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 14586bf464f9SBarry Smith ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 14598737fe31SLisandro Dalcin ts->vec_sol = x; 1460d763cef2SBarry Smith PetscFunctionReturn(0); 1461d763cef2SBarry Smith } 1462d763cef2SBarry Smith 1463e74ef692SMatthew Knepley #undef __FUNCT__ 1464e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep" 1465ac226902SBarry Smith /*@C 1466000e7ae3SMatthew Knepley TSSetPreStep - Sets the general-purpose function 14673f2090d5SJed Brown called once at the beginning of each time step. 1468000e7ae3SMatthew Knepley 14693f9fe445SBarry Smith Logically Collective on TS 1470000e7ae3SMatthew Knepley 1471000e7ae3SMatthew Knepley Input Parameters: 1472000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1473000e7ae3SMatthew Knepley - func - The function 1474000e7ae3SMatthew Knepley 1475000e7ae3SMatthew Knepley Calling sequence of func: 1476000e7ae3SMatthew Knepley . func (TS ts); 1477000e7ae3SMatthew Knepley 1478000e7ae3SMatthew Knepley Level: intermediate 1479000e7ae3SMatthew Knepley 1480000e7ae3SMatthew Knepley .keywords: TS, timestep 1481000e7ae3SMatthew Knepley @*/ 14827087cfbeSBarry Smith PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1483000e7ae3SMatthew Knepley { 1484000e7ae3SMatthew Knepley PetscFunctionBegin; 14850700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1486000e7ae3SMatthew Knepley ts->ops->prestep = func; 1487000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1488000e7ae3SMatthew Knepley } 1489000e7ae3SMatthew Knepley 1490e74ef692SMatthew Knepley #undef __FUNCT__ 14913f2090d5SJed Brown #define __FUNCT__ "TSPreStep" 14923f2090d5SJed Brown /*@C 14933f2090d5SJed Brown TSPreStep - Runs the user-defined pre-step function. 14943f2090d5SJed Brown 14953f2090d5SJed Brown Collective on TS 14963f2090d5SJed Brown 14973f2090d5SJed Brown Input Parameters: 14983f2090d5SJed Brown . ts - The TS context obtained from TSCreate() 14993f2090d5SJed Brown 15003f2090d5SJed Brown Notes: 15013f2090d5SJed Brown TSPreStep() is typically used within time stepping implementations, 15023f2090d5SJed Brown so most users would not generally call this routine themselves. 15033f2090d5SJed Brown 15043f2090d5SJed Brown Level: developer 15053f2090d5SJed Brown 15063f2090d5SJed Brown .keywords: TS, timestep 15073f2090d5SJed Brown @*/ 15087087cfbeSBarry Smith PetscErrorCode TSPreStep(TS ts) 15093f2090d5SJed Brown { 15103f2090d5SJed Brown PetscErrorCode ierr; 15113f2090d5SJed Brown 15123f2090d5SJed Brown PetscFunctionBegin; 15130700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 151472ac3e02SJed Brown if (ts->ops->prestep) { 15153f2090d5SJed Brown PetscStackPush("TS PreStep function"); 15163f2090d5SJed Brown ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 15173f2090d5SJed Brown PetscStackPop; 1518312ce896SJed Brown } 15193f2090d5SJed Brown PetscFunctionReturn(0); 15203f2090d5SJed Brown } 15213f2090d5SJed Brown 15223f2090d5SJed Brown #undef __FUNCT__ 1523e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep" 1524ac226902SBarry Smith /*@C 1525000e7ae3SMatthew Knepley TSSetPostStep - Sets the general-purpose function 15263f2090d5SJed Brown called once at the end of each time step. 1527000e7ae3SMatthew Knepley 15283f9fe445SBarry Smith Logically Collective on TS 1529000e7ae3SMatthew Knepley 1530000e7ae3SMatthew Knepley Input Parameters: 1531000e7ae3SMatthew Knepley + ts - The TS context obtained from TSCreate() 1532000e7ae3SMatthew Knepley - func - The function 1533000e7ae3SMatthew Knepley 1534000e7ae3SMatthew Knepley Calling sequence of func: 1535000e7ae3SMatthew Knepley . func (TS ts); 1536000e7ae3SMatthew Knepley 1537000e7ae3SMatthew Knepley Level: intermediate 1538000e7ae3SMatthew Knepley 1539000e7ae3SMatthew Knepley .keywords: TS, timestep 1540000e7ae3SMatthew Knepley @*/ 15417087cfbeSBarry Smith PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1542000e7ae3SMatthew Knepley { 1543000e7ae3SMatthew Knepley PetscFunctionBegin; 15440700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1545000e7ae3SMatthew Knepley ts->ops->poststep = func; 1546000e7ae3SMatthew Knepley PetscFunctionReturn(0); 1547000e7ae3SMatthew Knepley } 1548000e7ae3SMatthew Knepley 1549e74ef692SMatthew Knepley #undef __FUNCT__ 15503f2090d5SJed Brown #define __FUNCT__ "TSPostStep" 15513f2090d5SJed Brown /*@C 15523f2090d5SJed Brown TSPostStep - Runs the user-defined post-step function. 15533f2090d5SJed Brown 15543f2090d5SJed Brown Collective on TS 15553f2090d5SJed Brown 15563f2090d5SJed Brown Input Parameters: 15573f2090d5SJed Brown . ts - The TS context obtained from TSCreate() 15583f2090d5SJed Brown 15593f2090d5SJed Brown Notes: 15603f2090d5SJed Brown TSPostStep() is typically used within time stepping implementations, 15613f2090d5SJed Brown so most users would not generally call this routine themselves. 15623f2090d5SJed Brown 15633f2090d5SJed Brown Level: developer 15643f2090d5SJed Brown 15653f2090d5SJed Brown .keywords: TS, timestep 15663f2090d5SJed Brown @*/ 15677087cfbeSBarry Smith PetscErrorCode TSPostStep(TS ts) 15683f2090d5SJed Brown { 15693f2090d5SJed Brown PetscErrorCode ierr; 15703f2090d5SJed Brown 15713f2090d5SJed Brown PetscFunctionBegin; 15720700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 157372ac3e02SJed Brown if (ts->ops->poststep) { 15743f2090d5SJed Brown PetscStackPush("TS PostStep function"); 15753f2090d5SJed Brown ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 15763f2090d5SJed Brown PetscStackPop; 157772ac3e02SJed Brown } 15783f2090d5SJed Brown PetscFunctionReturn(0); 15793f2090d5SJed Brown } 15803f2090d5SJed Brown 1581d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 1582d763cef2SBarry Smith 15834a2ae208SSatish Balay #undef __FUNCT__ 1584a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet" 1585d763cef2SBarry Smith /*@C 1586a6570f20SBarry Smith TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1587d763cef2SBarry Smith timestep to display the iteration's progress. 1588d763cef2SBarry Smith 15893f9fe445SBarry Smith Logically Collective on TS 1590d763cef2SBarry Smith 1591d763cef2SBarry Smith Input Parameters: 1592d763cef2SBarry Smith + ts - the TS context obtained from TSCreate() 1593d763cef2SBarry Smith . func - monitoring routine 1594329f5518SBarry Smith . mctx - [optional] user-defined context for private data for the 1595b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desired) 1596b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1597b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 1598d763cef2SBarry Smith 1599d763cef2SBarry Smith Calling sequence of func: 1600a7cc72afSBarry Smith $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1601d763cef2SBarry Smith 1602d763cef2SBarry Smith + ts - the TS context 1603d763cef2SBarry Smith . steps - iteration number 16041f06c33eSBarry Smith . time - current time 1605d763cef2SBarry Smith . x - current iterate 1606d763cef2SBarry Smith - mctx - [optional] monitoring context 1607d763cef2SBarry Smith 1608d763cef2SBarry Smith Notes: 1609d763cef2SBarry Smith This routine adds an additional monitor to the list of monitors that 1610d763cef2SBarry Smith already has been loaded. 1611d763cef2SBarry Smith 1612025f1a04SBarry Smith Fortran notes: Only a single monitor function can be set for each TS object 1613025f1a04SBarry Smith 1614d763cef2SBarry Smith Level: intermediate 1615d763cef2SBarry Smith 1616d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1617d763cef2SBarry Smith 1618a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel() 1619d763cef2SBarry Smith @*/ 1620c2efdce3SBarry Smith PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1621d763cef2SBarry Smith { 1622d763cef2SBarry Smith PetscFunctionBegin; 16230700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 162417186662SBarry Smith if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1625d763cef2SBarry Smith ts->monitor[ts->numbermonitors] = monitor; 1626329f5518SBarry Smith ts->mdestroy[ts->numbermonitors] = mdestroy; 1627d763cef2SBarry Smith ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1628d763cef2SBarry Smith PetscFunctionReturn(0); 1629d763cef2SBarry Smith } 1630d763cef2SBarry Smith 16314a2ae208SSatish Balay #undef __FUNCT__ 1632a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel" 1633d763cef2SBarry Smith /*@C 1634a6570f20SBarry Smith TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1635d763cef2SBarry Smith 16363f9fe445SBarry Smith Logically Collective on TS 1637d763cef2SBarry Smith 1638d763cef2SBarry Smith Input Parameters: 1639d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1640d763cef2SBarry Smith 1641d763cef2SBarry Smith Notes: 1642d763cef2SBarry Smith There is no way to remove a single, specific monitor. 1643d763cef2SBarry Smith 1644d763cef2SBarry Smith Level: intermediate 1645d763cef2SBarry Smith 1646d763cef2SBarry Smith .keywords: TS, timestep, set, monitor 1647d763cef2SBarry Smith 1648a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet() 1649d763cef2SBarry Smith @*/ 16507087cfbeSBarry Smith PetscErrorCode TSMonitorCancel(TS ts) 1651d763cef2SBarry Smith { 1652d952e501SBarry Smith PetscErrorCode ierr; 1653d952e501SBarry Smith PetscInt i; 1654d952e501SBarry Smith 1655d763cef2SBarry Smith PetscFunctionBegin; 16560700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1657d952e501SBarry Smith for (i=0; i<ts->numbermonitors; i++) { 1658d952e501SBarry Smith if (ts->mdestroy[i]) { 16593c4aec1bSBarry Smith ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1660d952e501SBarry Smith } 1661d952e501SBarry Smith } 1662d763cef2SBarry Smith ts->numbermonitors = 0; 1663d763cef2SBarry Smith PetscFunctionReturn(0); 1664d763cef2SBarry Smith } 1665d763cef2SBarry Smith 16664a2ae208SSatish Balay #undef __FUNCT__ 1667a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault" 1668d8e5e3e6SSatish Balay /*@ 1669a6570f20SBarry Smith TSMonitorDefault - Sets the Default monitor 16705516499fSSatish Balay 16715516499fSSatish Balay Level: intermediate 167241251cbbSSatish Balay 16735516499fSSatish Balay .keywords: TS, set, monitor 16745516499fSSatish Balay 167541251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet() 167641251cbbSSatish Balay @*/ 1677649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1678d763cef2SBarry Smith { 1679dfbe8321SBarry Smith PetscErrorCode ierr; 1680649052a6SBarry Smith PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1681d132466eSBarry Smith 1682d763cef2SBarry Smith PetscFunctionBegin; 1683649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1684649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1685649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1686d763cef2SBarry Smith PetscFunctionReturn(0); 1687d763cef2SBarry Smith } 1688d763cef2SBarry Smith 16894a2ae208SSatish Balay #undef __FUNCT__ 16904a2ae208SSatish Balay #define __FUNCT__ "TSStep" 1691d763cef2SBarry Smith /*@ 1692d763cef2SBarry Smith TSStep - Steps the requested number of timesteps. 1693d763cef2SBarry Smith 1694d763cef2SBarry Smith Collective on TS 1695d763cef2SBarry Smith 1696d763cef2SBarry Smith Input Parameter: 1697d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1698d763cef2SBarry Smith 1699d763cef2SBarry Smith Output Parameters: 1700d763cef2SBarry Smith + steps - number of iterations until termination 1701142b95e3SSatish Balay - ptime - time until termination 1702d763cef2SBarry Smith 1703d763cef2SBarry Smith Level: beginner 1704d763cef2SBarry Smith 1705d763cef2SBarry Smith .keywords: TS, timestep, solve 1706d763cef2SBarry Smith 1707d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy() 1708d763cef2SBarry Smith @*/ 17097087cfbeSBarry Smith PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1710d763cef2SBarry Smith { 1711dfbe8321SBarry Smith PetscErrorCode ierr; 1712d763cef2SBarry Smith 1713d763cef2SBarry Smith PetscFunctionBegin; 17140700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1715277b19d0SLisandro Dalcin 1716d405a339SMatthew Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 1717d405a339SMatthew Knepley 1718d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1719000e7ae3SMatthew Knepley ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1720d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1721d405a339SMatthew Knepley 17224bb05414SBarry Smith if (!PetscPreLoadingOn) { 17237adad957SLisandro Dalcin ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1724d405a339SMatthew Knepley } 1725d763cef2SBarry Smith PetscFunctionReturn(0); 1726d763cef2SBarry Smith } 1727d763cef2SBarry Smith 17284a2ae208SSatish Balay #undef __FUNCT__ 17296a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve" 17306a4d4014SLisandro Dalcin /*@ 17316a4d4014SLisandro Dalcin TSSolve - Steps the requested number of timesteps. 17326a4d4014SLisandro Dalcin 17336a4d4014SLisandro Dalcin Collective on TS 17346a4d4014SLisandro Dalcin 17356a4d4014SLisandro Dalcin Input Parameter: 17366a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 17376a4d4014SLisandro Dalcin - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 17386a4d4014SLisandro Dalcin 17396a4d4014SLisandro Dalcin Level: beginner 17406a4d4014SLisandro Dalcin 17416a4d4014SLisandro Dalcin .keywords: TS, timestep, solve 17426a4d4014SLisandro Dalcin 17436a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep() 17446a4d4014SLisandro Dalcin @*/ 17457087cfbeSBarry Smith PetscErrorCode TSSolve(TS ts, Vec x) 17466a4d4014SLisandro Dalcin { 17476a4d4014SLisandro Dalcin PetscInt steps; 17486a4d4014SLisandro Dalcin PetscReal ptime; 17496a4d4014SLisandro Dalcin PetscErrorCode ierr; 1750f22f69f0SBarry Smith 17516a4d4014SLisandro Dalcin PetscFunctionBegin; 17520700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17536a4d4014SLisandro Dalcin /* set solution vector if provided */ 17546a4d4014SLisandro Dalcin if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 17556a4d4014SLisandro Dalcin /* reset time step and iteration counters */ 17566a4d4014SLisandro Dalcin ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 17576a4d4014SLisandro Dalcin /* steps the requested number of timesteps. */ 17586a4d4014SLisandro Dalcin ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 17596a4d4014SLisandro Dalcin PetscFunctionReturn(0); 17606a4d4014SLisandro Dalcin } 17616a4d4014SLisandro Dalcin 17626a4d4014SLisandro Dalcin #undef __FUNCT__ 17634a2ae208SSatish Balay #define __FUNCT__ "TSMonitor" 1764d763cef2SBarry Smith /* 1765d763cef2SBarry Smith Runs the user provided monitor routines, if they exists. 1766d763cef2SBarry Smith */ 1767a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1768d763cef2SBarry Smith { 17696849ba73SBarry Smith PetscErrorCode ierr; 1770a7cc72afSBarry Smith PetscInt i,n = ts->numbermonitors; 1771d763cef2SBarry Smith 1772d763cef2SBarry Smith PetscFunctionBegin; 1773d763cef2SBarry Smith for (i=0; i<n; i++) { 1774142b95e3SSatish Balay ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1775d763cef2SBarry Smith } 1776d763cef2SBarry Smith PetscFunctionReturn(0); 1777d763cef2SBarry Smith } 1778d763cef2SBarry Smith 1779d763cef2SBarry Smith /* ------------------------------------------------------------------------*/ 1780d763cef2SBarry Smith 17814a2ae208SSatish Balay #undef __FUNCT__ 1782a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate" 1783d763cef2SBarry Smith /*@C 1784a6570f20SBarry Smith TSMonitorLGCreate - Creates a line graph context for use with 1785d763cef2SBarry Smith TS to monitor convergence of preconditioned residual norms. 1786d763cef2SBarry Smith 1787d763cef2SBarry Smith Collective on TS 1788d763cef2SBarry Smith 1789d763cef2SBarry Smith Input Parameters: 1790d763cef2SBarry Smith + host - the X display to open, or null for the local machine 1791d763cef2SBarry Smith . label - the title to put in the title bar 17927c922b88SBarry Smith . x, y - the screen coordinates of the upper left coordinate of the window 1793d763cef2SBarry Smith - m, n - the screen width and height in pixels 1794d763cef2SBarry Smith 1795d763cef2SBarry Smith Output Parameter: 1796d763cef2SBarry Smith . draw - the drawing context 1797d763cef2SBarry Smith 1798d763cef2SBarry Smith Options Database Key: 1799a6570f20SBarry Smith . -ts_monitor_draw - automatically sets line graph monitor 1800d763cef2SBarry Smith 1801d763cef2SBarry Smith Notes: 1802a6570f20SBarry Smith Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1803d763cef2SBarry Smith 1804d763cef2SBarry Smith Level: intermediate 1805d763cef2SBarry Smith 18067c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso 1807d763cef2SBarry Smith 1808a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet() 18097c922b88SBarry Smith 1810d763cef2SBarry Smith @*/ 18117087cfbeSBarry Smith PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1812d763cef2SBarry Smith { 1813b0a32e0cSBarry Smith PetscDraw win; 1814dfbe8321SBarry Smith PetscErrorCode ierr; 1815d763cef2SBarry Smith 1816d763cef2SBarry Smith PetscFunctionBegin; 1817b0a32e0cSBarry Smith ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1818b0a32e0cSBarry Smith ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1819b0a32e0cSBarry Smith ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1820b0a32e0cSBarry Smith ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1821d763cef2SBarry Smith 182252e6d16bSBarry Smith ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1823d763cef2SBarry Smith PetscFunctionReturn(0); 1824d763cef2SBarry Smith } 1825d763cef2SBarry Smith 18264a2ae208SSatish Balay #undef __FUNCT__ 1827a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG" 1828a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1829d763cef2SBarry Smith { 1830b0a32e0cSBarry Smith PetscDrawLG lg = (PetscDrawLG) monctx; 183187828ca2SBarry Smith PetscReal x,y = ptime; 1832dfbe8321SBarry Smith PetscErrorCode ierr; 1833d763cef2SBarry Smith 1834d763cef2SBarry Smith PetscFunctionBegin; 18357c922b88SBarry Smith if (!monctx) { 18367c922b88SBarry Smith MPI_Comm comm; 1837b0a32e0cSBarry Smith PetscViewer viewer; 18387c922b88SBarry Smith 18397c922b88SBarry Smith ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1840b0a32e0cSBarry Smith viewer = PETSC_VIEWER_DRAW_(comm); 1841b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 18427c922b88SBarry Smith } 18437c922b88SBarry Smith 1844b0a32e0cSBarry Smith if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 184587828ca2SBarry Smith x = (PetscReal)n; 1846b0a32e0cSBarry Smith ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1847d763cef2SBarry Smith if (n < 20 || (n % 5)) { 1848b0a32e0cSBarry Smith ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1849d763cef2SBarry Smith } 1850d763cef2SBarry Smith PetscFunctionReturn(0); 1851d763cef2SBarry Smith } 1852d763cef2SBarry Smith 18534a2ae208SSatish Balay #undef __FUNCT__ 1854a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy" 1855d763cef2SBarry Smith /*@C 1856a6570f20SBarry Smith TSMonitorLGDestroy - Destroys a line graph context that was created 1857a6570f20SBarry Smith with TSMonitorLGCreate(). 1858d763cef2SBarry Smith 1859b0a32e0cSBarry Smith Collective on PetscDrawLG 1860d763cef2SBarry Smith 1861d763cef2SBarry Smith Input Parameter: 1862d763cef2SBarry Smith . draw - the drawing context 1863d763cef2SBarry Smith 1864d763cef2SBarry Smith Level: intermediate 1865d763cef2SBarry Smith 1866d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy 1867d763cef2SBarry Smith 1868a6570f20SBarry Smith .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1869d763cef2SBarry Smith @*/ 18706bf464f9SBarry Smith PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 1871d763cef2SBarry Smith { 1872b0a32e0cSBarry Smith PetscDraw draw; 1873dfbe8321SBarry Smith PetscErrorCode ierr; 1874d763cef2SBarry Smith 1875d763cef2SBarry Smith PetscFunctionBegin; 18766bf464f9SBarry Smith ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 18776bf464f9SBarry Smith ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 1878b0a32e0cSBarry Smith ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1879d763cef2SBarry Smith PetscFunctionReturn(0); 1880d763cef2SBarry Smith } 1881d763cef2SBarry Smith 18824a2ae208SSatish Balay #undef __FUNCT__ 18834a2ae208SSatish Balay #define __FUNCT__ "TSGetTime" 1884d763cef2SBarry Smith /*@ 1885d763cef2SBarry Smith TSGetTime - Gets the current time. 1886d763cef2SBarry Smith 1887d763cef2SBarry Smith Not Collective 1888d763cef2SBarry Smith 1889d763cef2SBarry Smith Input Parameter: 1890d763cef2SBarry Smith . ts - the TS context obtained from TSCreate() 1891d763cef2SBarry Smith 1892d763cef2SBarry Smith Output Parameter: 1893d763cef2SBarry Smith . t - the current time 1894d763cef2SBarry Smith 1895d763cef2SBarry Smith Level: beginner 1896d763cef2SBarry Smith 1897d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1898d763cef2SBarry Smith 1899d763cef2SBarry Smith .keywords: TS, get, time 1900d763cef2SBarry Smith @*/ 19017087cfbeSBarry Smith PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1902d763cef2SBarry Smith { 1903d763cef2SBarry Smith PetscFunctionBegin; 19040700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 19054482741eSBarry Smith PetscValidDoublePointer(t,2); 1906d763cef2SBarry Smith *t = ts->ptime; 1907d763cef2SBarry Smith PetscFunctionReturn(0); 1908d763cef2SBarry Smith } 1909d763cef2SBarry Smith 19104a2ae208SSatish Balay #undef __FUNCT__ 19116a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime" 19126a4d4014SLisandro Dalcin /*@ 19136a4d4014SLisandro Dalcin TSSetTime - Allows one to reset the time. 19146a4d4014SLisandro Dalcin 19153f9fe445SBarry Smith Logically Collective on TS 19166a4d4014SLisandro Dalcin 19176a4d4014SLisandro Dalcin Input Parameters: 19186a4d4014SLisandro Dalcin + ts - the TS context obtained from TSCreate() 19196a4d4014SLisandro Dalcin - time - the time 19206a4d4014SLisandro Dalcin 19216a4d4014SLisandro Dalcin Level: intermediate 19226a4d4014SLisandro Dalcin 19236a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration() 19246a4d4014SLisandro Dalcin 19256a4d4014SLisandro Dalcin .keywords: TS, set, time 19266a4d4014SLisandro Dalcin @*/ 19277087cfbeSBarry Smith PetscErrorCode TSSetTime(TS ts, PetscReal t) 19286a4d4014SLisandro Dalcin { 19296a4d4014SLisandro Dalcin PetscFunctionBegin; 19300700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1931c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,t,2); 19326a4d4014SLisandro Dalcin ts->ptime = t; 19336a4d4014SLisandro Dalcin PetscFunctionReturn(0); 19346a4d4014SLisandro Dalcin } 19356a4d4014SLisandro Dalcin 19366a4d4014SLisandro Dalcin #undef __FUNCT__ 19374a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix" 1938d763cef2SBarry Smith /*@C 1939d763cef2SBarry Smith TSSetOptionsPrefix - Sets the prefix used for searching for all 1940d763cef2SBarry Smith TS options in the database. 1941d763cef2SBarry Smith 19423f9fe445SBarry Smith Logically Collective on TS 1943d763cef2SBarry Smith 1944d763cef2SBarry Smith Input Parameter: 1945d763cef2SBarry Smith + ts - The TS context 1946d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1947d763cef2SBarry Smith 1948d763cef2SBarry Smith Notes: 1949d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1950d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1951d763cef2SBarry Smith hyphen. 1952d763cef2SBarry Smith 1953d763cef2SBarry Smith Level: advanced 1954d763cef2SBarry Smith 1955d763cef2SBarry Smith .keywords: TS, set, options, prefix, database 1956d763cef2SBarry Smith 1957d763cef2SBarry Smith .seealso: TSSetFromOptions() 1958d763cef2SBarry Smith 1959d763cef2SBarry Smith @*/ 19607087cfbeSBarry Smith PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1961d763cef2SBarry Smith { 1962dfbe8321SBarry Smith PetscErrorCode ierr; 1963d763cef2SBarry Smith 1964d763cef2SBarry Smith PetscFunctionBegin; 19650700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1966d763cef2SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1967d763cef2SBarry Smith switch(ts->problem_type) { 1968d763cef2SBarry Smith case TS_NONLINEAR: 196959580b9cSBarry Smith if (ts->snes) { 1970d763cef2SBarry Smith ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 197159580b9cSBarry Smith } 1972d763cef2SBarry Smith break; 1973d763cef2SBarry Smith case TS_LINEAR: 197459580b9cSBarry Smith if (ts->ksp) { 197594b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 197659580b9cSBarry Smith } 1977d763cef2SBarry Smith break; 1978d763cef2SBarry Smith } 1979d763cef2SBarry Smith PetscFunctionReturn(0); 1980d763cef2SBarry Smith } 1981d763cef2SBarry Smith 1982d763cef2SBarry Smith 19834a2ae208SSatish Balay #undef __FUNCT__ 19844a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix" 1985d763cef2SBarry Smith /*@C 1986d763cef2SBarry Smith TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1987d763cef2SBarry Smith TS options in the database. 1988d763cef2SBarry Smith 19893f9fe445SBarry Smith Logically Collective on TS 1990d763cef2SBarry Smith 1991d763cef2SBarry Smith Input Parameter: 1992d763cef2SBarry Smith + ts - The TS context 1993d763cef2SBarry Smith - prefix - The prefix to prepend to all option names 1994d763cef2SBarry Smith 1995d763cef2SBarry Smith Notes: 1996d763cef2SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1997d763cef2SBarry Smith The first character of all runtime options is AUTOMATICALLY the 1998d763cef2SBarry Smith hyphen. 1999d763cef2SBarry Smith 2000d763cef2SBarry Smith Level: advanced 2001d763cef2SBarry Smith 2002d763cef2SBarry Smith .keywords: TS, append, options, prefix, database 2003d763cef2SBarry Smith 2004d763cef2SBarry Smith .seealso: TSGetOptionsPrefix() 2005d763cef2SBarry Smith 2006d763cef2SBarry Smith @*/ 20077087cfbeSBarry Smith PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 2008d763cef2SBarry Smith { 2009dfbe8321SBarry Smith PetscErrorCode ierr; 2010d763cef2SBarry Smith 2011d763cef2SBarry Smith PetscFunctionBegin; 20120700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2013d763cef2SBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2014d763cef2SBarry Smith switch(ts->problem_type) { 2015d763cef2SBarry Smith case TS_NONLINEAR: 20161ac94b3bSBarry Smith if (ts->snes) { 2017d763cef2SBarry Smith ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 20181ac94b3bSBarry Smith } 2019d763cef2SBarry Smith break; 2020d763cef2SBarry Smith case TS_LINEAR: 20211ac94b3bSBarry Smith if (ts->ksp) { 202294b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 20231ac94b3bSBarry Smith } 2024d763cef2SBarry Smith break; 2025d763cef2SBarry Smith } 2026d763cef2SBarry Smith PetscFunctionReturn(0); 2027d763cef2SBarry Smith } 2028d763cef2SBarry Smith 20294a2ae208SSatish Balay #undef __FUNCT__ 20304a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix" 2031d763cef2SBarry Smith /*@C 2032d763cef2SBarry Smith TSGetOptionsPrefix - Sets the prefix used for searching for all 2033d763cef2SBarry Smith TS options in the database. 2034d763cef2SBarry Smith 2035d763cef2SBarry Smith Not Collective 2036d763cef2SBarry Smith 2037d763cef2SBarry Smith Input Parameter: 2038d763cef2SBarry Smith . ts - The TS context 2039d763cef2SBarry Smith 2040d763cef2SBarry Smith Output Parameter: 2041d763cef2SBarry Smith . prefix - A pointer to the prefix string used 2042d763cef2SBarry Smith 2043d763cef2SBarry Smith Notes: On the fortran side, the user should pass in a string 'prifix' of 2044d763cef2SBarry Smith sufficient length to hold the prefix. 2045d763cef2SBarry Smith 2046d763cef2SBarry Smith Level: intermediate 2047d763cef2SBarry Smith 2048d763cef2SBarry Smith .keywords: TS, get, options, prefix, database 2049d763cef2SBarry Smith 2050d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix() 2051d763cef2SBarry Smith @*/ 20527087cfbeSBarry Smith PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2053d763cef2SBarry Smith { 2054dfbe8321SBarry Smith PetscErrorCode ierr; 2055d763cef2SBarry Smith 2056d763cef2SBarry Smith PetscFunctionBegin; 20570700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 20584482741eSBarry Smith PetscValidPointer(prefix,2); 2059d763cef2SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2060d763cef2SBarry Smith PetscFunctionReturn(0); 2061d763cef2SBarry Smith } 2062d763cef2SBarry Smith 20634a2ae208SSatish Balay #undef __FUNCT__ 20644a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian" 2065d763cef2SBarry Smith /*@C 2066d763cef2SBarry Smith TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2067d763cef2SBarry Smith 2068d763cef2SBarry Smith Not Collective, but parallel objects are returned if TS is parallel 2069d763cef2SBarry Smith 2070d763cef2SBarry Smith Input Parameter: 2071d763cef2SBarry Smith . ts - The TS context obtained from TSCreate() 2072d763cef2SBarry Smith 2073d763cef2SBarry Smith Output Parameters: 2074d763cef2SBarry Smith + J - The Jacobian J of F, where U_t = F(U,t) 2075d763cef2SBarry Smith . M - The preconditioner matrix, usually the same as J 2076d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine 2077d763cef2SBarry Smith 2078d763cef2SBarry Smith Notes: You can pass in PETSC_NULL for any return argument you do not need. 2079d763cef2SBarry Smith 2080d763cef2SBarry Smith Level: intermediate 2081d763cef2SBarry Smith 208226d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2083d763cef2SBarry Smith 2084d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian 2085d763cef2SBarry Smith @*/ 20867087cfbeSBarry Smith PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 2087d763cef2SBarry Smith { 2088d763cef2SBarry Smith PetscFunctionBegin; 208926d46c62SHong Zhang if (J) *J = ts->Arhs; 209026d46c62SHong Zhang if (M) *M = ts->B; 209126d46c62SHong Zhang if (ctx) *ctx = ts->jacP; 2092d763cef2SBarry Smith PetscFunctionReturn(0); 2093d763cef2SBarry Smith } 2094d763cef2SBarry Smith 20951713a123SBarry Smith #undef __FUNCT__ 20962eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian" 20972eca1d9cSJed Brown /*@C 20982eca1d9cSJed Brown TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 20992eca1d9cSJed Brown 21002eca1d9cSJed Brown Not Collective, but parallel objects are returned if TS is parallel 21012eca1d9cSJed Brown 21022eca1d9cSJed Brown Input Parameter: 21032eca1d9cSJed Brown . ts - The TS context obtained from TSCreate() 21042eca1d9cSJed Brown 21052eca1d9cSJed Brown Output Parameters: 21062eca1d9cSJed Brown + A - The Jacobian of F(t,U,U_t) 21072eca1d9cSJed Brown . B - The preconditioner matrix, often the same as A 21082eca1d9cSJed Brown . f - The function to compute the matrices 21092eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine 21102eca1d9cSJed Brown 21112eca1d9cSJed Brown Notes: You can pass in PETSC_NULL for any return argument you do not need. 21122eca1d9cSJed Brown 21132eca1d9cSJed Brown Level: advanced 21142eca1d9cSJed Brown 21152eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 21162eca1d9cSJed Brown 21172eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian 21182eca1d9cSJed Brown @*/ 21197087cfbeSBarry Smith PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 21202eca1d9cSJed Brown { 21212eca1d9cSJed Brown PetscFunctionBegin; 21222eca1d9cSJed Brown if (A) *A = ts->A; 21232eca1d9cSJed Brown if (B) *B = ts->B; 21242eca1d9cSJed Brown if (f) *f = ts->ops->ijacobian; 21252eca1d9cSJed Brown if (ctx) *ctx = ts->jacP; 21262eca1d9cSJed Brown PetscFunctionReturn(0); 21272eca1d9cSJed Brown } 21282eca1d9cSJed Brown 21292eca1d9cSJed Brown #undef __FUNCT__ 2130a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution" 21311713a123SBarry Smith /*@C 2132a6570f20SBarry Smith TSMonitorSolution - Monitors progress of the TS solvers by calling 21331713a123SBarry Smith VecView() for the solution at each timestep 21341713a123SBarry Smith 21351713a123SBarry Smith Collective on TS 21361713a123SBarry Smith 21371713a123SBarry Smith Input Parameters: 21381713a123SBarry Smith + ts - the TS context 21391713a123SBarry Smith . step - current time-step 2140142b95e3SSatish Balay . ptime - current time 21411713a123SBarry Smith - dummy - either a viewer or PETSC_NULL 21421713a123SBarry Smith 21431713a123SBarry Smith Level: intermediate 21441713a123SBarry Smith 21451713a123SBarry Smith .keywords: TS, vector, monitor, view 21461713a123SBarry Smith 2147a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 21481713a123SBarry Smith @*/ 21497087cfbeSBarry Smith PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 21501713a123SBarry Smith { 2151dfbe8321SBarry Smith PetscErrorCode ierr; 21521713a123SBarry Smith PetscViewer viewer = (PetscViewer) dummy; 21531713a123SBarry Smith 21541713a123SBarry Smith PetscFunctionBegin; 2155a34d58ebSBarry Smith if (!dummy) { 21567adad957SLisandro Dalcin viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 21571713a123SBarry Smith } 21581713a123SBarry Smith ierr = VecView(x,viewer);CHKERRQ(ierr); 21591713a123SBarry Smith PetscFunctionReturn(0); 21601713a123SBarry Smith } 21611713a123SBarry Smith 21621713a123SBarry Smith 21636c699258SBarry Smith #undef __FUNCT__ 21646c699258SBarry Smith #define __FUNCT__ "TSSetDM" 21656c699258SBarry Smith /*@ 21666c699258SBarry Smith TSSetDM - Sets the DM that may be used by some preconditioners 21676c699258SBarry Smith 21683f9fe445SBarry Smith Logically Collective on TS and DM 21696c699258SBarry Smith 21706c699258SBarry Smith Input Parameters: 21716c699258SBarry Smith + ts - the preconditioner context 21726c699258SBarry Smith - dm - the dm 21736c699258SBarry Smith 21746c699258SBarry Smith Level: intermediate 21756c699258SBarry Smith 21766c699258SBarry Smith 21776c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 21786c699258SBarry Smith @*/ 21797087cfbeSBarry Smith PetscErrorCode TSSetDM(TS ts,DM dm) 21806c699258SBarry Smith { 21816c699258SBarry Smith PetscErrorCode ierr; 21826c699258SBarry Smith 21836c699258SBarry Smith PetscFunctionBegin; 21840700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 218570663e4aSLisandro Dalcin ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 21866bf464f9SBarry Smith ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 21876c699258SBarry Smith ts->dm = dm; 218870663e4aSLisandro Dalcin if (ts->snes) { 218970663e4aSLisandro Dalcin ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr); 219070663e4aSLisandro Dalcin } 219170663e4aSLisandro Dalcin if (ts->ksp) { 219270663e4aSLisandro Dalcin ierr = KSPSetDM(ts->ksp,dm);CHKERRQ(ierr); 219370663e4aSLisandro Dalcin ierr = KSPSetDMActive(ts->ksp,PETSC_FALSE);CHKERRQ(ierr); 219470663e4aSLisandro Dalcin } 21956c699258SBarry Smith PetscFunctionReturn(0); 21966c699258SBarry Smith } 21976c699258SBarry Smith 21986c699258SBarry Smith #undef __FUNCT__ 21996c699258SBarry Smith #define __FUNCT__ "TSGetDM" 22006c699258SBarry Smith /*@ 22016c699258SBarry Smith TSGetDM - Gets the DM that may be used by some preconditioners 22026c699258SBarry Smith 22033f9fe445SBarry Smith Not Collective 22046c699258SBarry Smith 22056c699258SBarry Smith Input Parameter: 22066c699258SBarry Smith . ts - the preconditioner context 22076c699258SBarry Smith 22086c699258SBarry Smith Output Parameter: 22096c699258SBarry Smith . dm - the dm 22106c699258SBarry Smith 22116c699258SBarry Smith Level: intermediate 22126c699258SBarry Smith 22136c699258SBarry Smith 22146c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 22156c699258SBarry Smith @*/ 22167087cfbeSBarry Smith PetscErrorCode TSGetDM(TS ts,DM *dm) 22176c699258SBarry Smith { 22186c699258SBarry Smith PetscFunctionBegin; 22190700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 22206c699258SBarry Smith *dm = ts->dm; 22216c699258SBarry Smith PetscFunctionReturn(0); 22226c699258SBarry Smith } 22231713a123SBarry Smith 22240f5c6efeSJed Brown #undef __FUNCT__ 22250f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction" 22260f5c6efeSJed Brown /*@ 22270f5c6efeSJed Brown SNESTSFormFunction - Function to evaluate nonlinear residual 22280f5c6efeSJed Brown 22293f9fe445SBarry Smith Logically Collective on SNES 22300f5c6efeSJed Brown 22310f5c6efeSJed Brown Input Parameter: 2232d42a1c89SJed Brown + snes - nonlinear solver 22330f5c6efeSJed Brown . X - the current state at which to evaluate the residual 2234d42a1c89SJed Brown - ctx - user context, must be a TS 22350f5c6efeSJed Brown 22360f5c6efeSJed Brown Output Parameter: 22370f5c6efeSJed Brown . F - the nonlinear residual 22380f5c6efeSJed Brown 22390f5c6efeSJed Brown Notes: 22400f5c6efeSJed Brown This function is not normally called by users and is automatically registered with the SNES used by TS. 22410f5c6efeSJed Brown It is most frequently passed to MatFDColoringSetFunction(). 22420f5c6efeSJed Brown 22430f5c6efeSJed Brown Level: advanced 22440f5c6efeSJed Brown 22450f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction() 22460f5c6efeSJed Brown @*/ 22477087cfbeSBarry Smith PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 22480f5c6efeSJed Brown { 22490f5c6efeSJed Brown TS ts = (TS)ctx; 22500f5c6efeSJed Brown PetscErrorCode ierr; 22510f5c6efeSJed Brown 22520f5c6efeSJed Brown PetscFunctionBegin; 22530f5c6efeSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 22540f5c6efeSJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 22550f5c6efeSJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 22560f5c6efeSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,4); 22570f5c6efeSJed Brown ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 22580f5c6efeSJed Brown PetscFunctionReturn(0); 22590f5c6efeSJed Brown } 22600f5c6efeSJed Brown 22610f5c6efeSJed Brown #undef __FUNCT__ 22620f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian" 22630f5c6efeSJed Brown /*@ 22640f5c6efeSJed Brown SNESTSFormJacobian - Function to evaluate the Jacobian 22650f5c6efeSJed Brown 22660f5c6efeSJed Brown Collective on SNES 22670f5c6efeSJed Brown 22680f5c6efeSJed Brown Input Parameter: 22690f5c6efeSJed Brown + snes - nonlinear solver 22700f5c6efeSJed Brown . X - the current state at which to evaluate the residual 22710f5c6efeSJed Brown - ctx - user context, must be a TS 22720f5c6efeSJed Brown 22730f5c6efeSJed Brown Output Parameter: 22740f5c6efeSJed Brown + A - the Jacobian 22750f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A) 22760f5c6efeSJed Brown - flag - indicates any structure change in the matrix 22770f5c6efeSJed Brown 22780f5c6efeSJed Brown Notes: 22790f5c6efeSJed Brown This function is not normally called by users and is automatically registered with the SNES used by TS. 22800f5c6efeSJed Brown 22810f5c6efeSJed Brown Level: developer 22820f5c6efeSJed Brown 22830f5c6efeSJed Brown .seealso: SNESSetJacobian() 22840f5c6efeSJed Brown @*/ 22857087cfbeSBarry Smith PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 22860f5c6efeSJed Brown { 22870f5c6efeSJed Brown TS ts = (TS)ctx; 22880f5c6efeSJed Brown PetscErrorCode ierr; 22890f5c6efeSJed Brown 22900f5c6efeSJed Brown PetscFunctionBegin; 22910f5c6efeSJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 22920f5c6efeSJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 22930f5c6efeSJed Brown PetscValidPointer(A,3); 22940f5c6efeSJed Brown PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 22950f5c6efeSJed Brown PetscValidPointer(B,4); 22960f5c6efeSJed Brown PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 22970f5c6efeSJed Brown PetscValidPointer(flag,5); 22980f5c6efeSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,6); 22990f5c6efeSJed Brown ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 23000f5c6efeSJed Brown PetscFunctionReturn(0); 23010f5c6efeSJed Brown } 2302325fc9f4SBarry Smith 2303325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 2304c6db04a5SJed Brown #include <mex.h> 2305325fc9f4SBarry Smith 2306325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2307325fc9f4SBarry Smith 2308325fc9f4SBarry Smith #undef __FUNCT__ 2309325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab" 2310325fc9f4SBarry Smith /* 2311325fc9f4SBarry Smith TSComputeFunction_Matlab - Calls the function that has been set with 2312325fc9f4SBarry Smith TSSetFunctionMatlab(). 2313325fc9f4SBarry Smith 2314325fc9f4SBarry Smith Collective on TS 2315325fc9f4SBarry Smith 2316325fc9f4SBarry Smith Input Parameters: 2317325fc9f4SBarry Smith + snes - the TS context 2318325fc9f4SBarry Smith - x - input vector 2319325fc9f4SBarry Smith 2320325fc9f4SBarry Smith Output Parameter: 2321325fc9f4SBarry Smith . y - function vector, as set by TSSetFunction() 2322325fc9f4SBarry Smith 2323325fc9f4SBarry Smith Notes: 2324325fc9f4SBarry Smith TSComputeFunction() is typically used within nonlinear solvers 2325325fc9f4SBarry Smith implementations, so most users would not generally call this routine 2326325fc9f4SBarry Smith themselves. 2327325fc9f4SBarry Smith 2328325fc9f4SBarry Smith Level: developer 2329325fc9f4SBarry Smith 2330325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function 2331325fc9f4SBarry Smith 2332325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction() 2333325fc9f4SBarry Smith */ 23347087cfbeSBarry Smith PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2335325fc9f4SBarry Smith { 2336325fc9f4SBarry Smith PetscErrorCode ierr; 2337325fc9f4SBarry Smith TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2338325fc9f4SBarry Smith int nlhs = 1,nrhs = 7; 2339325fc9f4SBarry Smith mxArray *plhs[1],*prhs[7]; 2340325fc9f4SBarry Smith long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2341325fc9f4SBarry Smith 2342325fc9f4SBarry Smith PetscFunctionBegin; 2343325fc9f4SBarry Smith PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2344325fc9f4SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2345325fc9f4SBarry Smith PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2346325fc9f4SBarry Smith PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2347325fc9f4SBarry Smith PetscCheckSameComm(snes,1,x,3); 2348325fc9f4SBarry Smith PetscCheckSameComm(snes,1,y,5); 2349325fc9f4SBarry Smith 2350325fc9f4SBarry Smith ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2351325fc9f4SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2352880f3077SBarry Smith ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2353325fc9f4SBarry Smith ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2354325fc9f4SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 2355325fc9f4SBarry Smith prhs[1] = mxCreateDoubleScalar(time); 2356325fc9f4SBarry Smith prhs[2] = mxCreateDoubleScalar((double)lx); 2357325fc9f4SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lxdot); 2358325fc9f4SBarry Smith prhs[4] = mxCreateDoubleScalar((double)ly); 2359325fc9f4SBarry Smith prhs[5] = mxCreateString(sctx->funcname); 2360325fc9f4SBarry Smith prhs[6] = sctx->ctx; 2361325fc9f4SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2362325fc9f4SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2363325fc9f4SBarry Smith mxDestroyArray(prhs[0]); 2364325fc9f4SBarry Smith mxDestroyArray(prhs[1]); 2365325fc9f4SBarry Smith mxDestroyArray(prhs[2]); 2366325fc9f4SBarry Smith mxDestroyArray(prhs[3]); 2367325fc9f4SBarry Smith mxDestroyArray(prhs[4]); 2368325fc9f4SBarry Smith mxDestroyArray(prhs[5]); 2369325fc9f4SBarry Smith mxDestroyArray(plhs[0]); 2370325fc9f4SBarry Smith PetscFunctionReturn(0); 2371325fc9f4SBarry Smith } 2372325fc9f4SBarry Smith 2373325fc9f4SBarry Smith 2374325fc9f4SBarry Smith #undef __FUNCT__ 2375325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab" 2376325fc9f4SBarry Smith /* 2377325fc9f4SBarry Smith TSSetFunctionMatlab - Sets the function evaluation routine and function 2378325fc9f4SBarry Smith vector for use by the TS routines in solving ODEs 2379e3c5b3baSBarry Smith equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2380325fc9f4SBarry Smith 2381325fc9f4SBarry Smith Logically Collective on TS 2382325fc9f4SBarry Smith 2383325fc9f4SBarry Smith Input Parameters: 2384325fc9f4SBarry Smith + ts - the TS context 2385325fc9f4SBarry Smith - func - function evaluation routine 2386325fc9f4SBarry Smith 2387325fc9f4SBarry Smith Calling sequence of func: 2388325fc9f4SBarry Smith $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2389325fc9f4SBarry Smith 2390325fc9f4SBarry Smith Level: beginner 2391325fc9f4SBarry Smith 2392325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function 2393325fc9f4SBarry Smith 2394325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2395325fc9f4SBarry Smith */ 23967087cfbeSBarry Smith PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2397325fc9f4SBarry Smith { 2398325fc9f4SBarry Smith PetscErrorCode ierr; 2399325fc9f4SBarry Smith TSMatlabContext *sctx; 2400325fc9f4SBarry Smith 2401325fc9f4SBarry Smith PetscFunctionBegin; 2402325fc9f4SBarry Smith /* currently sctx is memory bleed */ 2403325fc9f4SBarry Smith ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2404325fc9f4SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2405325fc9f4SBarry Smith /* 2406325fc9f4SBarry Smith This should work, but it doesn't 2407325fc9f4SBarry Smith sctx->ctx = ctx; 2408325fc9f4SBarry Smith mexMakeArrayPersistent(sctx->ctx); 2409325fc9f4SBarry Smith */ 2410325fc9f4SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 2411325fc9f4SBarry Smith ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2412325fc9f4SBarry Smith PetscFunctionReturn(0); 2413325fc9f4SBarry Smith } 2414325fc9f4SBarry Smith 2415325fc9f4SBarry Smith #undef __FUNCT__ 2416325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab" 2417325fc9f4SBarry Smith /* 2418325fc9f4SBarry Smith TSComputeJacobian_Matlab - Calls the function that has been set with 2419325fc9f4SBarry Smith TSSetJacobianMatlab(). 2420325fc9f4SBarry Smith 2421325fc9f4SBarry Smith Collective on TS 2422325fc9f4SBarry Smith 2423325fc9f4SBarry Smith Input Parameters: 2424325fc9f4SBarry Smith + snes - the TS context 2425325fc9f4SBarry Smith . x - input vector 2426325fc9f4SBarry Smith . A, B - the matrices 2427325fc9f4SBarry Smith - ctx - user context 2428325fc9f4SBarry Smith 2429325fc9f4SBarry Smith Output Parameter: 2430325fc9f4SBarry Smith . flag - structure of the matrix 2431325fc9f4SBarry Smith 2432325fc9f4SBarry Smith Level: developer 2433325fc9f4SBarry Smith 2434325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function 2435325fc9f4SBarry Smith 2436325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction() 2437325fc9f4SBarry Smith @*/ 24387087cfbeSBarry Smith PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2439325fc9f4SBarry Smith { 2440325fc9f4SBarry Smith PetscErrorCode ierr; 2441325fc9f4SBarry Smith TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2442325fc9f4SBarry Smith int nlhs = 2,nrhs = 9; 2443325fc9f4SBarry Smith mxArray *plhs[2],*prhs[9]; 2444325fc9f4SBarry Smith long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2445325fc9f4SBarry Smith 2446325fc9f4SBarry Smith PetscFunctionBegin; 2447325fc9f4SBarry Smith PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2448325fc9f4SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2449325fc9f4SBarry Smith 2450325fc9f4SBarry Smith /* call Matlab function in ctx with arguments x and y */ 2451325fc9f4SBarry Smith 2452325fc9f4SBarry Smith ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2453325fc9f4SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2454325fc9f4SBarry Smith ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2455325fc9f4SBarry Smith ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2456325fc9f4SBarry Smith ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2457325fc9f4SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 2458325fc9f4SBarry Smith prhs[1] = mxCreateDoubleScalar((double)time); 2459325fc9f4SBarry Smith prhs[2] = mxCreateDoubleScalar((double)lx); 2460325fc9f4SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lxdot); 2461325fc9f4SBarry Smith prhs[4] = mxCreateDoubleScalar((double)shift); 2462325fc9f4SBarry Smith prhs[5] = mxCreateDoubleScalar((double)lA); 2463325fc9f4SBarry Smith prhs[6] = mxCreateDoubleScalar((double)lB); 2464325fc9f4SBarry Smith prhs[7] = mxCreateString(sctx->funcname); 2465325fc9f4SBarry Smith prhs[8] = sctx->ctx; 2466325fc9f4SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2467325fc9f4SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2468325fc9f4SBarry Smith *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2469325fc9f4SBarry Smith mxDestroyArray(prhs[0]); 2470325fc9f4SBarry Smith mxDestroyArray(prhs[1]); 2471325fc9f4SBarry Smith mxDestroyArray(prhs[2]); 2472325fc9f4SBarry Smith mxDestroyArray(prhs[3]); 2473325fc9f4SBarry Smith mxDestroyArray(prhs[4]); 2474325fc9f4SBarry Smith mxDestroyArray(prhs[5]); 2475325fc9f4SBarry Smith mxDestroyArray(prhs[6]); 2476325fc9f4SBarry Smith mxDestroyArray(prhs[7]); 2477325fc9f4SBarry Smith mxDestroyArray(plhs[0]); 2478325fc9f4SBarry Smith mxDestroyArray(plhs[1]); 2479325fc9f4SBarry Smith PetscFunctionReturn(0); 2480325fc9f4SBarry Smith } 2481325fc9f4SBarry Smith 2482325fc9f4SBarry Smith 2483325fc9f4SBarry Smith #undef __FUNCT__ 2484325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab" 2485325fc9f4SBarry Smith /* 2486325fc9f4SBarry Smith TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2487e3c5b3baSBarry 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 2488325fc9f4SBarry Smith 2489325fc9f4SBarry Smith Logically Collective on TS 2490325fc9f4SBarry Smith 2491325fc9f4SBarry Smith Input Parameters: 2492325fc9f4SBarry Smith + snes - the TS context 2493325fc9f4SBarry Smith . A,B - Jacobian matrices 2494325fc9f4SBarry Smith . func - function evaluation routine 2495325fc9f4SBarry Smith - ctx - user context 2496325fc9f4SBarry Smith 2497325fc9f4SBarry Smith Calling sequence of func: 2498325fc9f4SBarry Smith $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2499325fc9f4SBarry Smith 2500325fc9f4SBarry Smith 2501325fc9f4SBarry Smith Level: developer 2502325fc9f4SBarry Smith 2503325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function 2504325fc9f4SBarry Smith 2505325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2506325fc9f4SBarry Smith */ 25077087cfbeSBarry Smith PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2508325fc9f4SBarry Smith { 2509325fc9f4SBarry Smith PetscErrorCode ierr; 2510325fc9f4SBarry Smith TSMatlabContext *sctx; 2511325fc9f4SBarry Smith 2512325fc9f4SBarry Smith PetscFunctionBegin; 2513325fc9f4SBarry Smith /* currently sctx is memory bleed */ 2514325fc9f4SBarry Smith ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2515325fc9f4SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2516325fc9f4SBarry Smith /* 2517325fc9f4SBarry Smith This should work, but it doesn't 2518325fc9f4SBarry Smith sctx->ctx = ctx; 2519325fc9f4SBarry Smith mexMakeArrayPersistent(sctx->ctx); 2520325fc9f4SBarry Smith */ 2521325fc9f4SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 2522325fc9f4SBarry Smith ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2523325fc9f4SBarry Smith PetscFunctionReturn(0); 2524325fc9f4SBarry Smith } 2525325fc9f4SBarry Smith 2526b5b1a830SBarry Smith #undef __FUNCT__ 2527b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab" 2528b5b1a830SBarry Smith /* 2529b5b1a830SBarry Smith TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2530b5b1a830SBarry Smith 2531b5b1a830SBarry Smith Collective on TS 2532b5b1a830SBarry Smith 2533b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction() 2534b5b1a830SBarry Smith @*/ 25357087cfbeSBarry Smith PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2536b5b1a830SBarry Smith { 2537b5b1a830SBarry Smith PetscErrorCode ierr; 2538b5b1a830SBarry Smith TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2539a530c242SBarry Smith int nlhs = 1,nrhs = 6; 2540b5b1a830SBarry Smith mxArray *plhs[1],*prhs[6]; 2541b5b1a830SBarry Smith long long int lx = 0,ls = 0; 2542b5b1a830SBarry Smith 2543b5b1a830SBarry Smith PetscFunctionBegin; 2544b5b1a830SBarry Smith PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2545b5b1a830SBarry Smith PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2546b5b1a830SBarry Smith 2547b5b1a830SBarry Smith ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2548b5b1a830SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2549b5b1a830SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 2550b5b1a830SBarry Smith prhs[1] = mxCreateDoubleScalar((double)it); 2551b5b1a830SBarry Smith prhs[2] = mxCreateDoubleScalar((double)time); 2552b5b1a830SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lx); 2553b5b1a830SBarry Smith prhs[4] = mxCreateString(sctx->funcname); 2554b5b1a830SBarry Smith prhs[5] = sctx->ctx; 2555b5b1a830SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2556b5b1a830SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2557b5b1a830SBarry Smith mxDestroyArray(prhs[0]); 2558b5b1a830SBarry Smith mxDestroyArray(prhs[1]); 2559b5b1a830SBarry Smith mxDestroyArray(prhs[2]); 2560b5b1a830SBarry Smith mxDestroyArray(prhs[3]); 2561b5b1a830SBarry Smith mxDestroyArray(prhs[4]); 2562b5b1a830SBarry Smith mxDestroyArray(plhs[0]); 2563b5b1a830SBarry Smith PetscFunctionReturn(0); 2564b5b1a830SBarry Smith } 2565b5b1a830SBarry Smith 2566b5b1a830SBarry Smith 2567b5b1a830SBarry Smith #undef __FUNCT__ 2568b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab" 2569b5b1a830SBarry Smith /* 2570b5b1a830SBarry Smith TSMonitorSetMatlab - Sets the monitor function from Matlab 2571b5b1a830SBarry Smith 2572b5b1a830SBarry Smith Level: developer 2573b5b1a830SBarry Smith 2574b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function 2575b5b1a830SBarry Smith 2576b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2577b5b1a830SBarry Smith */ 25787087cfbeSBarry Smith PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2579b5b1a830SBarry Smith { 2580b5b1a830SBarry Smith PetscErrorCode ierr; 2581b5b1a830SBarry Smith TSMatlabContext *sctx; 2582b5b1a830SBarry Smith 2583b5b1a830SBarry Smith PetscFunctionBegin; 2584b5b1a830SBarry Smith /* currently sctx is memory bleed */ 2585b5b1a830SBarry Smith ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2586b5b1a830SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2587b5b1a830SBarry Smith /* 2588b5b1a830SBarry Smith This should work, but it doesn't 2589b5b1a830SBarry Smith sctx->ctx = ctx; 2590b5b1a830SBarry Smith mexMakeArrayPersistent(sctx->ctx); 2591b5b1a830SBarry Smith */ 2592b5b1a830SBarry Smith sctx->ctx = mxDuplicateArray(ctx); 2593b5b1a830SBarry Smith ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2594b5b1a830SBarry Smith PetscFunctionReturn(0); 2595b5b1a830SBarry Smith } 2596b5b1a830SBarry Smith 2597325fc9f4SBarry Smith #endif 2598