Actual source code: ex9.c

  1: static char help[] = "Solves DAE with integrator only on non-algebraic terms \n";

  3: #include <petscts.h>

  5: /*
  6:         \dot{U} = f(U,V)
  7:         F(U,V)  = 0

  9:     Same as ex6.c and ex7.c except calls the ARKIMEX integrator on the entire DAE
 10: */

 12: /*
 13:    f(U,V) = U + V

 15: */
 16: PetscErrorCode f(PetscReal t,Vec U,Vec V,Vec F)
 17: {

 21:   VecWAXPY(F,1.0,U,V);
 22:   return(0);
 23: }

 25: /*
 26:    F(U,V) = U - V

 28: */
 29: PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F)
 30: {

 34:   VecWAXPY(F,-1.0,V,U);
 35:   return(0);
 36: }

 38: typedef struct {
 39:   Vec            U,V;
 40:   Vec            UF,VF;
 41:   VecScatter     scatterU,scatterV;
 42:   PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec);
 43:   PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec);
 44: } AppCtx;

 46: extern PetscErrorCode TSFunctionRHS(TS,PetscReal,Vec,Vec,void*);
 47: extern PetscErrorCode TSFunctionI(TS,PetscReal,Vec,Vec,Vec,void*);

 49: int main(int argc,char **argv)
 50: {
 52:   AppCtx         ctx;
 53:   TS             ts;
 54:   Vec            tsrhs,UV;
 55:   IS             is;
 56:   PetscInt       I;
 57:   PetscMPIInt    rank;

 59:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 60:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 61:   TSCreate(PETSC_COMM_WORLD,&ts);
 62:   TSSetProblemType(ts,TS_NONLINEAR);
 63:   TSSetType(ts,TSROSW);
 64:   TSSetFromOptions(ts);
 65:   VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&tsrhs);
 66:   VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&UV);
 67:   TSSetRHSFunction(ts,tsrhs,TSFunctionRHS,&ctx);
 68:   TSSetIFunction(ts,NULL,TSFunctionI,&ctx);
 69:   ctx.f = f;
 70:   ctx.F = F;

 72:   VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.U);
 73:   VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V);
 74:   VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.UF);
 75:   VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.VF);
 76:   I    = 2*rank;
 77:   ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is);
 78:   VecScatterCreate(ctx.U,NULL,UV,is,&ctx.scatterU);
 79:   ISDestroy(&is);
 80:   I    = 2*rank + 1;
 81:   ISCreateGeneral(PETSC_COMM_WORLD,1,&I,PETSC_COPY_VALUES,&is);
 82:   VecScatterCreate(ctx.V,NULL,UV,is,&ctx.scatterV);
 83:   ISDestroy(&is);

 85:   VecSet(UV,1.0);
 86:   TSSolve(ts,UV);
 87:   VecDestroy(&tsrhs);
 88:   VecDestroy(&UV);
 89:   VecDestroy(&ctx.U);
 90:   VecDestroy(&ctx.V);
 91:   VecDestroy(&ctx.UF);
 92:   VecDestroy(&ctx.VF);
 93:   VecScatterDestroy(&ctx.scatterU);
 94:   VecScatterDestroy(&ctx.scatterV);
 95:   TSDestroy(&ts);
 96:   PetscFinalize();
 97:   return ierr;
 98: }

100: /*
101:    Defines the RHS function that is passed to the time-integrator.

103: */
104: PetscErrorCode TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void *actx)
105: {
106:   AppCtx         *ctx = (AppCtx*)actx;

110:   VecSet(F,0.0);
111:   VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE);
112:   VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE);
113:   VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE);
114:   VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE);
115:   (*ctx->f)(t,ctx->U,ctx->V,ctx->UF);
116:   VecScatterBegin(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD);
117:   VecScatterEnd(ctx->scatterU,ctx->UF,F,INSERT_VALUES,SCATTER_FORWARD);
118:   return(0);
119: }

121: /*
122:    Defines the nonlinear function that is passed to the time-integrator

124: */
125: PetscErrorCode TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx)
126: {
127:   AppCtx         *ctx = (AppCtx*)actx;

131:   VecCopy(UVdot,F);
132:   VecScatterBegin(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE);
133:   VecScatterEnd(ctx->scatterU,UV,ctx->U,INSERT_VALUES,SCATTER_REVERSE);
134:   VecScatterBegin(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE);
135:   VecScatterEnd(ctx->scatterV,UV,ctx->V,INSERT_VALUES,SCATTER_REVERSE);
136:   (*ctx->F)(t,ctx->U,ctx->V,ctx->VF);
137:   VecScatterBegin(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD);
138:   VecScatterEnd(ctx->scatterV,ctx->VF,F,INSERT_VALUES,SCATTER_FORWARD);
139:   return(0);
140: }