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 */ 16d71ae5a4SJacob Faibussowitsch PetscErrorCode f(PetscReal t, Vec UV, Vec F) 17d71ae5a4SJacob Faibussowitsch { 18c4762a1bSJed Brown const PetscScalar *u, *v; 19c4762a1bSJed Brown PetscScalar *f; 20c4762a1bSJed Brown PetscInt n, i; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(UV, &n)); 24c4762a1bSJed Brown n = n / 2; 259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(UV, &u)); 26c4762a1bSJed Brown v = u + n; 279566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &f)); 28c4762a1bSJed Brown for (i = 0; i < n; i++) f[i] = u[i] + v[i]; 299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(UV, &u)); 309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &f)); 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown F(U,V) = U - V 36c4762a1bSJed Brown */ 37d71ae5a4SJacob Faibussowitsch PetscErrorCode F(PetscReal t, Vec UV, Vec F) 38d71ae5a4SJacob Faibussowitsch { 39c4762a1bSJed Brown const PetscScalar *u, *v; 40c4762a1bSJed Brown PetscScalar *f; 41c4762a1bSJed Brown PetscInt n, i; 42c4762a1bSJed Brown 43c4762a1bSJed Brown PetscFunctionBeginUser; 449566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(UV, &n)); 45c4762a1bSJed Brown n = n / 2; 469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(UV, &u)); 47c4762a1bSJed Brown v = u + n; 489566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &f)); 49c4762a1bSJed Brown for (i = 0; i < n; i++) f[i] = u[i] - v[i]; 509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(UV, &u)); 519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &f)); 523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 67d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 68d71ae5a4SJacob Faibussowitsch { 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 76327415f7SBarry Smith PetscFunctionBeginUser; 77*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 799566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 809566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 819566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSEULER)); 829566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 839566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD, 1, PETSC_DETERMINE, &tsrhs)); 849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tsrhs, &U)); 859566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, tsrhs, TSFunction, &ctx)); 869566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0)); 87c4762a1bSJed Brown ctx.f = f; 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &ctx.snes)); 909566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(ctx.snes)); 919566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(ctx.snes, NULL, SNESFunction, &ctx)); 929566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(ctx.snes, NULL, NULL, SNESComputeJacobianDefault, &ctx)); 93c4762a1bSJed Brown ctx.F = F; 949566063dSJacob Faibussowitsch PetscCall(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 */ 979566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD, 2, PETSC_DETERMINE, &ctx.UV)); 98303a5415SBarry Smith i = 2 * rank; 999566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 1, &i, PETSC_COPY_VALUES, &is)); 1009566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(U, NULL, ctx.UV, is, &ctx.scatterU)); 1019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 102303a5415SBarry Smith i = 2 * rank + 1; 1039566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 1, &i, PETSC_COPY_VALUES, &is)); 1049566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ctx.V, NULL, ctx.UV, is, &ctx.scatterV)); 1059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(VecSet(U, 1.0)); 1089566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.V)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.UV)); 1129566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx.scatterU)); 1139566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx.scatterV)); 1149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tsrhs)); 1159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1169566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&ctx.snes)); 1179566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 119b122ec5aSJacob Faibussowitsch return 0; 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 */ 127d71ae5a4SJacob Faibussowitsch PetscErrorCode TSFunction(TS ts, PetscReal t, Vec U, Vec F, void *actx) 128d71ae5a4SJacob Faibussowitsch { 129c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)actx; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBeginUser; 132c4762a1bSJed Brown ctx->t = t; 1339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterU, U, ctx->UV, INSERT_VALUES, SCATTER_FORWARD)); 1349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterU, U, ctx->UV, INSERT_VALUES, SCATTER_FORWARD)); 1359566063dSJacob Faibussowitsch PetscCall(SNESSolve(ctx->snes, NULL, ctx->V)); 1369566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV, ctx->V, ctx->UV, INSERT_VALUES, SCATTER_FORWARD)); 1379566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV, ctx->V, ctx->UV, INSERT_VALUES, SCATTER_FORWARD)); 1389566063dSJacob Faibussowitsch PetscCall((*ctx->f)(t, ctx->UV, F)); 1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* 143c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 144c4762a1bSJed Brown 145c4762a1bSJed Brown */ 146d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFunction(SNES snes, Vec V, Vec F, void *actx) 147d71ae5a4SJacob Faibussowitsch { 148c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)actx; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 1519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scatterV, V, ctx->UV, INSERT_VALUES, SCATTER_FORWARD)); 1529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scatterV, V, ctx->UV, INSERT_VALUES, SCATTER_FORWARD)); 1539566063dSJacob Faibussowitsch PetscCall((*ctx->F)(ctx->t, ctx->UV, F)); 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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