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 except the user provided functions take input values as a single vector instead of two vectors 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown f(U,V) = U + V 14c4762a1bSJed Brown 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec UV,Vec F) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown const PetscScalar *u,*v; 19c4762a1bSJed Brown PetscScalar *f; 20c4762a1bSJed Brown PetscInt n,i; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBeginUser; 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(UV,&n)); 24c4762a1bSJed Brown n = n/2; 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(UV,&u)); 26c4762a1bSJed Brown v = u + n; 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(F,&f)); 28c4762a1bSJed Brown for (i=0; i<n; i++) f[i] = u[i] + v[i]; 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(UV,&u)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(F,&f)); 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown F(U,V) = U - V 36c4762a1bSJed Brown */ 37c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec UV,Vec F) 38c4762a1bSJed Brown { 39c4762a1bSJed Brown const PetscScalar *u,*v; 40c4762a1bSJed Brown PetscScalar *f; 41c4762a1bSJed Brown PetscInt n,i; 42c4762a1bSJed Brown 43c4762a1bSJed Brown PetscFunctionBeginUser; 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(UV,&n)); 45c4762a1bSJed Brown n = n/2; 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(UV,&u)); 47c4762a1bSJed Brown v = u + n; 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(F,&f)); 49c4762a1bSJed Brown for (i=0; i<n; i++) f[i] = u[i] - v[i]; 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(UV,&u)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(F,&f)); 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown typedef struct { 56c4762a1bSJed Brown PetscReal t; 57c4762a1bSJed Brown SNES snes; 58c4762a1bSJed Brown Vec UV,V; 59c4762a1bSJed Brown VecScatter scatterU,scatterV; 60c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal,Vec,Vec); 61c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal,Vec,Vec); 62c4762a1bSJed Brown } AppCtx; 63c4762a1bSJed Brown 64c4762a1bSJed Brown extern PetscErrorCode TSFunction(TS,PetscReal,Vec,Vec,void*); 65c4762a1bSJed Brown extern PetscErrorCode SNESFunction(SNES,Vec,Vec,void*); 66c4762a1bSJed Brown 67c4762a1bSJed Brown int main(int argc,char **argv) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown PetscErrorCode ierr; 70c4762a1bSJed Brown AppCtx ctx; 71c4762a1bSJed Brown TS ts; 72c4762a1bSJed Brown Vec tsrhs,U; 73c4762a1bSJed Brown IS is; 74303a5415SBarry Smith PetscInt i; 75c4762a1bSJed Brown PetscMPIInt rank; 76c4762a1bSJed Brown 77c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 78*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSEULER)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&tsrhs)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tsrhs,&U)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,tsrhs,TSFunction,&ctx)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,1.0)); 87c4762a1bSJed Brown ctx.f = f; 88c4762a1bSJed Brown 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&ctx.snes)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(ctx.snes)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(ctx.snes,NULL,SNESFunction,&ctx)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(ctx.snes,NULL,NULL,SNESComputeJacobianDefault,&ctx)); 93c4762a1bSJed Brown ctx.F = F; 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Create scatters to move between separate U and V representation and UV representation of solution */ 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&ctx.UV)); 98303a5415SBarry Smith i = 2*rank; 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(U,NULL,ctx.UV,is,&ctx.scatterU)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 102303a5415SBarry Smith i = 2*rank + 1; 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(ctx.V,NULL,ctx.UV,is,&ctx.scatterV)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 106c4762a1bSJed Brown 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(U,1.0)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 109c4762a1bSJed Brown 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx.V)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx.UV)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&ctx.scatterU)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&ctx.scatterV)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tsrhs)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&ctx.snes)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 118c4762a1bSJed Brown ierr = PetscFinalize(); 119c4762a1bSJed Brown return ierr; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* 123c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 124c4762a1bSJed Brown 125c4762a1bSJed Brown Solves F(U,V) for V and then computes f(U,V) 126c4762a1bSJed Brown */ 127c4762a1bSJed Brown PetscErrorCode TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx) 128c4762a1bSJed Brown { 129c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBeginUser; 132c4762a1bSJed Brown ctx->t = t; 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(ctx->snes,NULL,ctx->V)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->f)(t,ctx->UV,F)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* 143c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 144c4762a1bSJed Brown 145c4762a1bSJed Brown */ 146c4762a1bSJed Brown PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx) 147c4762a1bSJed Brown { 148c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ((*ctx->F)(ctx->t,ctx->UV,F)); 154c4762a1bSJed Brown PetscFunctionReturn(0); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157303a5415SBarry Smith /*TEST 158303a5415SBarry Smith 159303a5415SBarry Smith test: 160303a5415SBarry Smith args: -ts_monitor -ts_view 161303a5415SBarry Smith 162303a5415SBarry Smith TEST*/ 163