1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 57445fe48SJed Brown #include <petscsnesfas.h> 6316643e7SJed Brown 7316643e7SJed Brown typedef struct { 8316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 9eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 10ace3abfcSBarry Smith PetscBool extrapolate; 11eb284becSJed Brown PetscBool endpoint; 12316643e7SJed Brown PetscReal Theta; 13316643e7SJed Brown PetscReal shift; 14316643e7SJed Brown PetscReal stage_time; 15316643e7SJed Brown } TS_Theta; 16316643e7SJed Brown 17316643e7SJed Brown #undef __FUNCT__ 187445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 197445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 207445fe48SJed Brown { 217445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 227445fe48SJed Brown PetscErrorCode ierr; 237445fe48SJed Brown 247445fe48SJed Brown PetscFunctionBegin; 257445fe48SJed Brown if (X0) { 267445fe48SJed Brown if (dm && dm != ts->dm) { 270d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 287445fe48SJed Brown } else *X0 = ts->vec_sol; 297445fe48SJed Brown } 307445fe48SJed Brown if (Xdot) { 317445fe48SJed Brown if (dm && dm != ts->dm) { 320d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 337445fe48SJed Brown } else *Xdot = th->Xdot; 347445fe48SJed Brown } 357445fe48SJed Brown PetscFunctionReturn(0); 367445fe48SJed Brown } 377445fe48SJed Brown 380d0b770aSPeter Brune 390d0b770aSPeter Brune #undef __FUNCT__ 400d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 410d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 420d0b770aSPeter Brune { 430d0b770aSPeter Brune PetscErrorCode ierr; 440d0b770aSPeter Brune 450d0b770aSPeter Brune PetscFunctionBegin; 460d0b770aSPeter Brune if (X0) { 470d0b770aSPeter Brune if (dm && dm != ts->dm) { 480d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 490d0b770aSPeter Brune } 500d0b770aSPeter Brune } 510d0b770aSPeter Brune if (Xdot) { 520d0b770aSPeter Brune if (dm && dm != ts->dm) { 530d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 540d0b770aSPeter Brune } 550d0b770aSPeter Brune } 560d0b770aSPeter Brune PetscFunctionReturn(0); 570d0b770aSPeter Brune } 580d0b770aSPeter Brune 597445fe48SJed Brown #undef __FUNCT__ 607445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 617445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 627445fe48SJed Brown { 637445fe48SJed Brown 647445fe48SJed Brown PetscFunctionBegin; 657445fe48SJed Brown PetscFunctionReturn(0); 667445fe48SJed Brown } 677445fe48SJed Brown 687445fe48SJed Brown #undef __FUNCT__ 697445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 707445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 717445fe48SJed Brown { 727445fe48SJed Brown TS ts = (TS)ctx; 737445fe48SJed Brown PetscErrorCode ierr; 747445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 757445fe48SJed Brown 767445fe48SJed Brown PetscFunctionBegin; 777445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 787445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 797445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 807445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 817445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 827445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 830d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 840d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 857445fe48SJed Brown PetscFunctionReturn(0); 867445fe48SJed Brown } 877445fe48SJed Brown 887445fe48SJed Brown #undef __FUNCT__ 89316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 90193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 91316643e7SJed Brown { 92316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 934e17fbf5SBarry Smith PetscInt its,lits,reject; 94*a3314f2cSMatthew G Knepley PetscReal next_time_step = 0.0; 95f1b97656SJed Brown SNESConvergedReason snesreason; 962b5a38e1SLisandro Dalcin PetscErrorCode ierr; 97316643e7SJed Brown 98316643e7SJed Brown PetscFunctionBegin; 9931748224SBarry Smith if (ts->time_steps_since_decrease > 3 && ts->time_step < ts->time_step_orig) { 10031748224SBarry Smith /* smaller time step has worked successfully for three time-steps, try increasing time step*/ 10131748224SBarry Smith ts->time_step = 2.0*ts->time_step; 10231748224SBarry Smith ts->time_steps_since_decrease = 0; /* don't want to increase time step two time steps in a row */ 10331748224SBarry Smith } 1044e17fbf5SBarry Smith for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 105cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 106eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 107316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 108b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 109b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 110316643e7SJed Brown 111eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 112eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 113eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 114eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 115eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 116eb284becSJed Brown } 117ace68cafSJed Brown if (th->extrapolate) { 1182b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 119ace68cafSJed Brown } else { 1202b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 121ace68cafSJed Brown } 122eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 123316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 124316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 125f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 1265ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 12731748224SBarry Smith if (its < 10) ts->time_steps_since_decrease++; 12831748224SBarry Smith else ts->time_steps_since_decrease = 0; 1294e17fbf5SBarry Smith if (snesreason > 0) break; 1304e17fbf5SBarry Smith ierr = PetscInfo3(ts,"Step=%D, Cutting time-step from %g to %g\n",ts->steps,(double)ts->time_step,(double).5*ts->time_step);CHKERRQ(ierr); 1314e17fbf5SBarry Smith ts->time_step = .5*ts->time_step; 13231748224SBarry Smith ts->time_steps_since_decrease = 0; 1334e17fbf5SBarry Smith } 134f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 135f1b97656SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 136f1b97656SJed Brown ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 137f1b97656SJed Brown PetscFunctionReturn(0); 138f1b97656SJed Brown } 139eb284becSJed Brown if (th->endpoint) { 140eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 141eb284becSJed Brown } else { 1422b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 1432b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 144eb284becSJed Brown } 1452b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 146cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 147316643e7SJed Brown ts->steps++; 148316643e7SJed Brown PetscFunctionReturn(0); 149316643e7SJed Brown } 150316643e7SJed Brown 151cd652676SJed Brown #undef __FUNCT__ 152cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 153cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 154cd652676SJed Brown { 155cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1565a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 157cd652676SJed Brown PetscErrorCode ierr; 158cd652676SJed Brown 159cd652676SJed Brown PetscFunctionBegin; 160a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 1615a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 1625a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 163cd652676SJed Brown PetscFunctionReturn(0); 164cd652676SJed Brown } 165cd652676SJed Brown 166316643e7SJed Brown /*------------------------------------------------------------*/ 167316643e7SJed Brown #undef __FUNCT__ 168277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 169277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 170316643e7SJed Brown { 171316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 172316643e7SJed Brown PetscErrorCode ierr; 173316643e7SJed Brown 174316643e7SJed Brown PetscFunctionBegin; 1756bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 1766bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 177eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 178277b19d0SLisandro Dalcin PetscFunctionReturn(0); 179277b19d0SLisandro Dalcin } 180277b19d0SLisandro Dalcin 181277b19d0SLisandro Dalcin #undef __FUNCT__ 182277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 183277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 184277b19d0SLisandro Dalcin { 185277b19d0SLisandro Dalcin PetscErrorCode ierr; 186277b19d0SLisandro Dalcin 187277b19d0SLisandro Dalcin PetscFunctionBegin; 188277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 189277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 190335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 191335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 19226f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 193eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 194316643e7SJed Brown PetscFunctionReturn(0); 195316643e7SJed Brown } 196316643e7SJed Brown 197316643e7SJed Brown /* 198316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1992b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 200316643e7SJed Brown */ 201316643e7SJed Brown #undef __FUNCT__ 2020f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 2030f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 204316643e7SJed Brown { 205316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 206316643e7SJed Brown PetscErrorCode ierr; 2077445fe48SJed Brown Vec X0,Xdot; 2087445fe48SJed Brown DM dm,dmsave; 209316643e7SJed Brown 210316643e7SJed Brown PetscFunctionBegin; 2117445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2125a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 2137445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 2147445fe48SJed Brown ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr); 2157445fe48SJed Brown 2167445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 2177445fe48SJed Brown dmsave = ts->dm; 2187445fe48SJed Brown ts->dm = dm; 2197445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 2207445fe48SJed Brown ts->dm = dmsave; 2210d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 222316643e7SJed Brown PetscFunctionReturn(0); 223316643e7SJed Brown } 224316643e7SJed Brown 225316643e7SJed Brown #undef __FUNCT__ 2260f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 2270f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 228316643e7SJed Brown { 229316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 230316643e7SJed Brown PetscErrorCode ierr; 2317445fe48SJed Brown Vec Xdot; 2327445fe48SJed Brown DM dm,dmsave; 233316643e7SJed Brown 234316643e7SJed Brown PetscFunctionBegin; 2357445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2367445fe48SJed Brown 2370f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 2387445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 2397445fe48SJed Brown 2407445fe48SJed Brown dmsave = ts->dm; 2417445fe48SJed Brown ts->dm = dm; 2427445fe48SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 2437445fe48SJed Brown ts->dm = dmsave; 2440d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 245316643e7SJed Brown PetscFunctionReturn(0); 246316643e7SJed Brown } 247316643e7SJed Brown 248316643e7SJed Brown #undef __FUNCT__ 249316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 250316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 251316643e7SJed Brown { 252316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 253316643e7SJed Brown PetscErrorCode ierr; 2547445fe48SJed Brown SNES snes; 2557445fe48SJed Brown DM dm; 256316643e7SJed Brown 257316643e7SJed Brown PetscFunctionBegin; 258316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 259316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 2607445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2617445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2627445fe48SJed Brown if (dm) { 2637445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 2647445fe48SJed Brown } 265316643e7SJed Brown PetscFunctionReturn(0); 266316643e7SJed Brown } 267316643e7SJed Brown /*------------------------------------------------------------*/ 268316643e7SJed Brown 269316643e7SJed Brown #undef __FUNCT__ 270316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 271316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 272316643e7SJed Brown { 273316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 274316643e7SJed Brown PetscErrorCode ierr; 275316643e7SJed Brown 276316643e7SJed Brown PetscFunctionBegin; 277d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 278316643e7SJed Brown { 279316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 280acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 281eb284becSJed Brown ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr); 282d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 283316643e7SJed Brown } 284316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 285316643e7SJed Brown PetscFunctionReturn(0); 286316643e7SJed Brown } 287316643e7SJed Brown 288316643e7SJed Brown #undef __FUNCT__ 289316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 290316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 291316643e7SJed Brown { 292316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 293ace3abfcSBarry Smith PetscBool iascii; 294316643e7SJed Brown PetscErrorCode ierr; 295316643e7SJed Brown 296316643e7SJed Brown PetscFunctionBegin; 297251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 298316643e7SJed Brown if (iascii) { 299316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 300ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 301316643e7SJed Brown } 302d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 303316643e7SJed Brown PetscFunctionReturn(0); 304316643e7SJed Brown } 305316643e7SJed Brown 3060de4c49aSJed Brown EXTERN_C_BEGIN 3070de4c49aSJed Brown #undef __FUNCT__ 3080de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 3097087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 3100de4c49aSJed Brown { 3110de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3120de4c49aSJed Brown 3130de4c49aSJed Brown PetscFunctionBegin; 3140de4c49aSJed Brown *theta = th->Theta; 3150de4c49aSJed Brown PetscFunctionReturn(0); 3160de4c49aSJed Brown } 3170de4c49aSJed Brown 3180de4c49aSJed Brown #undef __FUNCT__ 3190de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 3207087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 3210de4c49aSJed Brown { 3220de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3230de4c49aSJed Brown 3240de4c49aSJed Brown PetscFunctionBegin; 325e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 3260de4c49aSJed Brown th->Theta = theta; 3270de4c49aSJed Brown PetscFunctionReturn(0); 3280de4c49aSJed Brown } 329eb284becSJed Brown 330eb284becSJed Brown #undef __FUNCT__ 33178e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 33226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 33326f2ff8fSLisandro Dalcin { 33426f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 33526f2ff8fSLisandro Dalcin 33626f2ff8fSLisandro Dalcin PetscFunctionBegin; 33726f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 33826f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 33926f2ff8fSLisandro Dalcin } 34026f2ff8fSLisandro Dalcin 34126f2ff8fSLisandro Dalcin #undef __FUNCT__ 34226f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 343eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 344eb284becSJed Brown { 345eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 346eb284becSJed Brown 347eb284becSJed Brown PetscFunctionBegin; 348eb284becSJed Brown th->endpoint = flg; 349eb284becSJed Brown PetscFunctionReturn(0); 350eb284becSJed Brown } 3510de4c49aSJed Brown EXTERN_C_END 3520de4c49aSJed Brown 353316643e7SJed Brown /* ------------------------------------------------------------ */ 354316643e7SJed Brown /*MC 35596f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 356316643e7SJed Brown 357316643e7SJed Brown Level: beginner 358316643e7SJed Brown 3594eb428fdSBarry Smith Options Database: 3604eb428fdSBarry Smith -ts_theta_theta <Theta> - Location of stage (0<Theta<=1); Theta = 1.0 ( 3614eb428fdSBarry Smith -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 3624eb428fdSBarry Smith -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method 3634eb428fdSBarry Smith 364eb284becSJed Brown Notes: 3654eb428fdSBarry Smith $ -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER) 3664eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 3674eb428fdSBarry Smith 3684eb428fdSBarry Smith 3694eb428fdSBarry Smith 370eb284becSJed Brown This method can be applied to DAE. 371eb284becSJed Brown 372eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 373eb284becSJed Brown 374eb284becSJed Brown .vb 375eb284becSJed Brown Theta | Theta 376eb284becSJed Brown ------------- 377eb284becSJed Brown | 1 378eb284becSJed Brown .ve 379eb284becSJed Brown 380eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 381eb284becSJed Brown 382eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 383eb284becSJed Brown 384eb284becSJed Brown .vb 385eb284becSJed Brown 0 | 0 0 386eb284becSJed Brown 1 | 1-Theta Theta 387eb284becSJed Brown ------------------- 388eb284becSJed Brown | 1-Theta Theta 389eb284becSJed Brown .ve 390eb284becSJed Brown 391eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 392eb284becSJed Brown 393eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 394eb284becSJed Brown 395eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 396eb284becSJed Brown 3974eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 398eb284becSJed Brown 399eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 400316643e7SJed Brown 401316643e7SJed Brown M*/ 402316643e7SJed Brown EXTERN_C_BEGIN 403316643e7SJed Brown #undef __FUNCT__ 404316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 4057087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 406316643e7SJed Brown { 407316643e7SJed Brown TS_Theta *th; 408316643e7SJed Brown PetscErrorCode ierr; 409316643e7SJed Brown 410316643e7SJed Brown PetscFunctionBegin; 411277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 412316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 413316643e7SJed Brown ts->ops->view = TSView_Theta; 414316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 415316643e7SJed Brown ts->ops->step = TSStep_Theta; 416cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 417316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 4180f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 4190f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 420316643e7SJed Brown 421316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 422316643e7SJed Brown ts->data = (void*)th; 423316643e7SJed Brown 4246f700aefSJed Brown th->extrapolate = PETSC_FALSE; 425316643e7SJed Brown th->Theta = 0.5; 426316643e7SJed Brown 4270de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 4280de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 42926f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 430eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 431316643e7SJed Brown PetscFunctionReturn(0); 432316643e7SJed Brown } 433316643e7SJed Brown EXTERN_C_END 4340de4c49aSJed Brown 4350de4c49aSJed Brown #undef __FUNCT__ 4360de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 4370de4c49aSJed Brown /*@ 4380de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 4390de4c49aSJed Brown 4400de4c49aSJed Brown Not Collective 4410de4c49aSJed Brown 4420de4c49aSJed Brown Input Parameter: 4430de4c49aSJed Brown . ts - timestepping context 4440de4c49aSJed Brown 4450de4c49aSJed Brown Output Parameter: 4460de4c49aSJed Brown . theta - stage abscissa 4470de4c49aSJed Brown 4480de4c49aSJed Brown Note: 4490de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 4500de4c49aSJed Brown 4510de4c49aSJed Brown Level: Advanced 4520de4c49aSJed Brown 4530de4c49aSJed Brown .seealso: TSThetaSetTheta() 4540de4c49aSJed Brown @*/ 4557087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 4560de4c49aSJed Brown { 4574ac538c5SBarry Smith PetscErrorCode ierr; 4580de4c49aSJed Brown 4590de4c49aSJed Brown PetscFunctionBegin; 460afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4610de4c49aSJed Brown PetscValidPointer(theta,2); 4624ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 4630de4c49aSJed Brown PetscFunctionReturn(0); 4640de4c49aSJed Brown } 4650de4c49aSJed Brown 4660de4c49aSJed Brown #undef __FUNCT__ 4670de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 4680de4c49aSJed Brown /*@ 4690de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 4700de4c49aSJed Brown 4710de4c49aSJed Brown Not Collective 4720de4c49aSJed Brown 4730de4c49aSJed Brown Input Parameter: 4740de4c49aSJed Brown + ts - timestepping context 4750de4c49aSJed Brown - theta - stage abscissa 4760de4c49aSJed Brown 4770de4c49aSJed Brown Options Database: 4780de4c49aSJed Brown . -ts_theta_theta <theta> 4790de4c49aSJed Brown 4800de4c49aSJed Brown Level: Intermediate 4810de4c49aSJed Brown 4820de4c49aSJed Brown .seealso: TSThetaGetTheta() 4830de4c49aSJed Brown @*/ 4847087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 4850de4c49aSJed Brown { 4864ac538c5SBarry Smith PetscErrorCode ierr; 4870de4c49aSJed Brown 4880de4c49aSJed Brown PetscFunctionBegin; 489afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4904ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 4910de4c49aSJed Brown PetscFunctionReturn(0); 4920de4c49aSJed Brown } 493f33bbcb6SJed Brown 494eb284becSJed Brown #undef __FUNCT__ 49526f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 49626f2ff8fSLisandro Dalcin /*@ 49726f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 49826f2ff8fSLisandro Dalcin 49926f2ff8fSLisandro Dalcin Not Collective 50026f2ff8fSLisandro Dalcin 50126f2ff8fSLisandro Dalcin Input Parameter: 50226f2ff8fSLisandro Dalcin . ts - timestepping context 50326f2ff8fSLisandro Dalcin 50426f2ff8fSLisandro Dalcin Output Parameter: 50526f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 50626f2ff8fSLisandro Dalcin 50726f2ff8fSLisandro Dalcin Level: Advanced 50826f2ff8fSLisandro Dalcin 50926f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 51026f2ff8fSLisandro Dalcin @*/ 51126f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 51226f2ff8fSLisandro Dalcin { 51326f2ff8fSLisandro Dalcin PetscErrorCode ierr; 51426f2ff8fSLisandro Dalcin 51526f2ff8fSLisandro Dalcin PetscFunctionBegin; 51626f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 51726f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 51826f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 51926f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 52026f2ff8fSLisandro Dalcin } 52126f2ff8fSLisandro Dalcin 52226f2ff8fSLisandro Dalcin #undef __FUNCT__ 523eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 524eb284becSJed Brown /*@ 525eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 526eb284becSJed Brown 527eb284becSJed Brown Not Collective 528eb284becSJed Brown 529eb284becSJed Brown Input Parameter: 530eb284becSJed Brown + ts - timestepping context 531eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 532eb284becSJed Brown 533eb284becSJed Brown Options Database: 534eb284becSJed Brown . -ts_theta_endpoint <flg> 535eb284becSJed Brown 536eb284becSJed Brown Level: Intermediate 537eb284becSJed Brown 538eb284becSJed Brown .seealso: TSTHETA, TSCN 539eb284becSJed Brown @*/ 540eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 541eb284becSJed Brown { 542eb284becSJed Brown PetscErrorCode ierr; 543eb284becSJed Brown 544eb284becSJed Brown PetscFunctionBegin; 545eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 546eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 547eb284becSJed Brown PetscFunctionReturn(0); 548eb284becSJed Brown } 549eb284becSJed Brown 550f33bbcb6SJed Brown /* 551f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 552f33bbcb6SJed Brown * The creation functions for these specializations are below. 553f33bbcb6SJed Brown */ 554f33bbcb6SJed Brown 555f33bbcb6SJed Brown #undef __FUNCT__ 556f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 557f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 558f33bbcb6SJed Brown { 559d52bd9f3SBarry Smith PetscErrorCode ierr; 560d52bd9f3SBarry Smith 561f33bbcb6SJed Brown PetscFunctionBegin; 562d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 563f33bbcb6SJed Brown PetscFunctionReturn(0); 564f33bbcb6SJed Brown } 565f33bbcb6SJed Brown 566f33bbcb6SJed Brown /*MC 567f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 568f33bbcb6SJed Brown 569f33bbcb6SJed Brown Level: beginner 570f33bbcb6SJed Brown 5714eb428fdSBarry Smith Notes: 5724eb428fdSBarry Smith TSCN is equivalent to TSTHETA with Theta=1.0 5734eb428fdSBarry Smith 5744eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 1. 5754eb428fdSBarry Smith 576f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 577f33bbcb6SJed Brown 578f33bbcb6SJed Brown M*/ 579f33bbcb6SJed Brown EXTERN_C_BEGIN 580f33bbcb6SJed Brown #undef __FUNCT__ 581f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 582f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 583f33bbcb6SJed Brown { 584f33bbcb6SJed Brown PetscErrorCode ierr; 585f33bbcb6SJed Brown 586f33bbcb6SJed Brown PetscFunctionBegin; 587f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 588f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 589f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 590f33bbcb6SJed Brown PetscFunctionReturn(0); 591f33bbcb6SJed Brown } 592f33bbcb6SJed Brown EXTERN_C_END 593f33bbcb6SJed Brown 594f33bbcb6SJed Brown #undef __FUNCT__ 595f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 596f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 597f33bbcb6SJed Brown { 598d52bd9f3SBarry Smith PetscErrorCode ierr; 599d52bd9f3SBarry Smith 600f33bbcb6SJed Brown PetscFunctionBegin; 601d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 602f33bbcb6SJed Brown PetscFunctionReturn(0); 603f33bbcb6SJed Brown } 604f33bbcb6SJed Brown 605f33bbcb6SJed Brown /*MC 606f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 607f33bbcb6SJed Brown 608f33bbcb6SJed Brown Level: beginner 609f33bbcb6SJed Brown 610f33bbcb6SJed Brown Notes: 6117cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 6127cf5af47SJed Brown 6137cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 614f33bbcb6SJed Brown 615f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 616f33bbcb6SJed Brown 617f33bbcb6SJed Brown M*/ 618f33bbcb6SJed Brown EXTERN_C_BEGIN 619f33bbcb6SJed Brown #undef __FUNCT__ 620f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 621f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 622f33bbcb6SJed Brown { 623f33bbcb6SJed Brown PetscErrorCode ierr; 624f33bbcb6SJed Brown 625f33bbcb6SJed Brown PetscFunctionBegin; 626f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 627f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 628eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 629f33bbcb6SJed Brown ts->ops->view = TSView_CN; 630f33bbcb6SJed Brown PetscFunctionReturn(0); 631f33bbcb6SJed Brown } 632f33bbcb6SJed Brown EXTERN_C_END 633