1*c4762a1bSJed Brown static char help[] = "Solves DAE with integrator only on non-algebraic terms \n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscts.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown /* 6*c4762a1bSJed Brown \dot{U} = f(U,V) 7*c4762a1bSJed Brown F(U,V) = 0 8*c4762a1bSJed Brown */ 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown /* 12*c4762a1bSJed Brown f(U,V) = U + V 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown */ 15*c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec U,Vec V,Vec F) 16*c4762a1bSJed Brown { 17*c4762a1bSJed Brown PetscErrorCode ierr; 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown PetscFunctionBeginUser; 20*c4762a1bSJed Brown ierr = VecWAXPY(F,1.0,U,V);CHKERRQ(ierr); 21*c4762a1bSJed Brown PetscFunctionReturn(0); 22*c4762a1bSJed Brown } 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown /* 25*c4762a1bSJed Brown F(U,V) = U - V 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown */ 28*c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F) 29*c4762a1bSJed Brown { 30*c4762a1bSJed Brown PetscErrorCode ierr; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown PetscFunctionBeginUser; 33*c4762a1bSJed Brown ierr = VecWAXPY(F,-1.0,V,U);CHKERRQ(ierr); 34*c4762a1bSJed Brown PetscFunctionReturn(0); 35*c4762a1bSJed Brown } 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown typedef struct { 38*c4762a1bSJed Brown PetscReal t; 39*c4762a1bSJed Brown SNES snes; 40*c4762a1bSJed Brown Vec U,V; 41*c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec); 42*c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec); 43*c4762a1bSJed Brown } AppCtx; 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown extern PetscErrorCode TSFunction(TS,PetscReal,Vec,Vec,void*); 46*c4762a1bSJed Brown extern PetscErrorCode SNESFunction(SNES,Vec,Vec,void*); 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown int main(int argc,char **argv) 49*c4762a1bSJed Brown { 50*c4762a1bSJed Brown PetscErrorCode ierr; 51*c4762a1bSJed Brown AppCtx ctx; 52*c4762a1bSJed Brown TS ts; 53*c4762a1bSJed Brown Vec tsrhs,U; 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 56*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,1,&tsrhs);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,1,&U);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,tsrhs,TSFunction,&ctx);CHKERRQ(ierr); 63*c4762a1bSJed Brown ctx.f = f; 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&ctx.snes);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = SNESSetFromOptions(ctx.snes);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = SNESSetFunction(ctx.snes,NULL,SNESFunction,&ctx);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = SNESSetJacobian(ctx.snes,NULL,NULL,SNESComputeJacobianDefault,&ctx);CHKERRQ(ierr); 69*c4762a1bSJed Brown ctx.F = F; 70*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,1,&ctx.V);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown ierr = VecSet(U,1.0);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown ierr = VecDestroy(&ctx.V);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = VecDestroy(&tsrhs);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = SNESDestroy(&ctx.snes);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = PetscFinalize(); 81*c4762a1bSJed Brown return ierr; 82*c4762a1bSJed Brown } 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /* 85*c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown Solves F(U,V) for V and then computes f(U,V) 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown */ 90*c4762a1bSJed Brown PetscErrorCode TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx) 91*c4762a1bSJed Brown { 92*c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 93*c4762a1bSJed Brown PetscErrorCode ierr; 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown PetscFunctionBeginUser; 96*c4762a1bSJed Brown ctx->t = t; 97*c4762a1bSJed Brown ctx->U = U; 98*c4762a1bSJed Brown ierr = SNESSolve(ctx->snes,NULL,ctx->V);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = (*ctx->f)(t,U,ctx->V,F);CHKERRQ(ierr); 100*c4762a1bSJed Brown PetscFunctionReturn(0); 101*c4762a1bSJed Brown } 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown /* 104*c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown */ 107*c4762a1bSJed Brown PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx) 108*c4762a1bSJed Brown { 109*c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 110*c4762a1bSJed Brown PetscErrorCode ierr; 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown PetscFunctionBeginUser; 113*c4762a1bSJed Brown ierr = (*ctx->F)(ctx->t,ctx->U,V,F);CHKERRQ(ierr); 114*c4762a1bSJed Brown PetscFunctionReturn(0); 115*c4762a1bSJed Brown } 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown 118