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