1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4b45d2f2cSJed Brown #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 */ 132ca6e920SHong Zhang Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage*/ 142ca6e920SHong Zhang Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage*/ 152ca6e920SHong Zhang Vec *VecSensiTemp; /* 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 ? */ 25316643e7SJed Brown } TS_Theta; 26316643e7SJed Brown 27316643e7SJed Brown #undef __FUNCT__ 287445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 297445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 307445fe48SJed Brown { 317445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 327445fe48SJed Brown PetscErrorCode ierr; 337445fe48SJed Brown 347445fe48SJed Brown PetscFunctionBegin; 357445fe48SJed Brown if (X0) { 367445fe48SJed Brown if (dm && dm != ts->dm) { 370d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 387445fe48SJed Brown } else *X0 = ts->vec_sol; 397445fe48SJed Brown } 407445fe48SJed Brown if (Xdot) { 417445fe48SJed Brown if (dm && dm != ts->dm) { 420d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 437445fe48SJed Brown } else *Xdot = th->Xdot; 447445fe48SJed Brown } 457445fe48SJed Brown PetscFunctionReturn(0); 467445fe48SJed Brown } 477445fe48SJed Brown 480d0b770aSPeter Brune 490d0b770aSPeter Brune #undef __FUNCT__ 500d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 510d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 520d0b770aSPeter Brune { 530d0b770aSPeter Brune PetscErrorCode ierr; 540d0b770aSPeter Brune 550d0b770aSPeter Brune PetscFunctionBegin; 560d0b770aSPeter Brune if (X0) { 570d0b770aSPeter Brune if (dm && dm != ts->dm) { 580d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 590d0b770aSPeter Brune } 600d0b770aSPeter Brune } 610d0b770aSPeter Brune if (Xdot) { 620d0b770aSPeter Brune if (dm && dm != ts->dm) { 630d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 640d0b770aSPeter Brune } 650d0b770aSPeter Brune } 660d0b770aSPeter Brune PetscFunctionReturn(0); 670d0b770aSPeter Brune } 680d0b770aSPeter Brune 697445fe48SJed Brown #undef __FUNCT__ 707445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 717445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 727445fe48SJed Brown { 737445fe48SJed Brown 747445fe48SJed Brown PetscFunctionBegin; 757445fe48SJed Brown PetscFunctionReturn(0); 767445fe48SJed Brown } 777445fe48SJed Brown 787445fe48SJed Brown #undef __FUNCT__ 797445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 807445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 817445fe48SJed Brown { 827445fe48SJed Brown TS ts = (TS)ctx; 837445fe48SJed Brown PetscErrorCode ierr; 847445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 857445fe48SJed Brown 867445fe48SJed Brown PetscFunctionBegin; 877445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 887445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 897445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 907445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 917445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 927445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 930d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 940d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 957445fe48SJed Brown PetscFunctionReturn(0); 967445fe48SJed Brown } 977445fe48SJed Brown 987445fe48SJed Brown #undef __FUNCT__ 99258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta" 100258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 101258e1594SPeter Brune { 102258e1594SPeter Brune 103258e1594SPeter Brune PetscFunctionBegin; 104258e1594SPeter Brune PetscFunctionReturn(0); 105258e1594SPeter Brune } 106258e1594SPeter Brune 107258e1594SPeter Brune #undef __FUNCT__ 108258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta" 109258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 110258e1594SPeter Brune { 111258e1594SPeter Brune TS ts = (TS)ctx; 112258e1594SPeter Brune PetscErrorCode ierr; 113258e1594SPeter Brune Vec X0,Xdot,X0_sub,Xdot_sub; 114258e1594SPeter Brune 115258e1594SPeter Brune PetscFunctionBegin; 116258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 117258e1594SPeter Brune ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 118258e1594SPeter Brune 119258e1594SPeter Brune ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 120258e1594SPeter Brune ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121258e1594SPeter Brune 122258e1594SPeter Brune ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123258e1594SPeter Brune ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124258e1594SPeter Brune 125258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 126258e1594SPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 127258e1594SPeter Brune PetscFunctionReturn(0); 128258e1594SPeter Brune } 129258e1594SPeter Brune 1303b1890cdSShri Abhyankar #undef __FUNCT__ 1313b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta" 1323b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done) 1333b1890cdSShri Abhyankar { 1343b1890cdSShri Abhyankar PetscErrorCode ierr; 1353b1890cdSShri Abhyankar TS_Theta *th = (TS_Theta*)ts->data; 1363b1890cdSShri Abhyankar 1373b1890cdSShri Abhyankar PetscFunctionBegin; 138ce94432eSBarry 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"); 1393b1890cdSShri Abhyankar if (order == th->order) { 1403b1890cdSShri Abhyankar if (th->endpoint) { 1413b1890cdSShri Abhyankar ierr = VecCopy(th->X,U);CHKERRQ(ierr); 1423b1890cdSShri Abhyankar } else { 1433b1890cdSShri Abhyankar PetscReal shift = 1./(th->Theta*ts->time_step); 1443b1890cdSShri Abhyankar ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr); 1453b1890cdSShri Abhyankar ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr); 1463b1890cdSShri Abhyankar } 1473b1890cdSShri Abhyankar } else if (order == th->order-1 && order) { 1483b1890cdSShri Abhyankar ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr); 1493b1890cdSShri Abhyankar } 1503b1890cdSShri Abhyankar PetscFunctionReturn(0); 1513b1890cdSShri Abhyankar } 152258e1594SPeter Brune 153258e1594SPeter Brune #undef __FUNCT__ 15424655328SShri #define __FUNCT__ "TSRollBack_Theta" 15524655328SShri static PetscErrorCode TSRollBack_Theta(TS ts) 15624655328SShri { 15724655328SShri TS_Theta *th = (TS_Theta*)ts->data; 15824655328SShri PetscErrorCode ierr; 15924655328SShri 16024655328SShri PetscFunctionBegin; 16124655328SShri ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 16224655328SShri th->status = TS_STEP_INCOMPLETE; 16324655328SShri PetscFunctionReturn(0); 16424655328SShri } 16524655328SShri 16624655328SShri #undef __FUNCT__ 167316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 168193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 169316643e7SJed Brown { 170316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1713b1890cdSShri Abhyankar PetscInt its,lits,reject,next_scheme; 1723b1890cdSShri Abhyankar PetscReal next_time_step; 1733b1890cdSShri Abhyankar TSAdapt adapt; 1744957b756SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 175051f2191SLisandro Dalcin PetscErrorCode ierr; 176316643e7SJed Brown 177316643e7SJed Brown PetscFunctionBegin; 1783b1890cdSShri Abhyankar th->status = TS_STEP_INCOMPLETE; 1793b1890cdSShri Abhyankar ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 180051f2191SLisandro Dalcin for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) { 181b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 182eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 183b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 184b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 185316643e7SJed Brown 186eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 187eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 188eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 189eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 190eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 191eb284becSJed Brown } 192ace68cafSJed Brown if (th->extrapolate) { 193b296d7d5SJed Brown ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 194ace68cafSJed Brown } else { 1952b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 196ace68cafSJed Brown } 197eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 198316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 199316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 2005ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 201051f2191SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr); 202552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 2034957b756SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr); 204051f2191SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 205051f2191SLisandro Dalcin 2060298fd71SBarry Smith ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr); 207051f2191SLisandro Dalcin th->status = TS_STEP_PENDING; 2083b1890cdSShri Abhyankar /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 209552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 2103b1890cdSShri Abhyankar ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 2110298fd71SBarry Smith ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr); 2123b1890cdSShri Abhyankar ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 213051f2191SLisandro Dalcin if (!accept) { /* Roll back the current step */ 214051f2191SLisandro Dalcin ts->ptime += next_time_step; /* This will be undone in rollback */ 215051f2191SLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 216051f2191SLisandro Dalcin ierr = TSRollBack(ts);CHKERRQ(ierr); 217051f2191SLisandro Dalcin goto reject_step; 218051f2191SLisandro Dalcin } 2193b1890cdSShri Abhyankar 22017145e75SHong Zhang if (ts->vec_costintegral) { 22117145e75SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 22217145e75SHong Zhang if (th->endpoint) { 22317145e75SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 22417145e75SHong Zhang ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(1.-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 22517145e75SHong Zhang } 22617145e75SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 22717145e75SHong Zhang if (th->endpoint) { 22817145e75SHong Zhang ierr = VecAXPY(ts->vec_costintegral,ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 22917145e75SHong Zhang }else { 23017145e75SHong Zhang ierr = VecAXPY(ts->vec_costintegral,ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 23117145e75SHong Zhang } 23217145e75SHong Zhang } 23317145e75SHong Zhang 2343b1890cdSShri Abhyankar /* ignore next_scheme for now */ 2352b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 236cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 237316643e7SJed Brown ts->steps++; 2383b1890cdSShri Abhyankar th->status = TS_STEP_COMPLETE; 239051f2191SLisandro Dalcin break; 240051f2191SLisandro Dalcin 241051f2191SLisandro Dalcin reject_step: 242051f2191SLisandro Dalcin if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) { 243051f2191SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 244051f2191SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 2453b1890cdSShri Abhyankar } 246051f2191SLisandro Dalcin continue; 2473b1890cdSShri Abhyankar } 248316643e7SJed Brown PetscFunctionReturn(0); 249316643e7SJed Brown } 250316643e7SJed Brown 251cd652676SJed Brown #undef __FUNCT__ 25242f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_Theta" 25342f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts) 2542ca6e920SHong Zhang { 2552ca6e920SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 2563c54f92cSHong Zhang Vec *VecDeltaLam = th->VecDeltaLam,*VecDeltaMu = th->VecDeltaMu,*VecSensiTemp = th->VecSensiTemp; 2572ca6e920SHong Zhang PetscInt nadj; 2582ca6e920SHong Zhang PetscErrorCode ierr; 2592ca6e920SHong Zhang Mat J,Jp; 2602ca6e920SHong Zhang KSP ksp; 2612ca6e920SHong Zhang PetscReal shift; 2622ca6e920SHong Zhang 2632ca6e920SHong Zhang PetscFunctionBegin; 2642ca6e920SHong Zhang 2652ca6e920SHong Zhang th->status = TS_STEP_INCOMPLETE; 2662ca6e920SHong Zhang ierr = SNESGetKSP(ts->snes,&ksp); 2672ca6e920SHong Zhang ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 2683fd52205SHong Zhang th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/ 2692ca6e920SHong Zhang 2702ca6e920SHong Zhang ierr = TSPreStep(ts);CHKERRQ(ierr); 2712ca6e920SHong Zhang 272a4cab896SHong Zhang /* Build RHS */ 2732c39e106SBarry Smith if (ts->vec_costintegral) { /* Cost function has an integral term */ 27405755b9cSHong Zhang if (th->endpoint) { 275d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr); 27605755b9cSHong Zhang }else { 277d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 27805755b9cSHong Zhang } 27905755b9cSHong Zhang } 280abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 2812ca6e920SHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 2822ca6e920SHong Zhang ierr = VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr); 2832c39e106SBarry Smith if (ts->vec_costintegral) { 28436eaed60SHong Zhang ierr = VecAXPY(VecSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 28536eaed60SHong Zhang } 2862ca6e920SHong Zhang } 2873c54f92cSHong Zhang 2882ca6e920SHong Zhang /* Build LHS */ 2892ca6e920SHong Zhang shift = -1./(th->Theta*ts->time_step); 2903c54f92cSHong Zhang if (th->endpoint) { 2913c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2923c54f92cSHong Zhang }else { 2933c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 2943c54f92cSHong Zhang } 2952ca6e920SHong Zhang ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 2962ca6e920SHong Zhang 2972ca6e920SHong Zhang /* Solve LHS X = RHS */ 298abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 2992ca6e920SHong Zhang ierr = KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);CHKERRQ(ierr); 3002ca6e920SHong Zhang } 3013c54f92cSHong Zhang 30236eaed60SHong Zhang /* Update sensitivities, and evaluate integrals if there is any */ 3033fd52205SHong Zhang if(th->endpoint && th->Theta!=1.) { /* two-stage case */ 3043fd52205SHong Zhang shift = -1./((th->Theta-1.)*ts->time_step); 3053c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3062c39e106SBarry Smith if (ts->vec_costintegral) { 307d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 30817145e75SHong Zhang } 309*fb07d950SHong Zhang if (!ts->costintegralfwd) { 3102c39e106SBarry Smith /* Evolve ts->vec_costintegral to compute integrals */ 311d4aa0a58SBarry Smith ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 3122c39e106SBarry Smith ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 313d4aa0a58SBarry Smith ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 3142c39e106SBarry Smith ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);CHKERRQ(ierr); 31536eaed60SHong Zhang } 316abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 3173c54f92cSHong Zhang ierr = MatMultTranspose(J,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 3182c39e106SBarry Smith if (ts->vec_costintegral) { 31936eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 32036eaed60SHong Zhang } 3213c54f92cSHong Zhang ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 3222ca6e920SHong Zhang } 3233fd52205SHong Zhang 3243fd52205SHong Zhang if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 3255bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 326abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 3273fd52205SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr); 3283fd52205SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecDeltaMu[nadj]);CHKERRQ(ierr); 3293fd52205SHong Zhang } 3305bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 331abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 3323fd52205SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr); 3333fd52205SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecDeltaMu[nadj]);CHKERRQ(ierr); 3343fd52205SHong Zhang } 3352c39e106SBarry Smith if (ts->vec_costintegral) { 336d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 337abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 33836eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 33936eaed60SHong Zhang } 340d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 341abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 34236eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 34336eaed60SHong Zhang } 34436eaed60SHong Zhang } 3453fd52205SHong Zhang } 3463fd52205SHong Zhang }else { /* one-stage case */ 3473c54f92cSHong Zhang shift = 0.0; 348a2ae3dd2SHong Zhang ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 3492c39e106SBarry Smith if (ts->vec_costintegral) { 350d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr); 35117145e75SHong Zhang } 352*fb07d950SHong Zhang if (!ts->costintegralfwd) { 3532c39e106SBarry Smith /* Evolve ts->vec_costintegral to compute integrals */ 354d4aa0a58SBarry Smith ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 3552c39e106SBarry Smith ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 35636eaed60SHong Zhang } 3573fd52205SHong Zhang /* When th->endpoint is true and th->Theta==1 (beuler method), the Jacobian is supposed to be evaluated at ts->ptime like this: 3583c54f92cSHong Zhang if(th->endpoint) { 3593c54f92cSHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 3603c54f92cSHong Zhang } 3613fd52205SHong Zhang but ts->ptime and ts->vec_sol have the same values as th->stage_time and th->X in this case. So the code is simplified here. 3623fd52205SHong Zhang */ 363abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 3643c54f92cSHong Zhang ierr = MatMultTranspose(J,VecDeltaLam[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 3653c54f92cSHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr); 3662c39e106SBarry Smith if (ts->vec_costintegral) { 36736eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr); 36836eaed60SHong Zhang } 3692ca6e920SHong Zhang } 3703fd52205SHong Zhang if (ts->vecs_sensip) { 3715bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 372abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 3733fd52205SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr); 3743fd52205SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecDeltaMu[nadj]);CHKERRQ(ierr); 3753fd52205SHong Zhang } 3762c39e106SBarry Smith if (ts->vec_costintegral) { 377d4aa0a58SBarry Smith ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 378abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 37936eaed60SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 38036eaed60SHong Zhang } 38136eaed60SHong Zhang } 3823fd52205SHong Zhang } 3832ca6e920SHong Zhang } 3842ca6e920SHong Zhang 3852ca6e920SHong Zhang ts->ptime += ts->time_step; 3862ca6e920SHong Zhang ts->steps++; 3872ca6e920SHong Zhang th->status = TS_STEP_COMPLETE; 3882ca6e920SHong Zhang PetscFunctionReturn(0); 3892ca6e920SHong Zhang } 3902ca6e920SHong Zhang 3912ca6e920SHong Zhang #undef __FUNCT__ 392cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 393cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 394cd652676SJed Brown { 395cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3965a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 397cd652676SJed Brown PetscErrorCode ierr; 398cd652676SJed Brown 399cd652676SJed Brown PetscFunctionBegin; 400a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 4015a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 4025a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 403cd652676SJed Brown PetscFunctionReturn(0); 404cd652676SJed Brown } 405cd652676SJed Brown 406316643e7SJed Brown /*------------------------------------------------------------*/ 407316643e7SJed Brown #undef __FUNCT__ 408277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 409277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 410316643e7SJed Brown { 411316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 412316643e7SJed Brown PetscErrorCode ierr; 413316643e7SJed Brown 414316643e7SJed Brown PetscFunctionBegin; 4156bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 4166bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 4173b1890cdSShri Abhyankar ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 418eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 419abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&th->VecDeltaLam);CHKERRQ(ierr); 420abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&th->VecDeltaMu);CHKERRQ(ierr); 421abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&th->VecSensiTemp);CHKERRQ(ierr); 422277b19d0SLisandro Dalcin PetscFunctionReturn(0); 423277b19d0SLisandro Dalcin } 424277b19d0SLisandro Dalcin 425277b19d0SLisandro Dalcin #undef __FUNCT__ 426277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 427277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 428277b19d0SLisandro Dalcin { 429277b19d0SLisandro Dalcin PetscErrorCode ierr; 430277b19d0SLisandro Dalcin 431277b19d0SLisandro Dalcin PetscFunctionBegin; 432277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 433277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 434bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 435bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 436bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 437bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 438316643e7SJed Brown PetscFunctionReturn(0); 439316643e7SJed Brown } 440316643e7SJed Brown 441316643e7SJed Brown /* 442316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 4432b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 444316643e7SJed Brown */ 445316643e7SJed Brown #undef __FUNCT__ 4460f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 4470f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 448316643e7SJed Brown { 449316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 450316643e7SJed Brown PetscErrorCode ierr; 4517445fe48SJed Brown Vec X0,Xdot; 4527445fe48SJed Brown DM dm,dmsave; 453b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 454316643e7SJed Brown 455316643e7SJed Brown PetscFunctionBegin; 4567445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4575a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 4587445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 459b296d7d5SJed Brown ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 4607445fe48SJed Brown 4617445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 4627445fe48SJed Brown dmsave = ts->dm; 4637445fe48SJed Brown ts->dm = dm; 4647445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 4657445fe48SJed Brown ts->dm = dmsave; 4660d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 467316643e7SJed Brown PetscFunctionReturn(0); 468316643e7SJed Brown } 469316643e7SJed Brown 470316643e7SJed Brown #undef __FUNCT__ 4710f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 472d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 473316643e7SJed Brown { 474316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 475316643e7SJed Brown PetscErrorCode ierr; 4767445fe48SJed Brown Vec Xdot; 4777445fe48SJed Brown DM dm,dmsave; 478b296d7d5SJed Brown PetscReal shift = 1./(th->Theta*ts->time_step); 479316643e7SJed Brown 480316643e7SJed Brown PetscFunctionBegin; 4817445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4827445fe48SJed Brown 4830f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 4840298fd71SBarry Smith ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 4857445fe48SJed Brown 4867445fe48SJed Brown dmsave = ts->dm; 4877445fe48SJed Brown ts->dm = dm; 488d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 4897445fe48SJed Brown ts->dm = dmsave; 4900298fd71SBarry Smith ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 491316643e7SJed Brown PetscFunctionReturn(0); 492316643e7SJed Brown } 493316643e7SJed Brown 494316643e7SJed Brown #undef __FUNCT__ 495316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 496316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 497316643e7SJed Brown { 498316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 499316643e7SJed Brown PetscErrorCode ierr; 5007445fe48SJed Brown SNES snes; 501ef749922SLisandro Dalcin TSAdapt adapt; 5027445fe48SJed Brown DM dm; 503316643e7SJed Brown 504316643e7SJed Brown PetscFunctionBegin; 50539d13666SHong Zhang if (!th->X) { 506316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 50739d13666SHong Zhang } 50839d13666SHong Zhang if (!th->Xdot) { 509316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 51039d13666SHong Zhang } 51139d13666SHong Zhang if (!th->X0) { 5123b1890cdSShri Abhyankar ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 51339d13666SHong Zhang } 5147445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 5157445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 5167445fe48SJed Brown if (dm) { 5177445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 518258e1594SPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 5197445fe48SJed Brown } 5203b1890cdSShri Abhyankar if (th->Theta == 0.5 && th->endpoint) th->order = 2; 5213b1890cdSShri Abhyankar else th->order = 1; 5223b1890cdSShri Abhyankar 523552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 524ef749922SLisandro Dalcin if (!th->adapt) { 5253b1890cdSShri Abhyankar ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr); 5263b1890cdSShri Abhyankar } 527b4dd3bf9SBarry Smith PetscFunctionReturn(0); 528b4dd3bf9SBarry Smith } 529b4dd3bf9SBarry Smith /*------------------------------------------------------------*/ 530b4dd3bf9SBarry Smith 531b4dd3bf9SBarry Smith #undef __FUNCT__ 532b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta" 533b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 534b4dd3bf9SBarry Smith { 535b4dd3bf9SBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 536b4dd3bf9SBarry Smith PetscErrorCode ierr; 537b4dd3bf9SBarry Smith 538b4dd3bf9SBarry Smith PetscFunctionBegin; 539abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecDeltaLam);CHKERRQ(ierr); 5402ca6e920SHong Zhang if(ts->vecs_sensip) { 541abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecDeltaMu);CHKERRQ(ierr); 5422ca6e920SHong Zhang } 543abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecSensiTemp);CHKERRQ(ierr); 544316643e7SJed Brown PetscFunctionReturn(0); 545316643e7SJed Brown } 546316643e7SJed Brown /*------------------------------------------------------------*/ 547316643e7SJed Brown 548316643e7SJed Brown #undef __FUNCT__ 549316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 5508c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts) 551316643e7SJed Brown { 552316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 553316643e7SJed Brown PetscErrorCode ierr; 554316643e7SJed Brown 555316643e7SJed Brown PetscFunctionBegin; 556e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 557316643e7SJed Brown { 5580298fd71SBarry Smith ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 5590298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 5600298fd71SBarry 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); 5610298fd71SBarry Smith ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr); 562d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 563316643e7SJed Brown } 564316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 565316643e7SJed Brown PetscFunctionReturn(0); 566316643e7SJed Brown } 567316643e7SJed Brown 568316643e7SJed Brown #undef __FUNCT__ 569316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 570316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 571316643e7SJed Brown { 572316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 573ace3abfcSBarry Smith PetscBool iascii; 574316643e7SJed Brown PetscErrorCode ierr; 575316643e7SJed Brown 576316643e7SJed Brown PetscFunctionBegin; 577251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 578316643e7SJed Brown if (iascii) { 5797c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 580ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 581316643e7SJed Brown } 582ac75fa18SLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 583316643e7SJed Brown PetscFunctionReturn(0); 584316643e7SJed Brown } 585316643e7SJed Brown 5860de4c49aSJed Brown #undef __FUNCT__ 5870de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 5887087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 5890de4c49aSJed Brown { 5900de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 5910de4c49aSJed Brown 5920de4c49aSJed Brown PetscFunctionBegin; 5930de4c49aSJed Brown *theta = th->Theta; 5940de4c49aSJed Brown PetscFunctionReturn(0); 5950de4c49aSJed Brown } 5960de4c49aSJed Brown 5970de4c49aSJed Brown #undef __FUNCT__ 5980de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 5997087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 6000de4c49aSJed Brown { 6010de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 6020de4c49aSJed Brown 6030de4c49aSJed Brown PetscFunctionBegin; 6047c8652ddSBarry Smith if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 6050de4c49aSJed Brown th->Theta = theta; 6060de4c49aSJed Brown PetscFunctionReturn(0); 6070de4c49aSJed Brown } 608eb284becSJed Brown 609eb284becSJed Brown #undef __FUNCT__ 61078e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 61126f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 61226f2ff8fSLisandro Dalcin { 61326f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 61426f2ff8fSLisandro Dalcin 61526f2ff8fSLisandro Dalcin PetscFunctionBegin; 61626f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 61726f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 61826f2ff8fSLisandro Dalcin } 61926f2ff8fSLisandro Dalcin 62026f2ff8fSLisandro Dalcin #undef __FUNCT__ 62126f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 622eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 623eb284becSJed Brown { 624eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 625eb284becSJed Brown 626eb284becSJed Brown PetscFunctionBegin; 627eb284becSJed Brown th->endpoint = flg; 628eb284becSJed Brown PetscFunctionReturn(0); 629eb284becSJed Brown } 6300de4c49aSJed Brown 631f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 632f9c1d6abSBarry Smith #undef __FUNCT__ 633f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 634f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 635f9c1d6abSBarry Smith { 636f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 637f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 6383fd8ae06SJed Brown const PetscReal one = 1.0; 639f9c1d6abSBarry Smith 640f9c1d6abSBarry Smith PetscFunctionBegin; 6413fd8ae06SJed Brown f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 642f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 643f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 644f9c1d6abSBarry Smith PetscFunctionReturn(0); 645f9c1d6abSBarry Smith } 646f9c1d6abSBarry Smith #endif 647f9c1d6abSBarry Smith 64842682096SHong Zhang #undef __FUNCT__ 64942682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta" 65042682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 65142682096SHong Zhang { 65242682096SHong Zhang TS_Theta *th = (TS_Theta*)ts->data; 65342682096SHong Zhang 65442682096SHong Zhang PetscFunctionBegin; 6552ca6e920SHong Zhang *ns = 1; 6562ca6e920SHong Zhang if(Y) { 6572ca6e920SHong Zhang *Y = &(th->X); 6582ca6e920SHong Zhang } 65942682096SHong Zhang PetscFunctionReturn(0); 66042682096SHong Zhang } 661f9c1d6abSBarry Smith 662316643e7SJed Brown /* ------------------------------------------------------------ */ 663316643e7SJed Brown /*MC 66496f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 665316643e7SJed Brown 666316643e7SJed Brown Level: beginner 667316643e7SJed Brown 6684eb428fdSBarry Smith Options Database: 6690c3ba866SJed Brown -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 6704eb428fdSBarry Smith -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 6710c3ba866SJed Brown -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 6724eb428fdSBarry Smith 673eb284becSJed Brown Notes: 6740c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 6750c3ba866SJed Brown $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 6764eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 6774eb428fdSBarry Smith 6784eb428fdSBarry Smith 6794eb428fdSBarry Smith 680eb284becSJed Brown This method can be applied to DAE. 681eb284becSJed Brown 682eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 683eb284becSJed Brown 684eb284becSJed Brown .vb 685eb284becSJed Brown Theta | Theta 686eb284becSJed Brown ------------- 687eb284becSJed Brown | 1 688eb284becSJed Brown .ve 689eb284becSJed Brown 690eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 691eb284becSJed Brown 692eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 693eb284becSJed Brown 694eb284becSJed Brown .vb 695eb284becSJed Brown 0 | 0 0 696eb284becSJed Brown 1 | 1-Theta Theta 697eb284becSJed Brown ------------------- 698eb284becSJed Brown | 1-Theta Theta 699eb284becSJed Brown .ve 700eb284becSJed Brown 701eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 702eb284becSJed Brown 703eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 704eb284becSJed Brown 705eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 706eb284becSJed Brown 7074eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 708eb284becSJed Brown 709eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 710316643e7SJed Brown 711316643e7SJed Brown M*/ 712316643e7SJed Brown #undef __FUNCT__ 713316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 7148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 715316643e7SJed Brown { 716316643e7SJed Brown TS_Theta *th; 717316643e7SJed Brown PetscErrorCode ierr; 718316643e7SJed Brown 719316643e7SJed Brown PetscFunctionBegin; 720277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 721316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 722316643e7SJed Brown ts->ops->view = TSView_Theta; 723316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 72442f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_Theta; 725316643e7SJed Brown ts->ops->step = TSStep_Theta; 726cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 7273b1890cdSShri Abhyankar ts->ops->evaluatestep = TSEvaluateStep_Theta; 72824655328SShri ts->ops->rollback = TSRollBack_Theta; 729316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 7300f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 7310f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 732f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 733f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 734f9c1d6abSBarry Smith #endif 73542682096SHong Zhang ts->ops->getstages = TSGetStages_Theta; 73642f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_Theta; 737316643e7SJed Brown 738b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 739316643e7SJed Brown ts->data = (void*)th; 740316643e7SJed Brown 7416f700aefSJed Brown th->extrapolate = PETSC_FALSE; 742316643e7SJed Brown th->Theta = 0.5; 7433b1890cdSShri Abhyankar th->ccfl = 1.0; 7443b1890cdSShri Abhyankar th->adapt = PETSC_FALSE; 745bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 746bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 747bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 748bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 749316643e7SJed Brown PetscFunctionReturn(0); 750316643e7SJed Brown } 7510de4c49aSJed Brown 7520de4c49aSJed Brown #undef __FUNCT__ 7530de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 7540de4c49aSJed Brown /*@ 7550de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 7560de4c49aSJed Brown 7570de4c49aSJed Brown Not Collective 7580de4c49aSJed Brown 7590de4c49aSJed Brown Input Parameter: 7600de4c49aSJed Brown . ts - timestepping context 7610de4c49aSJed Brown 7620de4c49aSJed Brown Output Parameter: 7630de4c49aSJed Brown . theta - stage abscissa 7640de4c49aSJed Brown 7650de4c49aSJed Brown Note: 7660de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 7670de4c49aSJed Brown 7680de4c49aSJed Brown Level: Advanced 7690de4c49aSJed Brown 7700de4c49aSJed Brown .seealso: TSThetaSetTheta() 7710de4c49aSJed Brown @*/ 7727087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 7730de4c49aSJed Brown { 7744ac538c5SBarry Smith PetscErrorCode ierr; 7750de4c49aSJed Brown 7760de4c49aSJed Brown PetscFunctionBegin; 777afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7780de4c49aSJed Brown PetscValidPointer(theta,2); 7794ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 7800de4c49aSJed Brown PetscFunctionReturn(0); 7810de4c49aSJed Brown } 7820de4c49aSJed Brown 7830de4c49aSJed Brown #undef __FUNCT__ 7840de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 7850de4c49aSJed Brown /*@ 7860de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 7870de4c49aSJed Brown 7880de4c49aSJed Brown Not Collective 7890de4c49aSJed Brown 7900de4c49aSJed Brown Input Parameter: 7910de4c49aSJed Brown + ts - timestepping context 7920de4c49aSJed Brown - theta - stage abscissa 7930de4c49aSJed Brown 7940de4c49aSJed Brown Options Database: 7950de4c49aSJed Brown . -ts_theta_theta <theta> 7960de4c49aSJed Brown 7970de4c49aSJed Brown Level: Intermediate 7980de4c49aSJed Brown 7990de4c49aSJed Brown .seealso: TSThetaGetTheta() 8000de4c49aSJed Brown @*/ 8017087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 8020de4c49aSJed Brown { 8034ac538c5SBarry Smith PetscErrorCode ierr; 8040de4c49aSJed Brown 8050de4c49aSJed Brown PetscFunctionBegin; 806afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8074ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 8080de4c49aSJed Brown PetscFunctionReturn(0); 8090de4c49aSJed Brown } 810f33bbcb6SJed Brown 811eb284becSJed Brown #undef __FUNCT__ 81226f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 81326f2ff8fSLisandro Dalcin /*@ 81426f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 81526f2ff8fSLisandro Dalcin 81626f2ff8fSLisandro Dalcin Not Collective 81726f2ff8fSLisandro Dalcin 81826f2ff8fSLisandro Dalcin Input Parameter: 81926f2ff8fSLisandro Dalcin . ts - timestepping context 82026f2ff8fSLisandro Dalcin 82126f2ff8fSLisandro Dalcin Output Parameter: 82226f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 82326f2ff8fSLisandro Dalcin 82426f2ff8fSLisandro Dalcin Level: Advanced 82526f2ff8fSLisandro Dalcin 82626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 82726f2ff8fSLisandro Dalcin @*/ 82826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 82926f2ff8fSLisandro Dalcin { 83026f2ff8fSLisandro Dalcin PetscErrorCode ierr; 83126f2ff8fSLisandro Dalcin 83226f2ff8fSLisandro Dalcin PetscFunctionBegin; 83326f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 83426f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 83526f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 83626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 83726f2ff8fSLisandro Dalcin } 83826f2ff8fSLisandro Dalcin 83926f2ff8fSLisandro Dalcin #undef __FUNCT__ 840eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 841eb284becSJed Brown /*@ 842eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 843eb284becSJed Brown 844eb284becSJed Brown Not Collective 845eb284becSJed Brown 846eb284becSJed Brown Input Parameter: 847eb284becSJed Brown + ts - timestepping context 848eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 849eb284becSJed Brown 850eb284becSJed Brown Options Database: 851eb284becSJed Brown . -ts_theta_endpoint <flg> 852eb284becSJed Brown 853eb284becSJed Brown Level: Intermediate 854eb284becSJed Brown 855eb284becSJed Brown .seealso: TSTHETA, TSCN 856eb284becSJed Brown @*/ 857eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 858eb284becSJed Brown { 859eb284becSJed Brown PetscErrorCode ierr; 860eb284becSJed Brown 861eb284becSJed Brown PetscFunctionBegin; 862eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 863eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 864eb284becSJed Brown PetscFunctionReturn(0); 865eb284becSJed Brown } 866eb284becSJed Brown 867f33bbcb6SJed Brown /* 868f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 869f33bbcb6SJed Brown * The creation functions for these specializations are below. 870f33bbcb6SJed Brown */ 871f33bbcb6SJed Brown 872f33bbcb6SJed Brown #undef __FUNCT__ 873f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 874f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 875f33bbcb6SJed Brown { 876d52bd9f3SBarry Smith PetscErrorCode ierr; 877d52bd9f3SBarry Smith 878f33bbcb6SJed Brown PetscFunctionBegin; 879d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 880f33bbcb6SJed Brown PetscFunctionReturn(0); 881f33bbcb6SJed Brown } 882f33bbcb6SJed Brown 883f33bbcb6SJed Brown /*MC 884f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 885f33bbcb6SJed Brown 886f33bbcb6SJed Brown Level: beginner 887f33bbcb6SJed Brown 8884eb428fdSBarry Smith Notes: 889c7eb6c99SShri Abhyankar TSBEULER is equivalent to TSTHETA with Theta=1.0 8904eb428fdSBarry Smith 8914eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 1. 8924eb428fdSBarry Smith 893f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 894f33bbcb6SJed Brown 895f33bbcb6SJed Brown M*/ 896f33bbcb6SJed Brown #undef __FUNCT__ 897f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 8988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 899f33bbcb6SJed Brown { 900f33bbcb6SJed Brown PetscErrorCode ierr; 901f33bbcb6SJed Brown 902f33bbcb6SJed Brown PetscFunctionBegin; 903f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 904f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 905f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 906f33bbcb6SJed Brown PetscFunctionReturn(0); 907f33bbcb6SJed Brown } 908f33bbcb6SJed Brown 909f33bbcb6SJed Brown #undef __FUNCT__ 910f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 911f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 912f33bbcb6SJed Brown { 913d52bd9f3SBarry Smith PetscErrorCode ierr; 914d52bd9f3SBarry Smith 915f33bbcb6SJed Brown PetscFunctionBegin; 916d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 917f33bbcb6SJed Brown PetscFunctionReturn(0); 918f33bbcb6SJed Brown } 919f33bbcb6SJed Brown 920f33bbcb6SJed Brown /*MC 921f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 922f33bbcb6SJed Brown 923f33bbcb6SJed Brown Level: beginner 924f33bbcb6SJed Brown 925f33bbcb6SJed Brown Notes: 9267cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 9277cf5af47SJed Brown 9287cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 929f33bbcb6SJed Brown 930f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 931f33bbcb6SJed Brown 932f33bbcb6SJed Brown M*/ 933f33bbcb6SJed Brown #undef __FUNCT__ 934f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 9358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 936f33bbcb6SJed Brown { 937f33bbcb6SJed Brown PetscErrorCode ierr; 938f33bbcb6SJed Brown 939f33bbcb6SJed Brown PetscFunctionBegin; 940f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 941f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 942eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 943f33bbcb6SJed Brown ts->ops->view = TSView_CN; 944f33bbcb6SJed Brown PetscFunctionReturn(0); 945f33bbcb6SJed Brown } 946