1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Reaction Equation from Chemistry\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown 6c4762a1bSJed Brown Page 6, An example from Atomospheric Chemistry 7c4762a1bSJed Brown 8c4762a1bSJed Brown u_1_t = 9c4762a1bSJed Brown u_2_t = 10c4762a1bSJed Brown u_3_t = 11c4762a1bSJed Brown u_4_t = 12c4762a1bSJed Brown 13c4762a1bSJed Brown -ts_monitor_lg_error -ts_monitor_lg_solution -ts_view -ts_max_time 2.e4 14c4762a1bSJed Brown 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 19c4762a1bSJed Brown file automatically includes: 20c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 21c4762a1bSJed Brown petscmat.h - matrices 22c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 23c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 24c4762a1bSJed Brown petscksp.h - linear solvers 25c4762a1bSJed Brown */ 26c4762a1bSJed Brown 27c4762a1bSJed Brown #include <petscts.h> 28c4762a1bSJed Brown 29c4762a1bSJed Brown typedef struct { 30c4762a1bSJed Brown PetscScalar k1, k2, k3; 31c4762a1bSJed Brown PetscScalar sigma2; 32c4762a1bSJed Brown Vec initialsolution; 33c4762a1bSJed Brown } AppCtx; 34c4762a1bSJed Brown 35*d71ae5a4SJacob Faibussowitsch PetscScalar k1(AppCtx *ctx, PetscReal t) 36*d71ae5a4SJacob Faibussowitsch { 37c4762a1bSJed Brown PetscReal th = t / 3600.0; 38c4762a1bSJed Brown PetscReal barth = th - 24.0 * PetscFloorReal(th / 24.0); 39c4762a1bSJed Brown if (((((PetscInt)th) % 24) < 4) || ((((PetscInt)th) % 24) >= 20)) return (1.0e-40); 40c4762a1bSJed Brown else return (ctx->k1 * PetscExpReal(7.0 * PetscPowReal(PetscSinReal(.0625 * PETSC_PI * (barth - 4.0)), .2))); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43*d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) 44*d71ae5a4SJacob Faibussowitsch { 45c4762a1bSJed Brown PetscScalar *f; 46c4762a1bSJed Brown const PetscScalar *u, *udot; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 519566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &f)); 52c4762a1bSJed Brown f[0] = udot[0] - k1(ctx, t) * u[2] + ctx->k2 * u[0]; 53c4762a1bSJed Brown f[1] = udot[1] - k1(ctx, t) * u[2] + ctx->k3 * u[1] * u[3] - ctx->sigma2; 54c4762a1bSJed Brown f[2] = udot[2] - ctx->k3 * u[1] * u[3] + k1(ctx, t) * u[2]; 55c4762a1bSJed Brown f[3] = udot[3] - ctx->k2 * u[0] + ctx->k3 * u[1] * u[3]; 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &f)); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62*d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) 63*d71ae5a4SJacob Faibussowitsch { 64c4762a1bSJed Brown PetscInt rowcol[] = {0, 1, 2, 3}; 65c4762a1bSJed Brown PetscScalar J[4][4]; 66c4762a1bSJed Brown const PetscScalar *u, *udot; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBegin; 699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 719371c9d4SSatish Balay J[0][0] = a + ctx->k2; 729371c9d4SSatish Balay J[0][1] = 0.0; 739371c9d4SSatish Balay J[0][2] = -k1(ctx, t); 749371c9d4SSatish Balay J[0][3] = 0.0; 759371c9d4SSatish Balay J[1][0] = 0.0; 769371c9d4SSatish Balay J[1][1] = a + ctx->k3 * u[3]; 779371c9d4SSatish Balay J[1][2] = -k1(ctx, t); 789371c9d4SSatish Balay J[1][3] = ctx->k3 * u[1]; 799371c9d4SSatish Balay J[2][0] = 0.0; 809371c9d4SSatish Balay J[2][1] = -ctx->k3 * u[3]; 819371c9d4SSatish Balay J[2][2] = a + k1(ctx, t); 829371c9d4SSatish Balay J[2][3] = -ctx->k3 * u[1]; 839371c9d4SSatish Balay J[3][0] = -ctx->k2; 849371c9d4SSatish Balay J[3][1] = ctx->k3 * u[3]; 859371c9d4SSatish Balay J[3][2] = 0.0; 869371c9d4SSatish Balay J[3][3] = a + ctx->k3 * u[1]; 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 4, rowcol, 4, rowcol, &J[0][0], INSERT_VALUES)); 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 90c4762a1bSJed Brown 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 93c4762a1bSJed Brown if (A != B) { 949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown PetscFunctionReturn(0); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100*d71ae5a4SJacob Faibussowitsch static PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *ctx) 101*d71ae5a4SJacob Faibussowitsch { 102c4762a1bSJed Brown PetscFunctionBegin; 1039566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->initialsolution, U)); 1043c633725SBarry Smith PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Solution not given"); 105c4762a1bSJed Brown PetscFunctionReturn(0); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 109*d71ae5a4SJacob Faibussowitsch { 110c4762a1bSJed Brown TS ts; /* ODE integrator */ 111c4762a1bSJed Brown Vec U; /* solution */ 112c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 113c4762a1bSJed Brown PetscMPIInt size; 114c4762a1bSJed Brown PetscInt n = 4; 115c4762a1bSJed Brown AppCtx ctx; 116c4762a1bSJed Brown PetscScalar *u; 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 119c4762a1bSJed Brown Initialize program 120c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 121327415f7SBarry Smith PetscFunctionBeginUser; 1229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1243c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Create necessary matrix and vectors 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown ctx.k1 = 1.0e-5; 137c4762a1bSJed Brown ctx.k2 = 1.0e5; 138c4762a1bSJed Brown ctx.k3 = 1.0e-16; 139c4762a1bSJed Brown ctx.sigma2 = 1.0e6; 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &ctx.initialsolution)); 1429566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(ctx.initialsolution, &u)); 143c4762a1bSJed Brown u[0] = 0.0; 144c4762a1bSJed Brown u[1] = 1.3e8; 145c4762a1bSJed Brown u[2] = 5.0e11; 146c4762a1bSJed Brown u[3] = 8.0e11; 1479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(ctx.initialsolution, &u)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150c4762a1bSJed Brown Create timestepping solver context 151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1529566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1539566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1549566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW)); 1559566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx)); 1569566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159c4762a1bSJed Brown Set initial conditions 160c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1619566063dSJacob Faibussowitsch PetscCall(Solution(ts, 0, U, &ctx)); 1629566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 4.0 * 3600)); 1639566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1.0)); 1649566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 167c4762a1bSJed Brown Set solver options 168c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1699566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 518400.0)); 1709566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1719566063dSJacob Faibussowitsch PetscCall(TSSetMaxStepRejections(ts, 100)); 1729566063dSJacob Faibussowitsch PetscCall(TSSetMaxSNESFailures(ts, -1)); /* unlimited */ 1739566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176c4762a1bSJed Brown Solve nonlinear system 177c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1789566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 182c4762a1bSJed Brown are no longer needed. 183c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.initialsolution)); 1859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1879566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 190b122ec5aSJacob Faibussowitsch return 0; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown /*TEST 194c4762a1bSJed Brown 195c4762a1bSJed Brown test: 196c4762a1bSJed Brown args: -ts_view -ts_max_time 2.e4 197c4762a1bSJed Brown timeoutfactor: 15 198c4762a1bSJed Brown requires: !single 199c4762a1bSJed Brown 200c4762a1bSJed Brown TEST*/ 201