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 { 10316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 113b1890cdSShri Abhyankar Vec X0; /* work vector to store X0 */ 12eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 13b5e0532dSHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage*/ 14b5e0532dSHong Zhang Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage*/ 15b5e0532dSHong Zhang Vec *VecsSensiTemp; /* Vector to be timed with Jacobian transpose*/ 16ace3abfcSBarry Smith PetscBool extrapolate; 17eb284becSJed Brown PetscBool endpoint; 18316643e7SJed Brown PetscReal Theta; 19316643e7SJed Brown PetscReal stage_time; 203b1890cdSShri Abhyankar TSStepStatus status; 213b1890cdSShri Abhyankar char *name; 223b1890cdSShri Abhyankar PetscInt order; 233b1890cdSShri Abhyankar PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 243b1890cdSShri Abhyankar PetscBool adapt; /* use time-step adaptivity ? */ 25b5e0532dSHong Zhang PetscReal ptime; 26316643e7SJed Brown } TS_Theta; 27316643e7SJed Brown 28316643e7SJed Brown #undef __FUNCT__ 297445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 307445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 317445fe48SJed Brown { 327445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 337445fe48SJed Brown PetscErrorCode ierr; 347445fe48SJed Brown 357445fe48SJed Brown PetscFunctionBegin; 367445fe48SJed Brown if (X0) { 377445fe48SJed Brown if (dm && dm != ts->dm) { 380d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 397445fe48SJed Brown } else *X0 = ts->vec_sol; 407445fe48SJed Brown } 417445fe48SJed Brown if (Xdot) { 427445fe48SJed Brown if (dm && dm != ts->dm) { 430d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 447445fe48SJed Brown } else *Xdot = th->Xdot; 457445fe48SJed Brown } 467445fe48SJed Brown PetscFunctionReturn(0); 477445fe48SJed Brown } 487445fe48SJed Brown 490d0b770aSPeter Brune 500d0b770aSPeter Brune #undef __FUNCT__ 510d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 520d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 530d0b770aSPeter Brune { 540d0b770aSPeter Brune PetscErrorCode ierr; 550d0b770aSPeter Brune 560d0b770aSPeter Brune PetscFunctionBegin; 570d0b770aSPeter Brune if (X0) { 580d0b770aSPeter Brune if (dm && dm != ts->dm) { 590d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 600d0b770aSPeter Brune } 610d0b770aSPeter Brune } 620d0b770aSPeter Brune if (Xdot) { 630d0b770aSPeter Brune if (dm && dm != ts->dm) { 640d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 650d0b770aSPeter Brune } 660d0b770aSPeter Brune } 670d0b770aSPeter Brune PetscFunctionReturn(0); 680d0b770aSPeter Brune } 690d0b770aSPeter Brune 707445fe48SJed Brown #undef __FUNCT__ 717445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 727445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 737445fe48SJed Brown { 747445fe48SJed Brown 757445fe48SJed Brown PetscFunctionBegin; 767445fe48SJed Brown PetscFunctionReturn(0); 777445fe48SJed Brown } 787445fe48SJed Brown 797445fe48SJed Brown #undef __FUNCT__ 807445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 817445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 827445fe48SJed Brown { 837445fe48SJed Brown TS ts = (TS)ctx; 847445fe48SJed Brown PetscErrorCode ierr; 857445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 867445fe48SJed Brown 877445fe48SJed Brown PetscFunctionBegin; 887445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 897445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 907445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 917445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 927445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 937445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 940d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 950d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 967445fe48SJed Brown PetscFunctionReturn(0); 977445fe48SJed Brown } 987445fe48SJed Brown 997445fe48SJed Brown #undef __FUNCT__ 100258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta" 101258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 102258e1594SPeter Brune { 103258e1594SPeter Brune 104258e1594SPeter Brune PetscFunctionBegin; 105258e1594SPeter Brune PetscFunctionReturn(0); 106258e1594SPeter Brune } 107258e1594SPeter Brune 108258e1594SPeter Brune #undef __FUNCT__ 109258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta" 110258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 111258e1594SPeter Brune { 112258e1594SPeter Brune TS ts = (TS)ctx; 113258e1594SPeter Brune PetscErrorCode ierr; 114258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 115258e1594SPeter Brune 116258e1594SPeter Brune PetscFunctionBegin; 117258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 118258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 119258e1594SPeter Brune 120258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 122258e1594SPeter Brune 123258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125258e1594SPeter Brune 126258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 127258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 128258e1594SPeter Brune PetscFunctionReturn(0); 129258e1594SPeter Brune } 130258e1594SPeter Brune 1313b1890cdSShri Abhyankar #undef __FUNCT__ 1323b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta" 1333b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done) 1343b1890cdSShri Abhyankar { 1353b1890cdSShri Abhyankar PetscErrorCode ierr; 1363b1890cdSShri Abhyankar TS_Theta *th = (TS_Theta*)ts->data; 1373b1890cdSShri Abhyankar 1383b1890cdSShri Abhyankar PetscFunctionBegin; 139ce94432eSBarry Smith if (order == 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"No time-step adaptivity implemented for 1st order theta method; Run with -ts_adapt_type none"); 1403b1890cdSShri Abhyankar if (order == th->order) { 1413b1890cdSShri Abhyankar if (th->endpoint) { 1423b1890cdSShri Abhyankar ierr = VecCopy(th->X,U);CHKERRQ(ierr); 1433b1890cdSShri Abhyankar } else { 1443b1890cdSShri Abhyankar PetscReal shift = 1./(th->Theta*ts->time_step); 1453b1890cdSShri Abhyankar ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr); 1463b1890cdSShri Abhyankar ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr); 1473b1890cdSShri Abhyankar } 1483b1890cdSShri Abhyankar } else if (order == th->order-1 && order) { 1493b1890cdSShri Abhyankar ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr); 1503b1890cdSShri Abhyankar } 1513b1890cdSShri Abhyankar PetscFunctionReturn(0); 1523b1890cdSShri Abhyankar } 153258e1594SPeter Brune 154258e1594SPeter Brune #undef __FUNCT__ 15524655328SShri #define __FUNCT__ "TSRollBack_Theta" 15624655328SShri static PetscErrorCode TSRollBack_Theta(TS ts) 15724655328SShri { 15824655328SShri TS_Theta *th = (TS_Theta*)ts->data; 15924655328SShri PetscErrorCode ierr; 16024655328SShri 16124655328SShri PetscFunctionBegin; 16224655328SShri ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 16324655328SShri th->status = TS_STEP_INCOMPLETE; 16424655328SShri PetscFunctionReturn(0); 16524655328SShri } 16624655328SShri 16724655328SShri #undef __FUNCT__ 168316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 169193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 170316643e7SJed Brown { 171316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1723b1890cdSShri Abhyankar PetscInt its,lits,reject,next_scheme; 1733b1890cdSShri Abhyankar PetscReal next_time_step; 1743b1890cdSShri Abhyankar TSAdapt adapt; 1754957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 176051f2191SLisandro Dalcin PetscErrorCode ierr; 177316643e7SJed Brown 178316643e7SJed Brown PetscFunctionBegin; 1793b1890cdSShri Abhyankar th->status = TS_STEP_INCOMPLETE; 1803b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 181051f2191SLisandro Dalcin for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) { 182b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 183eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 184b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 185b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 186316643e7SJed Brown 187eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 188eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 189eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 190eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 191eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 192eb284becSJed Brown } 193ace68cafSJed Brown if (th->extrapolate) { 194b296d7d5SJed Brown ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 195ace68cafSJed Brown } else { 1962b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 197ace68cafSJed Brown } 198eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 199316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 200316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 2015ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 202051f2191SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr); 203552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 2044957b756SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr); 205051f2191SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 206051f2191SLisandro Dalcin 2070298fd71SBarry Smith ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr); 208051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2093b1890cdSShri Abhyankar /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 210552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 2113b1890cdSShri Abhyankar ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 2120298fd71SBarry Smith ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr); 2133b1890cdSShri Abhyankar ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 214051f2191SLisandro Dalcin if (!accept) { /* Roll back the current step */ 215051f2191SLisandro Dalcin ts->ptime += next_time_step; /* This will be undone in rollback */ 216051f2191SLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 217051f2191SLisandro Dalcin ierr = TSRollBack(ts);CHKERRQ(ierr); 218051f2191SLisandro Dalcin goto reject_step; 219051f2191SLisandro Dalcin } 2203b1890cdSShri Abhyankar 22117145e75SHong Zhang if (ts->vec_costintegral) { 22217145e75SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 22317145e75SHong Zhang if (th->endpoint) { 22417145e75SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 22517145e75SHong Zhang ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(1.-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 22617145e75SHong Zhang } 22717145e75SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 22817145e75SHong Zhang if (th->endpoint) { 22917145e75SHong Zhang ierr = VecAXPY(ts->vec_costintegral,ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 23017145e75SHong Zhang }else { 23117145e75SHong Zhang ierr = VecAXPY(ts->vec_costintegral,ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 23217145e75SHong Zhang } 23317145e75SHong Zhang } 23417145e75SHong Zhang 2353b1890cdSShri Abhyankar /* ignore next_scheme for now */ 2362b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 237cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 238316643e7SJed Brown ts->steps++; 2393b1890cdSShri Abhyankar th->status = TS_STEP_COMPLETE; 240051f2191SLisandro Dalcin break; 241051f2191SLisandro Dalcin 242051f2191SLisandro Dalcin reject_step: 243051f2191SLisandro Dalcin if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) { 244051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 245051f2191SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 2463b1890cdSShri Abhyankar } 247051f2191SLisandro Dalcin continue; 2483b1890cdSShri Abhyankar } 249316643e7SJed Brown PetscFunctionReturn(0); 250316643e7SJed Brown } 251316643e7SJed Brown 252cd652676SJed Brown #undef __FUNCT__ 25342f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_Theta" 25442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2552ca6e920SHong Zhang { 2562ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 257b5e0532dSHong Zhang Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 2582ca6e920SHong Zhang PetscInt nadj; 2592ca6e920SHong Zhang PetscErrorCode ierr; 2602ca6e920SHong Zhang Mat J,Jp; 2612ca6e920SHong Zhang KSP ksp; 2622ca6e920SHong Zhang PetscReal shift; 2632ca6e920SHong Zhang 2642ca6e920SHong Zhang PetscFunctionBegin; 2652ca6e920SHong Zhang 2662ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 267302440fdSBarry Smith ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 2682ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 269b5e0532dSHong Zhang 270b5e0532dSHong Zhang /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 2713fd52205SHong Zhang th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/ 272b5e0532dSHong Zhang th->ptime = ts->ptime + ts->time_step; 2732ca6e920SHong Zhang 2742ca6e920SHong Zhang ierr = TSPreStep(ts);CHKERRQ(ierr); 2752ca6e920SHong Zhang 276a4cab896SHong Zhang /* Build RHS */ 2772c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 27805755b9cSHong Zhang if (th->endpoint) { 279d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr); 28005755b9cSHong Zhang }else { 281d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 28205755b9cSHong Zhang } 28305755b9cSHong Zhang } 284abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 285b5e0532dSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 286b5e0532dSHong Zhang ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr); 2872c39e106SBarry Smith if (ts->vec_costintegral) { 288b5e0532dSHong Zhang ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 28936eaed60SHong Zhang } 2902ca6e920SHong Zhang } 2913c54f92cSHong Zhang 2922ca6e920SHong Zhang /* Build LHS */ 2932ca6e920SHong Zhang shift = -1./(th->Theta*ts->time_step); 2943c54f92cSHong Zhang if (th->endpoint) { 2953c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2963c54f92cSHong Zhang }else { 2973c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2983c54f92cSHong Zhang } 2992ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 3002ca6e920SHong Zhang 3012ca6e920SHong Zhang /* Solve LHS X = RHS */ 302abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 303b5e0532dSHong Zhang ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 3042ca6e920SHong Zhang } 3053c54f92cSHong Zhang 30636eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 307b5e0532dSHong Zhang if(th->endpoint) { /* two-stage case */ 308b5e0532dSHong Zhang if (th->Theta!=1.) { 3093fd52205SHong Zhang shift = -1./((th->Theta-1.)*ts->time_step); 310b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3112c39e106SBarry Smith if (ts->vec_costintegral) { 312b5e0532dSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr); 313fb07d950SHong Zhang if (!ts->costintegralfwd) { 3142c39e106SBarry Smith /* Evolve ts->vec_costintegral to compute integrals */ 315d4aa0a58SBarry Smith ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 3162c39e106SBarry Smith ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 317b5e0532dSHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 3182c39e106SBarry Smith ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);CHKERRQ(ierr); 31936eaed60SHong Zhang } 3208263dbbfSHong Zhang } 321abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 322b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3232c39e106SBarry Smith if (ts->vec_costintegral) { 32436eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 32536eaed60SHong Zhang } 3263c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 3272ca6e920SHong Zhang } 328b5e0532dSHong Zhang }else { /* backward Euler */ 329b5e0532dSHong Zhang shift = 0.0; 330b5e0532dSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 331b5e0532dSHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 332b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 333b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 334b5e0532dSHong Zhang if (ts->vec_costintegral) { 335b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 336b5e0532dSHong Zhang if (!ts->costintegralfwd) { 337b5e0532dSHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 338b5e0532dSHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 339b5e0532dSHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 340b5e0532dSHong Zhang } 341b5e0532dSHong Zhang } 342b5e0532dSHong Zhang } 343b5e0532dSHong Zhang } 3443fd52205SHong Zhang 3453fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 3465bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 347abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 348b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 349b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3503fd52205SHong Zhang } 351b5e0532dSHong Zhang if (th->Theta!=1.) { 352b5e0532dSHong Zhang ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,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*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 356b5e0532dSHong Zhang } 3573fd52205SHong Zhang } 3582c39e106SBarry Smith if (ts->vec_costintegral) { 359d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 360abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 36136eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 36236eaed60SHong Zhang } 363b5e0532dSHong Zhang if (th->Theta!=1.) { 364b5e0532dSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 365abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 36636eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 36736eaed60SHong Zhang } 36836eaed60SHong Zhang } 3693fd52205SHong Zhang } 370b5e0532dSHong Zhang } 3713fd52205SHong Zhang }else { /* one-stage case */ 3723c54f92cSHong Zhang shift = 0.0; 373a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 3742c39e106SBarry Smith if (ts->vec_costintegral) { 375d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 376fb07d950SHong Zhang if (!ts->costintegralfwd) { 3772c39e106SBarry Smith /* Evolve ts->vec_costintegral to compute integrals */ 378d4aa0a58SBarry Smith ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 3792c39e106SBarry Smith ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 38036eaed60SHong Zhang } 3818263dbbfSHong Zhang } 382abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 383b5e0532dSHong Zhang ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 384b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 3852c39e106SBarry Smith if (ts->vec_costintegral) { 38636eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 38736eaed60SHong Zhang } 3882ca6e920SHong Zhang } 3893fd52205SHong Zhang if (ts->vecs_sensip) { 3905bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 391abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 392b5e0532dSHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 393b5e0532dSHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 3943fd52205SHong Zhang } 3952c39e106SBarry Smith if (ts->vec_costintegral) { 396d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 397abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 39836eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 39936eaed60SHong Zhang } 40036eaed60SHong Zhang } 4013fd52205SHong Zhang } 4022ca6e920SHong Zhang } 4032ca6e920SHong Zhang 4042ca6e920SHong Zhang ts->ptime += ts->time_step; 4052ca6e920SHong Zhang ts->steps++; 4062ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 4072ca6e920SHong Zhang PetscFunctionReturn(0); 4082ca6e920SHong Zhang } 4092ca6e920SHong Zhang 4102ca6e920SHong Zhang #undef __FUNCT__ 411cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 412cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 413cd652676SJed Brown { 414cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 4155a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 416cd652676SJed Brown PetscErrorCode ierr; 417cd652676SJed Brown 418cd652676SJed Brown PetscFunctionBegin; 419a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 4205a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 4215a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 422cd652676SJed Brown PetscFunctionReturn(0); 423cd652676SJed Brown } 424cd652676SJed Brown 425316643e7SJed Brown /*------------------------------------------------------------*/ 426316643e7SJed Brown #undef __FUNCT__ 427277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 428277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 429316643e7SJed Brown { 430316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 431316643e7SJed Brown PetscErrorCode ierr; 432316643e7SJed Brown 433316643e7SJed Brown PetscFunctionBegin; 4346bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 4356bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 4363b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 437eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 438b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 439b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 440b5e0532dSHong Zhang ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 441277b19d0SLisandro Dalcin PetscFunctionReturn(0); 442277b19d0SLisandro Dalcin } 443277b19d0SLisandro Dalcin 444277b19d0SLisandro Dalcin #undef __FUNCT__ 445277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 446277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 447277b19d0SLisandro Dalcin { 448277b19d0SLisandro Dalcin PetscErrorCode ierr; 449277b19d0SLisandro Dalcin 450277b19d0SLisandro Dalcin PetscFunctionBegin; 451277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 452277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 453bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 454bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 455bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 456bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 457316643e7SJed Brown PetscFunctionReturn(0); 458316643e7SJed Brown } 459316643e7SJed Brown 460316643e7SJed Brown /* 461316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 4622b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 463316643e7SJed Brown */ 464316643e7SJed Brown #undef __FUNCT__ 4650f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 4660f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 467316643e7SJed Brown { 468316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 469316643e7SJed Brown PetscErrorCode ierr; 4707445fe48SJed Brown Vec X0,Xdot; 4717445fe48SJed Brown DM dm,dmsave; 472b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 473316643e7SJed Brown 474316643e7SJed Brown PetscFunctionBegin; 4757445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4765a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 4777445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 478b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 4797445fe48SJed Brown 4807445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 4817445fe48SJed Brown dmsave = ts->dm; 4827445fe48SJed Brown ts->dm = dm; 4837445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 4847445fe48SJed Brown ts->dm = dmsave; 4850d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 486316643e7SJed Brown PetscFunctionReturn(0); 487316643e7SJed Brown } 488316643e7SJed Brown 489316643e7SJed Brown #undef __FUNCT__ 4900f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 491d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 492316643e7SJed Brown { 493316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 494316643e7SJed Brown PetscErrorCode ierr; 4957445fe48SJed Brown Vec Xdot; 4967445fe48SJed Brown DM dm,dmsave; 497b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 498316643e7SJed Brown 499316643e7SJed Brown PetscFunctionBegin; 5007445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 5017445fe48SJed Brown 5020f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 5030298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 5047445fe48SJed Brown 5057445fe48SJed Brown dmsave = ts->dm; 5067445fe48SJed Brown ts->dm = dm; 507d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 5087445fe48SJed Brown ts->dm = dmsave; 5090298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 510316643e7SJed Brown PetscFunctionReturn(0); 511316643e7SJed Brown } 512316643e7SJed Brown 513316643e7SJed Brown #undef __FUNCT__ 514316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 515316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 516316643e7SJed Brown { 517316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 518316643e7SJed Brown PetscErrorCode ierr; 5197445fe48SJed Brown SNES snes; 520ef749922SLisandro Dalcin TSAdapt adapt; 5217445fe48SJed Brown DM dm; 522316643e7SJed Brown 523316643e7SJed Brown PetscFunctionBegin; 52439d13666SHong Zhang if (!th->X) { 525316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 52639d13666SHong Zhang } 52739d13666SHong Zhang if (!th->Xdot) { 528316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 52939d13666SHong Zhang } 53039d13666SHong Zhang if (!th->X0) { 5313b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 53239d13666SHong Zhang } 5337445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 5347445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 5357445fe48SJed Brown if (dm) { 5367445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 537258e1594SPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5387445fe48SJed Brown } 5393b1890cdSShri Abhyankar if (th->Theta == 0.5 && th->endpoint) th->order = 2; 5403b1890cdSShri Abhyankar else th->order = 1; 5413b1890cdSShri Abhyankar 542552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 543ef749922SLisandro Dalcin if (!th->adapt) { 5443b1890cdSShri Abhyankar ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr); 5453b1890cdSShri Abhyankar } 546b4dd3bf9SBarry Smith PetscFunctionReturn(0); 547b4dd3bf9SBarry Smith } 5480c7376e5SHong Zhang 5490c7376e5SHong Zhang #undef __FUNCT__ 5500c7376e5SHong Zhang #define __FUNCT__ "TSSetUp_BEuler" 5510c7376e5SHong Zhang static PetscErrorCode TSSetUp_BEuler(TS ts) 5520c7376e5SHong Zhang { 5530c7376e5SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 5540c7376e5SHong Zhang PetscErrorCode ierr; 5550c7376e5SHong Zhang 5560c7376e5SHong Zhang PetscFunctionBegin; 557*b9eb5ee8SHong Zhang 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"); 5580c7376e5SHong Zhang ierr = TSSetUp_Theta(ts); 5590c7376e5SHong Zhang PetscFunctionReturn(0); 5600c7376e5SHong Zhang } 5610c7376e5SHong Zhang 5620c7376e5SHong Zhang #undef __FUNCT__ 5630c7376e5SHong Zhang #define __FUNCT__ "TSSetUp_CN" 5640c7376e5SHong Zhang static PetscErrorCode TSSetUp_CN(TS ts) 5650c7376e5SHong Zhang { 5660c7376e5SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 5670c7376e5SHong Zhang PetscErrorCode ierr; 5680c7376e5SHong Zhang 5690c7376e5SHong Zhang PetscFunctionBegin; 570*b9eb5ee8SHong Zhang 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"); 571*b9eb5ee8SHong Zhang 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"); 5720c7376e5SHong Zhang ierr = TSSetUp_Theta(ts); 5730c7376e5SHong Zhang PetscFunctionReturn(0); 5740c7376e5SHong Zhang } 575b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 576b4dd3bf9SBarry Smith 577b4dd3bf9SBarry Smith #undef __FUNCT__ 578b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta" 579b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 580b4dd3bf9SBarry Smith { 581b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 582b4dd3bf9SBarry Smith PetscErrorCode ierr; 583b4dd3bf9SBarry Smith 584b4dd3bf9SBarry Smith PetscFunctionBegin; 585b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 5862ca6e920SHong Zhang if(ts->vecs_sensip) { 587b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 5882ca6e920SHong Zhang } 589b5e0532dSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 590316643e7SJed Brown PetscFunctionReturn(0); 591316643e7SJed Brown } 592316643e7SJed Brown /*------------------------------------------------------------*/ 593316643e7SJed Brown 594316643e7SJed Brown #undef __FUNCT__ 595316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 5968c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts) 597316643e7SJed Brown { 598316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 599316643e7SJed Brown PetscErrorCode ierr; 600316643e7SJed Brown 601316643e7SJed Brown PetscFunctionBegin; 602e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 603316643e7SJed Brown { 6040298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 6050298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 6060298fd71SBarry 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); 6070298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 608d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 609316643e7SJed Brown } 610316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 611316643e7SJed Brown PetscFunctionReturn(0); 612316643e7SJed Brown } 613316643e7SJed Brown 614316643e7SJed Brown #undef __FUNCT__ 615316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 616316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 617316643e7SJed Brown { 618316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 619ace3abfcSBarry Smith PetscBool iascii; 620316643e7SJed Brown PetscErrorCode ierr; 621316643e7SJed Brown 622316643e7SJed Brown PetscFunctionBegin; 623251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 624316643e7SJed Brown if (iascii) { 6257c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 626ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 627316643e7SJed Brown } 628ac75fa18SLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 629316643e7SJed Brown PetscFunctionReturn(0); 630316643e7SJed Brown } 631316643e7SJed Brown 6320de4c49aSJed Brown #undef __FUNCT__ 6330de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 6347087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 6350de4c49aSJed Brown { 6360de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6370de4c49aSJed Brown 6380de4c49aSJed Brown PetscFunctionBegin; 6390de4c49aSJed Brown *theta = th->Theta; 6400de4c49aSJed Brown PetscFunctionReturn(0); 6410de4c49aSJed Brown } 6420de4c49aSJed Brown 6430de4c49aSJed Brown #undef __FUNCT__ 6440de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 6457087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 6460de4c49aSJed Brown { 6470de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6480de4c49aSJed Brown 6490de4c49aSJed Brown PetscFunctionBegin; 6507c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 6510de4c49aSJed Brown th->Theta = theta; 6520de4c49aSJed Brown PetscFunctionReturn(0); 6530de4c49aSJed Brown } 654eb284becSJed Brown 655eb284becSJed Brown #undef __FUNCT__ 65678e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 65726f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 65826f2ff8fSLisandro Dalcin { 65926f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 66026f2ff8fSLisandro Dalcin 66126f2ff8fSLisandro Dalcin PetscFunctionBegin; 66226f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 66326f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 66426f2ff8fSLisandro Dalcin } 66526f2ff8fSLisandro Dalcin 66626f2ff8fSLisandro Dalcin #undef __FUNCT__ 66726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 668eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 669eb284becSJed Brown { 670eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 671eb284becSJed Brown 672eb284becSJed Brown PetscFunctionBegin; 673eb284becSJed Brown th->endpoint = flg; 674eb284becSJed Brown PetscFunctionReturn(0); 675eb284becSJed Brown } 6760de4c49aSJed Brown 677f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 678f9c1d6abSBarry Smith #undef __FUNCT__ 679f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 680f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 681f9c1d6abSBarry Smith { 682f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 683f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 6843fd8ae06SJed Brown const PetscReal one = 1.0; 685f9c1d6abSBarry Smith 686f9c1d6abSBarry Smith PetscFunctionBegin; 6873fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 688f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 689f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 690f9c1d6abSBarry Smith PetscFunctionReturn(0); 691f9c1d6abSBarry Smith } 692f9c1d6abSBarry Smith #endif 693f9c1d6abSBarry Smith 69442682096SHong Zhang #undef __FUNCT__ 69542682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta" 69642682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 69742682096SHong Zhang { 69842682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 69942682096SHong Zhang 70042682096SHong Zhang PetscFunctionBegin; 7012ca6e920SHong Zhang *ns = 1; 7022ca6e920SHong Zhang if(Y) { 703b5e0532dSHong Zhang *Y = (th->endpoint)?&(th->X0):&(th->X); 7042ca6e920SHong Zhang } 70542682096SHong Zhang PetscFunctionReturn(0); 70642682096SHong Zhang } 707f9c1d6abSBarry Smith 708316643e7SJed Brown /* ------------------------------------------------------------ */ 709316643e7SJed Brown /*MC 71096f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 711316643e7SJed Brown 712316643e7SJed Brown Level: beginner 713316643e7SJed Brown 7144eb428fdSBarry Smith Options Database: 7150c3ba866SJed Brown -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 7164eb428fdSBarry Smith -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 7170c3ba866SJed Brown -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 7184eb428fdSBarry Smith 719eb284becSJed Brown Notes: 7200c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 7210c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 7224eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 7234eb428fdSBarry Smith 7244eb428fdSBarry Smith 7254eb428fdSBarry Smith 726eb284becSJed Brown This method can be applied to DAE. 727eb284becSJed Brown 728eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 729eb284becSJed Brown 730eb284becSJed Brown .vb 731eb284becSJed Brown Theta | Theta 732eb284becSJed Brown ------------- 733eb284becSJed Brown | 1 734eb284becSJed Brown .ve 735eb284becSJed Brown 736eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 737eb284becSJed Brown 738eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 739eb284becSJed Brown 740eb284becSJed Brown .vb 741eb284becSJed Brown 0 | 0 0 742eb284becSJed Brown 1 | 1-Theta Theta 743eb284becSJed Brown ------------------- 744eb284becSJed Brown | 1-Theta Theta 745eb284becSJed Brown .ve 746eb284becSJed Brown 747eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 748eb284becSJed Brown 749eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 750eb284becSJed Brown 751eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 752eb284becSJed Brown 7534eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 754eb284becSJed Brown 755eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 756316643e7SJed Brown 757316643e7SJed Brown M*/ 758316643e7SJed Brown #undef __FUNCT__ 759316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 7608cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 761316643e7SJed Brown { 762316643e7SJed Brown TS_Theta *th; 763316643e7SJed Brown PetscErrorCode ierr; 764316643e7SJed Brown 765316643e7SJed Brown PetscFunctionBegin; 766277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 767316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 768316643e7SJed Brown ts->ops->view = TSView_Theta; 769316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 77042f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 771316643e7SJed Brown ts->ops->step = TSStep_Theta; 772cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 7733b1890cdSShri Abhyankar ts->ops->evaluatestep = TSEvaluateStep_Theta; 77424655328SShri ts->ops->rollback = TSRollBack_Theta; 775316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 7760f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 7770f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 778f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 779f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 780f9c1d6abSBarry Smith #endif 78142682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 78242f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 783316643e7SJed Brown 784b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 785316643e7SJed Brown ts->data = (void*)th; 786316643e7SJed Brown 7876f700aefSJed Brown th->extrapolate = PETSC_FALSE; 788316643e7SJed Brown th->Theta = 0.5; 7893b1890cdSShri Abhyankar th->ccfl = 1.0; 7903b1890cdSShri Abhyankar th->adapt = PETSC_FALSE; 791bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 792bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 793bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 794bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 795316643e7SJed Brown PetscFunctionReturn(0); 796316643e7SJed Brown } 7970de4c49aSJed Brown 7980de4c49aSJed Brown #undef __FUNCT__ 7990de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 8000de4c49aSJed Brown /*@ 8010de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 8020de4c49aSJed Brown 8030de4c49aSJed Brown Not Collective 8040de4c49aSJed Brown 8050de4c49aSJed Brown Input Parameter: 8060de4c49aSJed Brown . ts - timestepping context 8070de4c49aSJed Brown 8080de4c49aSJed Brown Output Parameter: 8090de4c49aSJed Brown . theta - stage abscissa 8100de4c49aSJed Brown 8110de4c49aSJed Brown Note: 8120de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 8130de4c49aSJed Brown 8140de4c49aSJed Brown Level: Advanced 8150de4c49aSJed Brown 8160de4c49aSJed Brown .seealso: TSThetaSetTheta() 8170de4c49aSJed Brown @*/ 8187087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 8190de4c49aSJed Brown { 8204ac538c5SBarry Smith PetscErrorCode ierr; 8210de4c49aSJed Brown 8220de4c49aSJed Brown PetscFunctionBegin; 823afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8240de4c49aSJed Brown PetscValidPointer(theta,2); 8254ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 8260de4c49aSJed Brown PetscFunctionReturn(0); 8270de4c49aSJed Brown } 8280de4c49aSJed Brown 8290de4c49aSJed Brown #undef __FUNCT__ 8300de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 8310de4c49aSJed Brown /*@ 8320de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 8330de4c49aSJed Brown 8340de4c49aSJed Brown Not Collective 8350de4c49aSJed Brown 8360de4c49aSJed Brown Input Parameter: 8370de4c49aSJed Brown + ts - timestepping context 8380de4c49aSJed Brown - theta - stage abscissa 8390de4c49aSJed Brown 8400de4c49aSJed Brown Options Database: 8410de4c49aSJed Brown . -ts_theta_theta <theta> 8420de4c49aSJed Brown 8430de4c49aSJed Brown Level: Intermediate 8440de4c49aSJed Brown 8450de4c49aSJed Brown .seealso: TSThetaGetTheta() 8460de4c49aSJed Brown @*/ 8477087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 8480de4c49aSJed Brown { 8494ac538c5SBarry Smith PetscErrorCode ierr; 8500de4c49aSJed Brown 8510de4c49aSJed Brown PetscFunctionBegin; 852afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8534ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 8540de4c49aSJed Brown PetscFunctionReturn(0); 8550de4c49aSJed Brown } 856f33bbcb6SJed Brown 857eb284becSJed Brown #undef __FUNCT__ 85826f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 85926f2ff8fSLisandro Dalcin /*@ 86026f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 86126f2ff8fSLisandro Dalcin 86226f2ff8fSLisandro Dalcin Not Collective 86326f2ff8fSLisandro Dalcin 86426f2ff8fSLisandro Dalcin Input Parameter: 86526f2ff8fSLisandro Dalcin . ts - timestepping context 86626f2ff8fSLisandro Dalcin 86726f2ff8fSLisandro Dalcin Output Parameter: 86826f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 86926f2ff8fSLisandro Dalcin 87026f2ff8fSLisandro Dalcin Level: Advanced 87126f2ff8fSLisandro Dalcin 87226f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 87326f2ff8fSLisandro Dalcin @*/ 87426f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 87526f2ff8fSLisandro Dalcin { 87626f2ff8fSLisandro Dalcin PetscErrorCode ierr; 87726f2ff8fSLisandro Dalcin 87826f2ff8fSLisandro Dalcin PetscFunctionBegin; 87926f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 88026f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 88126f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 88226f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 88326f2ff8fSLisandro Dalcin } 88426f2ff8fSLisandro Dalcin 88526f2ff8fSLisandro Dalcin #undef __FUNCT__ 886eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 887eb284becSJed Brown /*@ 888eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 889eb284becSJed Brown 890eb284becSJed Brown Not Collective 891eb284becSJed Brown 892eb284becSJed Brown Input Parameter: 893eb284becSJed Brown + ts - timestepping context 894eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 895eb284becSJed Brown 896eb284becSJed Brown Options Database: 897eb284becSJed Brown . -ts_theta_endpoint <flg> 898eb284becSJed Brown 899eb284becSJed Brown Level: Intermediate 900eb284becSJed Brown 901eb284becSJed Brown .seealso: TSTHETA, TSCN 902eb284becSJed Brown @*/ 903eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 904eb284becSJed Brown { 905eb284becSJed Brown PetscErrorCode ierr; 906eb284becSJed Brown 907eb284becSJed Brown PetscFunctionBegin; 908eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 909eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 910eb284becSJed Brown PetscFunctionReturn(0); 911eb284becSJed Brown } 912eb284becSJed Brown 913f33bbcb6SJed Brown /* 914f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 915f33bbcb6SJed Brown * The creation functions for these specializations are below. 916f33bbcb6SJed Brown */ 917f33bbcb6SJed Brown 918f33bbcb6SJed Brown #undef __FUNCT__ 919f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 920f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 921f33bbcb6SJed Brown { 922d52bd9f3SBarry Smith PetscErrorCode ierr; 923d52bd9f3SBarry Smith 924f33bbcb6SJed Brown PetscFunctionBegin; 925d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 926f33bbcb6SJed Brown PetscFunctionReturn(0); 927f33bbcb6SJed Brown } 928f33bbcb6SJed Brown 929f33bbcb6SJed Brown /*MC 930f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 931f33bbcb6SJed Brown 932f33bbcb6SJed Brown Level: beginner 933f33bbcb6SJed Brown 9344eb428fdSBarry Smith Notes: 935c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 9364eb428fdSBarry Smith 9374eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 1. 9384eb428fdSBarry Smith 939f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 940f33bbcb6SJed Brown 941f33bbcb6SJed Brown M*/ 942f33bbcb6SJed Brown #undef __FUNCT__ 943f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 9448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 945f33bbcb6SJed Brown { 946f33bbcb6SJed Brown PetscErrorCode ierr; 947f33bbcb6SJed Brown 948f33bbcb6SJed Brown PetscFunctionBegin; 949f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 950f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 9510c7376e5SHong Zhang ts->ops->setup = TSSetUp_BEuler; 952f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 953f33bbcb6SJed Brown PetscFunctionReturn(0); 954f33bbcb6SJed Brown } 955f33bbcb6SJed Brown 956f33bbcb6SJed Brown #undef __FUNCT__ 957f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 958f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 959f33bbcb6SJed Brown { 960d52bd9f3SBarry Smith PetscErrorCode ierr; 961d52bd9f3SBarry Smith 962f33bbcb6SJed Brown PetscFunctionBegin; 963d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 964f33bbcb6SJed Brown PetscFunctionReturn(0); 965f33bbcb6SJed Brown } 966f33bbcb6SJed Brown 967f33bbcb6SJed Brown /*MC 968f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 969f33bbcb6SJed Brown 970f33bbcb6SJed Brown Level: beginner 971f33bbcb6SJed Brown 972f33bbcb6SJed Brown Notes: 9737cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 9747cf5af47SJed Brown 9757cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 976f33bbcb6SJed Brown 977f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 978f33bbcb6SJed Brown 979f33bbcb6SJed Brown M*/ 980f33bbcb6SJed Brown #undef __FUNCT__ 981f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 9828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 983f33bbcb6SJed Brown { 984f33bbcb6SJed Brown PetscErrorCode ierr; 985f33bbcb6SJed Brown 986f33bbcb6SJed Brown PetscFunctionBegin; 987f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 988f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 989eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 9900c7376e5SHong Zhang ts->ops->setup = TSSetUp_CN; 991f33bbcb6SJed Brown ts->ops->view = TSView_CN; 992f33bbcb6SJed Brown PetscFunctionReturn(0); 993f33bbcb6SJed Brown } 994