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 */ 103b1890cdSShri Abhyankar Vec X0; /* work vector to store X0 */ 11eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 12ace3abfcSBarry Smith PetscBool extrapolate; 13eb284becSJed Brown PetscBool endpoint; 14316643e7SJed Brown PetscReal Theta; 15316643e7SJed Brown PetscReal stage_time; 163b1890cdSShri Abhyankar TSStepStatus status; 173b1890cdSShri Abhyankar char *name; 183b1890cdSShri Abhyankar PetscInt order; 193b1890cdSShri Abhyankar PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 203b1890cdSShri Abhyankar PetscBool adapt; /* use time-step adaptivity ? */ 21316643e7SJed Brown } TS_Theta; 22316643e7SJed Brown 23316643e7SJed Brown #undef __FUNCT__ 247445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 257445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 267445fe48SJed Brown { 277445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 287445fe48SJed Brown PetscErrorCode ierr; 297445fe48SJed Brown 307445fe48SJed Brown PetscFunctionBegin; 317445fe48SJed Brown if (X0) { 327445fe48SJed Brown if (dm && dm != ts->dm) { 330d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 347445fe48SJed Brown } else *X0 = ts->vec_sol; 357445fe48SJed Brown } 367445fe48SJed Brown if (Xdot) { 377445fe48SJed Brown if (dm && dm != ts->dm) { 380d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 397445fe48SJed Brown } else *Xdot = th->Xdot; 407445fe48SJed Brown } 417445fe48SJed Brown PetscFunctionReturn(0); 427445fe48SJed Brown } 437445fe48SJed Brown 440d0b770aSPeter Brune 450d0b770aSPeter Brune #undef __FUNCT__ 460d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 470d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 480d0b770aSPeter Brune { 490d0b770aSPeter Brune PetscErrorCode ierr; 500d0b770aSPeter Brune 510d0b770aSPeter Brune PetscFunctionBegin; 520d0b770aSPeter Brune if (X0) { 530d0b770aSPeter Brune if (dm && dm != ts->dm) { 540d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 550d0b770aSPeter Brune } 560d0b770aSPeter Brune } 570d0b770aSPeter Brune if (Xdot) { 580d0b770aSPeter Brune if (dm && dm != ts->dm) { 590d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 600d0b770aSPeter Brune } 610d0b770aSPeter Brune } 620d0b770aSPeter Brune PetscFunctionReturn(0); 630d0b770aSPeter Brune } 640d0b770aSPeter Brune 657445fe48SJed Brown #undef __FUNCT__ 667445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 677445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 687445fe48SJed Brown { 697445fe48SJed Brown 707445fe48SJed Brown PetscFunctionBegin; 717445fe48SJed Brown PetscFunctionReturn(0); 727445fe48SJed Brown } 737445fe48SJed Brown 747445fe48SJed Brown #undef __FUNCT__ 757445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 767445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 777445fe48SJed Brown { 787445fe48SJed Brown TS ts = (TS)ctx; 797445fe48SJed Brown PetscErrorCode ierr; 807445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 817445fe48SJed Brown 827445fe48SJed Brown PetscFunctionBegin; 837445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 847445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 857445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 867445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 877445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 887445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 890d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 900d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 917445fe48SJed Brown PetscFunctionReturn(0); 927445fe48SJed Brown } 937445fe48SJed Brown 947445fe48SJed Brown #undef __FUNCT__ 95258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta" 96258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 97258e1594SPeter Brune { 98258e1594SPeter Brune 99258e1594SPeter Brune PetscFunctionBegin; 100258e1594SPeter Brune PetscFunctionReturn(0); 101258e1594SPeter Brune } 102258e1594SPeter Brune 103258e1594SPeter Brune #undef __FUNCT__ 104258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta" 105258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 106258e1594SPeter Brune { 107258e1594SPeter Brune TS ts = (TS)ctx; 108258e1594SPeter Brune PetscErrorCode ierr; 109258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 110258e1594SPeter Brune 111258e1594SPeter Brune PetscFunctionBegin; 112258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 113258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 114258e1594SPeter Brune 115258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 116258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 117258e1594SPeter Brune 118258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 119258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 120258e1594SPeter Brune 121258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 122258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 123258e1594SPeter Brune PetscFunctionReturn(0); 124258e1594SPeter Brune } 125258e1594SPeter Brune 1263b1890cdSShri Abhyankar #undef __FUNCT__ 1273b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta" 1283b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done) 1293b1890cdSShri Abhyankar { 1303b1890cdSShri Abhyankar PetscErrorCode ierr; 1313b1890cdSShri Abhyankar TS_Theta *th = (TS_Theta*)ts->data; 1323b1890cdSShri Abhyankar 1333b1890cdSShri Abhyankar PetscFunctionBegin; 134*ce94432eSBarry Smith if (order == 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"No time-step adaptivity implemented for 1st order theta method; Run with -ts_adapt_type none"); 1353b1890cdSShri Abhyankar if (order == th->order) { 1363b1890cdSShri Abhyankar if (th->endpoint) { 1373b1890cdSShri Abhyankar ierr = VecCopy(th->X,U);CHKERRQ(ierr); 1383b1890cdSShri Abhyankar } else { 1393b1890cdSShri Abhyankar PetscReal shift = 1./(th->Theta*ts->time_step); 1403b1890cdSShri Abhyankar ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr); 1413b1890cdSShri Abhyankar ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr); 1423b1890cdSShri Abhyankar } 1433b1890cdSShri Abhyankar } else if (order == th->order-1 && order) { 1443b1890cdSShri Abhyankar ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr); 1453b1890cdSShri Abhyankar } 1463b1890cdSShri Abhyankar PetscFunctionReturn(0); 1473b1890cdSShri Abhyankar } 148258e1594SPeter Brune 149258e1594SPeter Brune #undef __FUNCT__ 150316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 151193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 152316643e7SJed Brown { 153316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1543b1890cdSShri Abhyankar PetscInt its,lits,reject,next_scheme; 1553b1890cdSShri Abhyankar PetscReal next_time_step; 156f1b97656SJed Brown SNESConvergedReason snesreason; 1572b5a38e1SLisandro Dalcin PetscErrorCode ierr; 1583b1890cdSShri Abhyankar TSAdapt adapt; 1593b1890cdSShri Abhyankar PetscBool accept; 160316643e7SJed Brown 161316643e7SJed Brown PetscFunctionBegin; 1623b1890cdSShri Abhyankar th->status = TS_STEP_INCOMPLETE; 1633b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 1643b1890cdSShri Abhyankar for (reject=0; reject<ts->max_reject && !ts->reason && th->status != TS_STEP_COMPLETE; reject++,ts->reject++) { 165b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 166cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 167eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 168b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 169b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 170316643e7SJed Brown 171eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 172eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 173eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 174eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 175eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 176eb284becSJed Brown } 177ace68cafSJed Brown if (th->extrapolate) { 178b296d7d5SJed Brown ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 179ace68cafSJed Brown } else { 1802b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 181ace68cafSJed Brown } 182eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 183316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 184316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 185f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 1865ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 1873b1890cdSShri Abhyankar ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr); 1883b1890cdSShri Abhyankar ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 1893b1890cdSShri Abhyankar if (!accept) continue; 1900298fd71SBarry Smith ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr); 1913b1890cdSShri Abhyankar /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 1923b1890cdSShri Abhyankar ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr); 1933b1890cdSShri Abhyankar ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 1940298fd71SBarry Smith ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr); 1953b1890cdSShri Abhyankar ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 1963b1890cdSShri Abhyankar 1973b1890cdSShri Abhyankar if (accept) { 1983b1890cdSShri Abhyankar /* ignore next_scheme for now */ 1992b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 200cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 201316643e7SJed Brown ts->steps++; 2023b1890cdSShri Abhyankar th->status = TS_STEP_COMPLETE; 2033b1890cdSShri Abhyankar } else { /* Roll back the current step */ 2043b1890cdSShri Abhyankar ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 2053b1890cdSShri Abhyankar ts->time_step = next_time_step; 2063b1890cdSShri Abhyankar th->status = TS_STEP_INCOMPLETE; 2073b1890cdSShri Abhyankar } 2083b1890cdSShri Abhyankar } 209316643e7SJed Brown PetscFunctionReturn(0); 210316643e7SJed Brown } 211316643e7SJed Brown 212cd652676SJed Brown #undef __FUNCT__ 213cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 214cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 215cd652676SJed Brown { 216cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2175a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 218cd652676SJed Brown PetscErrorCode ierr; 219cd652676SJed Brown 220cd652676SJed Brown PetscFunctionBegin; 221a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 2225a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 2235a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 224cd652676SJed Brown PetscFunctionReturn(0); 225cd652676SJed Brown } 226cd652676SJed Brown 227316643e7SJed Brown /*------------------------------------------------------------*/ 228316643e7SJed Brown #undef __FUNCT__ 229277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 230277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 231316643e7SJed Brown { 232316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 233316643e7SJed Brown PetscErrorCode ierr; 234316643e7SJed Brown 235316643e7SJed Brown PetscFunctionBegin; 2366bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 2376bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 2383b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 239eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 240277b19d0SLisandro Dalcin PetscFunctionReturn(0); 241277b19d0SLisandro Dalcin } 242277b19d0SLisandro Dalcin 243277b19d0SLisandro Dalcin #undef __FUNCT__ 244277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 245277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 246277b19d0SLisandro Dalcin { 247277b19d0SLisandro Dalcin PetscErrorCode ierr; 248277b19d0SLisandro Dalcin 249277b19d0SLisandro Dalcin PetscFunctionBegin; 250277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 251277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 2520298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",NULL);CHKERRQ(ierr); 2530298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",NULL);CHKERRQ(ierr); 2540298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",NULL);CHKERRQ(ierr); 2550298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",NULL);CHKERRQ(ierr); 256316643e7SJed Brown PetscFunctionReturn(0); 257316643e7SJed Brown } 258316643e7SJed Brown 259316643e7SJed Brown /* 260316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 2612b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 262316643e7SJed Brown */ 263316643e7SJed Brown #undef __FUNCT__ 2640f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 2650f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 266316643e7SJed Brown { 267316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 268316643e7SJed Brown PetscErrorCode ierr; 2697445fe48SJed Brown Vec X0,Xdot; 2707445fe48SJed Brown DM dm,dmsave; 271b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 272316643e7SJed Brown 273316643e7SJed Brown PetscFunctionBegin; 2747445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2755a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 2767445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 277b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 2787445fe48SJed Brown 2797445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 2807445fe48SJed Brown dmsave = ts->dm; 2817445fe48SJed Brown ts->dm = dm; 2827445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 2837445fe48SJed Brown ts->dm = dmsave; 2840d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 285316643e7SJed Brown PetscFunctionReturn(0); 286316643e7SJed Brown } 287316643e7SJed Brown 288316643e7SJed Brown #undef __FUNCT__ 2890f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 2900f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 291316643e7SJed Brown { 292316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 293316643e7SJed Brown PetscErrorCode ierr; 2947445fe48SJed Brown Vec Xdot; 2957445fe48SJed Brown DM dm,dmsave; 296b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 297316643e7SJed Brown 298316643e7SJed Brown PetscFunctionBegin; 2997445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 3007445fe48SJed Brown 3010f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 3020298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 3037445fe48SJed Brown 3047445fe48SJed Brown dmsave = ts->dm; 3057445fe48SJed Brown ts->dm = dm; 306b296d7d5SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 3077445fe48SJed Brown ts->dm = dmsave; 3080298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 309316643e7SJed Brown PetscFunctionReturn(0); 310316643e7SJed Brown } 311316643e7SJed Brown 312316643e7SJed Brown #undef __FUNCT__ 313316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 314316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 315316643e7SJed Brown { 316316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 317316643e7SJed Brown PetscErrorCode ierr; 3187445fe48SJed Brown SNES snes; 3197445fe48SJed Brown DM dm; 320316643e7SJed Brown 321316643e7SJed Brown PetscFunctionBegin; 322316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 323316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 3243b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 3257445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 3267445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 3277445fe48SJed Brown if (dm) { 3287445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 329258e1594SPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 3307445fe48SJed Brown } 3313b1890cdSShri Abhyankar if (th->Theta == 0.5 && th->endpoint) th->order = 2; 3323b1890cdSShri Abhyankar else th->order = 1; 3333b1890cdSShri Abhyankar 3343b1890cdSShri Abhyankar if (!th->adapt) { 3353b1890cdSShri Abhyankar TSAdapt adapt; 3363b1890cdSShri Abhyankar ierr = TSAdaptDestroy(&ts->adapt);CHKERRQ(ierr); 3373b1890cdSShri Abhyankar ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr); 3383b1890cdSShri Abhyankar ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr); 3393b1890cdSShri Abhyankar } 340316643e7SJed Brown PetscFunctionReturn(0); 341316643e7SJed Brown } 342316643e7SJed Brown /*------------------------------------------------------------*/ 343316643e7SJed Brown 344316643e7SJed Brown #undef __FUNCT__ 345316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 346316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 347316643e7SJed Brown { 348316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 349316643e7SJed Brown PetscErrorCode ierr; 350316643e7SJed Brown 351316643e7SJed Brown PetscFunctionBegin; 352d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 353316643e7SJed Brown { 3540298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 3550298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 3560298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr); 3570298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 358d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 359316643e7SJed Brown } 360316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 361316643e7SJed Brown PetscFunctionReturn(0); 362316643e7SJed Brown } 363316643e7SJed Brown 364316643e7SJed Brown #undef __FUNCT__ 365316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 366316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 367316643e7SJed Brown { 368316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 369ace3abfcSBarry Smith PetscBool iascii; 370316643e7SJed Brown PetscErrorCode ierr; 371316643e7SJed Brown 372316643e7SJed Brown PetscFunctionBegin; 373251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 374316643e7SJed Brown if (iascii) { 375316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 376ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 377316643e7SJed Brown } 378d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 379316643e7SJed Brown PetscFunctionReturn(0); 380316643e7SJed Brown } 381316643e7SJed Brown 3820de4c49aSJed Brown EXTERN_C_BEGIN 3830de4c49aSJed Brown #undef __FUNCT__ 3840de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 3857087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 3860de4c49aSJed Brown { 3870de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3880de4c49aSJed Brown 3890de4c49aSJed Brown PetscFunctionBegin; 3900de4c49aSJed Brown *theta = th->Theta; 3910de4c49aSJed Brown PetscFunctionReturn(0); 3920de4c49aSJed Brown } 3930de4c49aSJed Brown 3940de4c49aSJed Brown #undef __FUNCT__ 3950de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 3967087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 3970de4c49aSJed Brown { 3980de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3990de4c49aSJed Brown 4000de4c49aSJed Brown PetscFunctionBegin; 401*ce94432eSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 4020de4c49aSJed Brown th->Theta = theta; 4030de4c49aSJed Brown PetscFunctionReturn(0); 4040de4c49aSJed Brown } 405eb284becSJed Brown 406eb284becSJed Brown #undef __FUNCT__ 40778e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 40826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 40926f2ff8fSLisandro Dalcin { 41026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 41126f2ff8fSLisandro Dalcin 41226f2ff8fSLisandro Dalcin PetscFunctionBegin; 41326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 41426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 41526f2ff8fSLisandro Dalcin } 41626f2ff8fSLisandro Dalcin 41726f2ff8fSLisandro Dalcin #undef __FUNCT__ 41826f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 419eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 420eb284becSJed Brown { 421eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 422eb284becSJed Brown 423eb284becSJed Brown PetscFunctionBegin; 424eb284becSJed Brown th->endpoint = flg; 425eb284becSJed Brown PetscFunctionReturn(0); 426eb284becSJed Brown } 4270de4c49aSJed Brown EXTERN_C_END 4280de4c49aSJed Brown 429f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 430f9c1d6abSBarry Smith #undef __FUNCT__ 431f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 432f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 433f9c1d6abSBarry Smith { 434f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 435f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 4363fd8ae06SJed Brown const PetscReal one = 1.0; 437f9c1d6abSBarry Smith 438f9c1d6abSBarry Smith PetscFunctionBegin; 4393fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 440f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 441f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 442f9c1d6abSBarry Smith PetscFunctionReturn(0); 443f9c1d6abSBarry Smith } 444f9c1d6abSBarry Smith #endif 445f9c1d6abSBarry Smith 446f9c1d6abSBarry Smith 447316643e7SJed Brown /* ------------------------------------------------------------ */ 448316643e7SJed Brown /*MC 44996f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 450316643e7SJed Brown 451316643e7SJed Brown Level: beginner 452316643e7SJed Brown 4534eb428fdSBarry Smith Options Database: 4540c3ba866SJed Brown -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 4554eb428fdSBarry Smith -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 4560c3ba866SJed Brown -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 4574eb428fdSBarry Smith 458eb284becSJed Brown Notes: 4590c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 4600c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 4614eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 4624eb428fdSBarry Smith 4634eb428fdSBarry Smith 4644eb428fdSBarry Smith 465eb284becSJed Brown This method can be applied to DAE. 466eb284becSJed Brown 467eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 468eb284becSJed Brown 469eb284becSJed Brown .vb 470eb284becSJed Brown Theta | Theta 471eb284becSJed Brown ------------- 472eb284becSJed Brown | 1 473eb284becSJed Brown .ve 474eb284becSJed Brown 475eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 476eb284becSJed Brown 477eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 478eb284becSJed Brown 479eb284becSJed Brown .vb 480eb284becSJed Brown 0 | 0 0 481eb284becSJed Brown 1 | 1-Theta Theta 482eb284becSJed Brown ------------------- 483eb284becSJed Brown | 1-Theta Theta 484eb284becSJed Brown .ve 485eb284becSJed Brown 486eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 487eb284becSJed Brown 488eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 489eb284becSJed Brown 490eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 491eb284becSJed Brown 4924eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 493eb284becSJed Brown 494eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 495316643e7SJed Brown 496316643e7SJed Brown M*/ 497316643e7SJed Brown EXTERN_C_BEGIN 498316643e7SJed Brown #undef __FUNCT__ 499316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 5007087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 501316643e7SJed Brown { 502316643e7SJed Brown TS_Theta *th; 503316643e7SJed Brown PetscErrorCode ierr; 504316643e7SJed Brown 505316643e7SJed Brown PetscFunctionBegin; 506277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 507316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 508316643e7SJed Brown ts->ops->view = TSView_Theta; 509316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 510316643e7SJed Brown ts->ops->step = TSStep_Theta; 511cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 5123b1890cdSShri Abhyankar ts->ops->evaluatestep = TSEvaluateStep_Theta; 513316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 5140f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 5150f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 516f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 517f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 518f9c1d6abSBarry Smith #endif 519316643e7SJed Brown 520316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 521316643e7SJed Brown ts->data = (void*)th; 522316643e7SJed Brown 5236f700aefSJed Brown th->extrapolate = PETSC_FALSE; 524316643e7SJed Brown th->Theta = 0.5; 5253b1890cdSShri Abhyankar th->ccfl = 1.0; 5263b1890cdSShri Abhyankar th->adapt = PETSC_FALSE; 5270de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 5280de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 52926f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 530eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 531316643e7SJed Brown PetscFunctionReturn(0); 532316643e7SJed Brown } 533316643e7SJed Brown EXTERN_C_END 5340de4c49aSJed Brown 5350de4c49aSJed Brown #undef __FUNCT__ 5360de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 5370de4c49aSJed Brown /*@ 5380de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 5390de4c49aSJed Brown 5400de4c49aSJed Brown Not Collective 5410de4c49aSJed Brown 5420de4c49aSJed Brown Input Parameter: 5430de4c49aSJed Brown . ts - timestepping context 5440de4c49aSJed Brown 5450de4c49aSJed Brown Output Parameter: 5460de4c49aSJed Brown . theta - stage abscissa 5470de4c49aSJed Brown 5480de4c49aSJed Brown Note: 5490de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 5500de4c49aSJed Brown 5510de4c49aSJed Brown Level: Advanced 5520de4c49aSJed Brown 5530de4c49aSJed Brown .seealso: TSThetaSetTheta() 5540de4c49aSJed Brown @*/ 5557087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 5560de4c49aSJed Brown { 5574ac538c5SBarry Smith PetscErrorCode ierr; 5580de4c49aSJed Brown 5590de4c49aSJed Brown PetscFunctionBegin; 560afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5610de4c49aSJed Brown PetscValidPointer(theta,2); 5624ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 5630de4c49aSJed Brown PetscFunctionReturn(0); 5640de4c49aSJed Brown } 5650de4c49aSJed Brown 5660de4c49aSJed Brown #undef __FUNCT__ 5670de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 5680de4c49aSJed Brown /*@ 5690de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 5700de4c49aSJed Brown 5710de4c49aSJed Brown Not Collective 5720de4c49aSJed Brown 5730de4c49aSJed Brown Input Parameter: 5740de4c49aSJed Brown + ts - timestepping context 5750de4c49aSJed Brown - theta - stage abscissa 5760de4c49aSJed Brown 5770de4c49aSJed Brown Options Database: 5780de4c49aSJed Brown . -ts_theta_theta <theta> 5790de4c49aSJed Brown 5800de4c49aSJed Brown Level: Intermediate 5810de4c49aSJed Brown 5820de4c49aSJed Brown .seealso: TSThetaGetTheta() 5830de4c49aSJed Brown @*/ 5847087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 5850de4c49aSJed Brown { 5864ac538c5SBarry Smith PetscErrorCode ierr; 5870de4c49aSJed Brown 5880de4c49aSJed Brown PetscFunctionBegin; 589afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5904ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 5910de4c49aSJed Brown PetscFunctionReturn(0); 5920de4c49aSJed Brown } 593f33bbcb6SJed Brown 594eb284becSJed Brown #undef __FUNCT__ 59526f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 59626f2ff8fSLisandro Dalcin /*@ 59726f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 59826f2ff8fSLisandro Dalcin 59926f2ff8fSLisandro Dalcin Not Collective 60026f2ff8fSLisandro Dalcin 60126f2ff8fSLisandro Dalcin Input Parameter: 60226f2ff8fSLisandro Dalcin . ts - timestepping context 60326f2ff8fSLisandro Dalcin 60426f2ff8fSLisandro Dalcin Output Parameter: 60526f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 60626f2ff8fSLisandro Dalcin 60726f2ff8fSLisandro Dalcin Level: Advanced 60826f2ff8fSLisandro Dalcin 60926f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 61026f2ff8fSLisandro Dalcin @*/ 61126f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 61226f2ff8fSLisandro Dalcin { 61326f2ff8fSLisandro Dalcin PetscErrorCode ierr; 61426f2ff8fSLisandro Dalcin 61526f2ff8fSLisandro Dalcin PetscFunctionBegin; 61626f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 61726f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 61826f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 61926f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 62026f2ff8fSLisandro Dalcin } 62126f2ff8fSLisandro Dalcin 62226f2ff8fSLisandro Dalcin #undef __FUNCT__ 623eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 624eb284becSJed Brown /*@ 625eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 626eb284becSJed Brown 627eb284becSJed Brown Not Collective 628eb284becSJed Brown 629eb284becSJed Brown Input Parameter: 630eb284becSJed Brown + ts - timestepping context 631eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 632eb284becSJed Brown 633eb284becSJed Brown Options Database: 634eb284becSJed Brown . -ts_theta_endpoint <flg> 635eb284becSJed Brown 636eb284becSJed Brown Level: Intermediate 637eb284becSJed Brown 638eb284becSJed Brown .seealso: TSTHETA, TSCN 639eb284becSJed Brown @*/ 640eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 641eb284becSJed Brown { 642eb284becSJed Brown PetscErrorCode ierr; 643eb284becSJed Brown 644eb284becSJed Brown PetscFunctionBegin; 645eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 646eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 647eb284becSJed Brown PetscFunctionReturn(0); 648eb284becSJed Brown } 649eb284becSJed Brown 650f33bbcb6SJed Brown /* 651f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 652f33bbcb6SJed Brown * The creation functions for these specializations are below. 653f33bbcb6SJed Brown */ 654f33bbcb6SJed Brown 655f33bbcb6SJed Brown #undef __FUNCT__ 656f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 657f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 658f33bbcb6SJed Brown { 659d52bd9f3SBarry Smith PetscErrorCode ierr; 660d52bd9f3SBarry Smith 661f33bbcb6SJed Brown PetscFunctionBegin; 662d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 663f33bbcb6SJed Brown PetscFunctionReturn(0); 664f33bbcb6SJed Brown } 665f33bbcb6SJed Brown 666f33bbcb6SJed Brown /*MC 667f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 668f33bbcb6SJed Brown 669f33bbcb6SJed Brown Level: beginner 670f33bbcb6SJed Brown 6714eb428fdSBarry Smith Notes: 672c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 6734eb428fdSBarry Smith 6744eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 1. 6754eb428fdSBarry Smith 676f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 677f33bbcb6SJed Brown 678f33bbcb6SJed Brown M*/ 679f33bbcb6SJed Brown EXTERN_C_BEGIN 680f33bbcb6SJed Brown #undef __FUNCT__ 681f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 682f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 683f33bbcb6SJed Brown { 684f33bbcb6SJed Brown PetscErrorCode ierr; 685f33bbcb6SJed Brown 686f33bbcb6SJed Brown PetscFunctionBegin; 687f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 688f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 689f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 690f33bbcb6SJed Brown PetscFunctionReturn(0); 691f33bbcb6SJed Brown } 692f33bbcb6SJed Brown EXTERN_C_END 693f33bbcb6SJed Brown 694f33bbcb6SJed Brown #undef __FUNCT__ 695f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 696f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 697f33bbcb6SJed Brown { 698d52bd9f3SBarry Smith PetscErrorCode ierr; 699d52bd9f3SBarry Smith 700f33bbcb6SJed Brown PetscFunctionBegin; 701d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 702f33bbcb6SJed Brown PetscFunctionReturn(0); 703f33bbcb6SJed Brown } 704f33bbcb6SJed Brown 705f33bbcb6SJed Brown /*MC 706f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 707f33bbcb6SJed Brown 708f33bbcb6SJed Brown Level: beginner 709f33bbcb6SJed Brown 710f33bbcb6SJed Brown Notes: 7117cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 7127cf5af47SJed Brown 7137cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 714f33bbcb6SJed Brown 715f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 716f33bbcb6SJed Brown 717f33bbcb6SJed Brown M*/ 718f33bbcb6SJed Brown EXTERN_C_BEGIN 719f33bbcb6SJed Brown #undef __FUNCT__ 720f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 721f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 722f33bbcb6SJed Brown { 723f33bbcb6SJed Brown PetscErrorCode ierr; 724f33bbcb6SJed Brown 725f33bbcb6SJed Brown PetscFunctionBegin; 726f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 727f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 728eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 729f33bbcb6SJed Brown ts->ops->view = TSView_CN; 730f33bbcb6SJed Brown PetscFunctionReturn(0); 731f33bbcb6SJed Brown } 732f33bbcb6SJed Brown EXTERN_C_END 733