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 PetscReal time_step_prev; 281566a47fSLisandro Dalcin Vec vec_sol_prev; 291566a47fSLisandro Dalcin Vec vec_lte_work; 301566a47fSLisandro Dalcin 311566a47fSLisandro Dalcin TSStepStatus status; 32316643e7SJed Brown } TS_Theta; 33316643e7SJed Brown 34316643e7SJed Brown #undef __FUNCT__ 357445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 367445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 377445fe48SJed Brown { 387445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 397445fe48SJed Brown PetscErrorCode ierr; 407445fe48SJed Brown 417445fe48SJed Brown PetscFunctionBegin; 427445fe48SJed Brown if (X0) { 437445fe48SJed Brown if (dm && dm != ts->dm) { 440d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 457445fe48SJed Brown } else *X0 = ts->vec_sol; 467445fe48SJed Brown } 477445fe48SJed Brown if (Xdot) { 487445fe48SJed Brown if (dm && dm != ts->dm) { 490d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 507445fe48SJed Brown } else *Xdot = th->Xdot; 517445fe48SJed Brown } 527445fe48SJed Brown PetscFunctionReturn(0); 537445fe48SJed Brown } 547445fe48SJed Brown 550d0b770aSPeter Brune #undef __FUNCT__ 560d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 570d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 580d0b770aSPeter Brune { 590d0b770aSPeter Brune PetscErrorCode ierr; 600d0b770aSPeter Brune 610d0b770aSPeter Brune PetscFunctionBegin; 620d0b770aSPeter Brune if (X0) { 630d0b770aSPeter Brune if (dm && dm != ts->dm) { 640d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 650d0b770aSPeter Brune } 660d0b770aSPeter Brune } 670d0b770aSPeter Brune if (Xdot) { 680d0b770aSPeter Brune if (dm && dm != ts->dm) { 690d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 700d0b770aSPeter Brune } 710d0b770aSPeter Brune } 720d0b770aSPeter Brune PetscFunctionReturn(0); 730d0b770aSPeter Brune } 740d0b770aSPeter Brune 757445fe48SJed Brown #undef __FUNCT__ 767445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 777445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 787445fe48SJed Brown { 797445fe48SJed Brown 807445fe48SJed Brown PetscFunctionBegin; 817445fe48SJed Brown PetscFunctionReturn(0); 827445fe48SJed Brown } 837445fe48SJed Brown 847445fe48SJed Brown #undef __FUNCT__ 857445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 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 1047445fe48SJed Brown #undef __FUNCT__ 105258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta" 106258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 107258e1594SPeter Brune { 108258e1594SPeter Brune 109258e1594SPeter Brune PetscFunctionBegin; 110258e1594SPeter Brune PetscFunctionReturn(0); 111258e1594SPeter Brune } 112258e1594SPeter Brune 113258e1594SPeter Brune #undef __FUNCT__ 114258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta" 115258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 116258e1594SPeter Brune { 117258e1594SPeter Brune TS ts = (TS)ctx; 118258e1594SPeter Brune PetscErrorCode ierr; 119258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 120258e1594SPeter Brune 121258e1594SPeter Brune PetscFunctionBegin; 122258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 123258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 124258e1594SPeter Brune 125258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 127258e1594SPeter Brune 128258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 129258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 130258e1594SPeter Brune 131258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 132258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 133258e1594SPeter Brune PetscFunctionReturn(0); 134258e1594SPeter Brune } 135258e1594SPeter Brune 1363b1890cdSShri Abhyankar #undef __FUNCT__ 137b1cb36f3SHong Zhang #define __FUNCT__ "TSForwardCostIntegral_Theta" 138b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 139b1cb36f3SHong Zhang { 140b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 141b1cb36f3SHong Zhang PetscErrorCode ierr; 142b1cb36f3SHong Zhang 143b1cb36f3SHong Zhang PetscFunctionBegin; 144b1cb36f3SHong Zhang /* backup cost integral */ 145b1cb36f3SHong Zhang ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 146b1cb36f3SHong Zhang if (th->endpoint) { 147b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1481566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 149b1cb36f3SHong Zhang } 150b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 151b1cb36f3SHong Zhang if (th->endpoint) { 152b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 153b1cb36f3SHong Zhang } else { 154b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 155b1cb36f3SHong Zhang } 156b1cb36f3SHong Zhang PetscFunctionReturn(0); 157b1cb36f3SHong Zhang } 158b1cb36f3SHong Zhang 159b1cb36f3SHong Zhang #undef __FUNCT__ 160b1cb36f3SHong Zhang #define __FUNCT__ "TSAdjointCostIntegral_Theta" 161b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 162b1cb36f3SHong Zhang { 163b1cb36f3SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 164b1cb36f3SHong Zhang PetscErrorCode ierr; 165b1cb36f3SHong Zhang 166b1cb36f3SHong Zhang PetscFunctionBegin; 167b1cb36f3SHong Zhang if (th->endpoint) { 168b1cb36f3SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 169b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 170b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 171b1cb36f3SHong Zhang if (th->Theta!=1) { 172b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 1731566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr); 174b1cb36f3SHong Zhang } 175b1cb36f3SHong Zhang }else { 176b1cb36f3SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 177b1cb36f3SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 178b1cb36f3SHong Zhang } 17924655328SShri PetscFunctionReturn(0); 18024655328SShri } 18124655328SShri 18224655328SShri #undef __FUNCT__ 1831566a47fSLisandro Dalcin #define __FUNCT__ "TS_SNESSolve" 1841566a47fSLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 1851566a47fSLisandro Dalcin { 1861566a47fSLisandro Dalcin PetscInt nits,lits; 1871566a47fSLisandro Dalcin PetscErrorCode ierr; 1881566a47fSLisandro Dalcin 1891566a47fSLisandro Dalcin PetscFunctionBegin; 1901566a47fSLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 1911566a47fSLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 1921566a47fSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1931566a47fSLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1941566a47fSLisandro Dalcin PetscFunctionReturn(0); 1951566a47fSLisandro Dalcin } 1961566a47fSLisandro Dalcin 1971566a47fSLisandro Dalcin #undef __FUNCT__ 198316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 199193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 200316643e7SJed Brown { 201316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2021566a47fSLisandro Dalcin PetscInt rejections = 0; 2034957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 2041566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 205051f2191SLisandro Dalcin PetscErrorCode ierr; 206316643e7SJed Brown 207316643e7SJed Brown PetscFunctionBegin; 2081566a47fSLisandro Dalcin 2091566a47fSLisandro Dalcin if (!ts->steprollback) { 2101566a47fSLisandro Dalcin if (th->adapt) { th->time_step_prev = ts->time_step_prev; } 2111566a47fSLisandro Dalcin if (th->adapt) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 2123b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 2131566a47fSLisandro Dalcin } 2141566a47fSLisandro Dalcin 2151566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 2161566a47fSLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 2171566a47fSLisandro Dalcin 2181566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 2191566a47fSLisandro Dalcin th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 220b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 221316643e7SJed Brown 222eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 223eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 2241566a47fSLisandro Dalcin ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 2251566a47fSLisandro Dalcin ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 2261566a47fSLisandro Dalcin ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 2271566a47fSLisandro Dalcin } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 2281566a47fSLisandro Dalcin ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 229eb284becSJed Brown } 230ace68cafSJed Brown if (th->extrapolate) { 2311566a47fSLisandro Dalcin ierr = VecWAXPY(th->X,1/shift,th->Xdot,th->X0);CHKERRQ(ierr); 232ace68cafSJed Brown } else { 2331566a47fSLisandro Dalcin ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 234ace68cafSJed Brown } 2351566a47fSLisandro Dalcin 2361566a47fSLisandro Dalcin ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 237051f2191SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr); 2381566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 239051f2191SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 240051f2191SLisandro Dalcin 241051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2421566a47fSLisandro Dalcin if (th->endpoint) { 2431566a47fSLisandro Dalcin ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 2441566a47fSLisandro Dalcin } else { 245cb3a5882SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 2461566a47fSLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 2471566a47fSLisandro Dalcin } 2481566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2491566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2501566a47fSLisandro Dalcin if (!accept) { 2511566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 2521566a47fSLisandro Dalcin ts->time_step = next_time_step; goto reject_step; 253051f2191SLisandro Dalcin } 2543b1890cdSShri Abhyankar 255b1cb36f3SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 256b1cb36f3SHong Zhang th->ptime = ts->ptime; 257b1cb36f3SHong Zhang th->time_step = ts->time_step; 25817145e75SHong Zhang } 2591566a47fSLisandro Dalcin 2602b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 261cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 262316643e7SJed Brown ts->steps++; 263051f2191SLisandro Dalcin break; 264051f2191SLisandro Dalcin 265051f2191SLisandro Dalcin reject_step: 2661566a47fSLisandro Dalcin ts->reject++; 2671566a47fSLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 268051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2691566a47fSLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 2703b1890cdSShri Abhyankar } 2713b1890cdSShri Abhyankar } 272316643e7SJed Brown PetscFunctionReturn(0); 273316643e7SJed Brown } 274316643e7SJed Brown 275cd652676SJed Brown #undef __FUNCT__ 27642f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_Theta" 27742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2782ca6e920SHong Zhang { 2792ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 280b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 2812ca6e920SHong Zhang PetscInt nadj; 2822ca6e920SHong Zhang PetscErrorCode ierr; 2832ca6e920SHong Zhang Mat J,Jp; 2842ca6e920SHong Zhang KSP ksp; 2852ca6e920SHong Zhang PetscReal shift; 2862ca6e920SHong Zhang 2872ca6e920SHong Zhang PetscFunctionBegin; 2882ca6e920SHong Zhang 2892ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 290302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 2912ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 292b5e0532dSHong Zhang 293b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 2943fd52205SHong Zhang th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/ 295b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 2962ca6e920SHong Zhang 297a4cab896SHong Zhang /* Build RHS */ 2982c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 29905755b9cSHong Zhang if (th->endpoint) { 300d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr); 30105755b9cSHong Zhang }else { 302d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 30305755b9cSHong Zhang } 30405755b9cSHong Zhang } 305abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 306b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 307b5e0532dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr); 3082c39e106SBarry Smith if (ts->vec_costintegral) { 309b5e0532dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 31036eaed60SHong Zhang } 3112ca6e920SHong Zhang } 3123c54f92cSHong Zhang 3132ca6e920SHong Zhang /* Build LHS */ 3142ca6e920SHong Zhang shift = -1./(th->Theta*ts->time_step); 3153c54f92cSHong Zhang if (th->endpoint) { 3163c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3173c54f92cSHong Zhang }else { 3183c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3193c54f92cSHong Zhang } 3202ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 3212ca6e920SHong Zhang 3222ca6e920SHong Zhang /* Solve LHS X = RHS */ 323abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 324b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 3252ca6e920SHong Zhang } 3263c54f92cSHong Zhang 32736eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 328b5e0532dSHong Zhang if(th->endpoint) { /* two-stage case */ 329b5e0532dSHong Zhang if (th->Theta!=1.) { 3303fd52205SHong Zhang shift = -1./((th->Theta-1.)*ts->time_step); 331b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3322c39e106SBarry Smith if (ts->vec_costintegral) { 333b5e0532dSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 3348263dbbfSHong Zhang } 335abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 336b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3372c39e106SBarry Smith if (ts->vec_costintegral) { 33836eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 33936eaed60SHong Zhang } 3403c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 3412ca6e920SHong Zhang } 342b5e0532dSHong Zhang }else { /* backward Euler */ 343b5e0532dSHong Zhang shift = 0.0; 344b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 345b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 346b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 347b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 348b5e0532dSHong Zhang if (ts->vec_costintegral) { 349b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 350b5e0532dSHong Zhang } 351b5e0532dSHong Zhang } 352b5e0532dSHong Zhang } 3533fd52205SHong Zhang 3543fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 3555bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 356abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 357b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 358b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3593fd52205SHong Zhang } 360b5e0532dSHong Zhang if (th->Theta!=1.) { 361b5e0532dSHong Zhang ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 362abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 363b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 364b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 365b5e0532dSHong Zhang } 3663fd52205SHong Zhang } 3672c39e106SBarry Smith if (ts->vec_costintegral) { 368d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 369abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 37036eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 37136eaed60SHong Zhang } 372b5e0532dSHong Zhang if (th->Theta!=1.) { 373b5e0532dSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 374abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 37536eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 37636eaed60SHong Zhang } 37736eaed60SHong Zhang } 3783fd52205SHong Zhang } 379b5e0532dSHong Zhang } 3803fd52205SHong Zhang }else { /* one-stage case */ 3813c54f92cSHong Zhang shift = 0.0; 382a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 3832c39e106SBarry Smith if (ts->vec_costintegral) { 384d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 3853c54f92cSHong Zhang } 386abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 387b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 388b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 3892c39e106SBarry Smith if (ts->vec_costintegral) { 39036eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 39136eaed60SHong Zhang } 3922ca6e920SHong Zhang } 3933fd52205SHong Zhang if (ts->vecs_sensip) { 3945bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 395abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 396b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 397b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3983fd52205SHong Zhang } 3992c39e106SBarry Smith if (ts->vec_costintegral) { 400d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 401abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 40236eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 40336eaed60SHong Zhang } 40436eaed60SHong Zhang } 4053fd52205SHong Zhang } 4062ca6e920SHong Zhang } 4072ca6e920SHong Zhang 4082ca6e920SHong Zhang ts->steps++; 4092ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 4102ca6e920SHong Zhang PetscFunctionReturn(0); 4112ca6e920SHong Zhang } 4122ca6e920SHong Zhang 4132ca6e920SHong Zhang #undef __FUNCT__ 414cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 415cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 416cd652676SJed Brown { 417cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 4185a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 419cd652676SJed Brown PetscErrorCode ierr; 420cd652676SJed Brown 421cd652676SJed Brown PetscFunctionBegin; 422a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 4235a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 4245a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 425cd652676SJed Brown PetscFunctionReturn(0); 426cd652676SJed Brown } 427cd652676SJed Brown 4281566a47fSLisandro Dalcin #define TSEvent_Status(ts) (ts->event ? ts->event->status : TSEVENT_NONE) 4291566a47fSLisandro Dalcin 4301566a47fSLisandro Dalcin #undef __FUNCT__ 4311566a47fSLisandro Dalcin #define __FUNCT__ "TSEvaluateWLTE_Theta" 4321566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 4331566a47fSLisandro Dalcin { 4341566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4351566a47fSLisandro Dalcin Vec X = ts->vec_sol; /* X = solution */ 4361566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 4371566a47fSLisandro Dalcin PetscErrorCode ierr; 4381566a47fSLisandro Dalcin 4391566a47fSLisandro Dalcin PetscFunctionBegin; 4401566a47fSLisandro Dalcin if (!ts->steps || TSEvent_Status(ts) == TSEVENT_RESET_NEXTSTEP) { 4411566a47fSLisandro Dalcin /* Cannot compute LTE in first step or in restart after event */ 4421566a47fSLisandro Dalcin *wlte = -1; 4431566a47fSLisandro Dalcin } else { 4441566a47fSLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 4451566a47fSLisandro Dalcin PetscReal a = 1 + th->time_step_prev/ts->time_step; 4461566a47fSLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 4471566a47fSLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 4481566a47fSLisandro Dalcin vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 4491566a47fSLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 4501566a47fSLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 4511566a47fSLisandro Dalcin ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte);CHKERRQ(ierr); 4521566a47fSLisandro Dalcin } 4531566a47fSLisandro Dalcin if (order) *order = 2; 4541566a47fSLisandro Dalcin PetscFunctionReturn(0); 4551566a47fSLisandro Dalcin } 4561566a47fSLisandro Dalcin 4571566a47fSLisandro Dalcin #undef __FUNCT__ 4581566a47fSLisandro Dalcin #define __FUNCT__ "TSRollBack_Theta" 4591566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts) 4601566a47fSLisandro Dalcin { 4611566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 4621566a47fSLisandro Dalcin PetscErrorCode ierr; 4631566a47fSLisandro Dalcin 4641566a47fSLisandro Dalcin PetscFunctionBegin; 4651566a47fSLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 4661566a47fSLisandro Dalcin if (ts->vec_costintegral && ts->costintegralfwd) { 4671566a47fSLisandro Dalcin ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4681566a47fSLisandro Dalcin } 4691566a47fSLisandro Dalcin PetscFunctionReturn(0); 4701566a47fSLisandro Dalcin } 4711566a47fSLisandro Dalcin 472316643e7SJed Brown /*------------------------------------------------------------*/ 473316643e7SJed Brown #undef __FUNCT__ 474277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 475277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 476316643e7SJed Brown { 477316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 478316643e7SJed Brown PetscErrorCode ierr; 479316643e7SJed Brown 480316643e7SJed Brown PetscFunctionBegin; 4816bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 4826bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 4833b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 484eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 4851566a47fSLisandro Dalcin 4861566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 4871566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 4881566a47fSLisandro Dalcin 489b1cb36f3SHong Zhang ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 490b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 491b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 492b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 493277b19d0SLisandro Dalcin PetscFunctionReturn(0); 494277b19d0SLisandro Dalcin } 495277b19d0SLisandro Dalcin 496277b19d0SLisandro Dalcin #undef __FUNCT__ 497277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 498277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 499277b19d0SLisandro Dalcin { 500277b19d0SLisandro Dalcin PetscErrorCode ierr; 501277b19d0SLisandro Dalcin 502277b19d0SLisandro Dalcin PetscFunctionBegin; 503277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 504277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 505bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 506bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 507bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 508bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 509316643e7SJed Brown PetscFunctionReturn(0); 510316643e7SJed Brown } 511316643e7SJed Brown 512316643e7SJed Brown /* 513316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 5142b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 515316643e7SJed Brown */ 516316643e7SJed Brown #undef __FUNCT__ 5170f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 5180f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 519316643e7SJed Brown { 520316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 521316643e7SJed Brown PetscErrorCode ierr; 5227445fe48SJed Brown Vec X0,Xdot; 5237445fe48SJed Brown DM dm,dmsave; 5241566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 525316643e7SJed Brown 526316643e7SJed Brown PetscFunctionBegin; 5277445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 5285a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 5297445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 530b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 5317445fe48SJed Brown 5327445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 5337445fe48SJed Brown dmsave = ts->dm; 5347445fe48SJed Brown ts->dm = dm; 5357445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 5367445fe48SJed Brown ts->dm = dmsave; 5370d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 538316643e7SJed Brown PetscFunctionReturn(0); 539316643e7SJed Brown } 540316643e7SJed Brown 541316643e7SJed Brown #undef __FUNCT__ 5420f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 543d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 544316643e7SJed Brown { 545316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 546316643e7SJed Brown PetscErrorCode ierr; 5477445fe48SJed Brown Vec Xdot; 5487445fe48SJed Brown DM dm,dmsave; 5491566a47fSLisandro Dalcin PetscReal shift = 1/(th->Theta*ts->time_step); 550316643e7SJed Brown 551316643e7SJed Brown PetscFunctionBegin; 5527445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 5537445fe48SJed Brown 5540f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 5550298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 5567445fe48SJed Brown 5577445fe48SJed Brown dmsave = ts->dm; 5587445fe48SJed Brown ts->dm = dm; 559d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 5607445fe48SJed Brown ts->dm = dmsave; 5610298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 562316643e7SJed Brown PetscFunctionReturn(0); 563316643e7SJed Brown } 564316643e7SJed Brown 565316643e7SJed Brown #undef __FUNCT__ 566316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 567316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 568316643e7SJed Brown { 569316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 570316643e7SJed Brown PetscErrorCode ierr; 571316643e7SJed Brown 572316643e7SJed Brown PetscFunctionBegin; 573b1cb36f3SHong Zhang if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 574b1cb36f3SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 575b1cb36f3SHong Zhang } 57639d13666SHong Zhang if (!th->X) { 577316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 57839d13666SHong Zhang } 57939d13666SHong Zhang if (!th->Xdot) { 580316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 58139d13666SHong Zhang } 58239d13666SHong Zhang if (!th->X0) { 5833b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 58439d13666SHong Zhang } 5851566a47fSLisandro Dalcin if (th->endpoint) { 5861566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 5877445fe48SJed Brown } 5883b1890cdSShri Abhyankar 5891566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 5901566a47fSLisandro Dalcin 5911566a47fSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 5921566a47fSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5931566a47fSLisandro Dalcin ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5941566a47fSLisandro Dalcin 5951566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 5961566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 597ef749922SLisandro Dalcin if (!th->adapt) { 5981566a47fSLisandro Dalcin ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr); 5991566a47fSLisandro Dalcin } else { 6001566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 6011566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 6021566a47fSLisandro Dalcin if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) 6031566a47fSLisandro Dalcin ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP; 6043b1890cdSShri Abhyankar } 6051566a47fSLisandro Dalcin 6061566a47fSLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 607b4dd3bf9SBarry Smith PetscFunctionReturn(0); 608b4dd3bf9SBarry Smith } 6090c7376e5SHong Zhang 610b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 611b4dd3bf9SBarry Smith 612b4dd3bf9SBarry Smith #undef __FUNCT__ 613b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta" 614b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 615b4dd3bf9SBarry Smith { 616b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 617b4dd3bf9SBarry Smith PetscErrorCode ierr; 618b4dd3bf9SBarry Smith 619b4dd3bf9SBarry Smith PetscFunctionBegin; 620b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 6212ca6e920SHong Zhang if(ts->vecs_sensip) { 622b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 6232ca6e920SHong Zhang } 624b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 625316643e7SJed Brown PetscFunctionReturn(0); 626316643e7SJed Brown } 627316643e7SJed Brown /*------------------------------------------------------------*/ 628316643e7SJed Brown 629316643e7SJed Brown #undef __FUNCT__ 630316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 6314416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 632316643e7SJed Brown { 633316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 634316643e7SJed Brown PetscErrorCode ierr; 635316643e7SJed Brown 636316643e7SJed Brown PetscFunctionBegin; 637e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 638316643e7SJed Brown { 6390298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 6400298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 6410298fd71SBarry 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); 6420298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 643316643e7SJed Brown } 644316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 645316643e7SJed Brown PetscFunctionReturn(0); 646316643e7SJed Brown } 647316643e7SJed Brown 648316643e7SJed Brown #undef __FUNCT__ 649316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 650316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 651316643e7SJed Brown { 652316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 653ace3abfcSBarry Smith PetscBool iascii; 654316643e7SJed Brown PetscErrorCode ierr; 655316643e7SJed Brown 656316643e7SJed Brown PetscFunctionBegin; 657251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 658316643e7SJed Brown if (iascii) { 6597c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 660ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 661316643e7SJed Brown } 662ac75fa18SLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 663316643e7SJed Brown PetscFunctionReturn(0); 664316643e7SJed Brown } 665316643e7SJed Brown 6660de4c49aSJed Brown #undef __FUNCT__ 6670de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 6687087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 6690de4c49aSJed Brown { 6700de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6710de4c49aSJed Brown 6720de4c49aSJed Brown PetscFunctionBegin; 6730de4c49aSJed Brown *theta = th->Theta; 6740de4c49aSJed Brown PetscFunctionReturn(0); 6750de4c49aSJed Brown } 6760de4c49aSJed Brown 6770de4c49aSJed Brown #undef __FUNCT__ 6780de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 6797087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 6800de4c49aSJed Brown { 6810de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6820de4c49aSJed Brown 6830de4c49aSJed Brown PetscFunctionBegin; 6847c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 6850de4c49aSJed Brown th->Theta = theta; 6861566a47fSLisandro Dalcin th->order = (th->Theta == 0.5) ? 2 : 1; 6870de4c49aSJed Brown PetscFunctionReturn(0); 6880de4c49aSJed Brown } 689eb284becSJed Brown 690eb284becSJed Brown #undef __FUNCT__ 69178e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 69226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 69326f2ff8fSLisandro Dalcin { 69426f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 69526f2ff8fSLisandro Dalcin 69626f2ff8fSLisandro Dalcin PetscFunctionBegin; 69726f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 69826f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 69926f2ff8fSLisandro Dalcin } 70026f2ff8fSLisandro Dalcin 70126f2ff8fSLisandro Dalcin #undef __FUNCT__ 70226f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 703eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 704eb284becSJed Brown { 705eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 706eb284becSJed Brown 707eb284becSJed Brown PetscFunctionBegin; 708eb284becSJed Brown th->endpoint = flg; 709eb284becSJed Brown PetscFunctionReturn(0); 710eb284becSJed Brown } 7110de4c49aSJed Brown 712f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 713f9c1d6abSBarry Smith #undef __FUNCT__ 714f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 715f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 716f9c1d6abSBarry Smith { 717f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 718f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 7193fd8ae06SJed Brown const PetscReal one = 1.0; 720f9c1d6abSBarry Smith 721f9c1d6abSBarry Smith PetscFunctionBegin; 7223fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 723f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 724f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 725f9c1d6abSBarry Smith PetscFunctionReturn(0); 726f9c1d6abSBarry Smith } 727f9c1d6abSBarry Smith #endif 728f9c1d6abSBarry Smith 72942682096SHong Zhang #undef __FUNCT__ 73042682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta" 73142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 73242682096SHong Zhang { 73342682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 73442682096SHong Zhang 73542682096SHong Zhang PetscFunctionBegin; 7361566a47fSLisandro Dalcin if (ns) *ns = 1; 7371566a47fSLisandro Dalcin if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 73842682096SHong Zhang PetscFunctionReturn(0); 73942682096SHong Zhang } 740f9c1d6abSBarry Smith 741316643e7SJed Brown /* ------------------------------------------------------------ */ 742316643e7SJed Brown /*MC 74396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 744316643e7SJed Brown 745316643e7SJed Brown Level: beginner 746316643e7SJed Brown 7474eb428fdSBarry Smith Options Database: 748c82ed962SBarry Smith + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 749c82ed962SBarry Smith . -ts_theta_extrapolate <flg> - Extrapolate stage solution from previous solution (sometimes unstable) 750c82ed962SBarry Smith . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 751c82ed962SBarry Smith - -ts_theta_adapt <flg> - Use time-step adaptivity with the Theta method 7524eb428fdSBarry Smith 753eb284becSJed Brown Notes: 7540c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 7550c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 7564eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 7574eb428fdSBarry Smith 758eb284becSJed Brown This method can be applied to DAE. 759eb284becSJed Brown 760eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 761eb284becSJed Brown 762eb284becSJed Brown .vb 763eb284becSJed Brown Theta | Theta 764eb284becSJed Brown ------------- 765eb284becSJed Brown | 1 766eb284becSJed Brown .ve 767eb284becSJed Brown 768eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 769eb284becSJed Brown 770eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 771eb284becSJed Brown 772eb284becSJed Brown .vb 773eb284becSJed Brown 0 | 0 0 774eb284becSJed Brown 1 | 1-Theta Theta 775eb284becSJed Brown ------------------- 776eb284becSJed Brown | 1-Theta Theta 777eb284becSJed Brown .ve 778eb284becSJed Brown 779eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 780eb284becSJed Brown 781eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 782eb284becSJed Brown 783eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 784eb284becSJed Brown 7854eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 786eb284becSJed Brown 787eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 788316643e7SJed Brown 789316643e7SJed Brown M*/ 790316643e7SJed Brown #undef __FUNCT__ 791316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 7928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 793316643e7SJed Brown { 794316643e7SJed Brown TS_Theta *th; 795316643e7SJed Brown PetscErrorCode ierr; 796316643e7SJed Brown 797316643e7SJed Brown PetscFunctionBegin; 798277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 799316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 800316643e7SJed Brown ts->ops->view = TSView_Theta; 801316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 80242f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 803316643e7SJed Brown ts->ops->step = TSStep_Theta; 804cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 8051566a47fSLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 80624655328SShri ts->ops->rollback = TSRollBack_Theta; 807316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 8080f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 8090f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 810f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 811f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 812f9c1d6abSBarry Smith #endif 81342682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 81442f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 815b1cb36f3SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 816b1cb36f3SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 817316643e7SJed Brown 818b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 819316643e7SJed Brown ts->data = (void*)th; 820316643e7SJed Brown 8216f700aefSJed Brown th->extrapolate = PETSC_FALSE; 822316643e7SJed Brown th->Theta = 0.5; 8231566a47fSLisandro Dalcin th->order = 2; 8243b1890cdSShri Abhyankar th->adapt = PETSC_FALSE; 825bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 826bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 827bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 828bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 829316643e7SJed Brown PetscFunctionReturn(0); 830316643e7SJed Brown } 8310de4c49aSJed Brown 8320de4c49aSJed Brown #undef __FUNCT__ 8330de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 8340de4c49aSJed Brown /*@ 8350de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 8360de4c49aSJed Brown 8370de4c49aSJed Brown Not Collective 8380de4c49aSJed Brown 8390de4c49aSJed Brown Input Parameter: 8400de4c49aSJed Brown . ts - timestepping context 8410de4c49aSJed Brown 8420de4c49aSJed Brown Output Parameter: 8430de4c49aSJed Brown . theta - stage abscissa 8440de4c49aSJed Brown 8450de4c49aSJed Brown Note: 8460de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 8470de4c49aSJed Brown 8480de4c49aSJed Brown Level: Advanced 8490de4c49aSJed Brown 8500de4c49aSJed Brown .seealso: TSThetaSetTheta() 8510de4c49aSJed Brown @*/ 8527087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 8530de4c49aSJed Brown { 8544ac538c5SBarry Smith PetscErrorCode ierr; 8550de4c49aSJed Brown 8560de4c49aSJed Brown PetscFunctionBegin; 857afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8580de4c49aSJed Brown PetscValidPointer(theta,2); 8594ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 8600de4c49aSJed Brown PetscFunctionReturn(0); 8610de4c49aSJed Brown } 8620de4c49aSJed Brown 8630de4c49aSJed Brown #undef __FUNCT__ 8640de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 8650de4c49aSJed Brown /*@ 8660de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 8670de4c49aSJed Brown 8680de4c49aSJed Brown Not Collective 8690de4c49aSJed Brown 8700de4c49aSJed Brown Input Parameter: 8710de4c49aSJed Brown + ts - timestepping context 8720de4c49aSJed Brown - theta - stage abscissa 8730de4c49aSJed Brown 8740de4c49aSJed Brown Options Database: 8750de4c49aSJed Brown . -ts_theta_theta <theta> 8760de4c49aSJed Brown 8770de4c49aSJed Brown Level: Intermediate 8780de4c49aSJed Brown 8790de4c49aSJed Brown .seealso: TSThetaGetTheta() 8800de4c49aSJed Brown @*/ 8817087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 8820de4c49aSJed Brown { 8834ac538c5SBarry Smith PetscErrorCode ierr; 8840de4c49aSJed Brown 8850de4c49aSJed Brown PetscFunctionBegin; 886afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8874ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 8880de4c49aSJed Brown PetscFunctionReturn(0); 8890de4c49aSJed Brown } 890f33bbcb6SJed Brown 891eb284becSJed Brown #undef __FUNCT__ 89226f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 89326f2ff8fSLisandro Dalcin /*@ 89426f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 89526f2ff8fSLisandro Dalcin 89626f2ff8fSLisandro Dalcin Not Collective 89726f2ff8fSLisandro Dalcin 89826f2ff8fSLisandro Dalcin Input Parameter: 89926f2ff8fSLisandro Dalcin . ts - timestepping context 90026f2ff8fSLisandro Dalcin 90126f2ff8fSLisandro Dalcin Output Parameter: 90226f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 90326f2ff8fSLisandro Dalcin 90426f2ff8fSLisandro Dalcin Level: Advanced 90526f2ff8fSLisandro Dalcin 90626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 90726f2ff8fSLisandro Dalcin @*/ 90826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 90926f2ff8fSLisandro Dalcin { 91026f2ff8fSLisandro Dalcin PetscErrorCode ierr; 91126f2ff8fSLisandro Dalcin 91226f2ff8fSLisandro Dalcin PetscFunctionBegin; 91326f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 91426f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 915*163d334eSBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 91626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 91726f2ff8fSLisandro Dalcin } 91826f2ff8fSLisandro Dalcin 91926f2ff8fSLisandro Dalcin #undef __FUNCT__ 920eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 921eb284becSJed Brown /*@ 922eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 923eb284becSJed Brown 924eb284becSJed Brown Not Collective 925eb284becSJed Brown 926eb284becSJed Brown Input Parameter: 927eb284becSJed Brown + ts - timestepping context 928eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 929eb284becSJed Brown 930eb284becSJed Brown Options Database: 931eb284becSJed Brown . -ts_theta_endpoint <flg> 932eb284becSJed Brown 933eb284becSJed Brown Level: Intermediate 934eb284becSJed Brown 935eb284becSJed Brown .seealso: TSTHETA, TSCN 936eb284becSJed Brown @*/ 937eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 938eb284becSJed Brown { 939eb284becSJed Brown PetscErrorCode ierr; 940eb284becSJed Brown 941eb284becSJed Brown PetscFunctionBegin; 942eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 943eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 944eb284becSJed Brown PetscFunctionReturn(0); 945eb284becSJed Brown } 946eb284becSJed Brown 947f33bbcb6SJed Brown /* 948f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 949f33bbcb6SJed Brown * The creation functions for these specializations are below. 950f33bbcb6SJed Brown */ 951f33bbcb6SJed Brown 952f33bbcb6SJed Brown #undef __FUNCT__ 9531566a47fSLisandro Dalcin #define __FUNCT__ "TSSetUp_BEuler" 9541566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts) 9551566a47fSLisandro Dalcin { 9561566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 9571566a47fSLisandro Dalcin PetscErrorCode ierr; 9581566a47fSLisandro Dalcin 9591566a47fSLisandro Dalcin PetscFunctionBegin; 9601566a47fSLisandro 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"); 9611566a47fSLisandro 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"); 9621566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 9631566a47fSLisandro Dalcin PetscFunctionReturn(0); 9641566a47fSLisandro Dalcin } 9651566a47fSLisandro Dalcin 9661566a47fSLisandro Dalcin #undef __FUNCT__ 967f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 968f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 969f33bbcb6SJed Brown { 970d52bd9f3SBarry Smith PetscErrorCode ierr; 971d52bd9f3SBarry Smith 972f33bbcb6SJed Brown PetscFunctionBegin; 9731566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 974f33bbcb6SJed Brown PetscFunctionReturn(0); 975f33bbcb6SJed Brown } 976f33bbcb6SJed Brown 977f33bbcb6SJed Brown /*MC 978f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 979f33bbcb6SJed Brown 980f33bbcb6SJed Brown Level: beginner 981f33bbcb6SJed Brown 9824eb428fdSBarry Smith Notes: 983c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 9844eb428fdSBarry Smith 9851566a47fSLisandro Dalcin $ -ts_type theta -ts_theta_theta 1.0 9864eb428fdSBarry Smith 987f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 988f33bbcb6SJed Brown 989f33bbcb6SJed Brown M*/ 990f33bbcb6SJed Brown #undef __FUNCT__ 991f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 9928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 993f33bbcb6SJed Brown { 994f33bbcb6SJed Brown PetscErrorCode ierr; 995f33bbcb6SJed Brown 996f33bbcb6SJed Brown PetscFunctionBegin; 997f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 998f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 9991566a47fSLisandro Dalcin ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 10000c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 1001f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 1002f33bbcb6SJed Brown PetscFunctionReturn(0); 1003f33bbcb6SJed Brown } 1004f33bbcb6SJed Brown 1005f33bbcb6SJed Brown #undef __FUNCT__ 10061566a47fSLisandro Dalcin #define __FUNCT__ "TSSetUp_CN" 10071566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts) 10081566a47fSLisandro Dalcin { 10091566a47fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 10101566a47fSLisandro Dalcin PetscErrorCode ierr; 10111566a47fSLisandro Dalcin 10121566a47fSLisandro Dalcin PetscFunctionBegin; 10131566a47fSLisandro 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"); 10141566a47fSLisandro 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"); 10151566a47fSLisandro Dalcin ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 10161566a47fSLisandro Dalcin PetscFunctionReturn(0); 10171566a47fSLisandro Dalcin } 10181566a47fSLisandro Dalcin 10191566a47fSLisandro Dalcin #undef __FUNCT__ 1020f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 1021f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1022f33bbcb6SJed Brown { 1023d52bd9f3SBarry Smith PetscErrorCode ierr; 1024d52bd9f3SBarry Smith 1025f33bbcb6SJed Brown PetscFunctionBegin; 10261566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 1027f33bbcb6SJed Brown PetscFunctionReturn(0); 1028f33bbcb6SJed Brown } 1029f33bbcb6SJed Brown 1030f33bbcb6SJed Brown /*MC 1031f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 1032f33bbcb6SJed Brown 1033f33bbcb6SJed Brown Level: beginner 1034f33bbcb6SJed Brown 1035f33bbcb6SJed Brown Notes: 10367cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 10377cf5af47SJed Brown 10387cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1039f33bbcb6SJed Brown 1040f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1041f33bbcb6SJed Brown 1042f33bbcb6SJed Brown M*/ 1043f33bbcb6SJed Brown #undef __FUNCT__ 1044f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 10458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1046f33bbcb6SJed Brown { 1047f33bbcb6SJed Brown PetscErrorCode ierr; 1048f33bbcb6SJed Brown 1049f33bbcb6SJed Brown PetscFunctionBegin; 1050f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1051f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1052eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 10530c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 1054f33bbcb6SJed Brown ts->ops->view = TSView_CN; 1055f33bbcb6SJed Brown PetscFunctionReturn(0); 1056f33bbcb6SJed Brown } 1057