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) { 277445fe48SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"TSTheta_X0",(PetscObject*)X0);CHKERRQ(ierr); 287445fe48SJed Brown if (!*X0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"TSTheta_X0 has not been composed with DM from SNES"); 297445fe48SJed Brown } else *X0 = ts->vec_sol; 307445fe48SJed Brown } 317445fe48SJed Brown if (Xdot) { 327445fe48SJed Brown if (dm && dm != ts->dm) { 337445fe48SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"TSTheta_Xdot",(PetscObject*)Xdot);CHKERRQ(ierr); 347445fe48SJed Brown if (!*Xdot) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"TSTheta_Xdot has not been composed with DM from SNES"); 357445fe48SJed Brown } else *Xdot = th->Xdot; 367445fe48SJed Brown } 377445fe48SJed Brown PetscFunctionReturn(0); 387445fe48SJed Brown } 397445fe48SJed Brown 407445fe48SJed Brown #undef __FUNCT__ 417445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta" 427445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 437445fe48SJed Brown { 447445fe48SJed Brown Vec X0,Xdot; 457445fe48SJed Brown PetscErrorCode ierr; 467445fe48SJed Brown 477445fe48SJed Brown PetscFunctionBegin; 487445fe48SJed Brown ierr = DMCreateGlobalVector(coarse,&X0);CHKERRQ(ierr); 497445fe48SJed Brown ierr = DMCreateGlobalVector(coarse,&Xdot);CHKERRQ(ierr); 507445fe48SJed Brown /* Oh noes, this would create a loop because the Vec holds a reference to the DM. 517445fe48SJed Brown Making a PetscContainer to hold these Vecs would make the following call succeed, but would create a reference loop. 527445fe48SJed Brown Need to decide on a way to break the reference counting loop. 537445fe48SJed Brown */ 547445fe48SJed Brown ierr = PetscObjectCompose((PetscObject)coarse,"TSTheta_X0",(PetscObject)X0);CHKERRQ(ierr); 557445fe48SJed Brown ierr = PetscObjectCompose((PetscObject)coarse,"TSTheta_Xdot",(PetscObject)Xdot);CHKERRQ(ierr); 567445fe48SJed Brown ierr = VecDestroy(&X0);CHKERRQ(ierr); 577445fe48SJed Brown ierr = VecDestroy(&Xdot);CHKERRQ(ierr); 587445fe48SJed Brown PetscFunctionReturn(0); 597445fe48SJed Brown } 607445fe48SJed Brown 617445fe48SJed Brown #undef __FUNCT__ 627445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta" 637445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 647445fe48SJed Brown { 657445fe48SJed Brown TS ts = (TS)ctx; 667445fe48SJed Brown PetscErrorCode ierr; 677445fe48SJed Brown Vec X0,Xdot,X0_c,Xdot_c; 687445fe48SJed Brown 697445fe48SJed Brown PetscFunctionBegin; 707445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 717445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 727445fe48SJed Brown ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 737445fe48SJed Brown ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 747445fe48SJed Brown ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 757445fe48SJed Brown ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 767445fe48SJed Brown PetscFunctionReturn(0); 777445fe48SJed Brown } 787445fe48SJed Brown 797445fe48SJed Brown #undef __FUNCT__ 80316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 81193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 82316643e7SJed Brown { 83316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 84b70ae86eSJed Brown PetscInt its,lits; 85cdbf8f93SLisandro Dalcin PetscReal next_time_step; 86f1b97656SJed Brown SNESConvergedReason snesreason; 872b5a38e1SLisandro Dalcin PetscErrorCode ierr; 88316643e7SJed Brown 89316643e7SJed Brown PetscFunctionBegin; 90cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 91eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 92316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 93*b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 94*b8123daeSJed Brown ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 95316643e7SJed Brown 96eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 97eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 98eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 99eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 100eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 101eb284becSJed Brown } 102ace68cafSJed Brown if (th->extrapolate) { 1032b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 104ace68cafSJed Brown } else { 1052b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 106ace68cafSJed Brown } 107eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 108316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 109316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 110f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 1115ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 112f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 113f1b97656SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 114f1b97656SJed 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); 115f1b97656SJed Brown PetscFunctionReturn(0); 116f1b97656SJed Brown } 117eb284becSJed Brown if (th->endpoint) { 118eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 119eb284becSJed Brown } else { 1202b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 1212b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 122eb284becSJed Brown } 1232b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 124cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 125316643e7SJed Brown ts->steps++; 126316643e7SJed Brown PetscFunctionReturn(0); 127316643e7SJed Brown } 128316643e7SJed Brown 129cd652676SJed Brown #undef __FUNCT__ 130cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 131cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 132cd652676SJed Brown { 133cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1345a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 135cd652676SJed Brown PetscErrorCode ierr; 136cd652676SJed Brown 137cd652676SJed Brown PetscFunctionBegin; 138a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 1395a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 1405a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 141cd652676SJed Brown PetscFunctionReturn(0); 142cd652676SJed Brown } 143cd652676SJed Brown 144316643e7SJed Brown /*------------------------------------------------------------*/ 145316643e7SJed Brown #undef __FUNCT__ 146277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 147277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 148316643e7SJed Brown { 149316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 150316643e7SJed Brown PetscErrorCode ierr; 151316643e7SJed Brown 152316643e7SJed Brown PetscFunctionBegin; 1536bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 1546bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 155eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 156277b19d0SLisandro Dalcin PetscFunctionReturn(0); 157277b19d0SLisandro Dalcin } 158277b19d0SLisandro Dalcin 159277b19d0SLisandro Dalcin #undef __FUNCT__ 160277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 161277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 162277b19d0SLisandro Dalcin { 163277b19d0SLisandro Dalcin PetscErrorCode ierr; 164277b19d0SLisandro Dalcin 165277b19d0SLisandro Dalcin PetscFunctionBegin; 166277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 167277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 168335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 169335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 17026f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 171eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 172316643e7SJed Brown PetscFunctionReturn(0); 173316643e7SJed Brown } 174316643e7SJed Brown 175316643e7SJed Brown /* 176316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1772b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 178316643e7SJed Brown */ 179316643e7SJed Brown #undef __FUNCT__ 1800f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1810f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 182316643e7SJed Brown { 183316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 184316643e7SJed Brown PetscErrorCode ierr; 1857445fe48SJed Brown Vec X0,Xdot; 1867445fe48SJed Brown DM dm,dmsave; 187316643e7SJed Brown 188316643e7SJed Brown PetscFunctionBegin; 1897445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1905a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 1917445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 1927445fe48SJed Brown ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr); 1937445fe48SJed Brown 1947445fe48SJed Brown /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 1957445fe48SJed Brown dmsave = ts->dm; 1967445fe48SJed Brown ts->dm = dm; 1977445fe48SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 1987445fe48SJed Brown ts->dm = dmsave; 199316643e7SJed Brown PetscFunctionReturn(0); 200316643e7SJed Brown } 201316643e7SJed Brown 202316643e7SJed Brown #undef __FUNCT__ 2030f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 2040f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 205316643e7SJed Brown { 206316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 207316643e7SJed Brown PetscErrorCode ierr; 2087445fe48SJed Brown Vec Xdot; 2097445fe48SJed Brown DM dm,dmsave; 210316643e7SJed Brown 211316643e7SJed Brown PetscFunctionBegin; 2127445fe48SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2137445fe48SJed Brown 2140f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 2157445fe48SJed Brown ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr); 2167445fe48SJed Brown 2177445fe48SJed Brown dmsave = ts->dm; 2187445fe48SJed Brown ts->dm = dm; 2197445fe48SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 2207445fe48SJed Brown ts->dm = dmsave; 221316643e7SJed Brown PetscFunctionReturn(0); 222316643e7SJed Brown } 223316643e7SJed Brown 224316643e7SJed Brown #undef __FUNCT__ 225316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 226316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 227316643e7SJed Brown { 228316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 229316643e7SJed Brown PetscErrorCode ierr; 2307445fe48SJed Brown SNES snes; 2317445fe48SJed Brown DM dm; 232316643e7SJed Brown 233316643e7SJed Brown PetscFunctionBegin; 234316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 235316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 2367445fe48SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2377445fe48SJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2387445fe48SJed Brown if (dm) { 2397445fe48SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 2407445fe48SJed Brown } 241316643e7SJed Brown PetscFunctionReturn(0); 242316643e7SJed Brown } 243316643e7SJed Brown /*------------------------------------------------------------*/ 244316643e7SJed Brown 245316643e7SJed Brown #undef __FUNCT__ 246316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 247316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 248316643e7SJed Brown { 249316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 250316643e7SJed Brown PetscErrorCode ierr; 251316643e7SJed Brown 252316643e7SJed Brown PetscFunctionBegin; 253d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 254316643e7SJed Brown { 255316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 256acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 257eb284becSJed 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); 258d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 259316643e7SJed Brown } 260316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 261316643e7SJed Brown PetscFunctionReturn(0); 262316643e7SJed Brown } 263316643e7SJed Brown 264316643e7SJed Brown #undef __FUNCT__ 265316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 266316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 267316643e7SJed Brown { 268316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 269ace3abfcSBarry Smith PetscBool iascii; 270316643e7SJed Brown PetscErrorCode ierr; 271316643e7SJed Brown 272316643e7SJed Brown PetscFunctionBegin; 273251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 274316643e7SJed Brown if (iascii) { 275316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 276ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 277316643e7SJed Brown } 278d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 279316643e7SJed Brown PetscFunctionReturn(0); 280316643e7SJed Brown } 281316643e7SJed Brown 2820de4c49aSJed Brown EXTERN_C_BEGIN 2830de4c49aSJed Brown #undef __FUNCT__ 2840de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 2857087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 2860de4c49aSJed Brown { 2870de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2880de4c49aSJed Brown 2890de4c49aSJed Brown PetscFunctionBegin; 2900de4c49aSJed Brown *theta = th->Theta; 2910de4c49aSJed Brown PetscFunctionReturn(0); 2920de4c49aSJed Brown } 2930de4c49aSJed Brown 2940de4c49aSJed Brown #undef __FUNCT__ 2950de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 2967087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 2970de4c49aSJed Brown { 2980de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2990de4c49aSJed Brown 3000de4c49aSJed Brown PetscFunctionBegin; 301e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 3020de4c49aSJed Brown th->Theta = theta; 3030de4c49aSJed Brown PetscFunctionReturn(0); 3040de4c49aSJed Brown } 305eb284becSJed Brown 306eb284becSJed Brown #undef __FUNCT__ 30778e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 30826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 30926f2ff8fSLisandro Dalcin { 31026f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 31126f2ff8fSLisandro Dalcin 31226f2ff8fSLisandro Dalcin PetscFunctionBegin; 31326f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 31426f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 31526f2ff8fSLisandro Dalcin } 31626f2ff8fSLisandro Dalcin 31726f2ff8fSLisandro Dalcin #undef __FUNCT__ 31826f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 319eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 320eb284becSJed Brown { 321eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 322eb284becSJed Brown 323eb284becSJed Brown PetscFunctionBegin; 324eb284becSJed Brown th->endpoint = flg; 325eb284becSJed Brown PetscFunctionReturn(0); 326eb284becSJed Brown } 3270de4c49aSJed Brown EXTERN_C_END 3280de4c49aSJed Brown 329316643e7SJed Brown /* ------------------------------------------------------------ */ 330316643e7SJed Brown /*MC 33196f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 332316643e7SJed Brown 333316643e7SJed Brown Level: beginner 334316643e7SJed Brown 335eb284becSJed Brown Notes: 336eb284becSJed Brown This method can be applied to DAE. 337eb284becSJed Brown 338eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 339eb284becSJed Brown 340eb284becSJed Brown .vb 341eb284becSJed Brown Theta | Theta 342eb284becSJed Brown ------------- 343eb284becSJed Brown | 1 344eb284becSJed Brown .ve 345eb284becSJed Brown 346eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 347eb284becSJed Brown 348eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 349eb284becSJed Brown 350eb284becSJed Brown .vb 351eb284becSJed Brown 0 | 0 0 352eb284becSJed Brown 1 | 1-Theta Theta 353eb284becSJed Brown ------------------- 354eb284becSJed Brown | 1-Theta Theta 355eb284becSJed Brown .ve 356eb284becSJed Brown 357eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 358eb284becSJed Brown 359eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 360eb284becSJed Brown 361eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 362eb284becSJed Brown 363eb284becSJed Brown is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 364eb284becSJed Brown 365eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 366316643e7SJed Brown 367316643e7SJed Brown M*/ 368316643e7SJed Brown EXTERN_C_BEGIN 369316643e7SJed Brown #undef __FUNCT__ 370316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 3717087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 372316643e7SJed Brown { 373316643e7SJed Brown TS_Theta *th; 374316643e7SJed Brown PetscErrorCode ierr; 375316643e7SJed Brown 376316643e7SJed Brown PetscFunctionBegin; 377277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 378316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 379316643e7SJed Brown ts->ops->view = TSView_Theta; 380316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 381316643e7SJed Brown ts->ops->step = TSStep_Theta; 382cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 383316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 3840f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 3850f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 386316643e7SJed Brown 387316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 388316643e7SJed Brown ts->data = (void*)th; 389316643e7SJed Brown 3906f700aefSJed Brown th->extrapolate = PETSC_FALSE; 391316643e7SJed Brown th->Theta = 0.5; 392316643e7SJed Brown 3930de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 3940de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 39526f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 396eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 397316643e7SJed Brown PetscFunctionReturn(0); 398316643e7SJed Brown } 399316643e7SJed Brown EXTERN_C_END 4000de4c49aSJed Brown 4010de4c49aSJed Brown #undef __FUNCT__ 4020de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 4030de4c49aSJed Brown /*@ 4040de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 4050de4c49aSJed Brown 4060de4c49aSJed Brown Not Collective 4070de4c49aSJed Brown 4080de4c49aSJed Brown Input Parameter: 4090de4c49aSJed Brown . ts - timestepping context 4100de4c49aSJed Brown 4110de4c49aSJed Brown Output Parameter: 4120de4c49aSJed Brown . theta - stage abscissa 4130de4c49aSJed Brown 4140de4c49aSJed Brown Note: 4150de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 4160de4c49aSJed Brown 4170de4c49aSJed Brown Level: Advanced 4180de4c49aSJed Brown 4190de4c49aSJed Brown .seealso: TSThetaSetTheta() 4200de4c49aSJed Brown @*/ 4217087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 4220de4c49aSJed Brown { 4234ac538c5SBarry Smith PetscErrorCode ierr; 4240de4c49aSJed Brown 4250de4c49aSJed Brown PetscFunctionBegin; 426afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4270de4c49aSJed Brown PetscValidPointer(theta,2); 4284ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 4290de4c49aSJed Brown PetscFunctionReturn(0); 4300de4c49aSJed Brown } 4310de4c49aSJed Brown 4320de4c49aSJed Brown #undef __FUNCT__ 4330de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 4340de4c49aSJed Brown /*@ 4350de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 4360de4c49aSJed Brown 4370de4c49aSJed Brown Not Collective 4380de4c49aSJed Brown 4390de4c49aSJed Brown Input Parameter: 4400de4c49aSJed Brown + ts - timestepping context 4410de4c49aSJed Brown - theta - stage abscissa 4420de4c49aSJed Brown 4430de4c49aSJed Brown Options Database: 4440de4c49aSJed Brown . -ts_theta_theta <theta> 4450de4c49aSJed Brown 4460de4c49aSJed Brown Level: Intermediate 4470de4c49aSJed Brown 4480de4c49aSJed Brown .seealso: TSThetaGetTheta() 4490de4c49aSJed Brown @*/ 4507087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 4510de4c49aSJed Brown { 4524ac538c5SBarry Smith PetscErrorCode ierr; 4530de4c49aSJed Brown 4540de4c49aSJed Brown PetscFunctionBegin; 455afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4564ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 4570de4c49aSJed Brown PetscFunctionReturn(0); 4580de4c49aSJed Brown } 459f33bbcb6SJed Brown 460eb284becSJed Brown #undef __FUNCT__ 46126f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 46226f2ff8fSLisandro Dalcin /*@ 46326f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 46426f2ff8fSLisandro Dalcin 46526f2ff8fSLisandro Dalcin Not Collective 46626f2ff8fSLisandro Dalcin 46726f2ff8fSLisandro Dalcin Input Parameter: 46826f2ff8fSLisandro Dalcin . ts - timestepping context 46926f2ff8fSLisandro Dalcin 47026f2ff8fSLisandro Dalcin Output Parameter: 47126f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 47226f2ff8fSLisandro Dalcin 47326f2ff8fSLisandro Dalcin Level: Advanced 47426f2ff8fSLisandro Dalcin 47526f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 47626f2ff8fSLisandro Dalcin @*/ 47726f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 47826f2ff8fSLisandro Dalcin { 47926f2ff8fSLisandro Dalcin PetscErrorCode ierr; 48026f2ff8fSLisandro Dalcin 48126f2ff8fSLisandro Dalcin PetscFunctionBegin; 48226f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 48326f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 48426f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 48526f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 48626f2ff8fSLisandro Dalcin } 48726f2ff8fSLisandro Dalcin 48826f2ff8fSLisandro Dalcin #undef __FUNCT__ 489eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 490eb284becSJed Brown /*@ 491eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 492eb284becSJed Brown 493eb284becSJed Brown Not Collective 494eb284becSJed Brown 495eb284becSJed Brown Input Parameter: 496eb284becSJed Brown + ts - timestepping context 497eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 498eb284becSJed Brown 499eb284becSJed Brown Options Database: 500eb284becSJed Brown . -ts_theta_endpoint <flg> 501eb284becSJed Brown 502eb284becSJed Brown Level: Intermediate 503eb284becSJed Brown 504eb284becSJed Brown .seealso: TSTHETA, TSCN 505eb284becSJed Brown @*/ 506eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 507eb284becSJed Brown { 508eb284becSJed Brown PetscErrorCode ierr; 509eb284becSJed Brown 510eb284becSJed Brown PetscFunctionBegin; 511eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 512eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 513eb284becSJed Brown PetscFunctionReturn(0); 514eb284becSJed Brown } 515eb284becSJed Brown 516f33bbcb6SJed Brown /* 517f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 518f33bbcb6SJed Brown * The creation functions for these specializations are below. 519f33bbcb6SJed Brown */ 520f33bbcb6SJed Brown 521f33bbcb6SJed Brown #undef __FUNCT__ 522f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 523f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 524f33bbcb6SJed Brown { 525d52bd9f3SBarry Smith PetscErrorCode ierr; 526d52bd9f3SBarry Smith 527f33bbcb6SJed Brown PetscFunctionBegin; 528d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 529f33bbcb6SJed Brown PetscFunctionReturn(0); 530f33bbcb6SJed Brown } 531f33bbcb6SJed Brown 532f33bbcb6SJed Brown /*MC 533f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 534f33bbcb6SJed Brown 535f33bbcb6SJed Brown Level: beginner 536f33bbcb6SJed Brown 537f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 538f33bbcb6SJed Brown 539f33bbcb6SJed Brown M*/ 540f33bbcb6SJed Brown EXTERN_C_BEGIN 541f33bbcb6SJed Brown #undef __FUNCT__ 542f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 543f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 544f33bbcb6SJed Brown { 545f33bbcb6SJed Brown PetscErrorCode ierr; 546f33bbcb6SJed Brown 547f33bbcb6SJed Brown PetscFunctionBegin; 548f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 549f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 550f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 551f33bbcb6SJed Brown PetscFunctionReturn(0); 552f33bbcb6SJed Brown } 553f33bbcb6SJed Brown EXTERN_C_END 554f33bbcb6SJed Brown 555f33bbcb6SJed Brown #undef __FUNCT__ 556f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 557f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 558f33bbcb6SJed Brown { 559d52bd9f3SBarry Smith PetscErrorCode ierr; 560d52bd9f3SBarry Smith 561f33bbcb6SJed Brown PetscFunctionBegin; 562d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 563f33bbcb6SJed Brown PetscFunctionReturn(0); 564f33bbcb6SJed Brown } 565f33bbcb6SJed Brown 566f33bbcb6SJed Brown /*MC 567f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 568f33bbcb6SJed Brown 569f33bbcb6SJed Brown Level: beginner 570f33bbcb6SJed Brown 571f33bbcb6SJed Brown Notes: 5727cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 5737cf5af47SJed Brown 5747cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 575f33bbcb6SJed Brown 576f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 577f33bbcb6SJed Brown 578f33bbcb6SJed Brown M*/ 579f33bbcb6SJed Brown EXTERN_C_BEGIN 580f33bbcb6SJed Brown #undef __FUNCT__ 581f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 582f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 583f33bbcb6SJed Brown { 584f33bbcb6SJed Brown PetscErrorCode ierr; 585f33bbcb6SJed Brown 586f33bbcb6SJed Brown PetscFunctionBegin; 587f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 588f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 589eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 590f33bbcb6SJed Brown ts->ops->view = TSView_CN; 591f33bbcb6SJed Brown PetscFunctionReturn(0); 592f33bbcb6SJed Brown } 593f33bbcb6SJed Brown EXTERN_C_END 594