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; 120fd10508SMatthew Knepley Vec X0,X,Xdot; /* Storage for stage solution, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */ 131566a47fSLisandro Dalcin Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 141566a47fSLisandro Dalcin PetscReal Theta; 151a352fa8SHong Zhang PetscReal shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */ 161566a47fSLisandro Dalcin PetscInt order; 171566a47fSLisandro Dalcin PetscBool endpoint; 181566a47fSLisandro Dalcin PetscBool extrapolate; 19715f1b00SHong Zhang TSStepStatus status; 201a352fa8SHong Zhang Vec VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */ 211a352fa8SHong Zhang PetscReal ptime0; /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */ 221a352fa8SHong Zhang PetscReal time_step0; /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/ 231566a47fSLisandro Dalcin 24715f1b00SHong Zhang /* context for sensitivity analysis */ 25022c081aSHong Zhang PetscInt num_tlm; /* Total number of tangent linear equations */ 26b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 27b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 2813b1b0ffSHong Zhang Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 2913b1b0ffSHong Zhang Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 30b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 3113b1b0ffSHong Zhang Mat MatFwdSensip0; /* backup for roll-backs due to events */ 327207e0fdSHong Zhang Mat MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 337207e0fdSHong Zhang Mat MatIntegralSensip0; /* backup for roll-backs due to events */ 34b5b4017aSHong Zhang Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 35b5b4017aSHong Zhang Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 36b5b4017aSHong Zhang Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 37b5b4017aSHong Zhang Vec *VecsAffine; /* Working vectors to store residuals */ 38715f1b00SHong Zhang /* context for error estimation */ 391566a47fSLisandro Dalcin Vec vec_sol_prev; 401566a47fSLisandro Dalcin Vec vec_lte_work; 41316643e7SJed Brown } TS_Theta; 42316643e7SJed Brown 437445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 447445fe48SJed Brown { 457445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 467445fe48SJed Brown PetscErrorCode ierr; 477445fe48SJed Brown 487445fe48SJed Brown PetscFunctionBegin; 497445fe48SJed Brown if (X0) { 507445fe48SJed Brown if (dm && dm != ts->dm) { 510d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 527445fe48SJed Brown } else *X0 = ts->vec_sol; 537445fe48SJed Brown } 547445fe48SJed Brown if (Xdot) { 557445fe48SJed Brown if (dm && dm != ts->dm) { 560d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 577445fe48SJed Brown } else *Xdot = th->Xdot; 587445fe48SJed Brown } 597445fe48SJed Brown PetscFunctionReturn(0); 607445fe48SJed Brown } 617445fe48SJed Brown 620d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 630d0b770aSPeter Brune { 640d0b770aSPeter Brune PetscErrorCode ierr; 650d0b770aSPeter Brune 660d0b770aSPeter Brune PetscFunctionBegin; 670d0b770aSPeter Brune if (X0) { 680d0b770aSPeter Brune if (dm && dm != ts->dm) { 690d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 700d0b770aSPeter Brune } 710d0b770aSPeter Brune } 720d0b770aSPeter Brune if (Xdot) { 730d0b770aSPeter Brune if (dm && dm != ts->dm) { 740d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 750d0b770aSPeter Brune } 760d0b770aSPeter Brune } 770d0b770aSPeter Brune PetscFunctionReturn(0); 780d0b770aSPeter Brune } 790d0b770aSPeter Brune 807445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 817445fe48SJed Brown { 827445fe48SJed Brown PetscFunctionBegin; 837445fe48SJed Brown PetscFunctionReturn(0); 847445fe48SJed Brown } 857445fe48SJed Brown 867445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 877445fe48SJed Brown { 887445fe48SJed Brown TS ts = (TS)ctx; 897445fe48SJed Brown PetscErrorCode ierr; 907445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 917445fe48SJed Brown 927445fe48SJed Brown PetscFunctionBegin; 937445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 947445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 957445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 967445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 977445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 987445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 990d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 1000d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 1017445fe48SJed Brown PetscFunctionReturn(0); 1027445fe48SJed Brown } 1037445fe48SJed Brown 104258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 105258e1594SPeter Brune { 106258e1594SPeter Brune PetscFunctionBegin; 107258e1594SPeter Brune PetscFunctionReturn(0); 108258e1594SPeter Brune } 109258e1594SPeter Brune 110258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 111258e1594SPeter Brune { 112258e1594SPeter Brune TS ts = (TS)ctx; 113258e1594SPeter Brune PetscErrorCode ierr; 114258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 115258e1594SPeter Brune 116258e1594SPeter Brune PetscFunctionBegin; 117258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 118258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 119258e1594SPeter Brune 120258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122258e1594SPeter Brune 123258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125258e1594SPeter Brune 126258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 127258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 128258e1594SPeter Brune PetscFunctionReturn(0); 129258e1594SPeter Brune } 130258e1594SPeter Brune 131022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 132022c081aSHong Zhang { 133022c081aSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 134cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 135022c081aSHong Zhang PetscErrorCode ierr; 136022c081aSHong Zhang 137022c081aSHong Zhang PetscFunctionBegin; 138022c081aSHong Zhang if (th->endpoint) { 139022c081aSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 140022c081aSHong Zhang if (th->Theta!=1.0) { 1411a352fa8SHong Zhang ierr = TSComputeRHSFunction(quadts,th->ptime0,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1421a352fa8SHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step0*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 143022c081aSHong Zhang } 144cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 1451a352fa8SHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step0*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 146022c081aSHong Zhang } else { 147cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 1481a352fa8SHong Zhang ierr = VecAXPY(quadts->vec_sol,th->time_step0,ts->vec_costintegrand);CHKERRQ(ierr); 149022c081aSHong Zhang } 150022c081aSHong Zhang PetscFunctionReturn(0); 151022c081aSHong Zhang } 152022c081aSHong Zhang 153b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 154b1cb36f3SHong Zhang { 155b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 156cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 157b1cb36f3SHong Zhang PetscErrorCode ierr; 158b1cb36f3SHong Zhang 159b1cb36f3SHong Zhang PetscFunctionBegin; 160b1cb36f3SHong Zhang /* backup cost integral */ 161cd4cee2dSHong Zhang ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr); 162022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 163b1cb36f3SHong Zhang PetscFunctionReturn(0); 164b1cb36f3SHong Zhang } 165b1cb36f3SHong Zhang 166b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 167b1cb36f3SHong Zhang { 1681a352fa8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 169b1cb36f3SHong Zhang PetscErrorCode ierr; 170b1cb36f3SHong Zhang 171b1cb36f3SHong Zhang PetscFunctionBegin; 1721a352fa8SHong Zhang /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */ 1731a352fa8SHong Zhang th->ptime0 = ts->ptime + ts->time_step; 1741a352fa8SHong Zhang th->time_step0 = -ts->time_step; 175022c081aSHong Zhang ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 17624655328SShri PetscFunctionReturn(0); 17724655328SShri } 17824655328SShri 179470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 1801566a47fSLisandro Dalcin { 1811566a47fSLisandro Dalcin PetscInt nits,lits; 1821566a47fSLisandro Dalcin PetscErrorCode ierr; 1831566a47fSLisandro Dalcin 1841566a47fSLisandro Dalcin PetscFunctionBegin; 1851566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1861566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1871566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1881566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1891566a47fSLisandro Dalcin PetscFunctionReturn(0); 1901566a47fSLisandro Dalcin } 1911566a47fSLisandro Dalcin 192193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 193316643e7SJed Brown { 194316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1951566a47fSLisandro Dalcin PetscInt rejections = 0; 1964957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 1971566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 198051f2191SLisandro Dalcin PetscErrorCode ierr; 199316643e7SJed Brown 200316643e7SJed Brown PetscFunctionBegin; 2011566a47fSLisandro Dalcin if (!ts->steprollback) { 2022ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 2033b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 2041566a47fSLisandro Dalcin } 2051566a47fSLisandro Dalcin 2061566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 20799260352SHong Zhang while (!ts->reason && th->status != TS_STEP_COMPLETE) { 2083e05ceb1SHong Zhang th->shift = 1/(th->Theta*ts->time_step); 2091566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 210be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 211fecfb714SLisandro Dalcin if (th->extrapolate && !ts->steprestart) { 2123e05ceb1SHong Zhang ierr = VecAXPY(th->X,1/th->shift,th->Xdot);CHKERRQ(ierr); 213be5899b3SLisandro Dalcin } 214eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 215eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2161566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2171566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2181566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2191566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2201566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 221eb284becSJed Brown } 222be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 223470880abSPatrick Sanan ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 224be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2251566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 226fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 227051f2191SLisandro Dalcin 228051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2291566a47fSLisandro Dalcin if (th->endpoint) { 2301566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2311566a47fSLisandro Dalcin } else { 2323e05ceb1SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr); 2331566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2341566a47fSLisandro Dalcin } 2351566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2361566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2371566a47fSLisandro Dalcin if (!accept) { 2381566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 239be5899b3SLisandro Dalcin ts->time_step = next_time_step; 240be5899b3SLisandro Dalcin goto reject_step; 241051f2191SLisandro Dalcin } 2423b1890cdSShri Abhyankar 243715f1b00SHong Zhang if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 2441a352fa8SHong Zhang th->ptime0 = ts->ptime; 2451a352fa8SHong Zhang th->time_step0 = ts->time_step; 24617145e75SHong Zhang } 2472b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 248cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 249051f2191SLisandro Dalcin break; 250051f2191SLisandro Dalcin 251051f2191SLisandro Dalcin reject_step: 252fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 2531566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 254051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2551566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2563b1890cdSShri Abhyankar } 2573b1890cdSShri Abhyankar } 258316643e7SJed Brown PetscFunctionReturn(0); 259316643e7SJed Brown } 260316643e7SJed Brown 261c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) 262c9681e5dSHong Zhang { 263c9681e5dSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 264cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 265c9681e5dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 266c9681e5dSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 267c9681e5dSHong Zhang PetscInt nadj; 2687207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 269c9681e5dSHong Zhang KSP ksp; 270c9681e5dSHong Zhang PetscScalar *xarr; 271077c4feaSHong Zhang TSEquationType eqtype; 272077c4feaSHong Zhang PetscBool isexplicitode = PETSC_FALSE; 2731a352fa8SHong Zhang PetscReal adjoint_time_step; 274c9681e5dSHong Zhang PetscErrorCode ierr; 275c9681e5dSHong Zhang 276c9681e5dSHong Zhang PetscFunctionBegin; 277077c4feaSHong Zhang ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr); 278077c4feaSHong Zhang if (eqtype == TS_EQ_ODE_EXPLICIT) { 279077c4feaSHong Zhang isexplicitode = PETSC_TRUE; 280077c4feaSHong Zhang VecsDeltaLam = ts->vecs_sensi; 281077c4feaSHong Zhang VecsDeltaLam2 = ts->vecs_sensi2; 282077c4feaSHong Zhang } 283c9681e5dSHong Zhang th->status = TS_STEP_INCOMPLETE; 284c9681e5dSHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 285cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 2867207e0fdSHong Zhang if (quadts) { 287cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 288cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 2897207e0fdSHong Zhang } 290c9681e5dSHong Zhang 2911a352fa8SHong Zhang th->stage_time = ts->ptime; 2921a352fa8SHong Zhang adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 293c9681e5dSHong Zhang 29433f72c5dSHong Zhang /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */ 2957207e0fdSHong Zhang if (quadts) { 296*82ad9101SHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,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); 3011a352fa8SHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./adjoint_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 */ 312c87ba875SIvan Yashchuk th->shift = 1./adjoint_time_step; 313*82ad9101SHong Zhang ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr); 31404f7482eSHong 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; 32329b3c8a1SHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order 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 */ 332*82ad9101SHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 333c9681e5dSHong Zhang /* lambda_s^T F_UP w_2 */ 334*82ad9101SHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,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); 3371a352fa8SHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],1./adjoint_time_step);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; 35029b3c8a1SHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order 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 */ 356077c4feaSHong Zhang if (!isexplicitode) { 3573e05ceb1SHong Zhang th->shift = 0.0; 358*82ad9101SHong Zhang ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr); 35904f7482eSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 360c9681e5dSHong Zhang ierr = MatScale(J,-1.);CHKERRQ(ierr); 361c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 36233f72c5dSHong Zhang /* Add f_U \lambda_s to the original RHS */ 363c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 3641a352fa8SHong Zhang ierr = VecScale(VecsSensiTemp[nadj],adjoint_time_step);CHKERRQ(ierr); 365c9681e5dSHong Zhang ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 366c9681e5dSHong Zhang if (ts->vecs_sensi2) { 367c9681e5dSHong Zhang ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 3681a352fa8SHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],adjoint_time_step);CHKERRQ(ierr); 369c9681e5dSHong Zhang ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 370c9681e5dSHong Zhang } 371c9681e5dSHong Zhang } 372077c4feaSHong Zhang } 373c9681e5dSHong Zhang if (ts->vecs_sensip) { 374*82ad9101SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-1./adjoint_time_step,1.0/adjoint_time_step,0,th->X0,ts->vec_sol);CHKERRQ(ierr); 375*82ad9101SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,1./adjoint_time_step,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */ 3767207e0fdSHong Zhang if (quadts) { 377*82ad9101SHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 3787207e0fdSHong Zhang } 379c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 380c9681e5dSHong Zhang /* lambda_s^T F_PU w_1 */ 381*82ad9101SHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 38233f72c5dSHong Zhang /* lambda_s^T F_PP w_2 */ 383*82ad9101SHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 384c9681e5dSHong Zhang } 385cd4cee2dSHong Zhang 386c9681e5dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 387c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 3881a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 389cd4cee2dSHong Zhang if (quadJp) { 390cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 391cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 3921a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col);CHKERRQ(ierr); 393cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 394cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 39533f72c5dSHong Zhang } 396c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 397c9681e5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 3981a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 399c9681e5dSHong Zhang if (ts->vecs_fpu) { 4001a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 401c9681e5dSHong Zhang } 402c9681e5dSHong Zhang if (ts->vecs_fpp) { 4031a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 404c9681e5dSHong Zhang } 405c9681e5dSHong Zhang } 406c9681e5dSHong Zhang } 407c9681e5dSHong Zhang } 408c9681e5dSHong Zhang 409c9681e5dSHong Zhang if (ts->vecs_sensi2) { 410c9681e5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 411c9681e5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 412c9681e5dSHong Zhang } 413c9681e5dSHong Zhang th->status = TS_STEP_COMPLETE; 414c9681e5dSHong Zhang PetscFunctionReturn(0); 415c9681e5dSHong Zhang } 416c9681e5dSHong Zhang 41742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 4182ca6e920SHong Zhang { 4192ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 420cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 421b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 422b5b4017aSHong Zhang Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 4232ca6e920SHong Zhang PetscInt nadj; 4247207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 4252ca6e920SHong Zhang KSP ksp; 426b5b4017aSHong Zhang PetscScalar *xarr; 4271a352fa8SHong Zhang PetscReal adjoint_time_step; 4281a352fa8SHong Zhang PetscReal adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, ususally ts->ptime is larger than adjoint_ptime) */ 429b5b4017aSHong Zhang PetscErrorCode ierr; 4302ca6e920SHong Zhang 4312ca6e920SHong Zhang PetscFunctionBegin; 432c9681e5dSHong Zhang if (th->Theta == 1.) { 433c9681e5dSHong Zhang ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr); 434c9681e5dSHong Zhang PetscFunctionReturn(0); 435c9681e5dSHong Zhang } 4362ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 437302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 438cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 4397207e0fdSHong Zhang if (quadts) { 440cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 441cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 4427207e0fdSHong Zhang } 443b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 4441a352fa8SHong Zhang th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); 4451a352fa8SHong Zhang adjoint_ptime = ts->ptime + ts->time_step; 4461a352fa8SHong Zhang adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */ 4472ca6e920SHong Zhang 448*82ad9101SHong Zhang if (!th->endpoint) { 449*82ad9101SHong Zhang ierr = VecAXPBYPCZ(th->X,1.0-th->Theta,th->Theta,0,th->X0,ts->vec_sol);CHKERRQ(ierr); 450*82ad9101SHong Zhang } 451*82ad9101SHong Zhang 452b5b4017aSHong Zhang /* Build RHS for first-order adjoint */ 45333f72c5dSHong Zhang /* Cost function has an integral term */ 4547207e0fdSHong Zhang if (quadts) { 45505755b9cSHong Zhang if (th->endpoint) { 456cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 45705755b9cSHong Zhang } else { 458cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 45905755b9cSHong Zhang } 4607207e0fdSHong Zhang } 46133f72c5dSHong Zhang 462abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 463b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 4641a352fa8SHong Zhang ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*adjoint_time_step));CHKERRQ(ierr); 465cd4cee2dSHong Zhang if (quadJ) { 466cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 467cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 468cd4cee2dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr); 469cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 470cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 47136eaed60SHong Zhang } 4722ca6e920SHong Zhang } 4733c54f92cSHong Zhang 474b5b4017aSHong Zhang /* Build LHS for first-order adjoint */ 4751a352fa8SHong Zhang th->shift = 1./(th->Theta*adjoint_time_step); 4763c54f92cSHong Zhang if (th->endpoint) { 47704f7482eSHong Zhang ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr); 4783c54f92cSHong Zhang } else { 47904f7482eSHong Zhang ierr = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr); 4803c54f92cSHong Zhang } 481cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 4822ca6e920SHong Zhang 483b5b4017aSHong Zhang /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 484abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 4851d14d03bSHong Zhang KSPConvergedReason kspreason; 486b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 4871d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 4881d14d03bSHong Zhang if (kspreason < 0) { 4891d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 49029b3c8a1SHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 4911d14d03bSHong Zhang } 4922ca6e920SHong Zhang } 4933c54f92cSHong Zhang 4941d14d03bSHong Zhang /* Second-order adjoint */ 495b5b4017aSHong Zhang if (ts->vecs_sensi2) { /* U_{n+1} */ 496b5b4017aSHong Zhang if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 497b5b4017aSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 498b5b4017aSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 499b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 500b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 501b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 502b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 503b5b4017aSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 504b5b4017aSHong Zhang /* lambda_s^T F_UP w_2 */ 505b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 506b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 507b5b4017aSHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 5083e05ceb1SHong Zhang ierr = VecScale(VecsSensi2Temp[nadj],th->shift);CHKERRQ(ierr); 50933f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 510b5b4017aSHong Zhang if (ts->vecs_fup) { 51133f72c5dSHong Zhang ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 512b5b4017aSHong Zhang } 513b5b4017aSHong Zhang } 514b5b4017aSHong Zhang /* Solve stage equation LHS X = RHS for second-order adjoint */ 515b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 5161d14d03bSHong Zhang KSPConvergedReason kspreason; 5171d14d03bSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 5181d14d03bSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 5191d14d03bSHong Zhang if (kspreason < 0) { 5201d14d03bSHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 52129b3c8a1SHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr); 5221d14d03bSHong Zhang } 523b5b4017aSHong Zhang } 524b5b4017aSHong Zhang } 525b5b4017aSHong Zhang 52636eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 5271d14d03bSHong Zhang if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */ 5281a352fa8SHong Zhang th->shift = 1./((th->Theta-1.)*adjoint_time_step); 5291a352fa8SHong Zhang th->stage_time = adjoint_ptime; 53004f7482eSHong Zhang ierr = TSComputeSNESJacobian(ts,th->X0,J,Jpre);CHKERRQ(ierr); 53104f7482eSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 53233f72c5dSHong Zhang /* R_U at t_n */ 5337207e0fdSHong Zhang if (quadts) { 5341a352fa8SHong Zhang ierr = TSComputeRHSJacobian(quadts,adjoint_ptime,th->X0,quadJ,NULL);CHKERRQ(ierr); 5357207e0fdSHong Zhang } 536abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 537b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 538cd4cee2dSHong Zhang if (quadJ) { 539cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 540cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 541cd4cee2dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr); 542cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 543cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 544b5b4017aSHong Zhang } 5453e05ceb1SHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./th->shift);CHKERRQ(ierr); 546b5b4017aSHong Zhang } 5471d14d03bSHong Zhang 5481d14d03bSHong Zhang /* Second-order adjoint */ 5491d14d03bSHong Zhang if (ts->vecs_sensi2) { /* U_n */ 550b5b4017aSHong Zhang /* Get w1 at t_n from TLM matrix */ 551b5b4017aSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 552b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 553b5b4017aSHong Zhang /* lambda_s^T F_UU w_1 */ 5541a352fa8SHong Zhang ierr = TSComputeIHessianProductFunctionUU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 555b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 55633f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 557b5b4017aSHong Zhang /* lambda_s^T F_UU w_2 */ 5581a352fa8SHong Zhang ierr = TSComputeIHessianProductFunctionUP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 559b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 56033f72c5dSHong 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 */ 561b5b4017aSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 56233f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 563b5b4017aSHong Zhang if (ts->vecs_fup) { 56433f72c5dSHong Zhang ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 565b5b4017aSHong Zhang } 5663e05ceb1SHong Zhang ierr = VecScale(ts->vecs_sensi2[nadj],1./th->shift);CHKERRQ(ierr); 567b5b4017aSHong Zhang } 568b5e0532dSHong Zhang } 5693fd52205SHong Zhang 570c586f2e8SHong Zhang th->stage_time = ts->ptime; /* recover the old value */ 571c586f2e8SHong Zhang 5723fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 573b5b4017aSHong Zhang /* U_{n+1} */ 574*82ad9101SHong Zhang th->shift = 1.0/(adjoint_time_step*th->Theta); 575*82ad9101SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol);CHKERRQ(ierr); 5761a352fa8SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,-1./(th->Theta*adjoint_time_step),ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 5777207e0fdSHong Zhang if (quadts) { 578cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 5797207e0fdSHong Zhang } 580abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 581b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 5821a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 5833b711c3fSHong Zhang if (quadJp) { 5843b711c3fSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 5853b711c3fSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 5863b711c3fSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*th->Theta,ts->vec_drdp_col);CHKERRQ(ierr); 5873b711c3fSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 5883b711c3fSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 5893b711c3fSHong Zhang } 5903fd52205SHong Zhang } 59133f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 59233f72c5dSHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 59333f72c5dSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 59433f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 595b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 596b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 59733f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 59833f72c5dSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 59933f72c5dSHong Zhang 600b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 601b1e111ebSHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 602b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 60333f72c5dSHong 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) */ 604b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 6051a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 606b5b4017aSHong Zhang if (ts->vecs_fpu) { 6071a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 608b5b4017aSHong Zhang } 609b5b4017aSHong Zhang if (ts->vecs_fpp) { 6101a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 611b5b4017aSHong Zhang } 612b5b4017aSHong Zhang } 613b5b4017aSHong Zhang } 614b5b4017aSHong Zhang 615b5b4017aSHong Zhang /* U_s */ 616*82ad9101SHong Zhang ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 6171a352fa8SHong Zhang ierr = TSComputeIJacobianP(ts,adjoint_ptime,th->X0,th->Xdot,1./((th->Theta-1.0)*adjoint_time_step),ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 6187207e0fdSHong Zhang if (quadts) { 6191a352fa8SHong Zhang ierr = TSComputeRHSJacobianP(quadts,adjoint_ptime,th->X0,quadJp);CHKERRQ(ierr); 6207207e0fdSHong Zhang } 621abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 622b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 6231a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 6243b711c3fSHong Zhang if (quadJp) { 6253b711c3fSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 6263b711c3fSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 6273b711c3fSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*(1.0-th->Theta),ts->vec_drdp_col);CHKERRQ(ierr); 6283b711c3fSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 6293b711c3fSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 6303b711c3fSHong Zhang } 63133f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* second-order */ 63233f72c5dSHong Zhang /* Get w1 at t_n from TLM matrix */ 63333f72c5dSHong Zhang ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 63433f72c5dSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 635b5b4017aSHong Zhang /* lambda_s^T F_PU w_1 */ 6361a352fa8SHong Zhang ierr = TSComputeIHessianProductFunctionPU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 63733f72c5dSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 63833f72c5dSHong Zhang ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr); 639b5b4017aSHong Zhang /* lambda_s^T F_PP w_2 */ 6401a352fa8SHong Zhang ierr = TSComputeIHessianProductFunctionPP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 641b5b4017aSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 64233f72c5dSHong 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) */ 643b5b4017aSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 6441a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 645b5b4017aSHong Zhang if (ts->vecs_fpu) { 6461a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 647b5e0532dSHong Zhang } 648b5b4017aSHong Zhang if (ts->vecs_fpp) { 6491a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 650b5b4017aSHong Zhang } 65136eaed60SHong Zhang } 65236eaed60SHong Zhang } 6533fd52205SHong Zhang } 654b5e0532dSHong Zhang } 6553fd52205SHong Zhang } else { /* one-stage case */ 6563e05ceb1SHong Zhang th->shift = 0.0; 65704f7482eSHong Zhang ierr = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr); /* get -f_y */ 65804f7482eSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 6597207e0fdSHong Zhang if (quadts) { 660cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 6617207e0fdSHong Zhang } 662abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 663b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 6641a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-adjoint_time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 665cd4cee2dSHong Zhang if (quadJ) { 666cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr); 667cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr); 6681a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],adjoint_time_step,ts->vec_drdu_col);CHKERRQ(ierr); 669cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr); 670cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr); 67136eaed60SHong Zhang } 6722ca6e920SHong Zhang } 6733fd52205SHong Zhang if (ts->vecs_sensip) { 674*82ad9101SHong Zhang th->shift = 1.0/(adjoint_time_step*th->Theta); 675*82ad9101SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr); 6763e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 6777207e0fdSHong Zhang if (quadts) { 678cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 6797207e0fdSHong Zhang } 680abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 68133f72c5dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 6821a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 683cd4cee2dSHong Zhang if (quadJp) { 684cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr); 685cd4cee2dSHong Zhang ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr); 6861a352fa8SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col);CHKERRQ(ierr); 687cd4cee2dSHong Zhang ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr); 688cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr); 68936eaed60SHong Zhang } 69036eaed60SHong Zhang } 6913fd52205SHong Zhang } 6922ca6e920SHong Zhang } 6932ca6e920SHong Zhang 6942ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 6952ca6e920SHong Zhang PetscFunctionReturn(0); 6962ca6e920SHong Zhang } 6972ca6e920SHong Zhang 698cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 699cd652676SJed Brown { 700cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 701be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 702cd652676SJed Brown PetscErrorCode ierr; 703cd652676SJed Brown 704cd652676SJed Brown PetscFunctionBegin; 705a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 706be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 707be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 708cd652676SJed Brown PetscFunctionReturn(0); 709cd652676SJed Brown } 710cd652676SJed Brown 7111566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 7121566a47fSLisandro Dalcin { 7131566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 7141566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 7151566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 7167453f775SEmil Constantinescu PetscReal wltea,wlter; 7171566a47fSLisandro Dalcin PetscErrorCode ierr; 7181566a47fSLisandro Dalcin 7191566a47fSLisandro Dalcin PetscFunctionBegin; 7202ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 7211566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 722fecfb714SLisandro Dalcin if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 7231566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 724fecfb714SLisandro Dalcin { 725be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 726be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 7271566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 7281566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 7291566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 7301566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 7311566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 7327453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 7331566a47fSLisandro Dalcin } 7341566a47fSLisandro Dalcin if (order) *order = 2; 7351566a47fSLisandro Dalcin PetscFunctionReturn(0); 7361566a47fSLisandro Dalcin } 7371566a47fSLisandro Dalcin 7381566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 7391566a47fSLisandro Dalcin { 7401566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 7417207e0fdSHong Zhang TS quadts = ts->quadraturets; 7421566a47fSLisandro Dalcin PetscErrorCode ierr; 7431566a47fSLisandro Dalcin 7441566a47fSLisandro Dalcin PetscFunctionBegin; 7451566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 746cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 747cd4cee2dSHong Zhang ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr); 7481566a47fSLisandro Dalcin } 749715f1b00SHong Zhang th->status = TS_STEP_INCOMPLETE; 75013b1b0ffSHong Zhang if (ts->mat_sensip) { 75113b1b0ffSHong Zhang ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7526f92202bSHong Zhang } 7537207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7547207e0fdSHong Zhang ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7556f92202bSHong Zhang } 756715f1b00SHong Zhang PetscFunctionReturn(0); 757715f1b00SHong Zhang } 758715f1b00SHong Zhang 759715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts) 760715f1b00SHong Zhang { 761715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 762cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 76313b1b0ffSHong Zhang Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 764b5b4017aSHong Zhang Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 7657207e0fdSHong Zhang PetscInt ntlm; 766715f1b00SHong Zhang KSP ksp; 7677207e0fdSHong Zhang Mat J,Jpre,quadJ = NULL,quadJp = NULL; 76813b1b0ffSHong Zhang PetscScalar *barr,*xarr; 7691a352fa8SHong Zhang PetscReal previous_shift; 770715f1b00SHong Zhang PetscErrorCode ierr; 771715f1b00SHong Zhang 772715f1b00SHong Zhang PetscFunctionBegin; 7731a352fa8SHong Zhang previous_shift = th->shift; 77413b1b0ffSHong Zhang ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 77513b1b0ffSHong Zhang 7767207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 7777207e0fdSHong Zhang ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7786f92202bSHong Zhang } 779715f1b00SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 780cd4cee2dSHong Zhang ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 7817207e0fdSHong Zhang if (quadts) { 782cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr); 783cd4cee2dSHong Zhang ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr); 7847207e0fdSHong Zhang } 785715f1b00SHong Zhang 786715f1b00SHong Zhang /* Build RHS */ 787715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method*/ 7881a352fa8SHong Zhang th->shift = 1./((th->Theta-1.)*th->time_step0); 7891a352fa8SHong Zhang ierr = TSComputeIJacobian(ts,th->ptime0,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 79013b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 79113b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 792715f1b00SHong Zhang 793022c081aSHong Zhang /* Add the f_p forcing terms */ 7940e11c21fSHong Zhang if (ts->Jacp) { 795*82ad9101SHong Zhang ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 7961a352fa8SHong Zhang ierr = TSComputeIJacobianP(ts,th->ptime0,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 79733f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 798*82ad9101SHong Zhang th->shift = previous_shift; 799*82ad9101SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol);CHKERRQ(ierr); 8003e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 80133f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 8020e11c21fSHong Zhang } 803715f1b00SHong Zhang } else { /* 1-stage method */ 8043e05ceb1SHong Zhang th->shift = 0.0; 8053e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 80613b1b0ffSHong Zhang ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 80713b1b0ffSHong Zhang ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 80813b1b0ffSHong Zhang 809022c081aSHong Zhang /* Add the f_p forcing terms */ 8100e11c21fSHong Zhang if (ts->Jacp) { 811*82ad9101SHong Zhang th->shift = previous_shift; 812*82ad9101SHong Zhang ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr); 8133e05ceb1SHong Zhang ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); 81433f72c5dSHong Zhang ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 815715f1b00SHong Zhang } 8160e11c21fSHong Zhang } 817715f1b00SHong Zhang 818715f1b00SHong Zhang /* Build LHS */ 8191a352fa8SHong Zhang th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */ 820715f1b00SHong Zhang if (th->endpoint) { 8213e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 822715f1b00SHong Zhang } else { 8233e05ceb1SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); 824715f1b00SHong Zhang } 825cd4cee2dSHong Zhang ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr); 826715f1b00SHong Zhang 827715f1b00SHong Zhang /* 828715f1b00SHong Zhang Evaluate the first stage of integral gradients with the 2-stage method: 829c9ad9de0SHong Zhang drdu|t_n*S(t_n) + drdp|t_n 830715f1b00SHong Zhang This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 831715f1b00SHong Zhang */ 832715f1b00SHong Zhang if (th->endpoint) { /* 2-stage method only */ 8337207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 8341a352fa8SHong Zhang ierr = TSComputeRHSJacobian(quadts,th->ptime0,th->X0,quadJ,NULL);CHKERRQ(ierr); 8351a352fa8SHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->ptime0,th->X0,quadJp);CHKERRQ(ierr); 8367207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8377207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8381a352fa8SHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step0*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 839715f1b00SHong Zhang } 840715f1b00SHong Zhang } 841715f1b00SHong Zhang 842715f1b00SHong Zhang /* Solve the tangent linear equation for forward sensitivities to parameters */ 843022c081aSHong Zhang for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 844b5b4017aSHong Zhang KSPConvergedReason kspreason; 84513b1b0ffSHong Zhang ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 846b5b4017aSHong Zhang ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 847715f1b00SHong Zhang if (th->endpoint) { 84813b1b0ffSHong Zhang ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 849b5b4017aSHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 850b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 851b5b4017aSHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 85213b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 853715f1b00SHong Zhang } else { 854b5b4017aSHong Zhang ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 855715f1b00SHong Zhang } 856b5b4017aSHong Zhang ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 857b5b4017aSHong Zhang if (kspreason < 0) { 858b5b4017aSHong Zhang ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 859b5b4017aSHong Zhang ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 860b5b4017aSHong Zhang } 861b5b4017aSHong Zhang ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 86213b1b0ffSHong Zhang ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 863715f1b00SHong Zhang } 864715f1b00SHong Zhang 86513b1b0ffSHong Zhang /* 86613b1b0ffSHong Zhang Evaluate the second stage of integral gradients with the 2-stage method: 867c9ad9de0SHong Zhang drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 86813b1b0ffSHong Zhang */ 8697207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 87013b1b0ffSHong Zhang if (!th->endpoint) { 8717207e0fdSHong Zhang ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */ 872cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr); 873cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr); 8747207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8757207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8761a352fa8SHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step0,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 87713b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 87813b1b0ffSHong Zhang } else { 879cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr); 880cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr); 8817207e0fdSHong Zhang ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 8827207e0fdSHong Zhang ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8831a352fa8SHong Zhang ierr = MatAXPY(quadts->mat_sensip,th->time_step0*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 88413b1b0ffSHong Zhang } 88513b1b0ffSHong Zhang } else { 88613b1b0ffSHong Zhang if (!th->endpoint) { 88713b1b0ffSHong Zhang ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 888715f1b00SHong Zhang } 889715f1b00SHong Zhang } 8901566a47fSLisandro Dalcin PetscFunctionReturn(0); 8911566a47fSLisandro Dalcin } 8921566a47fSLisandro Dalcin 8936affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 8946affb6f8SHong Zhang { 8956affb6f8SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 8966affb6f8SHong Zhang 8976affb6f8SHong Zhang PetscFunctionBegin; 8986affb6f8SHong Zhang if (ns) *ns = 1; 8996affb6f8SHong Zhang if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 9006affb6f8SHong Zhang PetscFunctionReturn(0); 9016affb6f8SHong Zhang } 9026affb6f8SHong Zhang 903316643e7SJed Brown /*------------------------------------------------------------*/ 904277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 905316643e7SJed Brown { 906316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 907316643e7SJed Brown PetscErrorCode ierr; 908316643e7SJed Brown 909316643e7SJed Brown PetscFunctionBegin; 9106bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 9116bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 9123b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 913eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 9141566a47fSLisandro Dalcin 9151566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 9161566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 9171566a47fSLisandro Dalcin 918b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 919277b19d0SLisandro Dalcin PetscFunctionReturn(0); 920277b19d0SLisandro Dalcin } 921277b19d0SLisandro Dalcin 922ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts) 923ece46509SHong Zhang { 924ece46509SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 925ece46509SHong Zhang PetscErrorCode ierr; 926ece46509SHong Zhang 927ece46509SHong Zhang PetscFunctionBegin; 928ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 929ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 930ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 931ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 932ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 933ece46509SHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 934ece46509SHong Zhang PetscFunctionReturn(0); 935ece46509SHong Zhang } 936ece46509SHong Zhang 937277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 938277b19d0SLisandro Dalcin { 939277b19d0SLisandro Dalcin PetscErrorCode ierr; 940277b19d0SLisandro Dalcin 941277b19d0SLisandro Dalcin PetscFunctionBegin; 942277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 943b3a6b972SJed Brown if (ts->dm) { 944b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 945b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 946b3a6b972SJed Brown } 947277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 948bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 949bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 950bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 951bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 952316643e7SJed Brown PetscFunctionReturn(0); 953316643e7SJed Brown } 954316643e7SJed Brown 955316643e7SJed Brown /* 956316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 9572b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 9580fd10508SMatthew Knepley 9590fd10508SMatthew Knepley Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true, 9600fd10508SMatthew Knepley otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is 9610fd10508SMatthew Knepley U = (U_{n+1} + U0)/2 962316643e7SJed Brown */ 9630f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 964316643e7SJed Brown { 965316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 966316643e7SJed Brown PetscErrorCode ierr; 9677445fe48SJed Brown Vec X0,Xdot; 9687445fe48SJed Brown DM dm,dmsave; 9693e05ceb1SHong Zhang PetscReal shift = th->shift; 970316643e7SJed Brown 971316643e7SJed Brown PetscFunctionBegin; 9727445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 9735a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 9747445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 975494ce9fcSHong Zhang if (x != X0) { 976b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 977494ce9fcSHong Zhang } else { 978494ce9fcSHong Zhang ierr = VecZeroEntries(Xdot);CHKERRQ(ierr); 979494ce9fcSHong Zhang } 9807445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 9817445fe48SJed Brown dmsave = ts->dm; 9827445fe48SJed Brown ts->dm = dm; 9837445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 9847445fe48SJed Brown ts->dm = dmsave; 9850d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 986316643e7SJed Brown PetscFunctionReturn(0); 987316643e7SJed Brown } 988316643e7SJed Brown 989d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 990316643e7SJed Brown { 991316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 992316643e7SJed Brown PetscErrorCode ierr; 9937445fe48SJed Brown Vec Xdot; 9947445fe48SJed Brown DM dm,dmsave; 9953e05ceb1SHong Zhang PetscReal shift = th->shift; 996316643e7SJed Brown 997316643e7SJed Brown PetscFunctionBegin; 9987445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 999be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 10000298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 10017445fe48SJed Brown 10027445fe48SJed Brown dmsave = ts->dm; 10037445fe48SJed Brown ts->dm = dm; 1004d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 10057445fe48SJed Brown ts->dm = dmsave; 10060298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 1007316643e7SJed Brown PetscFunctionReturn(0); 1008316643e7SJed Brown } 1009316643e7SJed Brown 1010715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts) 1011715f1b00SHong Zhang { 1012715f1b00SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 10137207e0fdSHong Zhang TS quadts = ts->quadraturets; 1014715f1b00SHong Zhang PetscErrorCode ierr; 1015715f1b00SHong Zhang 1016715f1b00SHong Zhang PetscFunctionBegin; 1017022c081aSHong Zhang /* combine sensitivities to parameters and sensitivities to initial values into one array */ 101813b1b0ffSHong Zhang th->num_tlm = ts->num_parameters; 101913b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 10207207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 10217207e0fdSHong Zhang ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr); 10227207e0fdSHong Zhang ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr); 1023715f1b00SHong Zhang } 10246f92202bSHong Zhang /* backup sensitivity results for roll-backs */ 102513b1b0ffSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 102613b1b0ffSHong Zhang 1027b5b4017aSHong Zhang ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 1028715f1b00SHong Zhang PetscFunctionReturn(0); 1029715f1b00SHong Zhang } 1030715f1b00SHong Zhang 10317adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts) 10327adebddeSHong Zhang { 10337adebddeSHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 10347207e0fdSHong Zhang TS quadts = ts->quadraturets; 10357adebddeSHong Zhang PetscErrorCode ierr; 10367adebddeSHong Zhang 10377adebddeSHong Zhang PetscFunctionBegin; 10387207e0fdSHong Zhang if (quadts && quadts->mat_sensip) { 10397207e0fdSHong Zhang ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr); 10407207e0fdSHong Zhang ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr); 10417adebddeSHong Zhang } 10427adebddeSHong Zhang ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 10437adebddeSHong Zhang ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 10447adebddeSHong Zhang ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 10457adebddeSHong Zhang PetscFunctionReturn(0); 10467adebddeSHong Zhang } 10477adebddeSHong Zhang 1048316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 1049316643e7SJed Brown { 1050316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1051cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 10522ffb9264SLisandro Dalcin PetscBool match; 1053316643e7SJed Brown PetscErrorCode ierr; 1054316643e7SJed Brown 1055316643e7SJed Brown PetscFunctionBegin; 1056cd4cee2dSHong Zhang if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */ 1057cd4cee2dSHong Zhang ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr); 1058b1cb36f3SHong Zhang } 105939d13666SHong Zhang if (!th->X) { 1060316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 106139d13666SHong Zhang } 106239d13666SHong Zhang if (!th->Xdot) { 1063316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 106439d13666SHong Zhang } 106539d13666SHong Zhang if (!th->X0) { 10663b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 106739d13666SHong Zhang } 10681566a47fSLisandro Dalcin if (th->endpoint) { 10691566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 10707445fe48SJed Brown } 10713b1890cdSShri Abhyankar 10721566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 107320d49056SMatthew G. Knepley th->shift = 1/(th->Theta*ts->time_step); 10741566a47fSLisandro Dalcin 10751566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 10761566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10771566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 10781566a47fSLisandro Dalcin 10791566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 10801566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 10812ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 10822ffb9264SLisandro Dalcin if (!match) { 10831566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 10841566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 10853b1890cdSShri Abhyankar } 10861566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 1087b4dd3bf9SBarry Smith PetscFunctionReturn(0); 1088b4dd3bf9SBarry Smith } 10890c7376e5SHong Zhang 1090b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 1091b4dd3bf9SBarry Smith 1092b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 1093b4dd3bf9SBarry Smith { 1094b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 1095b4dd3bf9SBarry Smith PetscErrorCode ierr; 1096b4dd3bf9SBarry Smith 1097b4dd3bf9SBarry Smith PetscFunctionBegin; 1098b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 1099b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 11002ca6e920SHong Zhang if (ts->vecs_sensip) { 1101b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 11022ca6e920SHong Zhang } 1103b5b4017aSHong Zhang if (ts->vecs_sensi2) { 1104b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 1105b5b4017aSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 110667633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 110767633408SHong Zhang if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu; 110867633408SHong Zhang if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup; 1109b5b4017aSHong Zhang } 1110c9681e5dSHong Zhang if (ts->vecs_sensi2p) { 1111c9681e5dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 111267633408SHong Zhang /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */ 111367633408SHong Zhang if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu; 111467633408SHong Zhang if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp; 1115b5b4017aSHong Zhang } 1116316643e7SJed Brown PetscFunctionReturn(0); 1117316643e7SJed Brown } 1118316643e7SJed Brown 11194416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 1120316643e7SJed Brown { 1121316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1122316643e7SJed Brown PetscErrorCode ierr; 1123316643e7SJed Brown 1124316643e7SJed Brown PetscFunctionBegin; 1125e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 1126316643e7SJed Brown { 11270298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 11280298fd71SBarry 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); 112903842d09SLisandro 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); 1130316643e7SJed Brown } 1131316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1132316643e7SJed Brown PetscFunctionReturn(0); 1133316643e7SJed Brown } 1134316643e7SJed Brown 1135316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 1136316643e7SJed Brown { 1137316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1138ace3abfcSBarry Smith PetscBool iascii; 1139316643e7SJed Brown PetscErrorCode ierr; 1140316643e7SJed Brown 1141316643e7SJed Brown PetscFunctionBegin; 1142251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1143316643e7SJed Brown if (iascii) { 11447c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 1145ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1146316643e7SJed Brown } 1147316643e7SJed Brown PetscFunctionReturn(0); 1148316643e7SJed Brown } 1149316643e7SJed Brown 1150560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 11510de4c49aSJed Brown { 11520de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11530de4c49aSJed Brown 11540de4c49aSJed Brown PetscFunctionBegin; 11550de4c49aSJed Brown *theta = th->Theta; 11560de4c49aSJed Brown PetscFunctionReturn(0); 11570de4c49aSJed Brown } 11580de4c49aSJed Brown 1159560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 11600de4c49aSJed Brown { 11610de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 11620de4c49aSJed Brown 11630de4c49aSJed Brown PetscFunctionBegin; 11647c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 11650de4c49aSJed Brown th->Theta = theta; 11661566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 11670de4c49aSJed Brown PetscFunctionReturn(0); 11680de4c49aSJed Brown } 1169eb284becSJed Brown 1170560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 117126f2ff8fSLisandro Dalcin { 117226f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 117326f2ff8fSLisandro Dalcin 117426f2ff8fSLisandro Dalcin PetscFunctionBegin; 117526f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 117626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 117726f2ff8fSLisandro Dalcin } 117826f2ff8fSLisandro Dalcin 1179560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 1180eb284becSJed Brown { 1181eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1182eb284becSJed Brown 1183eb284becSJed Brown PetscFunctionBegin; 1184eb284becSJed Brown th->endpoint = flg; 1185eb284becSJed Brown PetscFunctionReturn(0); 1186eb284becSJed Brown } 11870de4c49aSJed Brown 1188f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1189f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 1190f9c1d6abSBarry Smith { 1191f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 1192f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 11933fd8ae06SJed Brown const PetscReal one = 1.0; 1194f9c1d6abSBarry Smith 1195f9c1d6abSBarry Smith PetscFunctionBegin; 11963fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 1197f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 1198f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 1199f9c1d6abSBarry Smith PetscFunctionReturn(0); 1200f9c1d6abSBarry Smith } 1201f9c1d6abSBarry Smith #endif 1202f9c1d6abSBarry Smith 120342682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 120442682096SHong Zhang { 120542682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 120642682096SHong Zhang 120742682096SHong Zhang PetscFunctionBegin; 12081566a47fSLisandro Dalcin if (ns) *ns = 1; 1209*82ad9101SHong Zhang if (Y) *Y = &(th->X0); 121042682096SHong Zhang PetscFunctionReturn(0); 121142682096SHong Zhang } 1212f9c1d6abSBarry Smith 1213316643e7SJed Brown /* ------------------------------------------------------------ */ 1214316643e7SJed Brown /*MC 121596f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 1216316643e7SJed Brown 1217316643e7SJed Brown Level: beginner 1218316643e7SJed Brown 12194eb428fdSBarry Smith Options Database: 1220c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 1221c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 122203842d09SLisandro Dalcin - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 12234eb428fdSBarry Smith 1224eb284becSJed Brown Notes: 12250c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 12260c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 12274eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 12284eb428fdSBarry Smith 1229eb284becSJed Brown This method can be applied to DAE. 1230eb284becSJed Brown 1231eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 1232eb284becSJed Brown 1233eb284becSJed Brown .vb 1234eb284becSJed Brown Theta | Theta 1235eb284becSJed Brown ------------- 1236eb284becSJed Brown | 1 1237eb284becSJed Brown .ve 1238eb284becSJed Brown 1239eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 1240eb284becSJed Brown 1241eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1242eb284becSJed Brown 1243eb284becSJed Brown .vb 1244eb284becSJed Brown 0 | 0 0 1245eb284becSJed Brown 1 | 1-Theta Theta 1246eb284becSJed Brown ------------------- 1247eb284becSJed Brown | 1-Theta Theta 1248eb284becSJed Brown .ve 1249eb284becSJed Brown 1250eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1251eb284becSJed Brown 1252eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 1253eb284becSJed Brown 1254eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 1255eb284becSJed Brown 12564eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1257eb284becSJed Brown 1258eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1259316643e7SJed Brown 1260316643e7SJed Brown M*/ 12618cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1262316643e7SJed Brown { 1263316643e7SJed Brown TS_Theta *th; 1264316643e7SJed Brown PetscErrorCode ierr; 1265316643e7SJed Brown 1266316643e7SJed Brown PetscFunctionBegin; 1267277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 1268ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1269316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 1270316643e7SJed Brown ts->ops->view = TSView_Theta; 1271316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 127242f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1273ece46509SHong Zhang ts->ops->adjointreset = TSAdjointReset_Theta; 1274316643e7SJed Brown ts->ops->step = TSStep_Theta; 1275cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 12761566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 127724655328SShri ts->ops->rollback = TSRollBack_Theta; 1278316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 12790f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 12800f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1281f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 1282f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 1283f9c1d6abSBarry Smith #endif 128442682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 128542f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 1286b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1287b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 12882ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 12896affb6f8SHong Zhang 1290715f1b00SHong Zhang ts->ops->forwardsetup = TSForwardSetUp_Theta; 12917adebddeSHong Zhang ts->ops->forwardreset = TSForwardReset_Theta; 1292715f1b00SHong Zhang ts->ops->forwardstep = TSForwardStep_Theta; 12936affb6f8SHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1294316643e7SJed Brown 1295efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1296efd4aadfSBarry Smith 1297b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1298316643e7SJed Brown ts->data = (void*)th; 1299316643e7SJed Brown 1300715f1b00SHong Zhang th->VecsDeltaLam = NULL; 1301715f1b00SHong Zhang th->VecsDeltaMu = NULL; 1302715f1b00SHong Zhang th->VecsSensiTemp = NULL; 1303b5b4017aSHong Zhang th->VecsSensi2Temp = NULL; 1304715f1b00SHong Zhang 13056f700aefSJed Brown th->extrapolate = PETSC_FALSE; 1306316643e7SJed Brown th->Theta = 0.5; 13071566a47fSLisandro Dalcin th->order = 2; 1308bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1309bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1310bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1311bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1312316643e7SJed Brown PetscFunctionReturn(0); 1313316643e7SJed Brown } 13140de4c49aSJed Brown 13150de4c49aSJed Brown /*@ 13160de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 13170de4c49aSJed Brown 13180de4c49aSJed Brown Not Collective 13190de4c49aSJed Brown 13200de4c49aSJed Brown Input Parameter: 13210de4c49aSJed Brown . ts - timestepping context 13220de4c49aSJed Brown 13230de4c49aSJed Brown Output Parameter: 13240de4c49aSJed Brown . theta - stage abscissa 13250de4c49aSJed Brown 13260de4c49aSJed Brown Note: 13270de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 13280de4c49aSJed Brown 13290de4c49aSJed Brown Level: Advanced 13300de4c49aSJed Brown 13310de4c49aSJed Brown .seealso: TSThetaSetTheta() 13320de4c49aSJed Brown @*/ 13337087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 13340de4c49aSJed Brown { 13354ac538c5SBarry Smith PetscErrorCode ierr; 13360de4c49aSJed Brown 13370de4c49aSJed Brown PetscFunctionBegin; 1338afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13390de4c49aSJed Brown PetscValidPointer(theta,2); 13404ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 13410de4c49aSJed Brown PetscFunctionReturn(0); 13420de4c49aSJed Brown } 13430de4c49aSJed Brown 13440de4c49aSJed Brown /*@ 13450de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 13460de4c49aSJed Brown 13470de4c49aSJed Brown Not Collective 13480de4c49aSJed Brown 13490de4c49aSJed Brown Input Parameter: 13500de4c49aSJed Brown + ts - timestepping context 13510de4c49aSJed Brown - theta - stage abscissa 13520de4c49aSJed Brown 13530de4c49aSJed Brown Options Database: 13540de4c49aSJed Brown . -ts_theta_theta <theta> 13550de4c49aSJed Brown 13560de4c49aSJed Brown Level: Intermediate 13570de4c49aSJed Brown 13580de4c49aSJed Brown .seealso: TSThetaGetTheta() 13590de4c49aSJed Brown @*/ 13607087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 13610de4c49aSJed Brown { 13624ac538c5SBarry Smith PetscErrorCode ierr; 13630de4c49aSJed Brown 13640de4c49aSJed Brown PetscFunctionBegin; 1365afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13664ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 13670de4c49aSJed Brown PetscFunctionReturn(0); 13680de4c49aSJed Brown } 1369f33bbcb6SJed Brown 137026f2ff8fSLisandro Dalcin /*@ 137126f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 137226f2ff8fSLisandro Dalcin 137326f2ff8fSLisandro Dalcin Not Collective 137426f2ff8fSLisandro Dalcin 137526f2ff8fSLisandro Dalcin Input Parameter: 137626f2ff8fSLisandro Dalcin . ts - timestepping context 137726f2ff8fSLisandro Dalcin 137826f2ff8fSLisandro Dalcin Output Parameter: 137926f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 138026f2ff8fSLisandro Dalcin 138126f2ff8fSLisandro Dalcin Level: Advanced 138226f2ff8fSLisandro Dalcin 138326f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 138426f2ff8fSLisandro Dalcin @*/ 138526f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 138626f2ff8fSLisandro Dalcin { 138726f2ff8fSLisandro Dalcin PetscErrorCode ierr; 138826f2ff8fSLisandro Dalcin 138926f2ff8fSLisandro Dalcin PetscFunctionBegin; 139026f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 139126f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 1392163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 139326f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 139426f2ff8fSLisandro Dalcin } 139526f2ff8fSLisandro Dalcin 1396eb284becSJed Brown /*@ 1397eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1398eb284becSJed Brown 1399eb284becSJed Brown Not Collective 1400eb284becSJed Brown 1401eb284becSJed Brown Input Parameter: 1402eb284becSJed Brown + ts - timestepping context 1403eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 1404eb284becSJed Brown 1405eb284becSJed Brown Options Database: 1406eb284becSJed Brown . -ts_theta_endpoint <flg> 1407eb284becSJed Brown 1408eb284becSJed Brown Level: Intermediate 1409eb284becSJed Brown 1410eb284becSJed Brown .seealso: TSTHETA, TSCN 1411eb284becSJed Brown @*/ 1412eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1413eb284becSJed Brown { 1414eb284becSJed Brown PetscErrorCode ierr; 1415eb284becSJed Brown 1416eb284becSJed Brown PetscFunctionBegin; 1417eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1418eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1419eb284becSJed Brown PetscFunctionReturn(0); 1420eb284becSJed Brown } 1421eb284becSJed Brown 1422f33bbcb6SJed Brown /* 1423f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1424f33bbcb6SJed Brown * The creation functions for these specializations are below. 1425f33bbcb6SJed Brown */ 1426f33bbcb6SJed Brown 14271566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 14281566a47fSLisandro Dalcin { 14291566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14301566a47fSLisandro Dalcin PetscErrorCode ierr; 14311566a47fSLisandro Dalcin 14321566a47fSLisandro Dalcin PetscFunctionBegin; 14331566a47fSLisandro 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"); 14341566a47fSLisandro 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"); 14351566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14361566a47fSLisandro Dalcin PetscFunctionReturn(0); 14371566a47fSLisandro Dalcin } 14381566a47fSLisandro Dalcin 1439f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1440f33bbcb6SJed Brown { 1441f33bbcb6SJed Brown PetscFunctionBegin; 1442f33bbcb6SJed Brown PetscFunctionReturn(0); 1443f33bbcb6SJed Brown } 1444f33bbcb6SJed Brown 1445f33bbcb6SJed Brown /*MC 1446f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 1447f33bbcb6SJed Brown 1448f33bbcb6SJed Brown Level: beginner 1449f33bbcb6SJed Brown 14504eb428fdSBarry Smith Notes: 1451c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 14524eb428fdSBarry Smith 14531566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 14544eb428fdSBarry Smith 1455f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1456f33bbcb6SJed Brown 1457f33bbcb6SJed Brown M*/ 14588cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1459f33bbcb6SJed Brown { 1460f33bbcb6SJed Brown PetscErrorCode ierr; 1461f33bbcb6SJed Brown 1462f33bbcb6SJed Brown PetscFunctionBegin; 1463f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1464f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 14651566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 14660c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1467f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1468f33bbcb6SJed Brown PetscFunctionReturn(0); 1469f33bbcb6SJed Brown } 1470f33bbcb6SJed Brown 14711566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 14721566a47fSLisandro Dalcin { 14731566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 14741566a47fSLisandro Dalcin PetscErrorCode ierr; 14751566a47fSLisandro Dalcin 14761566a47fSLisandro Dalcin PetscFunctionBegin; 14771566a47fSLisandro 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"); 14781566a47fSLisandro 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"); 14791566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 14801566a47fSLisandro Dalcin PetscFunctionReturn(0); 14811566a47fSLisandro Dalcin } 14821566a47fSLisandro Dalcin 1483f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1484f33bbcb6SJed Brown { 1485f33bbcb6SJed Brown PetscFunctionBegin; 1486f33bbcb6SJed Brown PetscFunctionReturn(0); 1487f33bbcb6SJed Brown } 1488f33bbcb6SJed Brown 1489f33bbcb6SJed Brown /*MC 1490f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1491f33bbcb6SJed Brown 1492f33bbcb6SJed Brown Level: beginner 1493f33bbcb6SJed Brown 1494f33bbcb6SJed Brown Notes: 14957cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 14967cf5af47SJed Brown 14977cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1498f33bbcb6SJed Brown 1499f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1500f33bbcb6SJed Brown 1501f33bbcb6SJed Brown M*/ 15028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1503f33bbcb6SJed Brown { 1504f33bbcb6SJed Brown PetscErrorCode ierr; 1505f33bbcb6SJed Brown 1506f33bbcb6SJed Brown PetscFunctionBegin; 1507f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1508f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1509eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 15100c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1511f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1512f33bbcb6SJed Brown PetscFunctionReturn(0); 1513f33bbcb6SJed Brown } 1514