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*9566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(UV,&n)); 24c4762a1bSJed Brown n = n/2; 25*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(UV,&u)); 26c4762a1bSJed Brown v = u + n; 27*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&f)); 28c4762a1bSJed Brown for (i=0; i<n; i++) f[i] = u[i] + v[i]; 29*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(UV,&u)); 30*9566063dSJacob Faibussowitsch PetscCall(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*9566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(UV,&n)); 45c4762a1bSJed Brown n = n/2; 46*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(UV,&u)); 47c4762a1bSJed Brown v = u + n; 48*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&f)); 49c4762a1bSJed Brown for (i=0; i<n; i++) f[i] = u[i] - v[i]; 50*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(UV,&u)); 51*9566063dSJacob Faibussowitsch PetscCall(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 AppCtx ctx; 70c4762a1bSJed Brown TS ts; 71c4762a1bSJed Brown Vec tsrhs,U; 72c4762a1bSJed Brown IS is; 73303a5415SBarry Smith PetscInt i; 74c4762a1bSJed Brown PetscMPIInt rank; 75c4762a1bSJed Brown 76*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 77*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 78*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 79*9566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 80*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSEULER)); 81*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 82*9566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&tsrhs)); 83*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tsrhs,&U)); 84*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,tsrhs,TSFunction,&ctx)); 85*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,1.0)); 86c4762a1bSJed Brown ctx.f = f; 87c4762a1bSJed Brown 88*9566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&ctx.snes)); 89*9566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(ctx.snes)); 90*9566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(ctx.snes,NULL,SNESFunction,&ctx)); 91*9566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(ctx.snes,NULL,NULL,SNESComputeJacobianDefault,&ctx)); 92c4762a1bSJed Brown ctx.F = F; 93*9566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&ctx.V)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* Create scatters to move between separate U and V representation and UV representation of solution */ 96*9566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD,2,PETSC_DETERMINE,&ctx.UV)); 97303a5415SBarry Smith i = 2*rank; 98*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is)); 99*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(U,NULL,ctx.UV,is,&ctx.scatterU)); 100*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 101303a5415SBarry Smith i = 2*rank + 1; 102*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,&i,PETSC_COPY_VALUES,&is)); 103*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ctx.V,NULL,ctx.UV,is,&ctx.scatterV)); 104*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 105c4762a1bSJed Brown 106*9566063dSJacob Faibussowitsch PetscCall(VecSet(U,1.0)); 107*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 108c4762a1bSJed Brown 109*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.V)); 110*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.UV)); 111*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx.scatterU)); 112*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx.scatterV)); 113*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tsrhs)); 114*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 115*9566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&ctx.snes)); 116*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 117*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 118b122ec5aSJacob Faibussowitsch return 0; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 123c4762a1bSJed Brown 124c4762a1bSJed Brown Solves F(U,V) for V and then computes f(U,V) 125c4762a1bSJed Brown */ 126c4762a1bSJed Brown PetscErrorCode TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscFunctionBeginUser; 131c4762a1bSJed Brown ctx->t = t; 132*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 133*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterU,U,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 134*9566063dSJacob Faibussowitsch PetscCall(SNESSolve(ctx->snes,NULL,ctx->V)); 135*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 136*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV,ctx->V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 137*9566063dSJacob Faibussowitsch PetscCall((*ctx->f)(t,ctx->UV,F)); 138c4762a1bSJed Brown PetscFunctionReturn(0); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* 142c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 143c4762a1bSJed Brown 144c4762a1bSJed Brown */ 145c4762a1bSJed Brown PetscErrorCode SNESFunction(SNES snes,Vec V,Vec F,void *actx) 146c4762a1bSJed Brown { 147c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)actx; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBeginUser; 150*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 151*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV,V,ctx->UV,INSERT_VALUES,SCATTER_FORWARD)); 152*9566063dSJacob Faibussowitsch PetscCall((*ctx->F)(ctx->t,ctx->UV,F)); 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156303a5415SBarry Smith /*TEST 157303a5415SBarry Smith 158303a5415SBarry Smith test: 159303a5415SBarry Smith args: -ts_monitor -ts_view 160303a5415SBarry Smith 161303a5415SBarry Smith TEST*/ 162