1*ca4445c7SIlya Fursov #include <petscts.h> 2*ca4445c7SIlya Fursov #include <stdio.h> 3*ca4445c7SIlya Fursov 4*ca4445c7SIlya Fursov #define NEW_VERSION // Applicable for the new features; avoid this for the old (current) TSEvent code 5*ca4445c7SIlya Fursov 6*ca4445c7SIlya Fursov static char help[] = "Simple linear problem with events\n" 7*ca4445c7SIlya Fursov "x_dot = 0.2*y\n" 8*ca4445c7SIlya Fursov "y_dot = -0.2*x\n" 9*ca4445c7SIlya Fursov "Using one or several event functions (on rank-0)\n" 10*ca4445c7SIlya Fursov "This program is mostly intended to test the Anderson-Bjorck iteration\n" 11*ca4445c7SIlya Fursov "Options:\n" 12*ca4445c7SIlya Fursov "-dir d : zero-crossing direction for events\n" 13*ca4445c7SIlya Fursov "-flg : additional output in Postevent\n" 14*ca4445c7SIlya Fursov "-restart : flag for TSRestartStep() in PostEvent\n" 15*ca4445c7SIlya Fursov "-dtpost x : if x > 0, then on even PostEvent calls dt_postevent = x is set, on odd PostEvent calls dt_postevent = 0 is set,\n" 16*ca4445c7SIlya Fursov " if x == 0, nothing happens\n" 17*ca4445c7SIlya Fursov "-func F : selects the event function [0, ..., 11], if F == -1 (default) is set, all event functions are taken\n"; 18*ca4445c7SIlya Fursov 19*ca4445c7SIlya Fursov #define MAX_NFUNC 100 // max event functions per rank 20*ca4445c7SIlya Fursov #define MAX_NEV 5000 // max zero crossings for each rank 21*ca4445c7SIlya Fursov 22*ca4445c7SIlya Fursov typedef struct { 23*ca4445c7SIlya Fursov PetscMPIInt rank, size; 24*ca4445c7SIlya Fursov PetscReal pi; 25*ca4445c7SIlya Fursov PetscReal fvals[MAX_NFUNC]; // helper array for reporting the residuals 26*ca4445c7SIlya Fursov PetscReal evres[MAX_NEV]; // times of found zero-crossings 27*ca4445c7SIlya Fursov PetscInt evnum[MAX_NEV]; // number of zero-crossings at each time 28*ca4445c7SIlya Fursov PetscInt cnt; // counter 29*ca4445c7SIlya Fursov PetscBool flg; // flag for additional print in PostEvent 30*ca4445c7SIlya Fursov PetscBool restart; // flag for TSRestartStep() in PostEvent 31*ca4445c7SIlya Fursov PetscReal dtpost; // post-event step 32*ca4445c7SIlya Fursov PetscInt postcnt; // counter for PostEvent calls 33*ca4445c7SIlya Fursov PetscInt F; // event-function index 34*ca4445c7SIlya Fursov PetscInt Fnum; // total available event functions 35*ca4445c7SIlya Fursov } AppCtx; 36*ca4445c7SIlya Fursov 37*ca4445c7SIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal gval[], void *ctx); 38*ca4445c7SIlya Fursov PetscErrorCode Postevent(TS ts, PetscInt nev_zero, PetscInt evs_zero[], PetscReal t, Vec U, PetscBool fwd, void *ctx); 39*ca4445c7SIlya Fursov 40*ca4445c7SIlya Fursov int main(int argc, char **argv) 41*ca4445c7SIlya Fursov { 42*ca4445c7SIlya Fursov TS ts; 43*ca4445c7SIlya Fursov Mat A; 44*ca4445c7SIlya Fursov Vec sol; 45*ca4445c7SIlya Fursov PetscInt n, dir0, m = 0; 46*ca4445c7SIlya Fursov PetscInt dir[MAX_NFUNC], inds[2]; 47*ca4445c7SIlya Fursov PetscBool term[MAX_NFUNC]; 48*ca4445c7SIlya Fursov PetscScalar *x, vals[4]; 49*ca4445c7SIlya Fursov AppCtx ctx; 50*ca4445c7SIlya Fursov 51*ca4445c7SIlya Fursov PetscFunctionBeginUser; 52*ca4445c7SIlya Fursov PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 53*ca4445c7SIlya Fursov setbuf(stdout, NULL); 54*ca4445c7SIlya Fursov PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &ctx.rank)); 55*ca4445c7SIlya Fursov PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &ctx.size)); 56*ca4445c7SIlya Fursov ctx.pi = PetscAcosReal(-1.0); 57*ca4445c7SIlya Fursov ctx.cnt = 0; 58*ca4445c7SIlya Fursov ctx.flg = PETSC_FALSE; 59*ca4445c7SIlya Fursov ctx.restart = PETSC_FALSE; 60*ca4445c7SIlya Fursov ctx.dtpost = 0; 61*ca4445c7SIlya Fursov ctx.postcnt = 0; 62*ca4445c7SIlya Fursov ctx.F = -1; 63*ca4445c7SIlya Fursov ctx.Fnum = 12; 64*ca4445c7SIlya Fursov 65*ca4445c7SIlya Fursov // The linear problem has a 2*2 matrix. The matrix is constant 66*ca4445c7SIlya Fursov if (ctx.rank == 0) m = 2; 67*ca4445c7SIlya Fursov inds[0] = 0; 68*ca4445c7SIlya Fursov inds[1] = 1; 69*ca4445c7SIlya Fursov vals[0] = 0; 70*ca4445c7SIlya Fursov vals[1] = 0.2; 71*ca4445c7SIlya Fursov vals[2] = -0.2; 72*ca4445c7SIlya Fursov vals[3] = 0; 73*ca4445c7SIlya Fursov PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, m, m, PETSC_DETERMINE, PETSC_DETERMINE, 2, NULL, 0, NULL, &A)); 74*ca4445c7SIlya Fursov PetscCall(MatSetValues(A, m, inds, m, inds, vals, INSERT_VALUES)); 75*ca4445c7SIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 76*ca4445c7SIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 77*ca4445c7SIlya Fursov PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 78*ca4445c7SIlya Fursov 79*ca4445c7SIlya Fursov PetscCall(MatCreateVecs(A, &sol, NULL)); 80*ca4445c7SIlya Fursov PetscCall(VecGetArray(sol, &x)); 81*ca4445c7SIlya Fursov if (ctx.rank == 0) { // initial conditions 82*ca4445c7SIlya Fursov x[0] = 0; // sin(0) 83*ca4445c7SIlya Fursov x[1] = 1; // cos(0) 84*ca4445c7SIlya Fursov } 85*ca4445c7SIlya Fursov PetscCall(VecRestoreArray(sol, &x)); 86*ca4445c7SIlya Fursov 87*ca4445c7SIlya Fursov PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 88*ca4445c7SIlya Fursov PetscCall(TSSetProblemType(ts, TS_LINEAR)); 89*ca4445c7SIlya Fursov 90*ca4445c7SIlya Fursov PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, NULL)); 91*ca4445c7SIlya Fursov PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, NULL)); 92*ca4445c7SIlya Fursov 93*ca4445c7SIlya Fursov PetscCall(TSSetTime(ts, 0.03)); 94*ca4445c7SIlya Fursov PetscCall(TSSetTimeStep(ts, 0.1)); 95*ca4445c7SIlya Fursov PetscCall(TSSetType(ts, TSBEULER)); 96*ca4445c7SIlya Fursov PetscCall(TSSetMaxSteps(ts, 10000)); 97*ca4445c7SIlya Fursov PetscCall(TSSetMaxTime(ts, 4.0)); 98*ca4445c7SIlya Fursov PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 99*ca4445c7SIlya Fursov PetscCall(TSSetFromOptions(ts)); 100*ca4445c7SIlya Fursov 101*ca4445c7SIlya Fursov // Set the event handling 102*ca4445c7SIlya Fursov dir0 = 0; 103*ca4445c7SIlya Fursov PetscCall(PetscOptionsGetInt(NULL, NULL, "-dir", &dir0, NULL)); // desired zero-crossing direction 104*ca4445c7SIlya Fursov PetscCall(PetscOptionsHasName(NULL, NULL, "-flg", &ctx.flg)); // flag for additional output 105*ca4445c7SIlya Fursov PetscCall(PetscOptionsGetBool(NULL, NULL, "-restart", &ctx.restart, NULL)); // flag for TSRestartStep() 106*ca4445c7SIlya Fursov PetscCall(PetscOptionsGetReal(NULL, NULL, "-dtpost", &ctx.dtpost, NULL)); // post-event step 107*ca4445c7SIlya Fursov PetscCall(PetscOptionsGetInt(NULL, NULL, "-F", &ctx.F, NULL)); // event-function index 108*ca4445c7SIlya Fursov PetscCheck(ctx.F >= -1 && ctx.F < ctx.Fnum, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Value of 'F' is out of range"); 109*ca4445c7SIlya Fursov 110*ca4445c7SIlya Fursov n = 0; // event counter 111*ca4445c7SIlya Fursov if (ctx.rank == 0) { // all events -- on rank-0 112*ca4445c7SIlya Fursov if (ctx.F == -1) 113*ca4445c7SIlya Fursov for (n = 0; n < ctx.Fnum; n++) { // all event-functions 114*ca4445c7SIlya Fursov dir[n] = dir0; 115*ca4445c7SIlya Fursov term[n] = PETSC_FALSE; 116*ca4445c7SIlya Fursov } 117*ca4445c7SIlya Fursov else { // single event-function 118*ca4445c7SIlya Fursov dir[n] = dir0; 119*ca4445c7SIlya Fursov term[n++] = PETSC_FALSE; 120*ca4445c7SIlya Fursov } 121*ca4445c7SIlya Fursov } 122*ca4445c7SIlya Fursov PetscCall(TSSetEventHandler(ts, n, dir, term, EventFunction, Postevent, &ctx)); 123*ca4445c7SIlya Fursov 124*ca4445c7SIlya Fursov // Solution 125*ca4445c7SIlya Fursov PetscCall(TSSolve(ts, sol)); 126*ca4445c7SIlya Fursov 127*ca4445c7SIlya Fursov // The 3 columns printed are: [RANK] [num. of events at the given time] [time of event] 128*ca4445c7SIlya Fursov for (PetscInt j = 0; j < ctx.cnt; j++) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%d\t%" PetscInt_FMT "\t%.5g\n", ctx.rank, ctx.evnum[j], (double)ctx.evres[j])); 129*ca4445c7SIlya Fursov PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 130*ca4445c7SIlya Fursov 131*ca4445c7SIlya Fursov PetscCall(MatDestroy(&A)); 132*ca4445c7SIlya Fursov PetscCall(TSDestroy(&ts)); 133*ca4445c7SIlya Fursov PetscCall(VecDestroy(&sol)); 134*ca4445c7SIlya Fursov 135*ca4445c7SIlya Fursov PetscCall(PetscFinalize()); 136*ca4445c7SIlya Fursov return 0; 137*ca4445c7SIlya Fursov } 138*ca4445c7SIlya Fursov 139*ca4445c7SIlya Fursov /* 140*ca4445c7SIlya Fursov User callback for defining the event-functions 141*ca4445c7SIlya Fursov */ 142*ca4445c7SIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal gval[], void *ctx) 143*ca4445c7SIlya Fursov { 144*ca4445c7SIlya Fursov PetscInt n = 0; 145*ca4445c7SIlya Fursov AppCtx *Ctx = (AppCtx *)ctx; 146*ca4445c7SIlya Fursov 147*ca4445c7SIlya Fursov PetscFunctionBeginUser; 148*ca4445c7SIlya Fursov // for the test purposes, event-functions are defined based on t 149*ca4445c7SIlya Fursov // all events -- on rank-0 150*ca4445c7SIlya Fursov if (Ctx->rank == 0) { 151*ca4445c7SIlya Fursov if (Ctx->F == 0 || Ctx->F == -1) gval[n++] = PetscSinReal(Ctx->pi * t) / Ctx->pi; // FUNC-0, roots 1, 2, 3, 4 152*ca4445c7SIlya Fursov if (Ctx->F == 1 || Ctx->F == -1) gval[n++] = PetscLogReal(t); // FUNC-2, root 1 153*ca4445c7SIlya Fursov if (Ctx->F == 2 || Ctx->F == -1) { // FUNC-3, root 1 154*ca4445c7SIlya Fursov if (t < 2) gval[n++] = (1 - PetscPowReal(t - 2, 12)) / 12.0; 155*ca4445c7SIlya Fursov else gval[n++] = 1 / 12.0; 156*ca4445c7SIlya Fursov } 157*ca4445c7SIlya Fursov if (Ctx->F == 3 || Ctx->F == -1) gval[n++] = t - PetscExpReal(PetscSinReal(t)) + 1; // FUNC-5, root 1.69681 158*ca4445c7SIlya Fursov if (Ctx->F == 4 || Ctx->F == -1) gval[n++] = (1e10 * PetscPowReal(t, 1 / t) - 1) / 100; // FUNC-6, root 0.1 159*ca4445c7SIlya Fursov if (Ctx->F == 5 || Ctx->F == -1) gval[n++] = PetscLogReal(t - 0.02) * PetscLogReal(t - 0.02) * PetscSignReal(t - 1.02) * 1e7; // FUNC-7, root 1.02 160*ca4445c7SIlya Fursov if (Ctx->F == 6 || Ctx->F == -1) gval[n++] = 4 * PetscCosReal(t) - PetscExpReal(t); // FUNC-14, root 0.904788 161*ca4445c7SIlya Fursov if (Ctx->F == 7 || Ctx->F == -1) gval[n++] = (20.0 * t - 1) / (19.0 * t) / 10; // FUNC-15, root 0.05 162*ca4445c7SIlya Fursov if (Ctx->F == 8 || Ctx->F == -1) gval[n++] = ((t - 1) * PetscExpReal(-20 * t) + PetscPowReal(t, 20)) * 1e4; // FUNC-16, root 0.552 163*ca4445c7SIlya Fursov if (Ctx->F == 9 || Ctx->F == -1) gval[n++] = (t * t * (t * t / 3.0 + PetscSqrtReal(2.0) * PetscSinReal(t)) - PetscSqrtReal(3.0) / 18) * 10; // FUNC-17, root 0.399 164*ca4445c7SIlya Fursov if (Ctx->F == 10 || Ctx->F == -1) gval[n++] = ((t * t + 1) * PetscSinReal(t) - PetscExpReal(PetscSqrtReal(t)) * (t - 1) * (t * t - 5)) / 10; // FUNC-18, roots 0.87, 2.388 165*ca4445c7SIlya Fursov if (Ctx->F == 11 || Ctx->F == -1) gval[n++] = 2 * t - 5; // FUNC-21, root 2.5 166*ca4445c7SIlya Fursov } 167*ca4445c7SIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 168*ca4445c7SIlya Fursov } 169*ca4445c7SIlya Fursov 170*ca4445c7SIlya Fursov /* 171*ca4445c7SIlya Fursov User callback for the post-event stuff 172*ca4445c7SIlya Fursov */ 173*ca4445c7SIlya Fursov PetscErrorCode Postevent(TS ts, PetscInt nev_zero, PetscInt evs_zero[], PetscReal t, Vec U, PetscBool fwd, void *ctx) 174*ca4445c7SIlya Fursov { 175*ca4445c7SIlya Fursov AppCtx *Ctx = (AppCtx *)ctx; 176*ca4445c7SIlya Fursov 177*ca4445c7SIlya Fursov PetscFunctionBeginUser; 178*ca4445c7SIlya Fursov if (Ctx->flg) { 179*ca4445c7SIlya Fursov PetscCallBack("EventFunction", EventFunction(ts, t, U, Ctx->fvals, ctx)); 180*ca4445c7SIlya Fursov PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] At t = %20.16g : %" PetscInt_FMT " events triggered, fvalues =", Ctx->rank, (double)t, nev_zero)); 181*ca4445c7SIlya Fursov for (PetscInt j = 0; j < nev_zero; j++) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\t%g", (double)Ctx->fvals[evs_zero[j]])); 182*ca4445c7SIlya Fursov PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\n")); 183*ca4445c7SIlya Fursov PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 184*ca4445c7SIlya Fursov } 185*ca4445c7SIlya Fursov 186*ca4445c7SIlya Fursov if (Ctx->cnt + nev_zero < MAX_NEV) { 187*ca4445c7SIlya Fursov for (PetscInt i = 0; i < nev_zero; i++) { // save the repeating zeros separately for easier/unified testing 188*ca4445c7SIlya Fursov Ctx->evres[Ctx->cnt] = t; 189*ca4445c7SIlya Fursov Ctx->evnum[Ctx->cnt++] = 1; 190*ca4445c7SIlya Fursov } 191*ca4445c7SIlya Fursov } 192*ca4445c7SIlya Fursov 193*ca4445c7SIlya Fursov #ifdef NEW_VERSION 194*ca4445c7SIlya Fursov Ctx->postcnt++; // sync 195*ca4445c7SIlya Fursov if (Ctx->dtpost > 0) { 196*ca4445c7SIlya Fursov if (Ctx->postcnt % 2 == 0) PetscCall(TSSetPostEventStep(ts, Ctx->dtpost)); 197*ca4445c7SIlya Fursov else PetscCall(TSSetPostEventStep(ts, 0)); 198*ca4445c7SIlya Fursov } 199*ca4445c7SIlya Fursov #endif 200*ca4445c7SIlya Fursov 201*ca4445c7SIlya Fursov if (Ctx->restart) PetscCall(TSRestartStep(ts)); 202*ca4445c7SIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 203*ca4445c7SIlya Fursov } 204*ca4445c7SIlya Fursov /*---------------------------------------------------------------------------------------------*/ 205*ca4445c7SIlya Fursov /*TEST 206*ca4445c7SIlya Fursov test: 207*ca4445c7SIlya Fursov suffix: 0 208*ca4445c7SIlya Fursov output_file: output/ex4_0.out 209*ca4445c7SIlya Fursov args: -dir 0 210*ca4445c7SIlya Fursov args: -ts_adapt_dt_min 1e-10 -ts_event_dt_min 1e-6 211*ca4445c7SIlya Fursov args: -ts_dt 0.4 212*ca4445c7SIlya Fursov args: -restart 1 213*ca4445c7SIlya Fursov args: -ts_event_tol {{1e-8 1e-15}} 214*ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 215*ca4445c7SIlya Fursov args: -dtpost {{0 0.25}} 216*ca4445c7SIlya Fursov args: -ts_event_post_event_step {{0 0.35}} 217*ca4445c7SIlya Fursov args: -ts_type {{beuler rk}} 218*ca4445c7SIlya Fursov nsize: {{1 4}} 219*ca4445c7SIlya Fursov filter: sort 220*ca4445c7SIlya Fursov filter_output: sort 221*ca4445c7SIlya Fursov 222*ca4445c7SIlya Fursov test: 223*ca4445c7SIlya Fursov suffix: F7 224*ca4445c7SIlya Fursov output_file: output/ex4_F7.out 225*ca4445c7SIlya Fursov args: -dir 0 226*ca4445c7SIlya Fursov args: -ts_adapt_dt_min 1e-10 -ts_event_dt_min 1e-6 227*ca4445c7SIlya Fursov args: -ts_dt 0.4 228*ca4445c7SIlya Fursov args: -F 7 229*ca4445c7SIlya Fursov args: -ts_event_tol {{1e-8 1e-15}} 230*ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 231*ca4445c7SIlya Fursov args: -ts_type {{beuler rk}} 232*ca4445c7SIlya Fursov nsize: 1 233*ca4445c7SIlya Fursov 234*ca4445c7SIlya Fursov test: 235*ca4445c7SIlya Fursov suffix: F7revisit 236*ca4445c7SIlya Fursov output_file: output/ex4_F7revisit.out 237*ca4445c7SIlya Fursov args: -ts_event_monitor -F 7 -ts_dt 0.04 -ts_event_dt_min 0.016 238*ca4445c7SIlya Fursov nsize: 1 239*ca4445c7SIlya Fursov 240*ca4445c7SIlya Fursov test: 241*ca4445c7SIlya Fursov suffix: pos 242*ca4445c7SIlya Fursov output_file: output/ex4_pos.out 243*ca4445c7SIlya Fursov args: -dir 1 244*ca4445c7SIlya Fursov args: -ts_adapt_dt_min 1e-10 -ts_event_dt_min 1e-6 245*ca4445c7SIlya Fursov args: -ts_dt 0.4 246*ca4445c7SIlya Fursov args: -restart 0 247*ca4445c7SIlya Fursov args: -ts_event_tol {{1e-8 1e-15}} 248*ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 249*ca4445c7SIlya Fursov args: -dtpost {{0 0.25}} 250*ca4445c7SIlya Fursov args: -ts_event_post_event_step {{0 0.35}} 251*ca4445c7SIlya Fursov args: -ts_type {{beuler rk}} 252*ca4445c7SIlya Fursov nsize: {{1 4}} 253*ca4445c7SIlya Fursov filter: sort 254*ca4445c7SIlya Fursov filter_output: sort 255*ca4445c7SIlya Fursov 256*ca4445c7SIlya Fursov test: 257*ca4445c7SIlya Fursov suffix: neg 258*ca4445c7SIlya Fursov output_file: output/ex4_neg.out 259*ca4445c7SIlya Fursov args: -dir -1 260*ca4445c7SIlya Fursov args: -ts_adapt_dt_min 1e-10 -ts_event_dt_min 1e-6 261*ca4445c7SIlya Fursov args: -ts_dt 0.4 262*ca4445c7SIlya Fursov args: -restart 1 263*ca4445c7SIlya Fursov args: -ts_event_tol {{1e-8 1e-15}} 264*ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 265*ca4445c7SIlya Fursov args: -dtpost {{0 0.25}} 266*ca4445c7SIlya Fursov args: -ts_event_post_event_step {{0 0.35}} 267*ca4445c7SIlya Fursov args: -ts_type {{beuler rk}} 268*ca4445c7SIlya Fursov nsize: {{1 4}} 269*ca4445c7SIlya Fursov filter: sort 270*ca4445c7SIlya Fursov filter_output: sort 271*ca4445c7SIlya Fursov TEST*/ 272