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 PetscErrorCode ierr; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBeginUser; 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(F,1.0,U,V)); 22c4762a1bSJed Brown PetscFunctionReturn(0); 23c4762a1bSJed Brown } 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown F(U,V) = U - V 27c4762a1bSJed Brown 28c4762a1bSJed Brown */ 29c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F) 30c4762a1bSJed Brown { 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscFunctionBeginUser; 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(F,-1.0,V,U)); 35c4762a1bSJed Brown PetscFunctionReturn(0); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown typedef struct { 39c4762a1bSJed Brown Vec U,V; 40c4762a1bSJed Brown Vec UF,VF; 41c4762a1bSJed Brown VecScatter scatterU,scatterV; 42c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec); 43c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec); 44c4762a1bSJed Brown } AppCtx; 45c4762a1bSJed Brown 46c4762a1bSJed Brown extern PetscErrorCode TSFunctionRHS(TS,PetscReal,Vec,Vec,void*); 47c4762a1bSJed Brown extern PetscErrorCode TSFunctionI(TS,PetscReal,Vec,Vec,Vec,void*); 48c4762a1bSJed Brown 49c4762a1bSJed Brown int main(int argc,char **argv) 50c4762a1bSJed Brown { 51c4762a1bSJed Brown PetscErrorCode ierr; 52c4762a1bSJed Brown AppCtx ctx; 53c4762a1bSJed Brown TS ts; 54c4762a1bSJed Brown Vec tsrhs,UV; 55c4762a1bSJed Brown IS is; 56c4762a1bSJed Brown PetscInt I; 57c4762a1bSJed Brown PetscMPIInt rank; 58c4762a1bSJed Brown 59c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 60*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSROSW)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&tsrhs)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&UV)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,tsrhs,TSFunctionRHS,&ctx)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,TSFunctionI,&ctx)); 69c4762a1bSJed Brown ctx.f = f; 70c4762a1bSJed Brown ctx.F = F; 71c4762a1bSJed Brown 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.U)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.UF)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.VF)); 76c4762a1bSJed Brown I = 2*rank; 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(ctx.U,NULL,UV,is,&ctx.scatterU)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 80c4762a1bSJed Brown I = 2*rank + 1; 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(ctx.V,NULL,UV,is,&ctx.scatterV)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 84c4762a1bSJed Brown 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(UV,1.0)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,UV)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tsrhs)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&UV)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx.U)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx.V)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx.UF)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx.VF)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&ctx.scatterU)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&ctx.scatterV)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 96c4762a1bSJed Brown ierr = PetscFinalize(); 97c4762a1bSJed Brown return ierr; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 102c4762a1bSJed Brown 103c4762a1bSJed Brown */ 104c4762a1bSJed Brown PetscErrorCode TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void *actx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 107c4762a1bSJed Brown PetscErrorCode ierr; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBeginUser; 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(F,0.0)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->f)(t,ctx->U,ctx->V,ctx->UF)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD)); 118c4762a1bSJed Brown PetscFunctionReturn(0); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown Defines the nonlinear function that is passed to the time-integrator 123c4762a1bSJed Brown 124c4762a1bSJed Brown */ 125c4762a1bSJed Brown PetscErrorCode TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx) 126c4762a1bSJed Brown { 127c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 128c4762a1bSJed Brown PetscErrorCode ierr; 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscFunctionBeginUser; 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(UVdot,F)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->F)(t,ctx->U,ctx->V,ctx->VF)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141