1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5a055c8caSBarry Smith #include <petscsnes.h> 61e25c274SJed Brown #include <petscdm.h> 73c54f92cSHong Zhang #include <petscmat.h> 8316643e7SJed Brown 9316643e7SJed Brown typedef struct { 101566a47fSLisandro Dalcin PetscReal stage_time; 111566a47fSLisandro Dalcin Vec X0,X,Xdot; /* Storage for stages and time derivative */ 121566a47fSLisandro Dalcin Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 131566a47fSLisandro Dalcin 141566a47fSLisandro Dalcin PetscReal Theta; 151566a47fSLisandro Dalcin PetscInt order; 161566a47fSLisandro Dalcin PetscBool endpoint; 171566a47fSLisandro Dalcin PetscBool extrapolate; 181566a47fSLisandro Dalcin 19b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 20b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 21b5e0532dSHong Zhang Vec *VecsSensiTemp; /* Vector to be timed with Jacobian transpose */ 221566a47fSLisandro Dalcin Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 23b5e0532dSHong Zhang PetscReal ptime; 24b1cb36f3SHong Zhang PetscReal time_step; 251566a47fSLisandro Dalcin 261566a47fSLisandro Dalcin Vec vec_sol_prev; 271566a47fSLisandro Dalcin Vec vec_lte_work; 281566a47fSLisandro Dalcin 291566a47fSLisandro Dalcin TSStepStatus status; 30316643e7SJed Brown } TS_Theta; 31316643e7SJed Brown 327445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 337445fe48SJed Brown { 347445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 357445fe48SJed Brown PetscErrorCode ierr; 367445fe48SJed Brown 377445fe48SJed Brown PetscFunctionBegin; 387445fe48SJed Brown if (X0) { 397445fe48SJed Brown if (dm && dm != ts->dm) { 400d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 417445fe48SJed Brown } else *X0 = ts->vec_sol; 427445fe48SJed Brown } 437445fe48SJed Brown if (Xdot) { 447445fe48SJed Brown if (dm && dm != ts->dm) { 450d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 467445fe48SJed Brown } else *Xdot = th->Xdot; 477445fe48SJed Brown } 487445fe48SJed Brown PetscFunctionReturn(0); 497445fe48SJed Brown } 507445fe48SJed Brown 510d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 520d0b770aSPeter Brune { 530d0b770aSPeter Brune PetscErrorCode ierr; 540d0b770aSPeter Brune 550d0b770aSPeter Brune PetscFunctionBegin; 560d0b770aSPeter Brune if (X0) { 570d0b770aSPeter Brune if (dm && dm != ts->dm) { 580d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 590d0b770aSPeter Brune } 600d0b770aSPeter Brune } 610d0b770aSPeter Brune if (Xdot) { 620d0b770aSPeter Brune if (dm && dm != ts->dm) { 630d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 640d0b770aSPeter Brune } 650d0b770aSPeter Brune } 660d0b770aSPeter Brune PetscFunctionReturn(0); 670d0b770aSPeter Brune } 680d0b770aSPeter Brune 697445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 707445fe48SJed Brown { 717445fe48SJed Brown 727445fe48SJed Brown PetscFunctionBegin; 737445fe48SJed Brown PetscFunctionReturn(0); 747445fe48SJed Brown } 757445fe48SJed Brown 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 94258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 95258e1594SPeter Brune { 96258e1594SPeter Brune 97258e1594SPeter Brune PetscFunctionBegin; 98258e1594SPeter Brune PetscFunctionReturn(0); 99258e1594SPeter Brune } 100258e1594SPeter Brune 101258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 102258e1594SPeter Brune { 103258e1594SPeter Brune TS ts = (TS)ctx; 104258e1594SPeter Brune PetscErrorCode ierr; 105258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 106258e1594SPeter Brune 107258e1594SPeter Brune PetscFunctionBegin; 108258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 109258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 110258e1594SPeter Brune 111258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 112258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 113258e1594SPeter Brune 114258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 115258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 116258e1594SPeter Brune 117258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 118258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 119258e1594SPeter Brune PetscFunctionReturn(0); 120258e1594SPeter Brune } 121258e1594SPeter Brune 122b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 123b1cb36f3SHong Zhang { 124b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 125b1cb36f3SHong Zhang PetscErrorCode ierr; 126b1cb36f3SHong Zhang 127b1cb36f3SHong Zhang PetscFunctionBegin; 128b1cb36f3SHong Zhang /* backup cost integral */ 129b1cb36f3SHong Zhang ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 130b1cb36f3SHong Zhang if (th->endpoint) { 131b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1321566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 133b1cb36f3SHong Zhang } 134b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 135b1cb36f3SHong Zhang if (th->endpoint) { 136b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 137b1cb36f3SHong Zhang } else { 138b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 139b1cb36f3SHong Zhang } 140b1cb36f3SHong Zhang PetscFunctionReturn(0); 141b1cb36f3SHong Zhang } 142b1cb36f3SHong Zhang 143b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 144b1cb36f3SHong Zhang { 145b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 146b1cb36f3SHong Zhang PetscErrorCode ierr; 147b1cb36f3SHong Zhang 148b1cb36f3SHong Zhang PetscFunctionBegin; 149b1cb36f3SHong Zhang if (th->endpoint) { 150b1cb36f3SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 151b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 152b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 153b1cb36f3SHong Zhang if (th->Theta!=1) { 154b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1551566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr); 156b1cb36f3SHong Zhang } 157b1cb36f3SHong Zhang }else { 158b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 159b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 160b1cb36f3SHong Zhang } 16124655328SShri PetscFunctionReturn(0); 16224655328SShri } 16324655328SShri 1641566a47fSLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 1651566a47fSLisandro Dalcin { 1661566a47fSLisandro Dalcin PetscInt nits,lits; 1671566a47fSLisandro Dalcin PetscErrorCode ierr; 1681566a47fSLisandro Dalcin 1691566a47fSLisandro Dalcin PetscFunctionBegin; 1701566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1711566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1721566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1731566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1741566a47fSLisandro Dalcin PetscFunctionReturn(0); 1751566a47fSLisandro Dalcin } 1761566a47fSLisandro Dalcin 177193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 178316643e7SJed Brown { 179316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1801566a47fSLisandro Dalcin PetscInt rejections = 0; 1814957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1821566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 183051f2191SLisandro Dalcin PetscErrorCode ierr; 184316643e7SJed Brown 185316643e7SJed Brown PetscFunctionBegin; 1861566a47fSLisandro Dalcin if (!ts->steprollback) { 1872ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 1883b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 1891566a47fSLisandro Dalcin } 1901566a47fSLisandro Dalcin 1911566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 1921566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 1931566a47fSLisandro Dalcin 1941566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 1951566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 196316643e7SJed Brown 197be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 198fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 199be5899b3SLisandro Dalcin ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 200be5899b3SLisandro Dalcin } 201eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 202eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2031566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2041566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2051566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2061566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2071566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 208eb284becSJed Brown } 209be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 2101566a47fSLisandro Dalcin ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 211be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2121566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 213fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 214051f2191SLisandro Dalcin 215051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2161566a47fSLisandro Dalcin if (th->endpoint) { 2171566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2181566a47fSLisandro Dalcin } else { 219cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2201566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2211566a47fSLisandro Dalcin } 2221566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2231566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2241566a47fSLisandro Dalcin if (!accept) { 2251566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 226be5899b3SLisandro Dalcin ts->time_step = next_time_step; 227be5899b3SLisandro Dalcin goto reject_step; 228051f2191SLisandro Dalcin } 2293b1890cdSShri Abhyankar 230b1cb36f3SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 231b1cb36f3SHong Zhang th->ptime = ts->ptime; 232b1cb36f3SHong Zhang th->time_step = ts->time_step; 23317145e75SHong Zhang } 2341566a47fSLisandro Dalcin 2352b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 236cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 237051f2191SLisandro Dalcin break; 238051f2191SLisandro Dalcin 239051f2191SLisandro Dalcin reject_step: 240fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2411566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 242051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2431566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2443b1890cdSShri Abhyankar } 2453b1890cdSShri Abhyankar } 246316643e7SJed Brown PetscFunctionReturn(0); 247316643e7SJed Brown } 248316643e7SJed Brown 24942f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2502ca6e920SHong Zhang { 2512ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 252b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 2532ca6e920SHong Zhang PetscInt nadj; 2542ca6e920SHong Zhang PetscErrorCode ierr; 2552ca6e920SHong Zhang Mat J,Jp; 2562ca6e920SHong Zhang KSP ksp; 2572ca6e920SHong Zhang PetscReal shift; 2582ca6e920SHong Zhang 2592ca6e920SHong Zhang PetscFunctionBegin; 2602ca6e920SHong Zhang 2612ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 262302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 2632ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 264b5e0532dSHong Zhang 265b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 2663fd52205SHong Zhang th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/ 267b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 2682ca6e920SHong Zhang 269a4cab896SHong Zhang /* Build RHS */ 2702c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 27105755b9cSHong Zhang if (th->endpoint) { 272d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr); 27305755b9cSHong Zhang }else { 274d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 27505755b9cSHong Zhang } 27605755b9cSHong Zhang } 277abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 278b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 279b5e0532dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr); 2802c39e106SBarry Smith if (ts->vec_costintegral) { 281b5e0532dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 28236eaed60SHong Zhang } 2832ca6e920SHong Zhang } 2843c54f92cSHong Zhang 2852ca6e920SHong Zhang /* Build LHS */ 2862ca6e920SHong Zhang shift = -1./(th->Theta*ts->time_step); 2873c54f92cSHong Zhang if (th->endpoint) { 2883c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2893c54f92cSHong Zhang }else { 2903c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2913c54f92cSHong Zhang } 2922ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 2932ca6e920SHong Zhang 2942ca6e920SHong Zhang /* Solve LHS X = RHS */ 295abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 296b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 2972ca6e920SHong Zhang } 2983c54f92cSHong Zhang 29936eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 300b5e0532dSHong Zhang if(th->endpoint) { /* two-stage case */ 301b5e0532dSHong Zhang if (th->Theta!=1.) { 3023fd52205SHong Zhang shift = -1./((th->Theta-1.)*ts->time_step); 303b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3042c39e106SBarry Smith if (ts->vec_costintegral) { 305b5e0532dSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 3068263dbbfSHong Zhang } 307abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 308b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3092c39e106SBarry Smith if (ts->vec_costintegral) { 31036eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 31136eaed60SHong Zhang } 3123c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 3132ca6e920SHong Zhang } 314b5e0532dSHong Zhang }else { /* backward Euler */ 315b5e0532dSHong Zhang shift = 0.0; 316b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 317b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 318b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 319b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 320b5e0532dSHong Zhang if (ts->vec_costintegral) { 321b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 322b5e0532dSHong Zhang } 323b5e0532dSHong Zhang } 324b5e0532dSHong Zhang } 3253fd52205SHong Zhang 3263fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 3275bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 328abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 329b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 330b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3313fd52205SHong Zhang } 332b5e0532dSHong Zhang if (th->Theta!=1.) { 333b5e0532dSHong Zhang ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 334abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 335b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 336b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 337b5e0532dSHong Zhang } 3383fd52205SHong Zhang } 3392c39e106SBarry Smith if (ts->vec_costintegral) { 340d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 341abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 34236eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 34336eaed60SHong Zhang } 344b5e0532dSHong Zhang if (th->Theta!=1.) { 345b5e0532dSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 346abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 34736eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 34836eaed60SHong Zhang } 34936eaed60SHong Zhang } 3503fd52205SHong Zhang } 351b5e0532dSHong Zhang } 3523fd52205SHong Zhang }else { /* one-stage case */ 3533c54f92cSHong Zhang shift = 0.0; 354a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 3552c39e106SBarry Smith if (ts->vec_costintegral) { 356d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 3573c54f92cSHong Zhang } 358abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 359b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 360b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 3612c39e106SBarry Smith if (ts->vec_costintegral) { 36236eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 36336eaed60SHong Zhang } 3642ca6e920SHong Zhang } 3653fd52205SHong Zhang if (ts->vecs_sensip) { 3665bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 367abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 368b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 369b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3703fd52205SHong Zhang } 3712c39e106SBarry Smith if (ts->vec_costintegral) { 372d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 373abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 37436eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 37536eaed60SHong Zhang } 37636eaed60SHong Zhang } 3773fd52205SHong Zhang } 3782ca6e920SHong Zhang } 3792ca6e920SHong Zhang 3802ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 3812ca6e920SHong Zhang PetscFunctionReturn(0); 3822ca6e920SHong Zhang } 3832ca6e920SHong Zhang 384cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 385cd652676SJed Brown { 386cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 387be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 388cd652676SJed Brown PetscErrorCode ierr; 389cd652676SJed Brown 390cd652676SJed Brown PetscFunctionBegin; 391a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 392be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 393be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 394cd652676SJed Brown PetscFunctionReturn(0); 395cd652676SJed Brown } 396cd652676SJed Brown 3971566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 3981566a47fSLisandro Dalcin { 3991566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4001566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 4011566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 4027453f775SEmil Constantinescu PetscReal wltea,wlter; 4031566a47fSLisandro Dalcin PetscErrorCode ierr; 4041566a47fSLisandro Dalcin 4051566a47fSLisandro Dalcin PetscFunctionBegin; 4062ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 4071566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 408fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 4091566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 410fecfb714SLisandro Dalcin { 411be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 412be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 4131566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 4141566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 4151566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 4161566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 4171566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 4187453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 4191566a47fSLisandro Dalcin } 4201566a47fSLisandro Dalcin if (order) *order = 2; 4211566a47fSLisandro Dalcin PetscFunctionReturn(0); 4221566a47fSLisandro Dalcin } 4231566a47fSLisandro Dalcin 4241566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 4251566a47fSLisandro Dalcin { 4261566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4271566a47fSLisandro Dalcin PetscErrorCode ierr; 4281566a47fSLisandro Dalcin 4291566a47fSLisandro Dalcin PetscFunctionBegin; 4301566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 4311566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 4321566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4331566a47fSLisandro Dalcin } 4341566a47fSLisandro Dalcin PetscFunctionReturn(0); 4351566a47fSLisandro Dalcin } 4361566a47fSLisandro Dalcin 437316643e7SJed Brown /*------------------------------------------------------------*/ 438277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 439316643e7SJed Brown { 440316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 441316643e7SJed Brown PetscErrorCode ierr; 442316643e7SJed Brown 443316643e7SJed Brown PetscFunctionBegin; 4446bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 4456bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 4463b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 447eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 4481566a47fSLisandro Dalcin 4491566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 4501566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 4511566a47fSLisandro Dalcin 452b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 453b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 454b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 455b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 456277b19d0SLisandro Dalcin PetscFunctionReturn(0); 457277b19d0SLisandro Dalcin } 458277b19d0SLisandro Dalcin 459277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 460277b19d0SLisandro Dalcin { 461277b19d0SLisandro Dalcin PetscErrorCode ierr; 462277b19d0SLisandro Dalcin 463277b19d0SLisandro Dalcin PetscFunctionBegin; 464277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 465277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 466bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 467bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 468bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 469bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 470316643e7SJed Brown PetscFunctionReturn(0); 471316643e7SJed Brown } 472316643e7SJed Brown 473316643e7SJed Brown /* 474316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 4752b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 476316643e7SJed Brown */ 4770f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 478316643e7SJed Brown { 479316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 480316643e7SJed Brown PetscErrorCode ierr; 4817445fe48SJed Brown Vec X0,Xdot; 4827445fe48SJed Brown DM dm,dmsave; 4831566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 484316643e7SJed Brown 485316643e7SJed Brown PetscFunctionBegin; 4867445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4875a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 4887445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 489b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 4907445fe48SJed Brown 4917445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 4927445fe48SJed Brown dmsave = ts->dm; 4937445fe48SJed Brown ts->dm = dm; 4947445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 4957445fe48SJed Brown ts->dm = dmsave; 4960d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 497316643e7SJed Brown PetscFunctionReturn(0); 498316643e7SJed Brown } 499316643e7SJed Brown 500d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 501316643e7SJed Brown { 502316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 503316643e7SJed Brown PetscErrorCode ierr; 5047445fe48SJed Brown Vec Xdot; 5057445fe48SJed Brown DM dm,dmsave; 5061566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 507316643e7SJed Brown 508316643e7SJed Brown PetscFunctionBegin; 5097445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 510be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 5110298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 5127445fe48SJed Brown 5137445fe48SJed Brown dmsave = ts->dm; 5147445fe48SJed Brown ts->dm = dm; 515d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 5167445fe48SJed Brown ts->dm = dmsave; 5170298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 518316643e7SJed Brown PetscFunctionReturn(0); 519316643e7SJed Brown } 520316643e7SJed Brown 521316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 522316643e7SJed Brown { 523316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 5242ffb9264SLisandro Dalcin PetscBool match; 525316643e7SJed Brown PetscErrorCode ierr; 526316643e7SJed Brown 527316643e7SJed Brown PetscFunctionBegin; 528b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 529b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 530b1cb36f3SHong Zhang } 53139d13666SHong Zhang if (!th->X) { 532316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 53339d13666SHong Zhang } 53439d13666SHong Zhang if (!th->Xdot) { 535316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 53639d13666SHong Zhang } 53739d13666SHong Zhang if (!th->X0) { 5383b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 53939d13666SHong Zhang } 5401566a47fSLisandro Dalcin if (th->endpoint) { 5411566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 5427445fe48SJed Brown } 5433b1890cdSShri Abhyankar 5441566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 5451566a47fSLisandro Dalcin 5461566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 5471566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5481566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5491566a47fSLisandro Dalcin 5501566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 5511566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 5522ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 5532ffb9264SLisandro Dalcin if (!match) { 5541566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 5551566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 5563b1890cdSShri Abhyankar } 5571566a47fSLisandro Dalcin 5581566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 559b4dd3bf9SBarry Smith PetscFunctionReturn(0); 560b4dd3bf9SBarry Smith } 5610c7376e5SHong Zhang 562b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 563b4dd3bf9SBarry Smith 564b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 565b4dd3bf9SBarry Smith { 566b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 567b4dd3bf9SBarry Smith PetscErrorCode ierr; 568b4dd3bf9SBarry Smith 569b4dd3bf9SBarry Smith PetscFunctionBegin; 570b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 5712ca6e920SHong Zhang if(ts->vecs_sensip) { 572b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 5732ca6e920SHong Zhang } 574b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 575316643e7SJed Brown PetscFunctionReturn(0); 576316643e7SJed Brown } 577316643e7SJed Brown /*------------------------------------------------------------*/ 578316643e7SJed Brown 5794416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 580316643e7SJed Brown { 581316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 582316643e7SJed Brown PetscErrorCode ierr; 583316643e7SJed Brown 584316643e7SJed Brown PetscFunctionBegin; 585e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 586316643e7SJed Brown { 5870298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 5880298fd71SBarry 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); 58903842d09SLisandro Dalcin ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 590316643e7SJed Brown } 591316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 592316643e7SJed Brown PetscFunctionReturn(0); 593316643e7SJed Brown } 594316643e7SJed Brown 595316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 596316643e7SJed Brown { 597316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 598ace3abfcSBarry Smith PetscBool iascii; 599316643e7SJed Brown PetscErrorCode ierr; 600316643e7SJed Brown 601316643e7SJed Brown PetscFunctionBegin; 602251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 603316643e7SJed Brown if (iascii) { 6047c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 605ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 606316643e7SJed Brown } 607316643e7SJed Brown PetscFunctionReturn(0); 608316643e7SJed Brown } 609316643e7SJed Brown 610560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 6110de4c49aSJed Brown { 6120de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6130de4c49aSJed Brown 6140de4c49aSJed Brown PetscFunctionBegin; 6150de4c49aSJed Brown *theta = th->Theta; 6160de4c49aSJed Brown PetscFunctionReturn(0); 6170de4c49aSJed Brown } 6180de4c49aSJed Brown 619560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 6200de4c49aSJed Brown { 6210de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6220de4c49aSJed Brown 6230de4c49aSJed Brown PetscFunctionBegin; 6247c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 6250de4c49aSJed Brown th->Theta = theta; 6261566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 6270de4c49aSJed Brown PetscFunctionReturn(0); 6280de4c49aSJed Brown } 629eb284becSJed Brown 630560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 63126f2ff8fSLisandro Dalcin { 63226f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 63326f2ff8fSLisandro Dalcin 63426f2ff8fSLisandro Dalcin PetscFunctionBegin; 63526f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 63626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 63726f2ff8fSLisandro Dalcin } 63826f2ff8fSLisandro Dalcin 639560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 640eb284becSJed Brown { 641eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 642eb284becSJed Brown 643eb284becSJed Brown PetscFunctionBegin; 644eb284becSJed Brown th->endpoint = flg; 645eb284becSJed Brown PetscFunctionReturn(0); 646eb284becSJed Brown } 6470de4c49aSJed Brown 648f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 649f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 650f9c1d6abSBarry Smith { 651f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 652f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 6533fd8ae06SJed Brown const PetscReal one = 1.0; 654f9c1d6abSBarry Smith 655f9c1d6abSBarry Smith PetscFunctionBegin; 6563fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 657f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 658f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 659f9c1d6abSBarry Smith PetscFunctionReturn(0); 660f9c1d6abSBarry Smith } 661f9c1d6abSBarry Smith #endif 662f9c1d6abSBarry Smith 66342682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 66442682096SHong Zhang { 66542682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 66642682096SHong Zhang 66742682096SHong Zhang PetscFunctionBegin; 6681566a47fSLisandro Dalcin if (ns) *ns = 1; 6691566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 67042682096SHong Zhang PetscFunctionReturn(0); 67142682096SHong Zhang } 672f9c1d6abSBarry Smith 673316643e7SJed Brown /* ------------------------------------------------------------ */ 674316643e7SJed Brown /*MC 67596f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 676316643e7SJed Brown 677316643e7SJed Brown Level: beginner 678316643e7SJed Brown 6794eb428fdSBarry Smith Options Database: 680c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 681c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 68203842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 6834eb428fdSBarry Smith 684eb284becSJed Brown Notes: 6850c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 6860c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 6874eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 6884eb428fdSBarry Smith 689eb284becSJed Brown This method can be applied to DAE. 690eb284becSJed Brown 691eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 692eb284becSJed Brown 693eb284becSJed Brown .vb 694eb284becSJed Brown Theta | Theta 695eb284becSJed Brown ------------- 696eb284becSJed Brown | 1 697eb284becSJed Brown .ve 698eb284becSJed Brown 699eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 700eb284becSJed Brown 701eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 702eb284becSJed Brown 703eb284becSJed Brown .vb 704eb284becSJed Brown 0 | 0 0 705eb284becSJed Brown 1 | 1-Theta Theta 706eb284becSJed Brown ------------------- 707eb284becSJed Brown | 1-Theta Theta 708eb284becSJed Brown .ve 709eb284becSJed Brown 710eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 711eb284becSJed Brown 712eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 713eb284becSJed Brown 714eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 715eb284becSJed Brown 7164eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 717eb284becSJed Brown 718eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 719316643e7SJed Brown 720316643e7SJed Brown M*/ 7218cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 722316643e7SJed Brown { 723316643e7SJed Brown TS_Theta *th; 724316643e7SJed Brown PetscErrorCode ierr; 725316643e7SJed Brown 726316643e7SJed Brown PetscFunctionBegin; 727277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 728316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 729316643e7SJed Brown ts->ops->view = TSView_Theta; 730316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 73142f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 732316643e7SJed Brown ts->ops->step = TSStep_Theta; 733cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 7341566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 73524655328SShri ts->ops->rollback = TSRollBack_Theta; 736316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 7370f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 7380f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 739f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 740f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 741f9c1d6abSBarry Smith #endif 74242682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 74342f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 744b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 745b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 7462ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 747316643e7SJed Brown 748*efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 749*efd4aadfSBarry Smith 750b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 751316643e7SJed Brown ts->data = (void*)th; 752316643e7SJed Brown 7536f700aefSJed Brown th->extrapolate = PETSC_FALSE; 754316643e7SJed Brown th->Theta = 0.5; 7551566a47fSLisandro Dalcin th->order = 2; 756bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 757bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 758bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 759bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 760316643e7SJed Brown PetscFunctionReturn(0); 761316643e7SJed Brown } 7620de4c49aSJed Brown 7630de4c49aSJed Brown /*@ 7640de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 7650de4c49aSJed Brown 7660de4c49aSJed Brown Not Collective 7670de4c49aSJed Brown 7680de4c49aSJed Brown Input Parameter: 7690de4c49aSJed Brown . ts - timestepping context 7700de4c49aSJed Brown 7710de4c49aSJed Brown Output Parameter: 7720de4c49aSJed Brown . theta - stage abscissa 7730de4c49aSJed Brown 7740de4c49aSJed Brown Note: 7750de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 7760de4c49aSJed Brown 7770de4c49aSJed Brown Level: Advanced 7780de4c49aSJed Brown 7790de4c49aSJed Brown .seealso: TSThetaSetTheta() 7800de4c49aSJed Brown @*/ 7817087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 7820de4c49aSJed Brown { 7834ac538c5SBarry Smith PetscErrorCode ierr; 7840de4c49aSJed Brown 7850de4c49aSJed Brown PetscFunctionBegin; 786afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7870de4c49aSJed Brown PetscValidPointer(theta,2); 7884ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 7890de4c49aSJed Brown PetscFunctionReturn(0); 7900de4c49aSJed Brown } 7910de4c49aSJed Brown 7920de4c49aSJed Brown /*@ 7930de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 7940de4c49aSJed Brown 7950de4c49aSJed Brown Not Collective 7960de4c49aSJed Brown 7970de4c49aSJed Brown Input Parameter: 7980de4c49aSJed Brown + ts - timestepping context 7990de4c49aSJed Brown - theta - stage abscissa 8000de4c49aSJed Brown 8010de4c49aSJed Brown Options Database: 8020de4c49aSJed Brown . -ts_theta_theta <theta> 8030de4c49aSJed Brown 8040de4c49aSJed Brown Level: Intermediate 8050de4c49aSJed Brown 8060de4c49aSJed Brown .seealso: TSThetaGetTheta() 8070de4c49aSJed Brown @*/ 8087087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 8090de4c49aSJed Brown { 8104ac538c5SBarry Smith PetscErrorCode ierr; 8110de4c49aSJed Brown 8120de4c49aSJed Brown PetscFunctionBegin; 813afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8144ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 8150de4c49aSJed Brown PetscFunctionReturn(0); 8160de4c49aSJed Brown } 817f33bbcb6SJed Brown 81826f2ff8fSLisandro Dalcin /*@ 81926f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 82026f2ff8fSLisandro Dalcin 82126f2ff8fSLisandro Dalcin Not Collective 82226f2ff8fSLisandro Dalcin 82326f2ff8fSLisandro Dalcin Input Parameter: 82426f2ff8fSLisandro Dalcin . ts - timestepping context 82526f2ff8fSLisandro Dalcin 82626f2ff8fSLisandro Dalcin Output Parameter: 82726f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 82826f2ff8fSLisandro Dalcin 82926f2ff8fSLisandro Dalcin Level: Advanced 83026f2ff8fSLisandro Dalcin 83126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 83226f2ff8fSLisandro Dalcin @*/ 83326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 83426f2ff8fSLisandro Dalcin { 83526f2ff8fSLisandro Dalcin PetscErrorCode ierr; 83626f2ff8fSLisandro Dalcin 83726f2ff8fSLisandro Dalcin PetscFunctionBegin; 83826f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 83926f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 840163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 84126f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 84226f2ff8fSLisandro Dalcin } 84326f2ff8fSLisandro Dalcin 844eb284becSJed Brown /*@ 845eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 846eb284becSJed Brown 847eb284becSJed Brown Not Collective 848eb284becSJed Brown 849eb284becSJed Brown Input Parameter: 850eb284becSJed Brown + ts - timestepping context 851eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 852eb284becSJed Brown 853eb284becSJed Brown Options Database: 854eb284becSJed Brown . -ts_theta_endpoint <flg> 855eb284becSJed Brown 856eb284becSJed Brown Level: Intermediate 857eb284becSJed Brown 858eb284becSJed Brown .seealso: TSTHETA, TSCN 859eb284becSJed Brown @*/ 860eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 861eb284becSJed Brown { 862eb284becSJed Brown PetscErrorCode ierr; 863eb284becSJed Brown 864eb284becSJed Brown PetscFunctionBegin; 865eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 866eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 867eb284becSJed Brown PetscFunctionReturn(0); 868eb284becSJed Brown } 869eb284becSJed Brown 870f33bbcb6SJed Brown /* 871f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 872f33bbcb6SJed Brown * The creation functions for these specializations are below. 873f33bbcb6SJed Brown */ 874f33bbcb6SJed Brown 8751566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 8761566a47fSLisandro Dalcin { 8771566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 8781566a47fSLisandro Dalcin PetscErrorCode ierr; 8791566a47fSLisandro Dalcin 8801566a47fSLisandro Dalcin PetscFunctionBegin; 8811566a47fSLisandro Dalcin if (th->Theta != 1.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler\n"); 8821566a47fSLisandro Dalcin if (th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler\n"); 8831566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 8841566a47fSLisandro Dalcin PetscFunctionReturn(0); 8851566a47fSLisandro Dalcin } 8861566a47fSLisandro Dalcin 887f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 888f33bbcb6SJed Brown { 889f33bbcb6SJed Brown PetscFunctionBegin; 890f33bbcb6SJed Brown PetscFunctionReturn(0); 891f33bbcb6SJed Brown } 892f33bbcb6SJed Brown 893f33bbcb6SJed Brown /*MC 894f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 895f33bbcb6SJed Brown 896f33bbcb6SJed Brown Level: beginner 897f33bbcb6SJed Brown 8984eb428fdSBarry Smith Notes: 899c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 9004eb428fdSBarry Smith 9011566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 9024eb428fdSBarry Smith 903f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 904f33bbcb6SJed Brown 905f33bbcb6SJed Brown M*/ 9068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 907f33bbcb6SJed Brown { 908f33bbcb6SJed Brown PetscErrorCode ierr; 909f33bbcb6SJed Brown 910f33bbcb6SJed Brown PetscFunctionBegin; 911f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 912f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 9131566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 9140c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 915f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 916f33bbcb6SJed Brown PetscFunctionReturn(0); 917f33bbcb6SJed Brown } 918f33bbcb6SJed Brown 9191566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 9201566a47fSLisandro Dalcin { 9211566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 9221566a47fSLisandro Dalcin PetscErrorCode ierr; 9231566a47fSLisandro Dalcin 9241566a47fSLisandro Dalcin PetscFunctionBegin; 9251566a47fSLisandro Dalcin if (th->Theta != 0.5) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson\n"); 9261566a47fSLisandro Dalcin if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson\n"); 9271566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 9281566a47fSLisandro Dalcin PetscFunctionReturn(0); 9291566a47fSLisandro Dalcin } 9301566a47fSLisandro Dalcin 931f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 932f33bbcb6SJed Brown { 933f33bbcb6SJed Brown PetscFunctionBegin; 934f33bbcb6SJed Brown PetscFunctionReturn(0); 935f33bbcb6SJed Brown } 936f33bbcb6SJed Brown 937f33bbcb6SJed Brown /*MC 938f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 939f33bbcb6SJed Brown 940f33bbcb6SJed Brown Level: beginner 941f33bbcb6SJed Brown 942f33bbcb6SJed Brown Notes: 9437cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 9447cf5af47SJed Brown 9457cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 946f33bbcb6SJed Brown 947f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 948f33bbcb6SJed Brown 949f33bbcb6SJed Brown M*/ 9508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 951f33bbcb6SJed Brown { 952f33bbcb6SJed Brown PetscErrorCode ierr; 953f33bbcb6SJed Brown 954f33bbcb6SJed Brown PetscFunctionBegin; 955f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 956f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 957eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 9580c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 959f33bbcb6SJed Brown ts->ops->view = TSView_CN; 960f33bbcb6SJed Brown PetscFunctionReturn(0); 961f33bbcb6SJed Brown } 962