1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 57445fe48SJed Brown #include <petscsnesfas.h> 6316643e7SJed Brown 7316643e7SJed Brown typedef struct { 8316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 9eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 10ace3abfcSBarry Smith PetscBool extrapolate; 11eb284becSJed Brown PetscBool endpoint; 12316643e7SJed Brown PetscReal Theta; 13316643e7SJed Brown PetscReal shift; 14316643e7SJed Brown PetscReal stage_time; 15316643e7SJed Brown } TS_Theta; 16316643e7SJed Brown 17316643e7SJed Brown #undef __FUNCT__ 187445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 197445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 207445fe48SJed Brown { 217445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 227445fe48SJed Brown PetscErrorCode ierr; 237445fe48SJed Brown 247445fe48SJed Brown PetscFunctionBegin; 257445fe48SJed Brown if (X0) { 267445fe48SJed Brown if (dm && dm != ts->dm) { 270d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 287445fe48SJed Brown } else *X0 = ts->vec_sol; 297445fe48SJed Brown } 307445fe48SJed Brown if (Xdot) { 317445fe48SJed Brown if (dm && dm != ts->dm) { 320d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 337445fe48SJed Brown } else *Xdot = th->Xdot; 347445fe48SJed Brown } 357445fe48SJed Brown PetscFunctionReturn(0); 367445fe48SJed Brown } 377445fe48SJed Brown 380d0b770aSPeter Brune 390d0b770aSPeter Brune #undef __FUNCT__ 400d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 410d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 420d0b770aSPeter Brune { 430d0b770aSPeter Brune TS_Theta *th = (TS_Theta*)ts->data; 440d0b770aSPeter Brune PetscErrorCode ierr; 450d0b770aSPeter Brune 460d0b770aSPeter Brune PetscFunctionBegin; 470d0b770aSPeter Brune if (X0) { 480d0b770aSPeter Brune if (dm && dm != ts->dm) { 490d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 500d0b770aSPeter Brune } 510d0b770aSPeter Brune } 520d0b770aSPeter Brune if (Xdot) { 530d0b770aSPeter Brune if (dm && dm != ts->dm) { 540d0b770aSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 550d0b770aSPeter Brune } 560d0b770aSPeter Brune } 570d0b770aSPeter Brune PetscFunctionReturn(0); 580d0b770aSPeter Brune } 590d0b770aSPeter Brune 607445fe48SJed Brown #undef __FUNCT__ 617445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 627445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 637445fe48SJed Brown { 647445fe48SJed Brown 657445fe48SJed Brown PetscFunctionBegin; 667445fe48SJed Brown PetscFunctionReturn(0); 677445fe48SJed Brown } 687445fe48SJed Brown 697445fe48SJed Brown #undef __FUNCT__ 707445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 717445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 727445fe48SJed Brown { 737445fe48SJed Brown TS ts = (TS)ctx; 747445fe48SJed Brown PetscErrorCode ierr; 757445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 767445fe48SJed Brown 777445fe48SJed Brown PetscFunctionBegin; 787445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 797445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 807445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 817445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 827445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 837445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 840d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 850d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 867445fe48SJed Brown PetscFunctionReturn(0); 877445fe48SJed Brown } 887445fe48SJed Brown 897445fe48SJed Brown #undef __FUNCT__ 90316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 91193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 92316643e7SJed Brown { 93316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 94b70ae86eSJed Brown PetscInt its,lits; 95cdbf8f93SLisandro Dalcin PetscReal next_time_step; 96f1b97656SJed Brown SNESConvergedReason snesreason; 972b5a38e1SLisandro Dalcin PetscErrorCode ierr; 98316643e7SJed Brown 99316643e7SJed Brown PetscFunctionBegin; 100cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 101eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 102316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 103316643e7SJed Brown 104eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 105eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 106eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 107eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 108eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 109eb284becSJed Brown } 110ace68cafSJed Brown if (th->extrapolate) { 1112b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 112ace68cafSJed Brown } else { 1132b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 114ace68cafSJed Brown } 115eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 116316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 117316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 118f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 119*5ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 120f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 121f1b97656SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 122f1b97656SJed Brown ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 123f1b97656SJed Brown PetscFunctionReturn(0); 124f1b97656SJed Brown } 125eb284becSJed Brown if (th->endpoint) { 126eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 127eb284becSJed Brown } else { 1282b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 1292b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 130eb284becSJed Brown } 1312b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 132cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 133316643e7SJed Brown ts->steps++; 134316643e7SJed Brown PetscFunctionReturn(0); 135316643e7SJed Brown } 136316643e7SJed Brown 137cd652676SJed Brown #undef __FUNCT__ 138cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 139cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 140cd652676SJed Brown { 141cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1425a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 143cd652676SJed Brown PetscErrorCode ierr; 144cd652676SJed Brown 145cd652676SJed Brown PetscFunctionBegin; 146a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 1475a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 1485a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 149cd652676SJed Brown PetscFunctionReturn(0); 150cd652676SJed Brown } 151cd652676SJed Brown 152316643e7SJed Brown /*------------------------------------------------------------*/ 153316643e7SJed Brown #undef __FUNCT__ 154277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 155277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 156316643e7SJed Brown { 157316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 158316643e7SJed Brown PetscErrorCode ierr; 159316643e7SJed Brown 160316643e7SJed Brown PetscFunctionBegin; 1616bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 1626bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 163eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 164277b19d0SLisandro Dalcin PetscFunctionReturn(0); 165277b19d0SLisandro Dalcin } 166277b19d0SLisandro Dalcin 167277b19d0SLisandro Dalcin #undef __FUNCT__ 168277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 169277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 170277b19d0SLisandro Dalcin { 171277b19d0SLisandro Dalcin PetscErrorCode ierr; 172277b19d0SLisandro Dalcin 173277b19d0SLisandro Dalcin PetscFunctionBegin; 174277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 175277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 176335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 177335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 17826f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 179eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 180316643e7SJed Brown PetscFunctionReturn(0); 181316643e7SJed Brown } 182316643e7SJed Brown 183316643e7SJed Brown /* 184316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1852b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 186316643e7SJed Brown */ 187316643e7SJed Brown #undef __FUNCT__ 1880f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1890f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 190316643e7SJed Brown { 191316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 192316643e7SJed Brown PetscErrorCode ierr; 1937445fe48SJed Brown Vec X0,Xdot; 1947445fe48SJed Brown DM dm,dmsave; 195316643e7SJed Brown 196316643e7SJed Brown PetscFunctionBegin; 1977445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1985a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 1997445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 2007445fe48SJed Brown ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr); 2017445fe48SJed Brown 2027445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 2037445fe48SJed Brown dmsave = ts->dm; 2047445fe48SJed Brown ts->dm = dm; 2057445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 2067445fe48SJed Brown ts->dm = dmsave; 2070d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 208316643e7SJed Brown PetscFunctionReturn(0); 209316643e7SJed Brown } 210316643e7SJed Brown 211316643e7SJed Brown #undef __FUNCT__ 2120f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 2130f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 214316643e7SJed Brown { 215316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 216316643e7SJed Brown PetscErrorCode ierr; 2177445fe48SJed Brown Vec Xdot; 2187445fe48SJed Brown DM dm,dmsave; 219316643e7SJed Brown 220316643e7SJed Brown PetscFunctionBegin; 2217445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2227445fe48SJed Brown 2230f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 2247445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 2257445fe48SJed Brown 2267445fe48SJed Brown dmsave = ts->dm; 2277445fe48SJed Brown ts->dm = dm; 2287445fe48SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 2297445fe48SJed Brown ts->dm = dmsave; 2300d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 231316643e7SJed Brown PetscFunctionReturn(0); 232316643e7SJed Brown } 233316643e7SJed Brown 234316643e7SJed Brown #undef __FUNCT__ 235316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 236316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 237316643e7SJed Brown { 238316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 239316643e7SJed Brown PetscErrorCode ierr; 2407445fe48SJed Brown SNES snes; 2417445fe48SJed Brown DM dm; 242316643e7SJed Brown 243316643e7SJed Brown PetscFunctionBegin; 244316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 245316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 2467445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2477445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2487445fe48SJed Brown if (dm) { 2497445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 2507445fe48SJed Brown } 251316643e7SJed Brown PetscFunctionReturn(0); 252316643e7SJed Brown } 253316643e7SJed Brown /*------------------------------------------------------------*/ 254316643e7SJed Brown 255316643e7SJed Brown #undef __FUNCT__ 256316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 257316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 258316643e7SJed Brown { 259316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 260316643e7SJed Brown PetscErrorCode ierr; 261316643e7SJed Brown 262316643e7SJed Brown PetscFunctionBegin; 263d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 264316643e7SJed Brown { 265316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 266acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 267eb284becSJed Brown ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr); 268d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 269316643e7SJed Brown } 270316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 271316643e7SJed Brown PetscFunctionReturn(0); 272316643e7SJed Brown } 273316643e7SJed Brown 274316643e7SJed Brown #undef __FUNCT__ 275316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 276316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 277316643e7SJed Brown { 278316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 279ace3abfcSBarry Smith PetscBool iascii; 280316643e7SJed Brown PetscErrorCode ierr; 281316643e7SJed Brown 282316643e7SJed Brown PetscFunctionBegin; 283251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 284316643e7SJed Brown if (iascii) { 285316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 286ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 287316643e7SJed Brown } 288d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 289316643e7SJed Brown PetscFunctionReturn(0); 290316643e7SJed Brown } 291316643e7SJed Brown 2920de4c49aSJed Brown EXTERN_C_BEGIN 2930de4c49aSJed Brown #undef __FUNCT__ 2940de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 2957087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 2960de4c49aSJed Brown { 2970de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2980de4c49aSJed Brown 2990de4c49aSJed Brown PetscFunctionBegin; 3000de4c49aSJed Brown *theta = th->Theta; 3010de4c49aSJed Brown PetscFunctionReturn(0); 3020de4c49aSJed Brown } 3030de4c49aSJed Brown 3040de4c49aSJed Brown #undef __FUNCT__ 3050de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 3067087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 3070de4c49aSJed Brown { 3080de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3090de4c49aSJed Brown 3100de4c49aSJed Brown PetscFunctionBegin; 311e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 3120de4c49aSJed Brown th->Theta = theta; 3130de4c49aSJed Brown PetscFunctionReturn(0); 3140de4c49aSJed Brown } 315eb284becSJed Brown 316eb284becSJed Brown #undef __FUNCT__ 31778e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 31826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 31926f2ff8fSLisandro Dalcin { 32026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 32126f2ff8fSLisandro Dalcin 32226f2ff8fSLisandro Dalcin PetscFunctionBegin; 32326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 32426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 32526f2ff8fSLisandro Dalcin } 32626f2ff8fSLisandro Dalcin 32726f2ff8fSLisandro Dalcin #undef __FUNCT__ 32826f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 329eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 330eb284becSJed Brown { 331eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 332eb284becSJed Brown 333eb284becSJed Brown PetscFunctionBegin; 334eb284becSJed Brown th->endpoint = flg; 335eb284becSJed Brown PetscFunctionReturn(0); 336eb284becSJed Brown } 3370de4c49aSJed Brown EXTERN_C_END 3380de4c49aSJed Brown 339316643e7SJed Brown /* ------------------------------------------------------------ */ 340316643e7SJed Brown /*MC 34196f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 342316643e7SJed Brown 343316643e7SJed Brown Level: beginner 344316643e7SJed Brown 345eb284becSJed Brown Notes: 346eb284becSJed Brown This method can be applied to DAE. 347eb284becSJed Brown 348eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 349eb284becSJed Brown 350eb284becSJed Brown .vb 351eb284becSJed Brown Theta | Theta 352eb284becSJed Brown ------------- 353eb284becSJed Brown | 1 354eb284becSJed Brown .ve 355eb284becSJed Brown 356eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 357eb284becSJed Brown 358eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 359eb284becSJed Brown 360eb284becSJed Brown .vb 361eb284becSJed Brown 0 | 0 0 362eb284becSJed Brown 1 | 1-Theta Theta 363eb284becSJed Brown ------------------- 364eb284becSJed Brown | 1-Theta Theta 365eb284becSJed Brown .ve 366eb284becSJed Brown 367eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 368eb284becSJed Brown 369eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 370eb284becSJed Brown 371eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 372eb284becSJed Brown 373eb284becSJed Brown is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 374eb284becSJed Brown 375eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 376316643e7SJed Brown 377316643e7SJed Brown M*/ 378316643e7SJed Brown EXTERN_C_BEGIN 379316643e7SJed Brown #undef __FUNCT__ 380316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 3817087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 382316643e7SJed Brown { 383316643e7SJed Brown TS_Theta *th; 384316643e7SJed Brown PetscErrorCode ierr; 385316643e7SJed Brown 386316643e7SJed Brown PetscFunctionBegin; 387277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 388316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 389316643e7SJed Brown ts->ops->view = TSView_Theta; 390316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 391316643e7SJed Brown ts->ops->step = TSStep_Theta; 392cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 393316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 3940f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 3950f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 396316643e7SJed Brown 397316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 398316643e7SJed Brown ts->data = (void*)th; 399316643e7SJed Brown 4006f700aefSJed Brown th->extrapolate = PETSC_FALSE; 401316643e7SJed Brown th->Theta = 0.5; 402316643e7SJed Brown 4030de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 4040de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 40526f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 406eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 407316643e7SJed Brown PetscFunctionReturn(0); 408316643e7SJed Brown } 409316643e7SJed Brown EXTERN_C_END 4100de4c49aSJed Brown 4110de4c49aSJed Brown #undef __FUNCT__ 4120de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 4130de4c49aSJed Brown /*@ 4140de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 4150de4c49aSJed Brown 4160de4c49aSJed Brown Not Collective 4170de4c49aSJed Brown 4180de4c49aSJed Brown Input Parameter: 4190de4c49aSJed Brown . ts - timestepping context 4200de4c49aSJed Brown 4210de4c49aSJed Brown Output Parameter: 4220de4c49aSJed Brown . theta - stage abscissa 4230de4c49aSJed Brown 4240de4c49aSJed Brown Note: 4250de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 4260de4c49aSJed Brown 4270de4c49aSJed Brown Level: Advanced 4280de4c49aSJed Brown 4290de4c49aSJed Brown .seealso: TSThetaSetTheta() 4300de4c49aSJed Brown @*/ 4317087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 4320de4c49aSJed Brown { 4334ac538c5SBarry Smith PetscErrorCode ierr; 4340de4c49aSJed Brown 4350de4c49aSJed Brown PetscFunctionBegin; 436afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4370de4c49aSJed Brown PetscValidPointer(theta,2); 4384ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 4390de4c49aSJed Brown PetscFunctionReturn(0); 4400de4c49aSJed Brown } 4410de4c49aSJed Brown 4420de4c49aSJed Brown #undef __FUNCT__ 4430de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 4440de4c49aSJed Brown /*@ 4450de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 4460de4c49aSJed Brown 4470de4c49aSJed Brown Not Collective 4480de4c49aSJed Brown 4490de4c49aSJed Brown Input Parameter: 4500de4c49aSJed Brown + ts - timestepping context 4510de4c49aSJed Brown - theta - stage abscissa 4520de4c49aSJed Brown 4530de4c49aSJed Brown Options Database: 4540de4c49aSJed Brown . -ts_theta_theta <theta> 4550de4c49aSJed Brown 4560de4c49aSJed Brown Level: Intermediate 4570de4c49aSJed Brown 4580de4c49aSJed Brown .seealso: TSThetaGetTheta() 4590de4c49aSJed Brown @*/ 4607087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 4610de4c49aSJed Brown { 4624ac538c5SBarry Smith PetscErrorCode ierr; 4630de4c49aSJed Brown 4640de4c49aSJed Brown PetscFunctionBegin; 465afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4664ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 4670de4c49aSJed Brown PetscFunctionReturn(0); 4680de4c49aSJed Brown } 469f33bbcb6SJed Brown 470eb284becSJed Brown #undef __FUNCT__ 47126f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 47226f2ff8fSLisandro Dalcin /*@ 47326f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 47426f2ff8fSLisandro Dalcin 47526f2ff8fSLisandro Dalcin Not Collective 47626f2ff8fSLisandro Dalcin 47726f2ff8fSLisandro Dalcin Input Parameter: 47826f2ff8fSLisandro Dalcin . ts - timestepping context 47926f2ff8fSLisandro Dalcin 48026f2ff8fSLisandro Dalcin Output Parameter: 48126f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 48226f2ff8fSLisandro Dalcin 48326f2ff8fSLisandro Dalcin Level: Advanced 48426f2ff8fSLisandro Dalcin 48526f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 48626f2ff8fSLisandro Dalcin @*/ 48726f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 48826f2ff8fSLisandro Dalcin { 48926f2ff8fSLisandro Dalcin PetscErrorCode ierr; 49026f2ff8fSLisandro Dalcin 49126f2ff8fSLisandro Dalcin PetscFunctionBegin; 49226f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 49326f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 49426f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 49526f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 49626f2ff8fSLisandro Dalcin } 49726f2ff8fSLisandro Dalcin 49826f2ff8fSLisandro Dalcin #undef __FUNCT__ 499eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 500eb284becSJed Brown /*@ 501eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 502eb284becSJed Brown 503eb284becSJed Brown Not Collective 504eb284becSJed Brown 505eb284becSJed Brown Input Parameter: 506eb284becSJed Brown + ts - timestepping context 507eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 508eb284becSJed Brown 509eb284becSJed Brown Options Database: 510eb284becSJed Brown . -ts_theta_endpoint <flg> 511eb284becSJed Brown 512eb284becSJed Brown Level: Intermediate 513eb284becSJed Brown 514eb284becSJed Brown .seealso: TSTHETA, TSCN 515eb284becSJed Brown @*/ 516eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 517eb284becSJed Brown { 518eb284becSJed Brown PetscErrorCode ierr; 519eb284becSJed Brown 520eb284becSJed Brown PetscFunctionBegin; 521eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 522eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 523eb284becSJed Brown PetscFunctionReturn(0); 524eb284becSJed Brown } 525eb284becSJed Brown 526f33bbcb6SJed Brown /* 527f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 528f33bbcb6SJed Brown * The creation functions for these specializations are below. 529f33bbcb6SJed Brown */ 530f33bbcb6SJed Brown 531f33bbcb6SJed Brown #undef __FUNCT__ 532f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 533f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 534f33bbcb6SJed Brown { 535d52bd9f3SBarry Smith PetscErrorCode ierr; 536d52bd9f3SBarry Smith 537f33bbcb6SJed Brown PetscFunctionBegin; 538d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 539f33bbcb6SJed Brown PetscFunctionReturn(0); 540f33bbcb6SJed Brown } 541f33bbcb6SJed Brown 542f33bbcb6SJed Brown /*MC 543f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 544f33bbcb6SJed Brown 545f33bbcb6SJed Brown Level: beginner 546f33bbcb6SJed Brown 547f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 548f33bbcb6SJed Brown 549f33bbcb6SJed Brown M*/ 550f33bbcb6SJed Brown EXTERN_C_BEGIN 551f33bbcb6SJed Brown #undef __FUNCT__ 552f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 553f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 554f33bbcb6SJed Brown { 555f33bbcb6SJed Brown PetscErrorCode ierr; 556f33bbcb6SJed Brown 557f33bbcb6SJed Brown PetscFunctionBegin; 558f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 559f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 560f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 561f33bbcb6SJed Brown PetscFunctionReturn(0); 562f33bbcb6SJed Brown } 563f33bbcb6SJed Brown EXTERN_C_END 564f33bbcb6SJed Brown 565f33bbcb6SJed Brown #undef __FUNCT__ 566f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 567f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 568f33bbcb6SJed Brown { 569d52bd9f3SBarry Smith PetscErrorCode ierr; 570d52bd9f3SBarry Smith 571f33bbcb6SJed Brown PetscFunctionBegin; 572d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 573f33bbcb6SJed Brown PetscFunctionReturn(0); 574f33bbcb6SJed Brown } 575f33bbcb6SJed Brown 576f33bbcb6SJed Brown /*MC 577f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 578f33bbcb6SJed Brown 579f33bbcb6SJed Brown Level: beginner 580f33bbcb6SJed Brown 581f33bbcb6SJed Brown Notes: 5827cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 5837cf5af47SJed Brown 5847cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 585f33bbcb6SJed Brown 586f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 587f33bbcb6SJed Brown 588f33bbcb6SJed Brown M*/ 589f33bbcb6SJed Brown EXTERN_C_BEGIN 590f33bbcb6SJed Brown #undef __FUNCT__ 591f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 592f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 593f33bbcb6SJed Brown { 594f33bbcb6SJed Brown PetscErrorCode ierr; 595f33bbcb6SJed Brown 596f33bbcb6SJed Brown PetscFunctionBegin; 597f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 598f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 599eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 600f33bbcb6SJed Brown ts->ops->view = TSView_CN; 601f33bbcb6SJed Brown PetscFunctionReturn(0); 602f33bbcb6SJed Brown } 603f33bbcb6SJed Brown EXTERN_C_END 604