1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5a055c8caSBarry Smith #include <petscsnes.h> 61e25c274SJed Brown #include <petscdm.h> 73c54f92cSHong Zhang #include <petscmat.h> 8316643e7SJed Brown 9316643e7SJed Brown typedef struct { 101566a47fSLisandro Dalcin PetscReal stage_time; 111566a47fSLisandro Dalcin Vec X0,X,Xdot; /* Storage for stages and time derivative */ 121566a47fSLisandro Dalcin Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 131566a47fSLisandro Dalcin 141566a47fSLisandro Dalcin PetscReal Theta; 151566a47fSLisandro Dalcin PetscInt order; 161566a47fSLisandro Dalcin PetscBool endpoint; 171566a47fSLisandro Dalcin PetscBool extrapolate; 181566a47fSLisandro Dalcin 19b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 20b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 21b5e0532dSHong Zhang Vec *VecsSensiTemp; /* Vector to be timed with Jacobian transpose */ 221566a47fSLisandro Dalcin Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 23b5e0532dSHong Zhang PetscReal ptime; 24b1cb36f3SHong Zhang PetscReal time_step; 251566a47fSLisandro Dalcin 261566a47fSLisandro Dalcin PetscBool adapt; /* Use time-step adaptivity ? */ 271566a47fSLisandro Dalcin Vec vec_sol_prev; 281566a47fSLisandro Dalcin Vec vec_lte_work; 291566a47fSLisandro Dalcin 301566a47fSLisandro Dalcin TSStepStatus status; 31316643e7SJed Brown } TS_Theta; 32316643e7SJed Brown 33316643e7SJed Brown #undef __FUNCT__ 347445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 357445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 367445fe48SJed Brown { 377445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 387445fe48SJed Brown PetscErrorCode ierr; 397445fe48SJed Brown 407445fe48SJed Brown PetscFunctionBegin; 417445fe48SJed Brown if (X0) { 427445fe48SJed Brown if (dm && dm != ts->dm) { 430d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 447445fe48SJed Brown } else *X0 = ts->vec_sol; 457445fe48SJed Brown } 467445fe48SJed Brown if (Xdot) { 477445fe48SJed Brown if (dm && dm != ts->dm) { 480d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 497445fe48SJed Brown } else *Xdot = th->Xdot; 507445fe48SJed Brown } 517445fe48SJed Brown PetscFunctionReturn(0); 527445fe48SJed Brown } 537445fe48SJed Brown 540d0b770aSPeter Brune #undef __FUNCT__ 550d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 560d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 570d0b770aSPeter Brune { 580d0b770aSPeter Brune PetscErrorCode ierr; 590d0b770aSPeter Brune 600d0b770aSPeter Brune PetscFunctionBegin; 610d0b770aSPeter Brune if (X0) { 620d0b770aSPeter Brune if (dm && dm != ts->dm) { 630d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 640d0b770aSPeter Brune } 650d0b770aSPeter Brune } 660d0b770aSPeter Brune if (Xdot) { 670d0b770aSPeter Brune if (dm && dm != ts->dm) { 680d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 690d0b770aSPeter Brune } 700d0b770aSPeter Brune } 710d0b770aSPeter Brune PetscFunctionReturn(0); 720d0b770aSPeter Brune } 730d0b770aSPeter Brune 747445fe48SJed Brown #undef __FUNCT__ 757445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 767445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 777445fe48SJed Brown { 787445fe48SJed Brown 797445fe48SJed Brown PetscFunctionBegin; 807445fe48SJed Brown PetscFunctionReturn(0); 817445fe48SJed Brown } 827445fe48SJed Brown 837445fe48SJed Brown #undef __FUNCT__ 847445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 857445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 867445fe48SJed Brown { 877445fe48SJed Brown TS ts = (TS)ctx; 887445fe48SJed Brown PetscErrorCode ierr; 897445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 907445fe48SJed Brown 917445fe48SJed Brown PetscFunctionBegin; 927445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 937445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 947445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 957445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 967445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 977445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 980d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 990d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 1007445fe48SJed Brown PetscFunctionReturn(0); 1017445fe48SJed Brown } 1027445fe48SJed Brown 1037445fe48SJed Brown #undef __FUNCT__ 104258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta" 105258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 106258e1594SPeter Brune { 107258e1594SPeter Brune 108258e1594SPeter Brune PetscFunctionBegin; 109258e1594SPeter Brune PetscFunctionReturn(0); 110258e1594SPeter Brune } 111258e1594SPeter Brune 112258e1594SPeter Brune #undef __FUNCT__ 113258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta" 114258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 115258e1594SPeter Brune { 116258e1594SPeter Brune TS ts = (TS)ctx; 117258e1594SPeter Brune PetscErrorCode ierr; 118258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 119258e1594SPeter Brune 120258e1594SPeter Brune PetscFunctionBegin; 121258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 122258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 123258e1594SPeter Brune 124258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126258e1594SPeter Brune 127258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 128258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 129258e1594SPeter Brune 130258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 131258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 132258e1594SPeter Brune PetscFunctionReturn(0); 133258e1594SPeter Brune } 134258e1594SPeter Brune 1353b1890cdSShri Abhyankar #undef __FUNCT__ 136b1cb36f3SHong Zhang #define __FUNCT__ "TSForwardCostIntegral_Theta" 137b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 138b1cb36f3SHong Zhang { 139b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 140b1cb36f3SHong Zhang PetscErrorCode ierr; 141b1cb36f3SHong Zhang 142b1cb36f3SHong Zhang PetscFunctionBegin; 143b1cb36f3SHong Zhang /* backup cost integral */ 144b1cb36f3SHong Zhang ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 145b1cb36f3SHong Zhang if (th->endpoint) { 146b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1471566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 148b1cb36f3SHong Zhang } 149b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 150b1cb36f3SHong Zhang if (th->endpoint) { 151b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 152b1cb36f3SHong Zhang } else { 153b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 154b1cb36f3SHong Zhang } 155b1cb36f3SHong Zhang PetscFunctionReturn(0); 156b1cb36f3SHong Zhang } 157b1cb36f3SHong Zhang 158b1cb36f3SHong Zhang #undef __FUNCT__ 159b1cb36f3SHong Zhang #define __FUNCT__ "TSAdjointCostIntegral_Theta" 160b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 161b1cb36f3SHong Zhang { 162b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 163b1cb36f3SHong Zhang PetscErrorCode ierr; 164b1cb36f3SHong Zhang 165b1cb36f3SHong Zhang PetscFunctionBegin; 166b1cb36f3SHong Zhang if (th->endpoint) { 167b1cb36f3SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 168b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 169b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 170b1cb36f3SHong Zhang if (th->Theta!=1) { 171b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1721566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr); 173b1cb36f3SHong Zhang } 174b1cb36f3SHong Zhang }else { 175b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 176b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 177b1cb36f3SHong Zhang } 17824655328SShri PetscFunctionReturn(0); 17924655328SShri } 18024655328SShri 181*be5899b3SLisandro Dalcin #define TSEvent_Status(ts) (ts->event ? ts->event->status : TSEVENT_NONE) 182*be5899b3SLisandro Dalcin 18324655328SShri #undef __FUNCT__ 1841566a47fSLisandro Dalcin #define __FUNCT__ "TS_SNESSolve" 1851566a47fSLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 1861566a47fSLisandro Dalcin { 1871566a47fSLisandro Dalcin PetscInt nits,lits; 1881566a47fSLisandro Dalcin PetscErrorCode ierr; 1891566a47fSLisandro Dalcin 1901566a47fSLisandro Dalcin PetscFunctionBegin; 1911566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1921566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1931566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1941566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1951566a47fSLisandro Dalcin PetscFunctionReturn(0); 1961566a47fSLisandro Dalcin } 1971566a47fSLisandro Dalcin 1981566a47fSLisandro Dalcin #undef __FUNCT__ 199316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 200193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 201316643e7SJed Brown { 202316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2031566a47fSLisandro Dalcin PetscInt rejections = 0; 2044957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 2051566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 206051f2191SLisandro Dalcin PetscErrorCode ierr; 207316643e7SJed Brown 208316643e7SJed Brown PetscFunctionBegin; 2091566a47fSLisandro Dalcin if (!ts->steprollback) { 2101566a47fSLisandro Dalcin if (th->adapt) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 2113b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 2121566a47fSLisandro Dalcin } 2131566a47fSLisandro Dalcin 2141566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 2151566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 2161566a47fSLisandro Dalcin 2171566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 2181566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 219316643e7SJed Brown 220*be5899b3SLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 221*be5899b3SLisandro Dalcin if (th->extrapolate && ts->steps > 0 && TSEvent_Status(ts) != TSEVENT_RESET_NEXTSTEP) { 222*be5899b3SLisandro Dalcin ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 223*be5899b3SLisandro Dalcin } 224eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 225eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2261566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2271566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2281566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2291566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2301566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 231eb284becSJed Brown } 232*be5899b3SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 2331566a47fSLisandro Dalcin ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 234*be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 2351566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 236051f2191SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 237051f2191SLisandro Dalcin 238051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2391566a47fSLisandro Dalcin if (th->endpoint) { 2401566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2411566a47fSLisandro Dalcin } else { 242cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2431566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2441566a47fSLisandro Dalcin } 2451566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2461566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2471566a47fSLisandro Dalcin if (!accept) { 2481566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 249*be5899b3SLisandro Dalcin ts->time_step = next_time_step; 250*be5899b3SLisandro Dalcin goto reject_step; 251051f2191SLisandro Dalcin } 2523b1890cdSShri Abhyankar 253b1cb36f3SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 254b1cb36f3SHong Zhang th->ptime = ts->ptime; 255b1cb36f3SHong Zhang th->time_step = ts->time_step; 25617145e75SHong Zhang } 2571566a47fSLisandro Dalcin 2582b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 259cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 260051f2191SLisandro Dalcin break; 261051f2191SLisandro Dalcin 262051f2191SLisandro Dalcin reject_step: 2631566a47fSLisandro Dalcin ts->reject++; 2641566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 265051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2661566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2673b1890cdSShri Abhyankar } 2683b1890cdSShri Abhyankar } 269316643e7SJed Brown PetscFunctionReturn(0); 270316643e7SJed Brown } 271316643e7SJed Brown 272cd652676SJed Brown #undef __FUNCT__ 27342f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_Theta" 27442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2752ca6e920SHong Zhang { 2762ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 277b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 2782ca6e920SHong Zhang PetscInt nadj; 2792ca6e920SHong Zhang PetscErrorCode ierr; 2802ca6e920SHong Zhang Mat J,Jp; 2812ca6e920SHong Zhang KSP ksp; 2822ca6e920SHong Zhang PetscReal shift; 2832ca6e920SHong Zhang 2842ca6e920SHong Zhang PetscFunctionBegin; 2852ca6e920SHong Zhang 2862ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 287302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 2882ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 289b5e0532dSHong Zhang 290b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 2913fd52205SHong Zhang th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/ 292b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 2932ca6e920SHong Zhang 294a4cab896SHong Zhang /* Build RHS */ 2952c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 29605755b9cSHong Zhang if (th->endpoint) { 297d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr); 29805755b9cSHong Zhang }else { 299d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 30005755b9cSHong Zhang } 30105755b9cSHong Zhang } 302abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 303b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 304b5e0532dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr); 3052c39e106SBarry Smith if (ts->vec_costintegral) { 306b5e0532dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 30736eaed60SHong Zhang } 3082ca6e920SHong Zhang } 3093c54f92cSHong Zhang 3102ca6e920SHong Zhang /* Build LHS */ 3112ca6e920SHong Zhang shift = -1./(th->Theta*ts->time_step); 3123c54f92cSHong Zhang if (th->endpoint) { 3133c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3143c54f92cSHong Zhang }else { 3153c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3163c54f92cSHong Zhang } 3172ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 3182ca6e920SHong Zhang 3192ca6e920SHong Zhang /* Solve LHS X = RHS */ 320abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 321b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 3222ca6e920SHong Zhang } 3233c54f92cSHong Zhang 32436eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 325b5e0532dSHong Zhang if(th->endpoint) { /* two-stage case */ 326b5e0532dSHong Zhang if (th->Theta!=1.) { 3273fd52205SHong Zhang shift = -1./((th->Theta-1.)*ts->time_step); 328b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3292c39e106SBarry Smith if (ts->vec_costintegral) { 330b5e0532dSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 3318263dbbfSHong Zhang } 332abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 333b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3342c39e106SBarry Smith if (ts->vec_costintegral) { 33536eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 33636eaed60SHong Zhang } 3373c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 3382ca6e920SHong Zhang } 339b5e0532dSHong Zhang }else { /* backward Euler */ 340b5e0532dSHong Zhang shift = 0.0; 341b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 342b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 343b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 344b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 345b5e0532dSHong Zhang if (ts->vec_costintegral) { 346b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 347b5e0532dSHong Zhang } 348b5e0532dSHong Zhang } 349b5e0532dSHong Zhang } 3503fd52205SHong Zhang 3513fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 3525bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 353abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 354b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 355b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3563fd52205SHong Zhang } 357b5e0532dSHong Zhang if (th->Theta!=1.) { 358b5e0532dSHong Zhang ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 359abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 360b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 361b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 362b5e0532dSHong Zhang } 3633fd52205SHong Zhang } 3642c39e106SBarry Smith if (ts->vec_costintegral) { 365d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 366abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 36736eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 36836eaed60SHong Zhang } 369b5e0532dSHong Zhang if (th->Theta!=1.) { 370b5e0532dSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 371abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 37236eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 37336eaed60SHong Zhang } 37436eaed60SHong Zhang } 3753fd52205SHong Zhang } 376b5e0532dSHong Zhang } 3773fd52205SHong Zhang }else { /* one-stage case */ 3783c54f92cSHong Zhang shift = 0.0; 379a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 3802c39e106SBarry Smith if (ts->vec_costintegral) { 381d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 3823c54f92cSHong Zhang } 383abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 384b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 385b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 3862c39e106SBarry Smith if (ts->vec_costintegral) { 38736eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 38836eaed60SHong Zhang } 3892ca6e920SHong Zhang } 3903fd52205SHong Zhang if (ts->vecs_sensip) { 3915bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 392abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 393b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 394b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3953fd52205SHong Zhang } 3962c39e106SBarry Smith if (ts->vec_costintegral) { 397d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 398abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 39936eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 40036eaed60SHong Zhang } 40136eaed60SHong Zhang } 4023fd52205SHong Zhang } 4032ca6e920SHong Zhang } 4042ca6e920SHong Zhang 4052ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 4062ca6e920SHong Zhang PetscFunctionReturn(0); 4072ca6e920SHong Zhang } 4082ca6e920SHong Zhang 4092ca6e920SHong Zhang #undef __FUNCT__ 410cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 411cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 412cd652676SJed Brown { 413cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 414*be5899b3SLisandro Dalcin PetscReal dt = t - ts->ptime; 415cd652676SJed Brown PetscErrorCode ierr; 416cd652676SJed Brown 417cd652676SJed Brown PetscFunctionBegin; 418a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 419*be5899b3SLisandro Dalcin if (th->endpoint) dt *= th->Theta; 420*be5899b3SLisandro Dalcin ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 421cd652676SJed Brown PetscFunctionReturn(0); 422cd652676SJed Brown } 423cd652676SJed Brown 4241566a47fSLisandro Dalcin #undef __FUNCT__ 4251566a47fSLisandro Dalcin #define __FUNCT__ "TSEvaluateWLTE_Theta" 4261566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 4271566a47fSLisandro Dalcin { 4281566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4291566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 4301566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 4311566a47fSLisandro Dalcin PetscErrorCode ierr; 4321566a47fSLisandro Dalcin 4331566a47fSLisandro Dalcin PetscFunctionBegin; 4341566a47fSLisandro Dalcin if (!ts->steps || TSEvent_Status(ts) == TSEVENT_RESET_NEXTSTEP) { 4351566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 4361566a47fSLisandro Dalcin *wlte = -1; 4371566a47fSLisandro Dalcin } else { 4381566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 439*be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 440*be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 4411566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 4421566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 4431566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 4441566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 4451566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 4461566a47fSLisandro Dalcin ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte);CHKERRQ(ierr); 4471566a47fSLisandro Dalcin } 4481566a47fSLisandro Dalcin if (order) *order = 2; 4491566a47fSLisandro Dalcin PetscFunctionReturn(0); 4501566a47fSLisandro Dalcin } 4511566a47fSLisandro Dalcin 4521566a47fSLisandro Dalcin #undef __FUNCT__ 4531566a47fSLisandro Dalcin #define __FUNCT__ "TSRollBack_Theta" 4541566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 4551566a47fSLisandro Dalcin { 4561566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4571566a47fSLisandro Dalcin PetscErrorCode ierr; 4581566a47fSLisandro Dalcin 4591566a47fSLisandro Dalcin PetscFunctionBegin; 4601566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 4611566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 4621566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4631566a47fSLisandro Dalcin } 4641566a47fSLisandro Dalcin PetscFunctionReturn(0); 4651566a47fSLisandro Dalcin } 4661566a47fSLisandro Dalcin 467316643e7SJed Brown /*------------------------------------------------------------*/ 468316643e7SJed Brown #undef __FUNCT__ 469277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 470277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 471316643e7SJed Brown { 472316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 473316643e7SJed Brown PetscErrorCode ierr; 474316643e7SJed Brown 475316643e7SJed Brown PetscFunctionBegin; 4766bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 4776bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 4783b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 479eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 4801566a47fSLisandro Dalcin 4811566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 4821566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 4831566a47fSLisandro Dalcin 484b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 485b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 486b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 487b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 488277b19d0SLisandro Dalcin PetscFunctionReturn(0); 489277b19d0SLisandro Dalcin } 490277b19d0SLisandro Dalcin 491277b19d0SLisandro Dalcin #undef __FUNCT__ 492277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 493277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 494277b19d0SLisandro Dalcin { 495277b19d0SLisandro Dalcin PetscErrorCode ierr; 496277b19d0SLisandro Dalcin 497277b19d0SLisandro Dalcin PetscFunctionBegin; 498277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 499277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 500bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 501bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 502bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 503bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 504316643e7SJed Brown PetscFunctionReturn(0); 505316643e7SJed Brown } 506316643e7SJed Brown 507316643e7SJed Brown /* 508316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 5092b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 510316643e7SJed Brown */ 511316643e7SJed Brown #undef __FUNCT__ 5120f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 5130f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 514316643e7SJed Brown { 515316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 516316643e7SJed Brown PetscErrorCode ierr; 5177445fe48SJed Brown Vec X0,Xdot; 5187445fe48SJed Brown DM dm,dmsave; 5191566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 520316643e7SJed Brown 521316643e7SJed Brown PetscFunctionBegin; 5227445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 5235a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 5247445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 525b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 5267445fe48SJed Brown 5277445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 5287445fe48SJed Brown dmsave = ts->dm; 5297445fe48SJed Brown ts->dm = dm; 5307445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 5317445fe48SJed Brown ts->dm = dmsave; 5320d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 533316643e7SJed Brown PetscFunctionReturn(0); 534316643e7SJed Brown } 535316643e7SJed Brown 536316643e7SJed Brown #undef __FUNCT__ 5370f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 538d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 539316643e7SJed Brown { 540316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 541316643e7SJed Brown PetscErrorCode ierr; 5427445fe48SJed Brown Vec Xdot; 5437445fe48SJed Brown DM dm,dmsave; 5441566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 545316643e7SJed Brown 546316643e7SJed Brown PetscFunctionBegin; 5477445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 548*be5899b3SLisandro Dalcin /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 5490298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 5507445fe48SJed Brown 5517445fe48SJed Brown dmsave = ts->dm; 5527445fe48SJed Brown ts->dm = dm; 553d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 5547445fe48SJed Brown ts->dm = dmsave; 5550298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 556316643e7SJed Brown PetscFunctionReturn(0); 557316643e7SJed Brown } 558316643e7SJed Brown 559316643e7SJed Brown #undef __FUNCT__ 560316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 561316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 562316643e7SJed Brown { 563316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 564316643e7SJed Brown PetscErrorCode ierr; 565316643e7SJed Brown 566316643e7SJed Brown PetscFunctionBegin; 567b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 568b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 569b1cb36f3SHong Zhang } 57039d13666SHong Zhang if (!th->X) { 571316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 57239d13666SHong Zhang } 57339d13666SHong Zhang if (!th->Xdot) { 574316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 57539d13666SHong Zhang } 57639d13666SHong Zhang if (!th->X0) { 5773b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 57839d13666SHong Zhang } 5791566a47fSLisandro Dalcin if (th->endpoint) { 5801566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 5817445fe48SJed Brown } 5823b1890cdSShri Abhyankar 5831566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 5841566a47fSLisandro Dalcin 5851566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 5861566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5871566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5881566a47fSLisandro Dalcin 5891566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 5901566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 591ef749922SLisandro Dalcin if (!th->adapt) { 5921566a47fSLisandro Dalcin ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr); 5931566a47fSLisandro Dalcin } else { 5941566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 5951566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 5961566a47fSLisandro Dalcin if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) 5971566a47fSLisandro Dalcin ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP; 5983b1890cdSShri Abhyankar } 5991566a47fSLisandro Dalcin 6001566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 601b4dd3bf9SBarry Smith PetscFunctionReturn(0); 602b4dd3bf9SBarry Smith } 6030c7376e5SHong Zhang 604b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 605b4dd3bf9SBarry Smith 606b4dd3bf9SBarry Smith #undef __FUNCT__ 607b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta" 608b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 609b4dd3bf9SBarry Smith { 610b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 611b4dd3bf9SBarry Smith PetscErrorCode ierr; 612b4dd3bf9SBarry Smith 613b4dd3bf9SBarry Smith PetscFunctionBegin; 614b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 6152ca6e920SHong Zhang if(ts->vecs_sensip) { 616b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 6172ca6e920SHong Zhang } 618b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 619316643e7SJed Brown PetscFunctionReturn(0); 620316643e7SJed Brown } 621316643e7SJed Brown /*------------------------------------------------------------*/ 622316643e7SJed Brown 623316643e7SJed Brown #undef __FUNCT__ 624316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 6254416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 626316643e7SJed Brown { 627316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 628316643e7SJed Brown PetscErrorCode ierr; 629316643e7SJed Brown 630316643e7SJed Brown PetscFunctionBegin; 631e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 632316643e7SJed Brown { 6330298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 6340298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 6350298fd71SBarry 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); 6360298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 637316643e7SJed Brown } 638316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 639316643e7SJed Brown PetscFunctionReturn(0); 640316643e7SJed Brown } 641316643e7SJed Brown 642316643e7SJed Brown #undef __FUNCT__ 643316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 644316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 645316643e7SJed Brown { 646316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 647ace3abfcSBarry Smith PetscBool iascii; 648316643e7SJed Brown PetscErrorCode ierr; 649316643e7SJed Brown 650316643e7SJed Brown PetscFunctionBegin; 651251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 652316643e7SJed Brown if (iascii) { 6537c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 654ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 655316643e7SJed Brown } 6569c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 657ac75fa18SLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 658316643e7SJed Brown PetscFunctionReturn(0); 659316643e7SJed Brown } 660316643e7SJed Brown 6610de4c49aSJed Brown #undef __FUNCT__ 6620de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 663560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 6640de4c49aSJed Brown { 6650de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6660de4c49aSJed Brown 6670de4c49aSJed Brown PetscFunctionBegin; 6680de4c49aSJed Brown *theta = th->Theta; 6690de4c49aSJed Brown PetscFunctionReturn(0); 6700de4c49aSJed Brown } 6710de4c49aSJed Brown 6720de4c49aSJed Brown #undef __FUNCT__ 6730de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 674560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 6750de4c49aSJed Brown { 6760de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6770de4c49aSJed Brown 6780de4c49aSJed Brown PetscFunctionBegin; 6797c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 6800de4c49aSJed Brown th->Theta = theta; 6811566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 6820de4c49aSJed Brown PetscFunctionReturn(0); 6830de4c49aSJed Brown } 684eb284becSJed Brown 685eb284becSJed Brown #undef __FUNCT__ 68678e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 687560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 68826f2ff8fSLisandro Dalcin { 68926f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 69026f2ff8fSLisandro Dalcin 69126f2ff8fSLisandro Dalcin PetscFunctionBegin; 69226f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 69326f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 69426f2ff8fSLisandro Dalcin } 69526f2ff8fSLisandro Dalcin 69626f2ff8fSLisandro Dalcin #undef __FUNCT__ 69726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 698560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 699eb284becSJed Brown { 700eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 701eb284becSJed Brown 702eb284becSJed Brown PetscFunctionBegin; 703eb284becSJed Brown th->endpoint = flg; 704eb284becSJed Brown PetscFunctionReturn(0); 705eb284becSJed Brown } 7060de4c49aSJed Brown 707f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 708f9c1d6abSBarry Smith #undef __FUNCT__ 709f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 710f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 711f9c1d6abSBarry Smith { 712f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 713f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 7143fd8ae06SJed Brown const PetscReal one = 1.0; 715f9c1d6abSBarry Smith 716f9c1d6abSBarry Smith PetscFunctionBegin; 7173fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 718f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 719f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 720f9c1d6abSBarry Smith PetscFunctionReturn(0); 721f9c1d6abSBarry Smith } 722f9c1d6abSBarry Smith #endif 723f9c1d6abSBarry Smith 72442682096SHong Zhang #undef __FUNCT__ 72542682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta" 72642682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 72742682096SHong Zhang { 72842682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 72942682096SHong Zhang 73042682096SHong Zhang PetscFunctionBegin; 7311566a47fSLisandro Dalcin if (ns) *ns = 1; 7321566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 73342682096SHong Zhang PetscFunctionReturn(0); 73442682096SHong Zhang } 735f9c1d6abSBarry Smith 736316643e7SJed Brown /* ------------------------------------------------------------ */ 737316643e7SJed Brown /*MC 73896f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 739316643e7SJed Brown 740316643e7SJed Brown Level: beginner 741316643e7SJed Brown 7424eb428fdSBarry Smith Options Database: 743c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 744c82ed962SBarry Smith . -ts_theta_extrapolate <flg> - Extrapolate stage solution from previous solution (sometimes unstable) 745c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 746c82ed962SBarry Smith - -ts_theta_adapt <flg> - Use time-step adaptivity with the Theta method 7474eb428fdSBarry Smith 748eb284becSJed Brown Notes: 7490c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 7500c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 7514eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 7524eb428fdSBarry Smith 753eb284becSJed Brown This method can be applied to DAE. 754eb284becSJed Brown 755eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 756eb284becSJed Brown 757eb284becSJed Brown .vb 758eb284becSJed Brown Theta | Theta 759eb284becSJed Brown ------------- 760eb284becSJed Brown | 1 761eb284becSJed Brown .ve 762eb284becSJed Brown 763eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 764eb284becSJed Brown 765eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 766eb284becSJed Brown 767eb284becSJed Brown .vb 768eb284becSJed Brown 0 | 0 0 769eb284becSJed Brown 1 | 1-Theta Theta 770eb284becSJed Brown ------------------- 771eb284becSJed Brown | 1-Theta Theta 772eb284becSJed Brown .ve 773eb284becSJed Brown 774eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 775eb284becSJed Brown 776eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 777eb284becSJed Brown 778eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 779eb284becSJed Brown 7804eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 781eb284becSJed Brown 782eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 783316643e7SJed Brown 784316643e7SJed Brown M*/ 785316643e7SJed Brown #undef __FUNCT__ 786316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 7878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 788316643e7SJed Brown { 789316643e7SJed Brown TS_Theta *th; 790316643e7SJed Brown PetscErrorCode ierr; 791316643e7SJed Brown 792316643e7SJed Brown PetscFunctionBegin; 793277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 794316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 795316643e7SJed Brown ts->ops->view = TSView_Theta; 796316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 79742f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 798316643e7SJed Brown ts->ops->step = TSStep_Theta; 799cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 8001566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 80124655328SShri ts->ops->rollback = TSRollBack_Theta; 802316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 8030f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 8040f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 805f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 806f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 807f9c1d6abSBarry Smith #endif 80842682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 80942f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 810b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 811b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 812316643e7SJed Brown 813b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 814316643e7SJed Brown ts->data = (void*)th; 815316643e7SJed Brown 8166f700aefSJed Brown th->extrapolate = PETSC_FALSE; 817316643e7SJed Brown th->Theta = 0.5; 8181566a47fSLisandro Dalcin th->order = 2; 8193b1890cdSShri Abhyankar th->adapt = PETSC_FALSE; 820bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 821bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 822bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 823bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 824316643e7SJed Brown PetscFunctionReturn(0); 825316643e7SJed Brown } 8260de4c49aSJed Brown 8270de4c49aSJed Brown #undef __FUNCT__ 8280de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 8290de4c49aSJed Brown /*@ 8300de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 8310de4c49aSJed Brown 8320de4c49aSJed Brown Not Collective 8330de4c49aSJed Brown 8340de4c49aSJed Brown Input Parameter: 8350de4c49aSJed Brown . ts - timestepping context 8360de4c49aSJed Brown 8370de4c49aSJed Brown Output Parameter: 8380de4c49aSJed Brown . theta - stage abscissa 8390de4c49aSJed Brown 8400de4c49aSJed Brown Note: 8410de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 8420de4c49aSJed Brown 8430de4c49aSJed Brown Level: Advanced 8440de4c49aSJed Brown 8450de4c49aSJed Brown .seealso: TSThetaSetTheta() 8460de4c49aSJed Brown @*/ 8477087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 8480de4c49aSJed Brown { 8494ac538c5SBarry Smith PetscErrorCode ierr; 8500de4c49aSJed Brown 8510de4c49aSJed Brown PetscFunctionBegin; 852afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8530de4c49aSJed Brown PetscValidPointer(theta,2); 8544ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 8550de4c49aSJed Brown PetscFunctionReturn(0); 8560de4c49aSJed Brown } 8570de4c49aSJed Brown 8580de4c49aSJed Brown #undef __FUNCT__ 8590de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 8600de4c49aSJed Brown /*@ 8610de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 8620de4c49aSJed Brown 8630de4c49aSJed Brown Not Collective 8640de4c49aSJed Brown 8650de4c49aSJed Brown Input Parameter: 8660de4c49aSJed Brown + ts - timestepping context 8670de4c49aSJed Brown - theta - stage abscissa 8680de4c49aSJed Brown 8690de4c49aSJed Brown Options Database: 8700de4c49aSJed Brown . -ts_theta_theta <theta> 8710de4c49aSJed Brown 8720de4c49aSJed Brown Level: Intermediate 8730de4c49aSJed Brown 8740de4c49aSJed Brown .seealso: TSThetaGetTheta() 8750de4c49aSJed Brown @*/ 8767087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 8770de4c49aSJed Brown { 8784ac538c5SBarry Smith PetscErrorCode ierr; 8790de4c49aSJed Brown 8800de4c49aSJed Brown PetscFunctionBegin; 881afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8824ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 8830de4c49aSJed Brown PetscFunctionReturn(0); 8840de4c49aSJed Brown } 885f33bbcb6SJed Brown 886eb284becSJed Brown #undef __FUNCT__ 88726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 88826f2ff8fSLisandro Dalcin /*@ 88926f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 89026f2ff8fSLisandro Dalcin 89126f2ff8fSLisandro Dalcin Not Collective 89226f2ff8fSLisandro Dalcin 89326f2ff8fSLisandro Dalcin Input Parameter: 89426f2ff8fSLisandro Dalcin . ts - timestepping context 89526f2ff8fSLisandro Dalcin 89626f2ff8fSLisandro Dalcin Output Parameter: 89726f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 89826f2ff8fSLisandro Dalcin 89926f2ff8fSLisandro Dalcin Level: Advanced 90026f2ff8fSLisandro Dalcin 90126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 90226f2ff8fSLisandro Dalcin @*/ 90326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 90426f2ff8fSLisandro Dalcin { 90526f2ff8fSLisandro Dalcin PetscErrorCode ierr; 90626f2ff8fSLisandro Dalcin 90726f2ff8fSLisandro Dalcin PetscFunctionBegin; 90826f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 90926f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 910163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 91126f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 91226f2ff8fSLisandro Dalcin } 91326f2ff8fSLisandro Dalcin 91426f2ff8fSLisandro Dalcin #undef __FUNCT__ 915eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 916eb284becSJed Brown /*@ 917eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 918eb284becSJed Brown 919eb284becSJed Brown Not Collective 920eb284becSJed Brown 921eb284becSJed Brown Input Parameter: 922eb284becSJed Brown + ts - timestepping context 923eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 924eb284becSJed Brown 925eb284becSJed Brown Options Database: 926eb284becSJed Brown . -ts_theta_endpoint <flg> 927eb284becSJed Brown 928eb284becSJed Brown Level: Intermediate 929eb284becSJed Brown 930eb284becSJed Brown .seealso: TSTHETA, TSCN 931eb284becSJed Brown @*/ 932eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 933eb284becSJed Brown { 934eb284becSJed Brown PetscErrorCode ierr; 935eb284becSJed Brown 936eb284becSJed Brown PetscFunctionBegin; 937eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 938eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 939eb284becSJed Brown PetscFunctionReturn(0); 940eb284becSJed Brown } 941eb284becSJed Brown 942f33bbcb6SJed Brown /* 943f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 944f33bbcb6SJed Brown * The creation functions for these specializations are below. 945f33bbcb6SJed Brown */ 946f33bbcb6SJed Brown 947f33bbcb6SJed Brown #undef __FUNCT__ 9481566a47fSLisandro Dalcin #define __FUNCT__ "TSSetUp_BEuler" 9491566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 9501566a47fSLisandro Dalcin { 9511566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 9521566a47fSLisandro Dalcin PetscErrorCode ierr; 9531566a47fSLisandro Dalcin 9541566a47fSLisandro Dalcin PetscFunctionBegin; 9551566a47fSLisandro 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"); 9561566a47fSLisandro 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"); 9571566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 9581566a47fSLisandro Dalcin PetscFunctionReturn(0); 9591566a47fSLisandro Dalcin } 9601566a47fSLisandro Dalcin 9611566a47fSLisandro Dalcin #undef __FUNCT__ 962f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 963f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 964f33bbcb6SJed Brown { 965d52bd9f3SBarry Smith PetscErrorCode ierr; 966d52bd9f3SBarry Smith 967f33bbcb6SJed Brown PetscFunctionBegin; 9689c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 9691566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 970f33bbcb6SJed Brown PetscFunctionReturn(0); 971f33bbcb6SJed Brown } 972f33bbcb6SJed Brown 973f33bbcb6SJed Brown /*MC 974f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 975f33bbcb6SJed Brown 976f33bbcb6SJed Brown Level: beginner 977f33bbcb6SJed Brown 9784eb428fdSBarry Smith Notes: 979c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 9804eb428fdSBarry Smith 9811566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 9824eb428fdSBarry Smith 983f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 984f33bbcb6SJed Brown 985f33bbcb6SJed Brown M*/ 986f33bbcb6SJed Brown #undef __FUNCT__ 987f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 9888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 989f33bbcb6SJed Brown { 990f33bbcb6SJed Brown PetscErrorCode ierr; 991f33bbcb6SJed Brown 992f33bbcb6SJed Brown PetscFunctionBegin; 993f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 994f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 9951566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 9960c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 997f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 998f33bbcb6SJed Brown PetscFunctionReturn(0); 999f33bbcb6SJed Brown } 1000f33bbcb6SJed Brown 1001f33bbcb6SJed Brown #undef __FUNCT__ 10021566a47fSLisandro Dalcin #define __FUNCT__ "TSSetUp_CN" 10031566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 10041566a47fSLisandro Dalcin { 10051566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 10061566a47fSLisandro Dalcin PetscErrorCode ierr; 10071566a47fSLisandro Dalcin 10081566a47fSLisandro Dalcin PetscFunctionBegin; 10091566a47fSLisandro 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"); 10101566a47fSLisandro 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"); 10111566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 10121566a47fSLisandro Dalcin PetscFunctionReturn(0); 10131566a47fSLisandro Dalcin } 10141566a47fSLisandro Dalcin 10151566a47fSLisandro Dalcin #undef __FUNCT__ 1016f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 1017f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1018f33bbcb6SJed Brown { 1019d52bd9f3SBarry Smith PetscErrorCode ierr; 1020d52bd9f3SBarry Smith 1021f33bbcb6SJed Brown PetscFunctionBegin; 10229c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 10231566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 1024f33bbcb6SJed Brown PetscFunctionReturn(0); 1025f33bbcb6SJed Brown } 1026f33bbcb6SJed Brown 1027f33bbcb6SJed Brown /*MC 1028f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1029f33bbcb6SJed Brown 1030f33bbcb6SJed Brown Level: beginner 1031f33bbcb6SJed Brown 1032f33bbcb6SJed Brown Notes: 10337cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 10347cf5af47SJed Brown 10357cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1036f33bbcb6SJed Brown 1037f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1038f33bbcb6SJed Brown 1039f33bbcb6SJed Brown M*/ 1040f33bbcb6SJed Brown #undef __FUNCT__ 1041f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 10428cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1043f33bbcb6SJed Brown { 1044f33bbcb6SJed Brown PetscErrorCode ierr; 1045f33bbcb6SJed Brown 1046f33bbcb6SJed Brown PetscFunctionBegin; 1047f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1048f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1049eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 10500c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1051f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1052f33bbcb6SJed Brown PetscFunctionReturn(0); 1053f33bbcb6SJed Brown } 1054