1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4f9c1d6abSBarry Smith #define PETSC_DESIRE_COMPLEX 5b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 67445fe48SJed Brown #include <petscsnesfas.h> 7316643e7SJed Brown 8316643e7SJed Brown typedef struct { 9316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 10eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 11ace3abfcSBarry Smith PetscBool extrapolate; 12eb284becSJed Brown PetscBool endpoint; 13316643e7SJed Brown PetscReal Theta; 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__ 89*258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta" 90*258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 91*258e1594SPeter Brune { 92*258e1594SPeter Brune 93*258e1594SPeter Brune PetscFunctionBegin; 94*258e1594SPeter Brune PetscFunctionReturn(0); 95*258e1594SPeter Brune } 96*258e1594SPeter Brune 97*258e1594SPeter Brune #undef __FUNCT__ 98*258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta" 99*258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 100*258e1594SPeter Brune { 101*258e1594SPeter Brune TS ts = (TS)ctx; 102*258e1594SPeter Brune PetscErrorCode ierr; 103*258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 104*258e1594SPeter Brune 105*258e1594SPeter Brune PetscFunctionBegin; 106*258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 107*258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 108*258e1594SPeter Brune 109*258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 110*258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 111*258e1594SPeter Brune 112*258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 113*258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 114*258e1594SPeter Brune 115*258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 116*258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 117*258e1594SPeter Brune PetscFunctionReturn(0); 118*258e1594SPeter Brune } 119*258e1594SPeter Brune 120*258e1594SPeter Brune 121*258e1594SPeter Brune #undef __FUNCT__ 122316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 123193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 124316643e7SJed Brown { 125316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1264e17fbf5SBarry Smith PetscInt its,lits,reject; 127a3314f2cSMatthew G Knepley PetscReal next_time_step = 0.0; 128f1b97656SJed Brown SNESConvergedReason snesreason; 1292b5a38e1SLisandro Dalcin PetscErrorCode ierr; 130316643e7SJed Brown 131316643e7SJed Brown PetscFunctionBegin; 13231748224SBarry Smith if (ts->time_steps_since_decrease > 3 && ts->time_step < ts->time_step_orig) { 13331748224SBarry Smith /* smaller time step has worked successfully for three time-steps, try increasing time step*/ 13431748224SBarry Smith ts->time_step = 2.0*ts->time_step; 13531748224SBarry Smith ts->time_steps_since_decrease = 0; /* don't want to increase time step two time steps in a row */ 13631748224SBarry Smith } 1374e17fbf5SBarry Smith for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 138b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 139cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 140eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 141b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 142b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 143316643e7SJed Brown 144eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 145eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 146eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 147eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 148eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 149eb284becSJed Brown } 150ace68cafSJed Brown if (th->extrapolate) { 151b296d7d5SJed Brown ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 152ace68cafSJed Brown } else { 1532b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 154ace68cafSJed Brown } 155eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 156316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 157316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 158f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 1595ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 16031748224SBarry Smith if (its < 10) ts->time_steps_since_decrease++; 16131748224SBarry Smith else ts->time_steps_since_decrease = 0; 1624e17fbf5SBarry Smith if (snesreason > 0) break; 1634e17fbf5SBarry 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); 1644e17fbf5SBarry Smith ts->time_step = .5*ts->time_step; 16531748224SBarry Smith ts->time_steps_since_decrease = 0; 1664e17fbf5SBarry Smith } 167f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 168f1b97656SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 169f1b97656SJed 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); 170f1b97656SJed Brown PetscFunctionReturn(0); 171f1b97656SJed Brown } 172eb284becSJed Brown if (th->endpoint) { 173eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 174eb284becSJed Brown } else { 175b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 176b296d7d5SJed Brown ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 1772b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 178eb284becSJed Brown } 1792b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 180cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 181316643e7SJed Brown ts->steps++; 182316643e7SJed Brown PetscFunctionReturn(0); 183316643e7SJed Brown } 184316643e7SJed Brown 185cd652676SJed Brown #undef __FUNCT__ 186cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 187cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 188cd652676SJed Brown { 189cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1905a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 191cd652676SJed Brown PetscErrorCode ierr; 192cd652676SJed Brown 193cd652676SJed Brown PetscFunctionBegin; 194a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 1955a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 1965a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 197cd652676SJed Brown PetscFunctionReturn(0); 198cd652676SJed Brown } 199cd652676SJed Brown 200316643e7SJed Brown /*------------------------------------------------------------*/ 201316643e7SJed Brown #undef __FUNCT__ 202277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 203277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 204316643e7SJed Brown { 205316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 206316643e7SJed Brown PetscErrorCode ierr; 207316643e7SJed Brown 208316643e7SJed Brown PetscFunctionBegin; 2096bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 2106bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 211eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 212277b19d0SLisandro Dalcin PetscFunctionReturn(0); 213277b19d0SLisandro Dalcin } 214277b19d0SLisandro Dalcin 215277b19d0SLisandro Dalcin #undef __FUNCT__ 216277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 217277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 218277b19d0SLisandro Dalcin { 219277b19d0SLisandro Dalcin PetscErrorCode ierr; 220277b19d0SLisandro Dalcin 221277b19d0SLisandro Dalcin PetscFunctionBegin; 222277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 223277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 224335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 225335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 22626f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 227eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 228316643e7SJed Brown PetscFunctionReturn(0); 229316643e7SJed Brown } 230316643e7SJed Brown 231316643e7SJed Brown /* 232316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 2332b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 234316643e7SJed Brown */ 235316643e7SJed Brown #undef __FUNCT__ 2360f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 2370f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 238316643e7SJed Brown { 239316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 240316643e7SJed Brown PetscErrorCode ierr; 2417445fe48SJed Brown Vec X0,Xdot; 2427445fe48SJed Brown DM dm,dmsave; 243b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 244316643e7SJed Brown 245316643e7SJed Brown PetscFunctionBegin; 2467445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2475a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 2487445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 249b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 2507445fe48SJed Brown 2517445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 2527445fe48SJed Brown dmsave = ts->dm; 2537445fe48SJed Brown ts->dm = dm; 2547445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 2557445fe48SJed Brown ts->dm = dmsave; 2560d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 257316643e7SJed Brown PetscFunctionReturn(0); 258316643e7SJed Brown } 259316643e7SJed Brown 260316643e7SJed Brown #undef __FUNCT__ 2610f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 2620f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 263316643e7SJed Brown { 264316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 265316643e7SJed Brown PetscErrorCode ierr; 2667445fe48SJed Brown Vec Xdot; 2677445fe48SJed Brown DM dm,dmsave; 268b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 269316643e7SJed Brown 270316643e7SJed Brown PetscFunctionBegin; 2717445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2727445fe48SJed Brown 2730f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 2747445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 2757445fe48SJed Brown 2767445fe48SJed Brown dmsave = ts->dm; 2777445fe48SJed Brown ts->dm = dm; 278b296d7d5SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 2797445fe48SJed Brown ts->dm = dmsave; 2800d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 281316643e7SJed Brown PetscFunctionReturn(0); 282316643e7SJed Brown } 283316643e7SJed Brown 284316643e7SJed Brown #undef __FUNCT__ 285316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 286316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 287316643e7SJed Brown { 288316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 289316643e7SJed Brown PetscErrorCode ierr; 2907445fe48SJed Brown SNES snes; 2917445fe48SJed Brown DM dm; 292316643e7SJed Brown 293316643e7SJed Brown PetscFunctionBegin; 294316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 295316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 2967445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2977445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2987445fe48SJed Brown if (dm) { 2997445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 300*258e1594SPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 3017445fe48SJed Brown } 302316643e7SJed Brown PetscFunctionReturn(0); 303316643e7SJed Brown } 304316643e7SJed Brown /*------------------------------------------------------------*/ 305316643e7SJed Brown 306316643e7SJed Brown #undef __FUNCT__ 307316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 308316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 309316643e7SJed Brown { 310316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 311316643e7SJed Brown PetscErrorCode ierr; 312316643e7SJed Brown 313316643e7SJed Brown PetscFunctionBegin; 314d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 315316643e7SJed Brown { 316316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 317acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 318eb284becSJed 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); 319d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 320316643e7SJed Brown } 321316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 322316643e7SJed Brown PetscFunctionReturn(0); 323316643e7SJed Brown } 324316643e7SJed Brown 325316643e7SJed Brown #undef __FUNCT__ 326316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 327316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 328316643e7SJed Brown { 329316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 330ace3abfcSBarry Smith PetscBool iascii; 331316643e7SJed Brown PetscErrorCode ierr; 332316643e7SJed Brown 333316643e7SJed Brown PetscFunctionBegin; 334251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 335316643e7SJed Brown if (iascii) { 336316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 337ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 338316643e7SJed Brown } 339d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 340316643e7SJed Brown PetscFunctionReturn(0); 341316643e7SJed Brown } 342316643e7SJed Brown 3430de4c49aSJed Brown EXTERN_C_BEGIN 3440de4c49aSJed Brown #undef __FUNCT__ 3450de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 3467087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 3470de4c49aSJed Brown { 3480de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3490de4c49aSJed Brown 3500de4c49aSJed Brown PetscFunctionBegin; 3510de4c49aSJed Brown *theta = th->Theta; 3520de4c49aSJed Brown PetscFunctionReturn(0); 3530de4c49aSJed Brown } 3540de4c49aSJed Brown 3550de4c49aSJed Brown #undef __FUNCT__ 3560de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 3577087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 3580de4c49aSJed Brown { 3590de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3600de4c49aSJed Brown 3610de4c49aSJed Brown PetscFunctionBegin; 362e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 3630de4c49aSJed Brown th->Theta = theta; 3640de4c49aSJed Brown PetscFunctionReturn(0); 3650de4c49aSJed Brown } 366eb284becSJed Brown 367eb284becSJed Brown #undef __FUNCT__ 36878e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 36926f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 37026f2ff8fSLisandro Dalcin { 37126f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 37226f2ff8fSLisandro Dalcin 37326f2ff8fSLisandro Dalcin PetscFunctionBegin; 37426f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 37526f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 37626f2ff8fSLisandro Dalcin } 37726f2ff8fSLisandro Dalcin 37826f2ff8fSLisandro Dalcin #undef __FUNCT__ 37926f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 380eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 381eb284becSJed Brown { 382eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 383eb284becSJed Brown 384eb284becSJed Brown PetscFunctionBegin; 385eb284becSJed Brown th->endpoint = flg; 386eb284becSJed Brown PetscFunctionReturn(0); 387eb284becSJed Brown } 3880de4c49aSJed Brown EXTERN_C_END 3890de4c49aSJed Brown 390f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 391f9c1d6abSBarry Smith #undef __FUNCT__ 392f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 393f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 394f9c1d6abSBarry Smith { 395f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 396f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 397f9c1d6abSBarry Smith 398f9c1d6abSBarry Smith PetscFunctionBegin; 399f9c1d6abSBarry Smith f = (1.0 + (1.0 - th->Theta)*z)/(1.0 - th->Theta*z); 400f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 401f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 402f9c1d6abSBarry Smith PetscFunctionReturn(0); 403f9c1d6abSBarry Smith } 404f9c1d6abSBarry Smith #endif 405f9c1d6abSBarry Smith 406f9c1d6abSBarry Smith 407316643e7SJed Brown /* ------------------------------------------------------------ */ 408316643e7SJed Brown /*MC 40996f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 410316643e7SJed Brown 411316643e7SJed Brown Level: beginner 412316643e7SJed Brown 4134eb428fdSBarry Smith Options Database: 4144eb428fdSBarry Smith -ts_theta_theta <Theta> - Location of stage (0<Theta<=1); Theta = 1.0 ( 4154eb428fdSBarry Smith -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 4164eb428fdSBarry Smith -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method 4174eb428fdSBarry Smith 418eb284becSJed Brown Notes: 4194eb428fdSBarry Smith $ -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER) 4204eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 4214eb428fdSBarry Smith 4224eb428fdSBarry Smith 4234eb428fdSBarry Smith 424eb284becSJed Brown This method can be applied to DAE. 425eb284becSJed Brown 426eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 427eb284becSJed Brown 428eb284becSJed Brown .vb 429eb284becSJed Brown Theta | Theta 430eb284becSJed Brown ------------- 431eb284becSJed Brown | 1 432eb284becSJed Brown .ve 433eb284becSJed Brown 434eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 435eb284becSJed Brown 436eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 437eb284becSJed Brown 438eb284becSJed Brown .vb 439eb284becSJed Brown 0 | 0 0 440eb284becSJed Brown 1 | 1-Theta Theta 441eb284becSJed Brown ------------------- 442eb284becSJed Brown | 1-Theta Theta 443eb284becSJed Brown .ve 444eb284becSJed Brown 445eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 446eb284becSJed Brown 447eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 448eb284becSJed Brown 449eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 450eb284becSJed Brown 4514eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 452eb284becSJed Brown 453eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 454316643e7SJed Brown 455316643e7SJed Brown M*/ 456316643e7SJed Brown EXTERN_C_BEGIN 457316643e7SJed Brown #undef __FUNCT__ 458316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 4597087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 460316643e7SJed Brown { 461316643e7SJed Brown TS_Theta *th; 462316643e7SJed Brown PetscErrorCode ierr; 463316643e7SJed Brown 464316643e7SJed Brown PetscFunctionBegin; 465277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 466316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 467316643e7SJed Brown ts->ops->view = TSView_Theta; 468316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 469316643e7SJed Brown ts->ops->step = TSStep_Theta; 470cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 471316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 4720f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 4730f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 474f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 475f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 476f9c1d6abSBarry Smith #endif 477316643e7SJed Brown 478316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 479316643e7SJed Brown ts->data = (void*)th; 480316643e7SJed Brown 4816f700aefSJed Brown th->extrapolate = PETSC_FALSE; 482316643e7SJed Brown th->Theta = 0.5; 483316643e7SJed Brown 4840de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 4850de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 48626f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 487eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 488316643e7SJed Brown PetscFunctionReturn(0); 489316643e7SJed Brown } 490316643e7SJed Brown EXTERN_C_END 4910de4c49aSJed Brown 4920de4c49aSJed Brown #undef __FUNCT__ 4930de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 4940de4c49aSJed Brown /*@ 4950de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 4960de4c49aSJed Brown 4970de4c49aSJed Brown Not Collective 4980de4c49aSJed Brown 4990de4c49aSJed Brown Input Parameter: 5000de4c49aSJed Brown . ts - timestepping context 5010de4c49aSJed Brown 5020de4c49aSJed Brown Output Parameter: 5030de4c49aSJed Brown . theta - stage abscissa 5040de4c49aSJed Brown 5050de4c49aSJed Brown Note: 5060de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 5070de4c49aSJed Brown 5080de4c49aSJed Brown Level: Advanced 5090de4c49aSJed Brown 5100de4c49aSJed Brown .seealso: TSThetaSetTheta() 5110de4c49aSJed Brown @*/ 5127087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 5130de4c49aSJed Brown { 5144ac538c5SBarry Smith PetscErrorCode ierr; 5150de4c49aSJed Brown 5160de4c49aSJed Brown PetscFunctionBegin; 517afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5180de4c49aSJed Brown PetscValidPointer(theta,2); 5194ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 5200de4c49aSJed Brown PetscFunctionReturn(0); 5210de4c49aSJed Brown } 5220de4c49aSJed Brown 5230de4c49aSJed Brown #undef __FUNCT__ 5240de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 5250de4c49aSJed Brown /*@ 5260de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 5270de4c49aSJed Brown 5280de4c49aSJed Brown Not Collective 5290de4c49aSJed Brown 5300de4c49aSJed Brown Input Parameter: 5310de4c49aSJed Brown + ts - timestepping context 5320de4c49aSJed Brown - theta - stage abscissa 5330de4c49aSJed Brown 5340de4c49aSJed Brown Options Database: 5350de4c49aSJed Brown . -ts_theta_theta <theta> 5360de4c49aSJed Brown 5370de4c49aSJed Brown Level: Intermediate 5380de4c49aSJed Brown 5390de4c49aSJed Brown .seealso: TSThetaGetTheta() 5400de4c49aSJed Brown @*/ 5417087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 5420de4c49aSJed Brown { 5434ac538c5SBarry Smith PetscErrorCode ierr; 5440de4c49aSJed Brown 5450de4c49aSJed Brown PetscFunctionBegin; 546afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5474ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 5480de4c49aSJed Brown PetscFunctionReturn(0); 5490de4c49aSJed Brown } 550f33bbcb6SJed Brown 551eb284becSJed Brown #undef __FUNCT__ 55226f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 55326f2ff8fSLisandro Dalcin /*@ 55426f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 55526f2ff8fSLisandro Dalcin 55626f2ff8fSLisandro Dalcin Not Collective 55726f2ff8fSLisandro Dalcin 55826f2ff8fSLisandro Dalcin Input Parameter: 55926f2ff8fSLisandro Dalcin . ts - timestepping context 56026f2ff8fSLisandro Dalcin 56126f2ff8fSLisandro Dalcin Output Parameter: 56226f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 56326f2ff8fSLisandro Dalcin 56426f2ff8fSLisandro Dalcin Level: Advanced 56526f2ff8fSLisandro Dalcin 56626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 56726f2ff8fSLisandro Dalcin @*/ 56826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 56926f2ff8fSLisandro Dalcin { 57026f2ff8fSLisandro Dalcin PetscErrorCode ierr; 57126f2ff8fSLisandro Dalcin 57226f2ff8fSLisandro Dalcin PetscFunctionBegin; 57326f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 57426f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 57526f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 57626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 57726f2ff8fSLisandro Dalcin } 57826f2ff8fSLisandro Dalcin 57926f2ff8fSLisandro Dalcin #undef __FUNCT__ 580eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 581eb284becSJed Brown /*@ 582eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 583eb284becSJed Brown 584eb284becSJed Brown Not Collective 585eb284becSJed Brown 586eb284becSJed Brown Input Parameter: 587eb284becSJed Brown + ts - timestepping context 588eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 589eb284becSJed Brown 590eb284becSJed Brown Options Database: 591eb284becSJed Brown . -ts_theta_endpoint <flg> 592eb284becSJed Brown 593eb284becSJed Brown Level: Intermediate 594eb284becSJed Brown 595eb284becSJed Brown .seealso: TSTHETA, TSCN 596eb284becSJed Brown @*/ 597eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 598eb284becSJed Brown { 599eb284becSJed Brown PetscErrorCode ierr; 600eb284becSJed Brown 601eb284becSJed Brown PetscFunctionBegin; 602eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 603eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 604eb284becSJed Brown PetscFunctionReturn(0); 605eb284becSJed Brown } 606eb284becSJed Brown 607f33bbcb6SJed Brown /* 608f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 609f33bbcb6SJed Brown * The creation functions for these specializations are below. 610f33bbcb6SJed Brown */ 611f33bbcb6SJed Brown 612f33bbcb6SJed Brown #undef __FUNCT__ 613f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 614f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 615f33bbcb6SJed Brown { 616d52bd9f3SBarry Smith PetscErrorCode ierr; 617d52bd9f3SBarry Smith 618f33bbcb6SJed Brown PetscFunctionBegin; 619d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 620f33bbcb6SJed Brown PetscFunctionReturn(0); 621f33bbcb6SJed Brown } 622f33bbcb6SJed Brown 623f33bbcb6SJed Brown /*MC 624f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 625f33bbcb6SJed Brown 626f33bbcb6SJed Brown Level: beginner 627f33bbcb6SJed Brown 6284eb428fdSBarry Smith Notes: 6294eb428fdSBarry Smith TSCN is equivalent to TSTHETA with Theta=1.0 6304eb428fdSBarry Smith 6314eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 1. 6324eb428fdSBarry Smith 633f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 634f33bbcb6SJed Brown 635f33bbcb6SJed Brown M*/ 636f33bbcb6SJed Brown EXTERN_C_BEGIN 637f33bbcb6SJed Brown #undef __FUNCT__ 638f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 639f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 640f33bbcb6SJed Brown { 641f33bbcb6SJed Brown PetscErrorCode ierr; 642f33bbcb6SJed Brown 643f33bbcb6SJed Brown PetscFunctionBegin; 644f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 645f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 646f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 647f33bbcb6SJed Brown PetscFunctionReturn(0); 648f33bbcb6SJed Brown } 649f33bbcb6SJed Brown EXTERN_C_END 650f33bbcb6SJed Brown 651f33bbcb6SJed Brown #undef __FUNCT__ 652f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 653f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 654f33bbcb6SJed Brown { 655d52bd9f3SBarry Smith PetscErrorCode ierr; 656d52bd9f3SBarry Smith 657f33bbcb6SJed Brown PetscFunctionBegin; 658d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 659f33bbcb6SJed Brown PetscFunctionReturn(0); 660f33bbcb6SJed Brown } 661f33bbcb6SJed Brown 662f33bbcb6SJed Brown /*MC 663f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 664f33bbcb6SJed Brown 665f33bbcb6SJed Brown Level: beginner 666f33bbcb6SJed Brown 667f33bbcb6SJed Brown Notes: 6687cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 6697cf5af47SJed Brown 6707cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 671f33bbcb6SJed Brown 672f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 673f33bbcb6SJed Brown 674f33bbcb6SJed Brown M*/ 675f33bbcb6SJed Brown EXTERN_C_BEGIN 676f33bbcb6SJed Brown #undef __FUNCT__ 677f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 678f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 679f33bbcb6SJed Brown { 680f33bbcb6SJed Brown PetscErrorCode ierr; 681f33bbcb6SJed Brown 682f33bbcb6SJed Brown PetscFunctionBegin; 683f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 684f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 685eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 686f33bbcb6SJed Brown ts->ops->view = TSView_CN; 687f33bbcb6SJed Brown PetscFunctionReturn(0); 688f33bbcb6SJed Brown } 689f33bbcb6SJed Brown EXTERN_C_END 690