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