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 Same as ex6.c and ex7.c except calls the ARKIMEX integrator on the entire DAE 10*c4762a1bSJed Brown */ 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown /* 14*c4762a1bSJed Brown f(U,V) = U + V 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown */ 17*c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec UV,Vec F) 18*c4762a1bSJed Brown { 19*c4762a1bSJed Brown PetscErrorCode ierr; 20*c4762a1bSJed Brown const PetscScalar *u,*v; 21*c4762a1bSJed Brown PetscScalar *f; 22*c4762a1bSJed Brown PetscInt n,i; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown PetscFunctionBeginUser; 25*c4762a1bSJed Brown ierr = VecGetLocalSize(UV,&n);CHKERRQ(ierr); 26*c4762a1bSJed Brown n = n/2; 27*c4762a1bSJed Brown ierr = VecGetArrayRead(UV,&u);CHKERRQ(ierr); 28*c4762a1bSJed Brown v = u + n; 29*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 30*c4762a1bSJed Brown for (i=0; i<n; i++) f[i] = u[i] + v[i]; 31*c4762a1bSJed Brown ierr = VecRestoreArrayRead(UV,&u);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 33*c4762a1bSJed Brown PetscFunctionReturn(0); 34*c4762a1bSJed Brown } 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown /* 37*c4762a1bSJed Brown F(U,V) = U - V 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown */ 40*c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec UV,Vec F) 41*c4762a1bSJed Brown { 42*c4762a1bSJed Brown PetscErrorCode ierr; 43*c4762a1bSJed Brown const PetscScalar *u,*v; 44*c4762a1bSJed Brown PetscScalar *f; 45*c4762a1bSJed Brown PetscInt n,i; 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown PetscFunctionBeginUser; 48*c4762a1bSJed Brown ierr = VecGetLocalSize(UV,&n);CHKERRQ(ierr); 49*c4762a1bSJed Brown n = n/2; 50*c4762a1bSJed Brown ierr = VecGetArrayRead(UV,&u);CHKERRQ(ierr); 51*c4762a1bSJed Brown v = u + n; 52*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 53*c4762a1bSJed Brown f = f + n; 54*c4762a1bSJed Brown for (i=0; i<n; i++) f[i] = u[i] - v[i]; 55*c4762a1bSJed Brown f = f - n; 56*c4762a1bSJed Brown ierr = VecRestoreArrayRead(UV,&u);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 58*c4762a1bSJed Brown PetscFunctionReturn(0); 59*c4762a1bSJed Brown } 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown typedef struct { 62*c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal,Vec,Vec); 63*c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal,Vec,Vec); 64*c4762a1bSJed Brown } AppCtx; 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown extern PetscErrorCode TSFunctionRHS(TS,PetscReal,Vec,Vec,void*); 67*c4762a1bSJed Brown extern PetscErrorCode TSFunctionI(TS,PetscReal,Vec,Vec,Vec,void*); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown int main(int argc,char **argv) 70*c4762a1bSJed Brown { 71*c4762a1bSJed Brown PetscErrorCode ierr; 72*c4762a1bSJed Brown AppCtx ctx; 73*c4762a1bSJed Brown TS ts; 74*c4762a1bSJed Brown Vec tsrhs,UV; 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 77*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&tsrhs);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&UV);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,tsrhs,TSFunctionRHS,&ctx);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,TSFunctionI,&ctx);CHKERRQ(ierr); 85*c4762a1bSJed Brown ctx.f = f; 86*c4762a1bSJed Brown ctx.F = F; 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown ierr = VecSet(UV,1.0);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = TSSolve(ts,UV);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = VecDestroy(&tsrhs);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = VecDestroy(&UV);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = PetscFinalize(); 94*c4762a1bSJed Brown return ierr; 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /* 98*c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown */ 103*c4762a1bSJed Brown PetscErrorCode TSFunctionRHS(TS ts,PetscReal t,Vec UV,Vec F,void *actx) 104*c4762a1bSJed Brown { 105*c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 106*c4762a1bSJed Brown PetscErrorCode ierr; 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown PetscFunctionBeginUser; 109*c4762a1bSJed Brown ierr = VecSet(F,0.0);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = (*ctx->f)(t,UV,F);CHKERRQ(ierr); 111*c4762a1bSJed Brown PetscFunctionReturn(0); 112*c4762a1bSJed Brown } 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown /* 115*c4762a1bSJed Brown Defines the nonlinear function that is passed to the time-integrator 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown */ 118*c4762a1bSJed Brown PetscErrorCode TSFunctionI(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx) 119*c4762a1bSJed Brown { 120*c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 121*c4762a1bSJed Brown PetscErrorCode ierr; 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown PetscFunctionBeginUser; 124*c4762a1bSJed Brown ierr = VecCopy(UVdot,F);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = (*ctx->F)(t,UV,F);CHKERRQ(ierr); 126*c4762a1bSJed Brown PetscFunctionReturn(0); 127*c4762a1bSJed Brown } 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown 130