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 { 10715f1b00SHong Zhang /* context for time stepping */ 111566a47fSLisandro Dalcin PetscReal stage_time; 121566a47fSLisandro Dalcin Vec X0,X,Xdot; /* Storage for stages and time derivative */ 131566a47fSLisandro Dalcin Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 141566a47fSLisandro Dalcin PetscReal Theta; 15715f1b00SHong Zhang PetscReal ptime; 16715f1b00SHong Zhang PetscReal time_step; 171566a47fSLisandro Dalcin PetscInt order; 181566a47fSLisandro Dalcin PetscBool endpoint; 191566a47fSLisandro Dalcin PetscBool extrapolate; 20715f1b00SHong Zhang TSStepStatus status; 21715f1b00SHong Zhang Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 221566a47fSLisandro Dalcin 23715f1b00SHong Zhang /* context for sensitivity analysis */ 24022c081aSHong Zhang PetscInt num_tlm; /* Total number of tangent linear equations */ 25b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 26b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 2713b1b0ffSHong Zhang Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 2813b1b0ffSHong Zhang Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 29b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 3013b1b0ffSHong Zhang Mat MatFwdSensip0; /* backup for roll-backs due to events */ 317207e0fdSHong Zhang Mat MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 327207e0fdSHong Zhang Mat MatIntegralSensip0; /* backup for roll-backs due to events */ 33b5b4017aSHong Zhang Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 34b5b4017aSHong Zhang Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 35b5b4017aSHong Zhang Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 36b5b4017aSHong Zhang Vec *VecsAffine; /* Working vectors to store residuals */ 37715f1b00SHong Zhang /* context for error estimation */ 381566a47fSLisandro Dalcin Vec vec_sol_prev; 391566a47fSLisandro Dalcin Vec vec_lte_work; 40316643e7SJed Brown } TS_Theta; 41316643e7SJed Brown 427445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 437445fe48SJed Brown { 447445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 457445fe48SJed Brown PetscErrorCode ierr; 467445fe48SJed Brown 477445fe48SJed Brown PetscFunctionBegin; 487445fe48SJed Brown if (X0) { 497445fe48SJed Brown if (dm && dm != ts->dm) { 500d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 517445fe48SJed Brown } else *X0 = ts->vec_sol; 527445fe48SJed Brown } 537445fe48SJed Brown if (Xdot) { 547445fe48SJed Brown if (dm && dm != ts->dm) { 550d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 567445fe48SJed Brown } else *Xdot = th->Xdot; 577445fe48SJed Brown } 587445fe48SJed Brown PetscFunctionReturn(0); 597445fe48SJed Brown } 607445fe48SJed Brown 610d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 620d0b770aSPeter Brune { 630d0b770aSPeter Brune PetscErrorCode ierr; 640d0b770aSPeter Brune 650d0b770aSPeter Brune PetscFunctionBegin; 660d0b770aSPeter Brune if (X0) { 670d0b770aSPeter Brune if (dm && dm != ts->dm) { 680d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 690d0b770aSPeter Brune } 700d0b770aSPeter Brune } 710d0b770aSPeter Brune if (Xdot) { 720d0b770aSPeter Brune if (dm && dm != ts->dm) { 730d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 740d0b770aSPeter Brune } 750d0b770aSPeter Brune } 760d0b770aSPeter Brune PetscFunctionReturn(0); 770d0b770aSPeter Brune } 780d0b770aSPeter Brune 797445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 807445fe48SJed Brown { 817445fe48SJed Brown PetscFunctionBegin; 827445fe48SJed Brown PetscFunctionReturn(0); 837445fe48SJed Brown } 847445fe48SJed Brown 857445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 867445fe48SJed Brown { 877445fe48SJed Brown TS ts = (TS)ctx; 887445fe48SJed Brown PetscErrorCode ierr; 897445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 907445fe48SJed Brown 917445fe48SJed Brown PetscFunctionBegin; 927445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 937445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 947445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 957445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 967445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 977445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 980d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 990d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 1007445fe48SJed Brown PetscFunctionReturn(0); 1017445fe48SJed Brown } 1027445fe48SJed Brown 103258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 104258e1594SPeter Brune { 105258e1594SPeter Brune PetscFunctionBegin; 106258e1594SPeter Brune PetscFunctionReturn(0); 107258e1594SPeter Brune } 108258e1594SPeter Brune 109258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 110258e1594SPeter Brune { 111258e1594SPeter Brune TS ts = (TS)ctx; 112258e1594SPeter Brune PetscErrorCode ierr; 113258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 114258e1594SPeter Brune 115258e1594SPeter Brune PetscFunctionBegin; 116258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 117258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 118258e1594SPeter Brune 119258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 120258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune 122258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124258e1594SPeter Brune 125258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 126258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 127258e1594SPeter Brune PetscFunctionReturn(0); 128258e1594SPeter Brune } 129258e1594SPeter Brune 130022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 131022c081aSHong Zhang { 132022c081aSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 133cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 134022c081aSHong Zhang PetscErrorCode ierr; 135022c081aSHong Zhang 136022c081aSHong Zhang PetscFunctionBegin; 137022c081aSHong Zhang if (th->endpoint) { 138022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 139022c081aSHong Zhang if (th->Theta!=1.0) { 140cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 141cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 142022c081aSHong Zhang } 143cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 144cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 145022c081aSHong Zhang } else { 146cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 147cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 148022c081aSHong Zhang } 149022c081aSHong Zhang PetscFunctionReturn(0); 150022c081aSHong Zhang } 151022c081aSHong Zhang 152b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 153b1cb36f3SHong Zhang { 154b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 155cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 156b1cb36f3SHong Zhang PetscErrorCode ierr; 157b1cb36f3SHong Zhang 158b1cb36f3SHong Zhang PetscFunctionBegin; 159b1cb36f3SHong Zhang /* backup cost integral */ 160cd4cee2dSHong Zhang ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr); 161022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 162b1cb36f3SHong Zhang PetscFunctionReturn(0); 163b1cb36f3SHong Zhang } 164b1cb36f3SHong Zhang 165b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 166b1cb36f3SHong Zhang { 167b1cb36f3SHong Zhang PetscErrorCode ierr; 168b1cb36f3SHong Zhang 169b1cb36f3SHong Zhang PetscFunctionBegin; 170022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 17124655328SShri PetscFunctionReturn(0); 17224655328SShri } 17324655328SShri 174470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 1751566a47fSLisandro Dalcin { 1761566a47fSLisandro Dalcin PetscInt nits,lits; 1771566a47fSLisandro Dalcin PetscErrorCode ierr; 1781566a47fSLisandro Dalcin 1791566a47fSLisandro Dalcin PetscFunctionBegin; 1801566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1811566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1821566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1831566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1841566a47fSLisandro Dalcin PetscFunctionReturn(0); 1851566a47fSLisandro Dalcin } 1861566a47fSLisandro Dalcin 187193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 188316643e7SJed Brown { 189316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1901566a47fSLisandro Dalcin PetscInt rejections = 0; 1914957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1921566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 193051f2191SLisandro Dalcin PetscErrorCode ierr; 194316643e7SJed Brown 195316643e7SJed Brown PetscFunctionBegin; 1961566a47fSLisandro Dalcin if (!ts->steprollback) { 1972ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 1983b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 1991566a47fSLisandro Dalcin } 2001566a47fSLisandro Dalcin 2011566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 2021566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 203715f1b00SHong Zhang 2041566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 2051566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 206316643e7SJed Brown 207be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 208fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 209be5899b3SLisandro Dalcin ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 210be5899b3SLisandro Dalcin } 211eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 212eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2131566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2141566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2151566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2161566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2171566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 218eb284becSJed Brown } 219be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 220470880abSPatrick Sanan ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 221be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2221566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 223fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 224051f2191SLisandro Dalcin 225051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2261566a47fSLisandro Dalcin if (th->endpoint) { 2271566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2281566a47fSLisandro Dalcin } else { 229cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2301566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2311566a47fSLisandro Dalcin } 2321566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2331566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2341566a47fSLisandro Dalcin if (!accept) { 2351566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 236be5899b3SLisandro Dalcin ts->time_step = next_time_step; 237be5899b3SLisandro Dalcin goto reject_step; 238051f2191SLisandro Dalcin } 2393b1890cdSShri Abhyankar 240715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 241b1cb36f3SHong Zhang th->ptime = ts->ptime; 242b1cb36f3SHong Zhang th->time_step = ts->time_step; 24317145e75SHong Zhang } 2441566a47fSLisandro Dalcin 2452b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 246cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 247051f2191SLisandro Dalcin break; 248051f2191SLisandro Dalcin 249051f2191SLisandro Dalcin reject_step: 250fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2511566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 252051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2531566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2543b1890cdSShri Abhyankar } 2553b1890cdSShri Abhyankar } 256316643e7SJed Brown PetscFunctionReturn(0); 257316643e7SJed Brown } 258316643e7SJed Brown 259c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 260c9681e5dSHong Zhang { 261c9681e5dSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 262cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 263c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 264c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 265c9681e5dSHong Zhang PetscInt nadj; 2667207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 267c9681e5dSHong Zhang KSP ksp; 268c9681e5dSHong Zhang PetscReal shift; 269c9681e5dSHong Zhang PetscScalar *xarr; 270*077c4feaSHong Zhang TSEquationType eqtype; 271*077c4feaSHong Zhang PetscBool isexplicitode = PETSC_FALSE; 272c9681e5dSHong Zhang PetscErrorCode ierr; 273c9681e5dSHong Zhang 274c9681e5dSHong Zhang PetscFunctionBegin; 275*077c4feaSHong Zhang ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr); 276*077c4feaSHong Zhang if (eqtype == TS_EQ_ODE_EXPLICIT) { 277*077c4feaSHong Zhang isexplicitode = PETSC_TRUE; 278*077c4feaSHong Zhang VecsDeltaLam = ts->vecs_sensi; 279*077c4feaSHong Zhang VecsDeltaLam2 = ts->vecs_sensi2; 280*077c4feaSHong Zhang } 281c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 282c9681e5dSHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 283cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 2847207e0fdSHong Zhang if (quadts) { 285cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 286cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 2877207e0fdSHong Zhang } 288c9681e5dSHong Zhang 289c9681e5dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 290c9681e5dSHong Zhang th->stage_time = ts->ptime; /* time_step is negative*/ 291c9681e5dSHong Zhang th->ptime = ts->ptime + ts->time_step; 292c9681e5dSHong Zhang th->time_step = -ts->time_step; 293c9681e5dSHong Zhang 29433f72c5dSHong Zhang /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 2957207e0fdSHong Zhang if (quadts) { 296cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 2977207e0fdSHong Zhang } 298cd4cee2dSHong Zhang 299c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 300c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 301c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */ 302cd4cee2dSHong Zhang if (quadJ) { 303cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 304cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 305cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 306cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 307cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 308c9681e5dSHong Zhang } 309c9681e5dSHong Zhang } 310c9681e5dSHong Zhang 311c9681e5dSHong Zhang /* Build LHS for first-order adjoint */ 312c9681e5dSHong Zhang shift = 1./th->time_step; 313cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 314cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 315c9681e5dSHong Zhang 316c9681e5dSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 317c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 318c9681e5dSHong Zhang KSPConvergedReason kspreason; 319c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 320c9681e5dSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 321c9681e5dSHong Zhang if (kspreason < 0) { 322c9681e5dSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 323c9681e5dSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 324c9681e5dSHong Zhang } 325c9681e5dSHong Zhang } 326c9681e5dSHong Zhang 327c9681e5dSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 328c9681e5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 329c9681e5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 330c9681e5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 331c9681e5dSHong Zhang /* lambda_s^T F_UU w_1 */ 332b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 333c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 334b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 335c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 336c9681e5dSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 337c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 33833f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 339c9681e5dSHong Zhang if (ts->vecs_fup) { 34033f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 341c9681e5dSHong Zhang } 342c9681e5dSHong Zhang } 343c9681e5dSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 344c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 34505c9950eSHong Zhang KSPConvergedReason kspreason; 346c9681e5dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 34705c9950eSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 34805c9950eSHong Zhang if (kspreason < 0) { 34905c9950eSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 35005c9950eSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 35105c9950eSHong Zhang } 352c9681e5dSHong Zhang } 353c9681e5dSHong Zhang } 354c9681e5dSHong Zhang 355c9681e5dSHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 356*077c4feaSHong Zhang if (!isexplicitode) { 357c9681e5dSHong Zhang shift = 0.0; 358cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_U */ 359c9681e5dSHong Zhang ierr = MatScale(J,-1.);CHKERRQ(ierr); 360c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 36133f72c5dSHong Zhang /* Add f_U \lambda_s to the original RHS */ 362c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 363c9681e5dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr); 364c9681e5dSHong Zhang ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 365c9681e5dSHong Zhang if (ts->vecs_sensi2) { 366c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 367c9681e5dSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr); 368c9681e5dSHong Zhang ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 369c9681e5dSHong Zhang } 370c9681e5dSHong Zhang } 371*077c4feaSHong Zhang } 372c9681e5dSHong Zhang if (ts->vecs_sensip) { 373c9681e5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 3747207e0fdSHong Zhang if (quadts) { 375cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 3767207e0fdSHong Zhang } 377c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 378c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 379b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 38033f72c5dSHong Zhang /* lambda_s^T F_PP w_2 */ 381b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 382c9681e5dSHong Zhang } 383cd4cee2dSHong Zhang 384c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 385c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 386c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 387cd4cee2dSHong Zhang if (quadJp) { 388cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 389cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 390cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr); 391cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 392cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 39333f72c5dSHong Zhang } 394c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 395c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 396c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 397c9681e5dSHong Zhang if (ts->vecs_fpu) { 398c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 399c9681e5dSHong Zhang } 400c9681e5dSHong Zhang if (ts->vecs_fpp) { 401c9681e5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 402c9681e5dSHong Zhang } 403c9681e5dSHong Zhang } 404c9681e5dSHong Zhang } 405c9681e5dSHong Zhang } 406c9681e5dSHong Zhang 407c9681e5dSHong Zhang if (ts->vecs_sensi2) { 408c9681e5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 409c9681e5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 410c9681e5dSHong Zhang } 411c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 412c9681e5dSHong Zhang PetscFunctionReturn(0); 413c9681e5dSHong Zhang } 414c9681e5dSHong Zhang 41542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 4162ca6e920SHong Zhang { 4172ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 418cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 419b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 420b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 4212ca6e920SHong Zhang PetscInt nadj; 4227207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 4232ca6e920SHong Zhang KSP ksp; 4242ca6e920SHong Zhang PetscReal shift; 425b5b4017aSHong Zhang PetscScalar *xarr; 426b5b4017aSHong Zhang PetscErrorCode ierr; 4272ca6e920SHong Zhang 4282ca6e920SHong Zhang PetscFunctionBegin; 429c9681e5dSHong Zhang if (th->Theta == 1.) { 430c9681e5dSHong Zhang ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 431c9681e5dSHong Zhang PetscFunctionReturn(0); 432c9681e5dSHong Zhang } 4332ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 434302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 435cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 4367207e0fdSHong Zhang if (quadts) { 437cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 438cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 4397207e0fdSHong Zhang } 440b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 441022c081aSHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 442b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 443022c081aSHong Zhang th->time_step = -ts->time_step; 4442ca6e920SHong Zhang 445b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 44633f72c5dSHong Zhang /* Cost function has an integral term */ 4477207e0fdSHong Zhang if (quadts) { 44805755b9cSHong Zhang if (th->endpoint) { 449cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 45005755b9cSHong Zhang } else { 451cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 45205755b9cSHong Zhang } 4537207e0fdSHong Zhang } 45433f72c5dSHong Zhang 455abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 456b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 457022c081aSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 458cd4cee2dSHong Zhang if (quadJ) { 459cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 460cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 461cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 462cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 463cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 46436eaed60SHong Zhang } 4652ca6e920SHong Zhang } 4663c54f92cSHong Zhang 467b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 468022c081aSHong Zhang shift = 1./(th->Theta*th->time_step); 4693c54f92cSHong Zhang if (th->endpoint) { 470cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 4713c54f92cSHong Zhang } else { 472cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 4733c54f92cSHong Zhang } 474cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 4752ca6e920SHong Zhang 476b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 477abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 4781d14d03bSHong Zhang KSPConvergedReason kspreason; 479b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 4801d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 4811d14d03bSHong Zhang if (kspreason < 0) { 4821d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 4831d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 4841d14d03bSHong Zhang } 4852ca6e920SHong Zhang } 4863c54f92cSHong Zhang 4871d14d03bSHong Zhang /* Second-order adjoint */ 488b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 489b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 490b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 491b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 492b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 493b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 494b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 495b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 496b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 497b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 498b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 499b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 500b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 5011d14d03bSHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr); 50233f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 503b5b4017aSHong Zhang if (ts->vecs_fup) { 50433f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 505b5b4017aSHong Zhang } 506b5b4017aSHong Zhang } 507b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 508b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 5091d14d03bSHong Zhang KSPConvergedReason kspreason; 5101d14d03bSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 5111d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 5121d14d03bSHong Zhang if (kspreason < 0) { 5131d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 5141d14d03bSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 5151d14d03bSHong Zhang } 516b5b4017aSHong Zhang } 517b5b4017aSHong Zhang } 518b5b4017aSHong Zhang 51936eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 5201d14d03bSHong Zhang if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 521022c081aSHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 522cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 52333f72c5dSHong Zhang /* R_U at t_n */ 5247207e0fdSHong Zhang if (quadts) { 525cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 5267207e0fdSHong Zhang } 527abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 528b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 529cd4cee2dSHong Zhang if (quadJ) { 530cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 531cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 532cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr); 533cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 534cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 535b5b4017aSHong Zhang } 5361d14d03bSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 537b5b4017aSHong Zhang } 5381d14d03bSHong Zhang 5391d14d03bSHong Zhang /* Second-order adjoint */ 5401d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 541b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 542b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 543b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 544b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 545b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 546b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 54733f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 548b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 549b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 550b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 55133f72c5dSHong Zhang /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2 */ 552b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 55333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 554b5b4017aSHong Zhang if (ts->vecs_fup) { 55533f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 556b5b4017aSHong Zhang } 5571d14d03bSHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 558b5b4017aSHong Zhang } 559b5e0532dSHong Zhang } 5603fd52205SHong Zhang 5613fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 562b5b4017aSHong Zhang /* U_{n+1} */ 56333f72c5dSHong Zhang shift = -1./(th->Theta*th->time_step); 56433f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 5657207e0fdSHong Zhang if (quadts) { 566cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 5677207e0fdSHong Zhang } 568abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 569b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 57033f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5713fd52205SHong Zhang } 57233f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 57333f72c5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 57433f72c5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 57533f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 576b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 577b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 57833f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 57933f72c5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 58033f72c5dSHong Zhang 581b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 582b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 583b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 58433f72c5dSHong Zhang /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */ 585b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 58633f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 587b5b4017aSHong Zhang if (ts->vecs_fpu) { 58833f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 589b5b4017aSHong Zhang } 590b5b4017aSHong Zhang if (ts->vecs_fpp) { 59133f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 592b5b4017aSHong Zhang } 593b5b4017aSHong Zhang } 594b5b4017aSHong Zhang } 595b5b4017aSHong Zhang 596b5b4017aSHong Zhang /* U_s */ 59733f72c5dSHong Zhang shift = 1./((th->Theta-1.0)*th->time_step); 59833f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 5997207e0fdSHong Zhang if (quadts) { 600cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr); 6017207e0fdSHong Zhang } 602abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 603b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 60433f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 60533f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 60633f72c5dSHong Zhang /* Get w1 at t_n from TLM matrix */ 60733f72c5dSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 60833f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 609b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 610b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 61133f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 61233f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 613b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 614b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 615b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 61633f72c5dSHong Zhang /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */ 617b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 61833f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 619b5b4017aSHong Zhang if (ts->vecs_fpu) { 62033f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 621b5e0532dSHong Zhang } 622b5b4017aSHong Zhang if (ts->vecs_fpp) { 62333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 624b5b4017aSHong Zhang } 62536eaed60SHong Zhang } 62636eaed60SHong Zhang } 6273fd52205SHong Zhang } 628b5e0532dSHong Zhang } 6293fd52205SHong Zhang } else { /* one-stage case */ 6303c54f92cSHong Zhang shift = 0.0; 631cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 6327207e0fdSHong Zhang if (quadts) { 633cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 6347207e0fdSHong Zhang } 635abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 636b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 637022c081aSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 638cd4cee2dSHong Zhang if (quadJ) { 639cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 640cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 641cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vec_drdu_col);CHKERRQ(ierr); 642cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 643cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 64436eaed60SHong Zhang } 6452ca6e920SHong Zhang } 6463fd52205SHong Zhang if (ts->vecs_sensip) { 64733f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 6487207e0fdSHong Zhang if (quadts) { 649cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 6507207e0fdSHong Zhang } 651abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 65233f72c5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 65333f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 654cd4cee2dSHong Zhang if (quadJp) { 655cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 656cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 657cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr); 658cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 659cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 66036eaed60SHong Zhang } 66136eaed60SHong Zhang } 6623fd52205SHong Zhang } 6632ca6e920SHong Zhang } 6642ca6e920SHong Zhang 6652ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6662ca6e920SHong Zhang PetscFunctionReturn(0); 6672ca6e920SHong Zhang } 6682ca6e920SHong Zhang 669cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 670cd652676SJed Brown { 671cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 672be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 673cd652676SJed Brown PetscErrorCode ierr; 674cd652676SJed Brown 675cd652676SJed Brown PetscFunctionBegin; 676a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 677be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 678be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 679cd652676SJed Brown PetscFunctionReturn(0); 680cd652676SJed Brown } 681cd652676SJed Brown 6821566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 6831566a47fSLisandro Dalcin { 6841566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 6851566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 6861566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 6877453f775SEmil Constantinescu PetscReal wltea,wlter; 6881566a47fSLisandro Dalcin PetscErrorCode ierr; 6891566a47fSLisandro Dalcin 6901566a47fSLisandro Dalcin PetscFunctionBegin; 6912ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 6921566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 693fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 6941566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 695fecfb714SLisandro Dalcin { 696be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 697be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 6981566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 6991566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 7001566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 7011566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 7021566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 7037453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 7041566a47fSLisandro Dalcin } 7051566a47fSLisandro Dalcin if (order) *order = 2; 7061566a47fSLisandro Dalcin PetscFunctionReturn(0); 7071566a47fSLisandro Dalcin } 7081566a47fSLisandro Dalcin 7091566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 7101566a47fSLisandro Dalcin { 7111566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 7127207e0fdSHong Zhang TS quadts = ts->quadraturets; 7131566a47fSLisandro Dalcin PetscErrorCode ierr; 7141566a47fSLisandro Dalcin 7151566a47fSLisandro Dalcin PetscFunctionBegin; 7161566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 717cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 718cd4cee2dSHong Zhang ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr); 7191566a47fSLisandro Dalcin } 720715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 72113b1b0ffSHong Zhang if (ts->mat_sensip) { 72213b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7236f92202bSHong Zhang } 7247207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7257207e0fdSHong Zhang ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7266f92202bSHong Zhang } 727715f1b00SHong Zhang PetscFunctionReturn(0); 728715f1b00SHong Zhang } 729715f1b00SHong Zhang 730715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 731715f1b00SHong Zhang { 732715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 733cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 73413b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 735b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 7367207e0fdSHong Zhang PetscInt ntlm; 737715f1b00SHong Zhang KSP ksp; 7387207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 739715f1b00SHong Zhang PetscReal shift; 74013b1b0ffSHong Zhang PetscScalar *barr,*xarr; 741715f1b00SHong Zhang PetscErrorCode ierr; 742715f1b00SHong Zhang 743715f1b00SHong Zhang PetscFunctionBegin; 74413b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 74513b1b0ffSHong Zhang 7467207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7477207e0fdSHong Zhang ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7486f92202bSHong Zhang } 749715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 750cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 7517207e0fdSHong Zhang if (quadts) { 752cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 753cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 7547207e0fdSHong Zhang } 755715f1b00SHong Zhang 756715f1b00SHong Zhang /* Build RHS */ 757715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 758715f1b00SHong Zhang shift = 1./((th->Theta-1.)*th->time_step); 759cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 76013b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 76113b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 762715f1b00SHong Zhang 763022c081aSHong Zhang /* Add the f_p forcing terms */ 7640e11c21fSHong Zhang if (ts->Jacp) { 76533f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 76633f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 76733f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 76833f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 7690e11c21fSHong Zhang } 770715f1b00SHong Zhang } else { /* 1-stage method */ 771715f1b00SHong Zhang shift = 0.0; 772cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 77313b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 77413b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 77513b1b0ffSHong Zhang 776022c081aSHong Zhang /* Add the f_p forcing terms */ 7770e11c21fSHong Zhang if (ts->Jacp) { 77833f72c5dSHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 77933f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 780715f1b00SHong Zhang } 7810e11c21fSHong Zhang } 782715f1b00SHong Zhang 783715f1b00SHong Zhang /* Build LHS */ 784715f1b00SHong Zhang shift = 1/(th->Theta*th->time_step); 785715f1b00SHong Zhang if (th->endpoint) { 786cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 787715f1b00SHong Zhang } else { 788cd4cee2dSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 789715f1b00SHong Zhang } 790cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 791715f1b00SHong Zhang 792715f1b00SHong Zhang /* 793715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 794c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 795715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 796715f1b00SHong Zhang */ 797715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 7987207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 799cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 800cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr); 8017207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8027207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8037207e0fdSHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 804715f1b00SHong Zhang } 805715f1b00SHong Zhang } 806715f1b00SHong Zhang 807715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 808022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 809b5b4017aSHong Zhang KSPConvergedReason kspreason; 81013b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 811b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 812715f1b00SHong Zhang if (th->endpoint) { 81313b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 814b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 815b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 816b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 81713b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 818715f1b00SHong Zhang } else { 819b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 820715f1b00SHong Zhang } 821b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 822b5b4017aSHong Zhang if (kspreason < 0) { 823b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 824b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 825b5b4017aSHong Zhang } 826b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 82713b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 828715f1b00SHong Zhang } 829715f1b00SHong Zhang 830b5b4017aSHong Zhang 83113b1b0ffSHong Zhang /* 83213b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 833c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 83413b1b0ffSHong Zhang */ 8357207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 83613b1b0ffSHong Zhang if (!th->endpoint) { 8377207e0fdSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */ 838cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 839cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 8407207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8417207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8427207e0fdSHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 84313b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 84413b1b0ffSHong Zhang } else { 845cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 846cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 8477207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8487207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8497207e0fdSHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 85013b1b0ffSHong Zhang } 85113b1b0ffSHong Zhang } else { 85213b1b0ffSHong Zhang if (!th->endpoint) { 85313b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 854715f1b00SHong Zhang } 855715f1b00SHong Zhang } 8561566a47fSLisandro Dalcin PetscFunctionReturn(0); 8571566a47fSLisandro Dalcin } 8581566a47fSLisandro Dalcin 8596affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8606affb6f8SHong Zhang { 8616affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8626affb6f8SHong Zhang 8636affb6f8SHong Zhang PetscFunctionBegin; 8646affb6f8SHong Zhang if (ns) *ns = 1; 8656affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 8666affb6f8SHong Zhang PetscFunctionReturn(0); 8676affb6f8SHong Zhang } 8686affb6f8SHong Zhang 869316643e7SJed Brown /*------------------------------------------------------------*/ 870277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 871316643e7SJed Brown { 872316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 873316643e7SJed Brown PetscErrorCode ierr; 874316643e7SJed Brown 875316643e7SJed Brown PetscFunctionBegin; 8766bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 8776bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 8783b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 879eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 8801566a47fSLisandro Dalcin 8811566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 8821566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 8831566a47fSLisandro Dalcin 884b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 885277b19d0SLisandro Dalcin PetscFunctionReturn(0); 886277b19d0SLisandro Dalcin } 887277b19d0SLisandro Dalcin 888ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 889ece46509SHong Zhang { 890ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 891ece46509SHong Zhang PetscErrorCode ierr; 892ece46509SHong Zhang 893ece46509SHong Zhang PetscFunctionBegin; 894ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 895ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 896ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 897ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 898ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 899ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 900ece46509SHong Zhang PetscFunctionReturn(0); 901ece46509SHong Zhang } 902ece46509SHong Zhang 903277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 904277b19d0SLisandro Dalcin { 905277b19d0SLisandro Dalcin PetscErrorCode ierr; 906277b19d0SLisandro Dalcin 907277b19d0SLisandro Dalcin PetscFunctionBegin; 908277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 909b3a6b972SJed Brown if (ts->dm) { 910b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 911b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 912b3a6b972SJed Brown } 913277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 914bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 915bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 916bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 917bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 918316643e7SJed Brown PetscFunctionReturn(0); 919316643e7SJed Brown } 920316643e7SJed Brown 921316643e7SJed Brown /* 922316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9232b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 924316643e7SJed Brown */ 9250f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 926316643e7SJed Brown { 927316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 928316643e7SJed Brown PetscErrorCode ierr; 9297445fe48SJed Brown Vec X0,Xdot; 9307445fe48SJed Brown DM dm,dmsave; 9311566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 932316643e7SJed Brown 933316643e7SJed Brown PetscFunctionBegin; 9347445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9355a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9367445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 937b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 9387445fe48SJed Brown 9397445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9407445fe48SJed Brown dmsave = ts->dm; 9417445fe48SJed Brown ts->dm = dm; 9427445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9437445fe48SJed Brown ts->dm = dmsave; 9440d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 945316643e7SJed Brown PetscFunctionReturn(0); 946316643e7SJed Brown } 947316643e7SJed Brown 948d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 949316643e7SJed Brown { 950316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 951316643e7SJed Brown PetscErrorCode ierr; 9527445fe48SJed Brown Vec Xdot; 9537445fe48SJed Brown DM dm,dmsave; 9541566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 955316643e7SJed Brown 956316643e7SJed Brown PetscFunctionBegin; 9577445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 958be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 9590298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 9607445fe48SJed Brown 9617445fe48SJed Brown dmsave = ts->dm; 9627445fe48SJed Brown ts->dm = dm; 963d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 9647445fe48SJed Brown ts->dm = dmsave; 9650298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 966316643e7SJed Brown PetscFunctionReturn(0); 967316643e7SJed Brown } 968316643e7SJed Brown 969715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 970715f1b00SHong Zhang { 971715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 9727207e0fdSHong Zhang TS quadts = ts->quadraturets; 973715f1b00SHong Zhang PetscErrorCode ierr; 974715f1b00SHong Zhang 975715f1b00SHong Zhang PetscFunctionBegin; 976022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 97713b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 97813b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 9797207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 9807207e0fdSHong Zhang ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 9817207e0fdSHong Zhang ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr); 982715f1b00SHong Zhang } 9836f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 98413b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 98513b1b0ffSHong Zhang 986b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 987715f1b00SHong Zhang PetscFunctionReturn(0); 988715f1b00SHong Zhang } 989715f1b00SHong Zhang 9907adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 9917adebddeSHong Zhang { 9927adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 9937207e0fdSHong Zhang TS quadts = ts->quadraturets; 9947adebddeSHong Zhang PetscErrorCode ierr; 9957adebddeSHong Zhang 9967adebddeSHong Zhang PetscFunctionBegin; 9977207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 9987207e0fdSHong Zhang ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr); 9997207e0fdSHong Zhang ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr); 10007adebddeSHong Zhang } 10017adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 10027adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 10037adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 10047adebddeSHong Zhang PetscFunctionReturn(0); 10057adebddeSHong Zhang } 10067adebddeSHong Zhang 1007316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 1008316643e7SJed Brown { 1009316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1010cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 10112ffb9264SLisandro Dalcin PetscBool match; 1012316643e7SJed Brown PetscErrorCode ierr; 1013316643e7SJed Brown 1014316643e7SJed Brown PetscFunctionBegin; 1015cd4cee2dSHong Zhang if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 1016cd4cee2dSHong Zhang ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr); 1017b1cb36f3SHong Zhang } 101839d13666SHong Zhang if (!th->X) { 1019316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 102039d13666SHong Zhang } 102139d13666SHong Zhang if (!th->Xdot) { 1022316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 102339d13666SHong Zhang } 102439d13666SHong Zhang if (!th->X0) { 10253b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 102639d13666SHong Zhang } 10271566a47fSLisandro Dalcin if (th->endpoint) { 10281566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10297445fe48SJed Brown } 10303b1890cdSShri Abhyankar 10311566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 10321566a47fSLisandro Dalcin 10331566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10341566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10351566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10361566a47fSLisandro Dalcin 10371566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10381566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10392ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10402ffb9264SLisandro Dalcin if (!match) { 10411566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10421566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10433b1890cdSShri Abhyankar } 10441566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1045b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1046b4dd3bf9SBarry Smith } 10470c7376e5SHong Zhang 1048b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1049b4dd3bf9SBarry Smith 1050b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1051b4dd3bf9SBarry Smith { 1052b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1053b4dd3bf9SBarry Smith PetscErrorCode ierr; 1054b4dd3bf9SBarry Smith 1055b4dd3bf9SBarry Smith PetscFunctionBegin; 1056b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1057b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 10582ca6e920SHong Zhang if (ts->vecs_sensip) { 1059b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 10602ca6e920SHong Zhang } 1061b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1062b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1063b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 106467633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 106567633408SHong Zhang if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 106667633408SHong Zhang if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1067b5b4017aSHong Zhang } 1068c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1069c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 107067633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 107167633408SHong Zhang if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 107267633408SHong Zhang if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1073b5b4017aSHong Zhang } 1074316643e7SJed Brown PetscFunctionReturn(0); 1075316643e7SJed Brown } 1076316643e7SJed Brown 10774416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1078316643e7SJed Brown { 1079316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1080316643e7SJed Brown PetscErrorCode ierr; 1081316643e7SJed Brown 1082316643e7SJed Brown PetscFunctionBegin; 1083e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1084316643e7SJed Brown { 10850298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 10860298fd71SBarry 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); 108703842d09SLisandro 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); 1088316643e7SJed Brown } 1089316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1090316643e7SJed Brown PetscFunctionReturn(0); 1091316643e7SJed Brown } 1092316643e7SJed Brown 1093316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1094316643e7SJed Brown { 1095316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1096ace3abfcSBarry Smith PetscBool iascii; 1097316643e7SJed Brown PetscErrorCode ierr; 1098316643e7SJed Brown 1099316643e7SJed Brown PetscFunctionBegin; 1100251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1101316643e7SJed Brown if (iascii) { 11027c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1103ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1104316643e7SJed Brown } 1105316643e7SJed Brown PetscFunctionReturn(0); 1106316643e7SJed Brown } 1107316643e7SJed Brown 1108560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 11090de4c49aSJed Brown { 11100de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11110de4c49aSJed Brown 11120de4c49aSJed Brown PetscFunctionBegin; 11130de4c49aSJed Brown *theta = th->Theta; 11140de4c49aSJed Brown PetscFunctionReturn(0); 11150de4c49aSJed Brown } 11160de4c49aSJed Brown 1117560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 11180de4c49aSJed Brown { 11190de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11200de4c49aSJed Brown 11210de4c49aSJed Brown PetscFunctionBegin; 11227c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 11230de4c49aSJed Brown th->Theta = theta; 11241566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11250de4c49aSJed Brown PetscFunctionReturn(0); 11260de4c49aSJed Brown } 1127eb284becSJed Brown 1128560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 112926f2ff8fSLisandro Dalcin { 113026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 113126f2ff8fSLisandro Dalcin 113226f2ff8fSLisandro Dalcin PetscFunctionBegin; 113326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 113426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 113526f2ff8fSLisandro Dalcin } 113626f2ff8fSLisandro Dalcin 1137560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1138eb284becSJed Brown { 1139eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1140eb284becSJed Brown 1141eb284becSJed Brown PetscFunctionBegin; 1142eb284becSJed Brown th->endpoint = flg; 1143eb284becSJed Brown PetscFunctionReturn(0); 1144eb284becSJed Brown } 11450de4c49aSJed Brown 1146f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1147f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1148f9c1d6abSBarry Smith { 1149f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1150f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11513fd8ae06SJed Brown const PetscReal one = 1.0; 1152f9c1d6abSBarry Smith 1153f9c1d6abSBarry Smith PetscFunctionBegin; 11543fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1155f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1156f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1157f9c1d6abSBarry Smith PetscFunctionReturn(0); 1158f9c1d6abSBarry Smith } 1159f9c1d6abSBarry Smith #endif 1160f9c1d6abSBarry Smith 116142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 116242682096SHong Zhang { 116342682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 116442682096SHong Zhang 116542682096SHong Zhang PetscFunctionBegin; 11661566a47fSLisandro Dalcin if (ns) *ns = 1; 11671566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 116842682096SHong Zhang PetscFunctionReturn(0); 116942682096SHong Zhang } 1170f9c1d6abSBarry Smith 1171316643e7SJed Brown /* ------------------------------------------------------------ */ 1172316643e7SJed Brown /*MC 117396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1174316643e7SJed Brown 1175316643e7SJed Brown Level: beginner 1176316643e7SJed Brown 11774eb428fdSBarry Smith Options Database: 1178c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1179c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 118003842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 11814eb428fdSBarry Smith 1182eb284becSJed Brown Notes: 11830c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 11840c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 11854eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 11864eb428fdSBarry Smith 1187eb284becSJed Brown This method can be applied to DAE. 1188eb284becSJed Brown 1189eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1190eb284becSJed Brown 1191eb284becSJed Brown .vb 1192eb284becSJed Brown Theta | Theta 1193eb284becSJed Brown ------------- 1194eb284becSJed Brown | 1 1195eb284becSJed Brown .ve 1196eb284becSJed Brown 1197eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1198eb284becSJed Brown 1199eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1200eb284becSJed Brown 1201eb284becSJed Brown .vb 1202eb284becSJed Brown 0 | 0 0 1203eb284becSJed Brown 1 | 1-Theta Theta 1204eb284becSJed Brown ------------------- 1205eb284becSJed Brown | 1-Theta Theta 1206eb284becSJed Brown .ve 1207eb284becSJed Brown 1208eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1209eb284becSJed Brown 1210eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1211eb284becSJed Brown 1212eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1213eb284becSJed Brown 12144eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1215eb284becSJed Brown 1216eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1217316643e7SJed Brown 1218316643e7SJed Brown M*/ 12198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1220316643e7SJed Brown { 1221316643e7SJed Brown TS_Theta *th; 1222316643e7SJed Brown PetscErrorCode ierr; 1223316643e7SJed Brown 1224316643e7SJed Brown PetscFunctionBegin; 1225277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1226ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1227316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1228316643e7SJed Brown ts->ops->view = TSView_Theta; 1229316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 123042f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1231ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1232316643e7SJed Brown ts->ops->step = TSStep_Theta; 1233cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12341566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 123524655328SShri ts->ops->rollback = TSRollBack_Theta; 1236316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12370f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12380f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1239f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1240f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1241f9c1d6abSBarry Smith #endif 124242682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 124342f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1244b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1245b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12462ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12476affb6f8SHong Zhang 1248715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12497adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1250715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12516affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1252316643e7SJed Brown 1253efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1254efd4aadfSBarry Smith 1255b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1256316643e7SJed Brown ts->data = (void*)th; 1257316643e7SJed Brown 1258715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1259715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1260715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1261b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1262715f1b00SHong Zhang 12636f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1264316643e7SJed Brown th->Theta = 0.5; 12651566a47fSLisandro Dalcin th->order = 2; 1266bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1267bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1268bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1269bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1270316643e7SJed Brown PetscFunctionReturn(0); 1271316643e7SJed Brown } 12720de4c49aSJed Brown 12730de4c49aSJed Brown /*@ 12740de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 12750de4c49aSJed Brown 12760de4c49aSJed Brown Not Collective 12770de4c49aSJed Brown 12780de4c49aSJed Brown Input Parameter: 12790de4c49aSJed Brown . ts - timestepping context 12800de4c49aSJed Brown 12810de4c49aSJed Brown Output Parameter: 12820de4c49aSJed Brown . theta - stage abscissa 12830de4c49aSJed Brown 12840de4c49aSJed Brown Note: 12850de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 12860de4c49aSJed Brown 12870de4c49aSJed Brown Level: Advanced 12880de4c49aSJed Brown 12890de4c49aSJed Brown .seealso: TSThetaSetTheta() 12900de4c49aSJed Brown @*/ 12917087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 12920de4c49aSJed Brown { 12934ac538c5SBarry Smith PetscErrorCode ierr; 12940de4c49aSJed Brown 12950de4c49aSJed Brown PetscFunctionBegin; 1296afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12970de4c49aSJed Brown PetscValidPointer(theta,2); 12984ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 12990de4c49aSJed Brown PetscFunctionReturn(0); 13000de4c49aSJed Brown } 13010de4c49aSJed Brown 13020de4c49aSJed Brown /*@ 13030de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 13040de4c49aSJed Brown 13050de4c49aSJed Brown Not Collective 13060de4c49aSJed Brown 13070de4c49aSJed Brown Input Parameter: 13080de4c49aSJed Brown + ts - timestepping context 13090de4c49aSJed Brown - theta - stage abscissa 13100de4c49aSJed Brown 13110de4c49aSJed Brown Options Database: 13120de4c49aSJed Brown . -ts_theta_theta <theta> 13130de4c49aSJed Brown 13140de4c49aSJed Brown Level: Intermediate 13150de4c49aSJed Brown 13160de4c49aSJed Brown .seealso: TSThetaGetTheta() 13170de4c49aSJed Brown @*/ 13187087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 13190de4c49aSJed Brown { 13204ac538c5SBarry Smith PetscErrorCode ierr; 13210de4c49aSJed Brown 13220de4c49aSJed Brown PetscFunctionBegin; 1323afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13244ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 13250de4c49aSJed Brown PetscFunctionReturn(0); 13260de4c49aSJed Brown } 1327f33bbcb6SJed Brown 132826f2ff8fSLisandro Dalcin /*@ 132926f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 133026f2ff8fSLisandro Dalcin 133126f2ff8fSLisandro Dalcin Not Collective 133226f2ff8fSLisandro Dalcin 133326f2ff8fSLisandro Dalcin Input Parameter: 133426f2ff8fSLisandro Dalcin . ts - timestepping context 133526f2ff8fSLisandro Dalcin 133626f2ff8fSLisandro Dalcin Output Parameter: 133726f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 133826f2ff8fSLisandro Dalcin 133926f2ff8fSLisandro Dalcin Level: Advanced 134026f2ff8fSLisandro Dalcin 134126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 134226f2ff8fSLisandro Dalcin @*/ 134326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 134426f2ff8fSLisandro Dalcin { 134526f2ff8fSLisandro Dalcin PetscErrorCode ierr; 134626f2ff8fSLisandro Dalcin 134726f2ff8fSLisandro Dalcin PetscFunctionBegin; 134826f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 134926f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1350163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 135126f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 135226f2ff8fSLisandro Dalcin } 135326f2ff8fSLisandro Dalcin 1354eb284becSJed Brown /*@ 1355eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1356eb284becSJed Brown 1357eb284becSJed Brown Not Collective 1358eb284becSJed Brown 1359eb284becSJed Brown Input Parameter: 1360eb284becSJed Brown + ts - timestepping context 1361eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1362eb284becSJed Brown 1363eb284becSJed Brown Options Database: 1364eb284becSJed Brown . -ts_theta_endpoint <flg> 1365eb284becSJed Brown 1366eb284becSJed Brown Level: Intermediate 1367eb284becSJed Brown 1368eb284becSJed Brown .seealso: TSTHETA, TSCN 1369eb284becSJed Brown @*/ 1370eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1371eb284becSJed Brown { 1372eb284becSJed Brown PetscErrorCode ierr; 1373eb284becSJed Brown 1374eb284becSJed Brown PetscFunctionBegin; 1375eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1376eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1377eb284becSJed Brown PetscFunctionReturn(0); 1378eb284becSJed Brown } 1379eb284becSJed Brown 1380f33bbcb6SJed Brown /* 1381f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1382f33bbcb6SJed Brown * The creation functions for these specializations are below. 1383f33bbcb6SJed Brown */ 1384f33bbcb6SJed Brown 13851566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 13861566a47fSLisandro Dalcin { 13871566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 13881566a47fSLisandro Dalcin PetscErrorCode ierr; 13891566a47fSLisandro Dalcin 13901566a47fSLisandro Dalcin PetscFunctionBegin; 13911566a47fSLisandro 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"); 13921566a47fSLisandro 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"); 13931566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 13941566a47fSLisandro Dalcin PetscFunctionReturn(0); 13951566a47fSLisandro Dalcin } 13961566a47fSLisandro Dalcin 1397f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1398f33bbcb6SJed Brown { 1399f33bbcb6SJed Brown PetscFunctionBegin; 1400f33bbcb6SJed Brown PetscFunctionReturn(0); 1401f33bbcb6SJed Brown } 1402f33bbcb6SJed Brown 1403f33bbcb6SJed Brown /*MC 1404f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1405f33bbcb6SJed Brown 1406f33bbcb6SJed Brown Level: beginner 1407f33bbcb6SJed Brown 14084eb428fdSBarry Smith Notes: 1409c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 14104eb428fdSBarry Smith 14111566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 14124eb428fdSBarry Smith 1413f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1414f33bbcb6SJed Brown 1415f33bbcb6SJed Brown M*/ 14168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1417f33bbcb6SJed Brown { 1418f33bbcb6SJed Brown PetscErrorCode ierr; 1419f33bbcb6SJed Brown 1420f33bbcb6SJed Brown PetscFunctionBegin; 1421f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1422f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 14231566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 14240c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1425f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1426f33bbcb6SJed Brown PetscFunctionReturn(0); 1427f33bbcb6SJed Brown } 1428f33bbcb6SJed Brown 14291566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14301566a47fSLisandro Dalcin { 14311566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14321566a47fSLisandro Dalcin PetscErrorCode ierr; 14331566a47fSLisandro Dalcin 14341566a47fSLisandro Dalcin PetscFunctionBegin; 14351566a47fSLisandro 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"); 14361566a47fSLisandro 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"); 14371566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14381566a47fSLisandro Dalcin PetscFunctionReturn(0); 14391566a47fSLisandro Dalcin } 14401566a47fSLisandro Dalcin 1441f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1442f33bbcb6SJed Brown { 1443f33bbcb6SJed Brown PetscFunctionBegin; 1444f33bbcb6SJed Brown PetscFunctionReturn(0); 1445f33bbcb6SJed Brown } 1446f33bbcb6SJed Brown 1447f33bbcb6SJed Brown /*MC 1448f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1449f33bbcb6SJed Brown 1450f33bbcb6SJed Brown Level: beginner 1451f33bbcb6SJed Brown 1452f33bbcb6SJed Brown Notes: 14537cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14547cf5af47SJed Brown 14557cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1456f33bbcb6SJed Brown 1457f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1458f33bbcb6SJed Brown 1459f33bbcb6SJed Brown M*/ 14608cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1461f33bbcb6SJed Brown { 1462f33bbcb6SJed Brown PetscErrorCode ierr; 1463f33bbcb6SJed Brown 1464f33bbcb6SJed Brown PetscFunctionBegin; 1465f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1466f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1467eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 14680c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1469f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1470f33bbcb6SJed Brown PetscFunctionReturn(0); 1471f33bbcb6SJed Brown } 1472