1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4*f9c1d6abSBarry Smith #define PETSC_DESIRE_COMPLEX 5b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 67445fe48SJed Brown #include <petscsnesfas.h> 7316643e7SJed Brown 8316643e7SJed Brown typedef struct { 9316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 10eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 11ace3abfcSBarry Smith PetscBool extrapolate; 12eb284becSJed Brown PetscBool endpoint; 13316643e7SJed Brown PetscReal Theta; 14316643e7SJed Brown PetscReal shift; 15316643e7SJed Brown PetscReal stage_time; 16316643e7SJed Brown } TS_Theta; 17316643e7SJed Brown 18316643e7SJed Brown #undef __FUNCT__ 197445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot" 207445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 217445fe48SJed Brown { 227445fe48SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 237445fe48SJed Brown PetscErrorCode ierr; 247445fe48SJed Brown 257445fe48SJed Brown PetscFunctionBegin; 267445fe48SJed Brown if (X0) { 277445fe48SJed Brown if (dm && dm != ts->dm) { 280d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 297445fe48SJed Brown } else *X0 = ts->vec_sol; 307445fe48SJed Brown } 317445fe48SJed Brown if (Xdot) { 327445fe48SJed Brown if (dm && dm != ts->dm) { 330d0b770aSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 347445fe48SJed Brown } else *Xdot = th->Xdot; 357445fe48SJed Brown } 367445fe48SJed Brown PetscFunctionReturn(0); 377445fe48SJed Brown } 387445fe48SJed Brown 390d0b770aSPeter Brune 400d0b770aSPeter Brune #undef __FUNCT__ 410d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot" 420d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 430d0b770aSPeter Brune { 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; 944e17fbf5SBarry Smith PetscInt its,lits,reject; 95a3314f2cSMatthew G Knepley PetscReal next_time_step = 0.0; 96f1b97656SJed Brown SNESConvergedReason snesreason; 972b5a38e1SLisandro Dalcin PetscErrorCode ierr; 98316643e7SJed Brown 99316643e7SJed Brown PetscFunctionBegin; 10031748224SBarry Smith if (ts->time_steps_since_decrease > 3 && ts->time_step < ts->time_step_orig) { 10131748224SBarry Smith /* smaller time step has worked successfully for three time-steps, try increasing time step*/ 10231748224SBarry Smith ts->time_step = 2.0*ts->time_step; 10331748224SBarry Smith ts->time_steps_since_decrease = 0; /* don't want to increase time step two time steps in a row */ 10431748224SBarry Smith } 1054e17fbf5SBarry Smith for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 106cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 107eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 108316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 109b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 110b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 111316643e7SJed Brown 112eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 113eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 114eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 115eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 116eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 117eb284becSJed Brown } 118ace68cafSJed Brown if (th->extrapolate) { 1192b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 120ace68cafSJed Brown } else { 1212b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 122ace68cafSJed Brown } 123eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 124316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 125316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 126f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 1275ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 12831748224SBarry Smith if (its < 10) ts->time_steps_since_decrease++; 12931748224SBarry Smith else ts->time_steps_since_decrease = 0; 1304e17fbf5SBarry Smith if (snesreason > 0) break; 1314e17fbf5SBarry Smith ierr = PetscInfo3(ts,"Step=%D, Cutting time-step from %g to %g\n",ts->steps,(double)ts->time_step,(double).5*ts->time_step);CHKERRQ(ierr); 1324e17fbf5SBarry Smith ts->time_step = .5*ts->time_step; 13331748224SBarry Smith ts->time_steps_since_decrease = 0; 1344e17fbf5SBarry Smith } 135f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 136f1b97656SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 137f1b97656SJed 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); 138f1b97656SJed Brown PetscFunctionReturn(0); 139f1b97656SJed Brown } 140eb284becSJed Brown if (th->endpoint) { 141eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 142eb284becSJed Brown } else { 1432b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 1442b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 145eb284becSJed Brown } 1462b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 147cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 148316643e7SJed Brown ts->steps++; 149316643e7SJed Brown PetscFunctionReturn(0); 150316643e7SJed Brown } 151316643e7SJed Brown 152cd652676SJed Brown #undef __FUNCT__ 153cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 154cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 155cd652676SJed Brown { 156cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1575a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 158cd652676SJed Brown PetscErrorCode ierr; 159cd652676SJed Brown 160cd652676SJed Brown PetscFunctionBegin; 161a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 1625a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 1635a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 164cd652676SJed Brown PetscFunctionReturn(0); 165cd652676SJed Brown } 166cd652676SJed Brown 167316643e7SJed Brown /*------------------------------------------------------------*/ 168316643e7SJed Brown #undef __FUNCT__ 169277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 170277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 171316643e7SJed Brown { 172316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 173316643e7SJed Brown PetscErrorCode ierr; 174316643e7SJed Brown 175316643e7SJed Brown PetscFunctionBegin; 1766bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 1776bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 178eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 179277b19d0SLisandro Dalcin PetscFunctionReturn(0); 180277b19d0SLisandro Dalcin } 181277b19d0SLisandro Dalcin 182277b19d0SLisandro Dalcin #undef __FUNCT__ 183277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 184277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 185277b19d0SLisandro Dalcin { 186277b19d0SLisandro Dalcin PetscErrorCode ierr; 187277b19d0SLisandro Dalcin 188277b19d0SLisandro Dalcin PetscFunctionBegin; 189277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 190277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 191335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 192335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 19326f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 194eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 195316643e7SJed Brown PetscFunctionReturn(0); 196316643e7SJed Brown } 197316643e7SJed Brown 198316643e7SJed Brown /* 199316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 2002b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 201316643e7SJed Brown */ 202316643e7SJed Brown #undef __FUNCT__ 2030f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 2040f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 205316643e7SJed Brown { 206316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 207316643e7SJed Brown PetscErrorCode ierr; 2087445fe48SJed Brown Vec X0,Xdot; 2097445fe48SJed Brown DM dm,dmsave; 210316643e7SJed Brown 211316643e7SJed Brown PetscFunctionBegin; 2127445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2135a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 2147445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 2157445fe48SJed Brown ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr); 2167445fe48SJed Brown 2177445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 2187445fe48SJed Brown dmsave = ts->dm; 2197445fe48SJed Brown ts->dm = dm; 2207445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 2217445fe48SJed Brown ts->dm = dmsave; 2220d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 223316643e7SJed Brown PetscFunctionReturn(0); 224316643e7SJed Brown } 225316643e7SJed Brown 226316643e7SJed Brown #undef __FUNCT__ 2270f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 2280f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 229316643e7SJed Brown { 230316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 231316643e7SJed Brown PetscErrorCode ierr; 2327445fe48SJed Brown Vec Xdot; 2337445fe48SJed Brown DM dm,dmsave; 234316643e7SJed Brown 235316643e7SJed Brown PetscFunctionBegin; 2367445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2377445fe48SJed Brown 2380f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 2397445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 2407445fe48SJed Brown 2417445fe48SJed Brown dmsave = ts->dm; 2427445fe48SJed Brown ts->dm = dm; 2437445fe48SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 2447445fe48SJed Brown ts->dm = dmsave; 2450d0b770aSPeter Brune ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 246316643e7SJed Brown PetscFunctionReturn(0); 247316643e7SJed Brown } 248316643e7SJed Brown 249316643e7SJed Brown #undef __FUNCT__ 250316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 251316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 252316643e7SJed Brown { 253316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 254316643e7SJed Brown PetscErrorCode ierr; 2557445fe48SJed Brown SNES snes; 2567445fe48SJed Brown DM dm; 257316643e7SJed Brown 258316643e7SJed Brown PetscFunctionBegin; 259316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 260316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 2617445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2627445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2637445fe48SJed Brown if (dm) { 2647445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 2657445fe48SJed Brown } 266316643e7SJed Brown PetscFunctionReturn(0); 267316643e7SJed Brown } 268316643e7SJed Brown /*------------------------------------------------------------*/ 269316643e7SJed Brown 270316643e7SJed Brown #undef __FUNCT__ 271316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 272316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 273316643e7SJed Brown { 274316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 275316643e7SJed Brown PetscErrorCode ierr; 276316643e7SJed Brown 277316643e7SJed Brown PetscFunctionBegin; 278d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 279316643e7SJed Brown { 280316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 281acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 282eb284becSJed 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); 283d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 284316643e7SJed Brown } 285316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 286316643e7SJed Brown PetscFunctionReturn(0); 287316643e7SJed Brown } 288316643e7SJed Brown 289316643e7SJed Brown #undef __FUNCT__ 290316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 291316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 292316643e7SJed Brown { 293316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 294ace3abfcSBarry Smith PetscBool iascii; 295316643e7SJed Brown PetscErrorCode ierr; 296316643e7SJed Brown 297316643e7SJed Brown PetscFunctionBegin; 298251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 299316643e7SJed Brown if (iascii) { 300316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 301ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 302316643e7SJed Brown } 303d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 304316643e7SJed Brown PetscFunctionReturn(0); 305316643e7SJed Brown } 306316643e7SJed Brown 3070de4c49aSJed Brown EXTERN_C_BEGIN 3080de4c49aSJed Brown #undef __FUNCT__ 3090de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 3107087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 3110de4c49aSJed Brown { 3120de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3130de4c49aSJed Brown 3140de4c49aSJed Brown PetscFunctionBegin; 3150de4c49aSJed Brown *theta = th->Theta; 3160de4c49aSJed Brown PetscFunctionReturn(0); 3170de4c49aSJed Brown } 3180de4c49aSJed Brown 3190de4c49aSJed Brown #undef __FUNCT__ 3200de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 3217087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 3220de4c49aSJed Brown { 3230de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 3240de4c49aSJed Brown 3250de4c49aSJed Brown PetscFunctionBegin; 326e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 3270de4c49aSJed Brown th->Theta = theta; 3280de4c49aSJed Brown PetscFunctionReturn(0); 3290de4c49aSJed Brown } 330eb284becSJed Brown 331eb284becSJed Brown #undef __FUNCT__ 33278e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 33326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 33426f2ff8fSLisandro Dalcin { 33526f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 33626f2ff8fSLisandro Dalcin 33726f2ff8fSLisandro Dalcin PetscFunctionBegin; 33826f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 33926f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 34026f2ff8fSLisandro Dalcin } 34126f2ff8fSLisandro Dalcin 34226f2ff8fSLisandro Dalcin #undef __FUNCT__ 34326f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 344eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 345eb284becSJed Brown { 346eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 347eb284becSJed Brown 348eb284becSJed Brown PetscFunctionBegin; 349eb284becSJed Brown th->endpoint = flg; 350eb284becSJed Brown PetscFunctionReturn(0); 351eb284becSJed Brown } 3520de4c49aSJed Brown EXTERN_C_END 3530de4c49aSJed Brown 354*f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 355*f9c1d6abSBarry Smith #undef __FUNCT__ 356*f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta" 357*f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 358*f9c1d6abSBarry Smith { 359*f9c1d6abSBarry Smith PetscComplex z = xr + xi*PETSC_i,f; 360*f9c1d6abSBarry Smith TS_Theta *th = (TS_Theta*)ts->data; 361*f9c1d6abSBarry Smith 362*f9c1d6abSBarry Smith PetscFunctionBegin; 363*f9c1d6abSBarry Smith f = (1.0 + (1.0 - th->Theta)*z)/(1.0 - th->Theta*z); 364*f9c1d6abSBarry Smith *yr = PetscRealPartComplex(f); 365*f9c1d6abSBarry Smith *yi = PetscImaginaryPartComplex(f); 366*f9c1d6abSBarry Smith PetscFunctionReturn(0); 367*f9c1d6abSBarry Smith } 368*f9c1d6abSBarry Smith #endif 369*f9c1d6abSBarry Smith 370*f9c1d6abSBarry Smith 371316643e7SJed Brown /* ------------------------------------------------------------ */ 372316643e7SJed Brown /*MC 37396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 374316643e7SJed Brown 375316643e7SJed Brown Level: beginner 376316643e7SJed Brown 3774eb428fdSBarry Smith Options Database: 3784eb428fdSBarry Smith -ts_theta_theta <Theta> - Location of stage (0<Theta<=1); Theta = 1.0 ( 3794eb428fdSBarry Smith -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable) 3804eb428fdSBarry Smith -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method 3814eb428fdSBarry Smith 382eb284becSJed Brown Notes: 3834eb428fdSBarry Smith $ -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER) 3844eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 3854eb428fdSBarry Smith 3864eb428fdSBarry Smith 3874eb428fdSBarry Smith 388eb284becSJed Brown This method can be applied to DAE. 389eb284becSJed Brown 390eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 391eb284becSJed Brown 392eb284becSJed Brown .vb 393eb284becSJed Brown Theta | Theta 394eb284becSJed Brown ------------- 395eb284becSJed Brown | 1 396eb284becSJed Brown .ve 397eb284becSJed Brown 398eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 399eb284becSJed Brown 400eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 401eb284becSJed Brown 402eb284becSJed Brown .vb 403eb284becSJed Brown 0 | 0 0 404eb284becSJed Brown 1 | 1-Theta Theta 405eb284becSJed Brown ------------------- 406eb284becSJed Brown | 1-Theta Theta 407eb284becSJed Brown .ve 408eb284becSJed Brown 409eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 410eb284becSJed Brown 411eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 412eb284becSJed Brown 413eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 414eb284becSJed Brown 4154eb428fdSBarry Smith is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 416eb284becSJed Brown 417eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 418316643e7SJed Brown 419316643e7SJed Brown M*/ 420316643e7SJed Brown EXTERN_C_BEGIN 421316643e7SJed Brown #undef __FUNCT__ 422316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 4237087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 424316643e7SJed Brown { 425316643e7SJed Brown TS_Theta *th; 426316643e7SJed Brown PetscErrorCode ierr; 427316643e7SJed Brown 428316643e7SJed Brown PetscFunctionBegin; 429277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 430316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 431316643e7SJed Brown ts->ops->view = TSView_Theta; 432316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 433316643e7SJed Brown ts->ops->step = TSStep_Theta; 434cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 435316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 4360f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 4370f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 438*f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX) 439*f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Theta; 440*f9c1d6abSBarry Smith #endif 441316643e7SJed Brown 442316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 443316643e7SJed Brown ts->data = (void*)th; 444316643e7SJed Brown 4456f700aefSJed Brown th->extrapolate = PETSC_FALSE; 446316643e7SJed Brown th->Theta = 0.5; 447316643e7SJed Brown 4480de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 4490de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 45026f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 451eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 452316643e7SJed Brown PetscFunctionReturn(0); 453316643e7SJed Brown } 454316643e7SJed Brown EXTERN_C_END 4550de4c49aSJed Brown 4560de4c49aSJed Brown #undef __FUNCT__ 4570de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 4580de4c49aSJed Brown /*@ 4590de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 4600de4c49aSJed Brown 4610de4c49aSJed Brown Not Collective 4620de4c49aSJed Brown 4630de4c49aSJed Brown Input Parameter: 4640de4c49aSJed Brown . ts - timestepping context 4650de4c49aSJed Brown 4660de4c49aSJed Brown Output Parameter: 4670de4c49aSJed Brown . theta - stage abscissa 4680de4c49aSJed Brown 4690de4c49aSJed Brown Note: 4700de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 4710de4c49aSJed Brown 4720de4c49aSJed Brown Level: Advanced 4730de4c49aSJed Brown 4740de4c49aSJed Brown .seealso: TSThetaSetTheta() 4750de4c49aSJed Brown @*/ 4767087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 4770de4c49aSJed Brown { 4784ac538c5SBarry Smith PetscErrorCode ierr; 4790de4c49aSJed Brown 4800de4c49aSJed Brown PetscFunctionBegin; 481afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4820de4c49aSJed Brown PetscValidPointer(theta,2); 4834ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 4840de4c49aSJed Brown PetscFunctionReturn(0); 4850de4c49aSJed Brown } 4860de4c49aSJed Brown 4870de4c49aSJed Brown #undef __FUNCT__ 4880de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 4890de4c49aSJed Brown /*@ 4900de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 4910de4c49aSJed Brown 4920de4c49aSJed Brown Not Collective 4930de4c49aSJed Brown 4940de4c49aSJed Brown Input Parameter: 4950de4c49aSJed Brown + ts - timestepping context 4960de4c49aSJed Brown - theta - stage abscissa 4970de4c49aSJed Brown 4980de4c49aSJed Brown Options Database: 4990de4c49aSJed Brown . -ts_theta_theta <theta> 5000de4c49aSJed Brown 5010de4c49aSJed Brown Level: Intermediate 5020de4c49aSJed Brown 5030de4c49aSJed Brown .seealso: TSThetaGetTheta() 5040de4c49aSJed Brown @*/ 5057087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 5060de4c49aSJed Brown { 5074ac538c5SBarry Smith PetscErrorCode ierr; 5080de4c49aSJed Brown 5090de4c49aSJed Brown PetscFunctionBegin; 510afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5114ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 5120de4c49aSJed Brown PetscFunctionReturn(0); 5130de4c49aSJed Brown } 514f33bbcb6SJed Brown 515eb284becSJed Brown #undef __FUNCT__ 51626f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 51726f2ff8fSLisandro Dalcin /*@ 51826f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 51926f2ff8fSLisandro Dalcin 52026f2ff8fSLisandro Dalcin Not Collective 52126f2ff8fSLisandro Dalcin 52226f2ff8fSLisandro Dalcin Input Parameter: 52326f2ff8fSLisandro Dalcin . ts - timestepping context 52426f2ff8fSLisandro Dalcin 52526f2ff8fSLisandro Dalcin Output Parameter: 52626f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 52726f2ff8fSLisandro Dalcin 52826f2ff8fSLisandro Dalcin Level: Advanced 52926f2ff8fSLisandro Dalcin 53026f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 53126f2ff8fSLisandro Dalcin @*/ 53226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 53326f2ff8fSLisandro Dalcin { 53426f2ff8fSLisandro Dalcin PetscErrorCode ierr; 53526f2ff8fSLisandro Dalcin 53626f2ff8fSLisandro Dalcin PetscFunctionBegin; 53726f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 53826f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 53926f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 54026f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 54126f2ff8fSLisandro Dalcin } 54226f2ff8fSLisandro Dalcin 54326f2ff8fSLisandro Dalcin #undef __FUNCT__ 544eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 545eb284becSJed Brown /*@ 546eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 547eb284becSJed Brown 548eb284becSJed Brown Not Collective 549eb284becSJed Brown 550eb284becSJed Brown Input Parameter: 551eb284becSJed Brown + ts - timestepping context 552eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 553eb284becSJed Brown 554eb284becSJed Brown Options Database: 555eb284becSJed Brown . -ts_theta_endpoint <flg> 556eb284becSJed Brown 557eb284becSJed Brown Level: Intermediate 558eb284becSJed Brown 559eb284becSJed Brown .seealso: TSTHETA, TSCN 560eb284becSJed Brown @*/ 561eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 562eb284becSJed Brown { 563eb284becSJed Brown PetscErrorCode ierr; 564eb284becSJed Brown 565eb284becSJed Brown PetscFunctionBegin; 566eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 567eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 568eb284becSJed Brown PetscFunctionReturn(0); 569eb284becSJed Brown } 570eb284becSJed Brown 571f33bbcb6SJed Brown /* 572f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 573f33bbcb6SJed Brown * The creation functions for these specializations are below. 574f33bbcb6SJed Brown */ 575f33bbcb6SJed Brown 576f33bbcb6SJed Brown #undef __FUNCT__ 577f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 578f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 579f33bbcb6SJed Brown { 580d52bd9f3SBarry Smith PetscErrorCode ierr; 581d52bd9f3SBarry Smith 582f33bbcb6SJed Brown PetscFunctionBegin; 583d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 584f33bbcb6SJed Brown PetscFunctionReturn(0); 585f33bbcb6SJed Brown } 586f33bbcb6SJed Brown 587f33bbcb6SJed Brown /*MC 588f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 589f33bbcb6SJed Brown 590f33bbcb6SJed Brown Level: beginner 591f33bbcb6SJed Brown 5924eb428fdSBarry Smith Notes: 5934eb428fdSBarry Smith TSCN is equivalent to TSTHETA with Theta=1.0 5944eb428fdSBarry Smith 5954eb428fdSBarry Smith $ -ts_type theta -ts_theta_theta 1. 5964eb428fdSBarry Smith 597f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 598f33bbcb6SJed Brown 599f33bbcb6SJed Brown M*/ 600f33bbcb6SJed Brown EXTERN_C_BEGIN 601f33bbcb6SJed Brown #undef __FUNCT__ 602f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 603f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 604f33bbcb6SJed Brown { 605f33bbcb6SJed Brown PetscErrorCode ierr; 606f33bbcb6SJed Brown 607f33bbcb6SJed Brown PetscFunctionBegin; 608f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 609f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 610f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 611f33bbcb6SJed Brown PetscFunctionReturn(0); 612f33bbcb6SJed Brown } 613f33bbcb6SJed Brown EXTERN_C_END 614f33bbcb6SJed Brown 615f33bbcb6SJed Brown #undef __FUNCT__ 616f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 617f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 618f33bbcb6SJed Brown { 619d52bd9f3SBarry Smith PetscErrorCode ierr; 620d52bd9f3SBarry Smith 621f33bbcb6SJed Brown PetscFunctionBegin; 622d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 623f33bbcb6SJed Brown PetscFunctionReturn(0); 624f33bbcb6SJed Brown } 625f33bbcb6SJed Brown 626f33bbcb6SJed Brown /*MC 627f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 628f33bbcb6SJed Brown 629f33bbcb6SJed Brown Level: beginner 630f33bbcb6SJed Brown 631f33bbcb6SJed Brown Notes: 6327cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 6337cf5af47SJed Brown 6347cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 635f33bbcb6SJed Brown 636f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 637f33bbcb6SJed Brown 638f33bbcb6SJed Brown M*/ 639f33bbcb6SJed Brown EXTERN_C_BEGIN 640f33bbcb6SJed Brown #undef __FUNCT__ 641f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 642f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 643f33bbcb6SJed Brown { 644f33bbcb6SJed Brown PetscErrorCode ierr; 645f33bbcb6SJed Brown 646f33bbcb6SJed Brown PetscFunctionBegin; 647f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 648f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 649eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 650f33bbcb6SJed Brown ts->ops->view = TSView_CN; 651f33bbcb6SJed Brown PetscFunctionReturn(0); 652f33bbcb6SJed Brown } 653f33bbcb6SJed Brown EXTERN_C_END 654