1c4762a1bSJed Brown static char help[] = "Solves DAE with integrator only on non-algebraic terms \n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscts.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* 6c4762a1bSJed Brown \dot{U} = f(U,V) 7c4762a1bSJed Brown F(U,V) = 0 8c4762a1bSJed Brown 9c4762a1bSJed Brown Same as ex6.c and ex7.c except calls the ARKIMEX integrator on the entire DAE 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown f(U,V) = U + V 14c4762a1bSJed Brown 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec U,Vec V,Vec F) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(VecWAXPY(F,1.0,U,V)); 20c4762a1bSJed Brown PetscFunctionReturn(0); 21c4762a1bSJed Brown } 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* 24c4762a1bSJed Brown F(U,V) = U - V 25c4762a1bSJed Brown 26c4762a1bSJed Brown */ 27c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown PetscFunctionBeginUser; 309566063dSJacob Faibussowitsch PetscCall(VecWAXPY(F,-1.0,V,U)); 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown typedef struct { 35c4762a1bSJed Brown Vec U,V; 36c4762a1bSJed Brown Vec UF,VF; 37c4762a1bSJed Brown VecScatter scatterU,scatterV; 38c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec); 39c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec); 40c4762a1bSJed Brown } AppCtx; 41c4762a1bSJed Brown 42c4762a1bSJed Brown extern PetscErrorCode TSFunctionRHS(TS,PetscReal,Vec,Vec,void*); 43c4762a1bSJed Brown extern PetscErrorCode TSFunctionI(TS,PetscReal,Vec,Vec,Vec,void*); 44c4762a1bSJed Brown 45c4762a1bSJed Brown int main(int argc,char **argv) 46c4762a1bSJed Brown { 47c4762a1bSJed Brown AppCtx ctx; 48c4762a1bSJed Brown TS ts; 49c4762a1bSJed Brown Vec tsrhs,UV; 50c4762a1bSJed Brown IS is; 51c4762a1bSJed Brown PetscInt I; 52c4762a1bSJed Brown PetscMPIInt rank; 53c4762a1bSJed Brown 54*327415f7SBarry Smith PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 579566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 589566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 599566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSROSW)); 609566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 619566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&tsrhs)); 629566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&UV)); 639566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,tsrhs,TSFunctionRHS,&ctx)); 649566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,TSFunctionI,&ctx)); 65c4762a1bSJed Brown ctx.f = f; 66c4762a1bSJed Brown ctx.F = F; 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.U)); 699566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V)); 709566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.UF)); 719566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.VF)); 72c4762a1bSJed Brown I = 2*rank; 739566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is)); 749566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ctx.U,NULL,UV,is,&ctx.scatterU)); 759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 76c4762a1bSJed Brown I = 2*rank + 1; 779566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is)); 789566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ctx.V,NULL,UV,is,&ctx.scatterV)); 799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(VecSet(UV,1.0)); 829566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,UV)); 839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tsrhs)); 849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&UV)); 859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.U)); 869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.V)); 879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.UF)); 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.VF)); 899566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx.scatterU)); 909566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx.scatterV)); 919566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 929566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 93b122ec5aSJacob Faibussowitsch return 0; 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 98c4762a1bSJed Brown 99c4762a1bSJed Brown */ 100c4762a1bSJed Brown PetscErrorCode TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void *actx) 101c4762a1bSJed Brown { 102c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBeginUser; 1059566063dSJacob Faibussowitsch PetscCall(VecSet(F,0.0)); 1069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 1079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 1089566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 1099566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 1109566063dSJacob Faibussowitsch PetscCall((*ctx->f)(t,ctx->U,ctx->V,ctx->UF)); 1119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD)); 1129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD)); 113c4762a1bSJed Brown PetscFunctionReturn(0); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* 117c4762a1bSJed Brown Defines the nonlinear function that is passed to the time-integrator 118c4762a1bSJed Brown 119c4762a1bSJed Brown */ 120c4762a1bSJed Brown PetscErrorCode TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx) 121c4762a1bSJed Brown { 122c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscFunctionBeginUser; 1259566063dSJacob Faibussowitsch PetscCall(VecCopy(UVdot,F)); 1269566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 1279566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 1289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 1299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 1309566063dSJacob Faibussowitsch PetscCall((*ctx->F)(t,ctx->U,ctx->V,ctx->VF)); 1319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD)); 1329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD)); 133c4762a1bSJed Brown PetscFunctionReturn(0); 134c4762a1bSJed Brown } 135