1*fa9584fbSIlya Fursov static char help[] = "Solves a time-dependent linear PDE with discontinuous right hand side.\n"; 2*fa9584fbSIlya Fursov 3*fa9584fbSIlya Fursov /* ------------------------------------------------------------------------ 4*fa9584fbSIlya Fursov 5*fa9584fbSIlya Fursov This program solves the one-dimensional quench front problem modeling a cooled 6*fa9584fbSIlya Fursov liquid rising on a hot metal rod 7*fa9584fbSIlya Fursov u_t = u_xx + g(u), 8*fa9584fbSIlya Fursov with 9*fa9584fbSIlya Fursov g(u) = -Au if u <= u_c, 10*fa9584fbSIlya Fursov = 0 if u > u_c 11*fa9584fbSIlya Fursov on the domain 0 <= x <= 1, with the boundary conditions 12*fa9584fbSIlya Fursov u(t,0) = 0, u_x(t,1) = 0, 13*fa9584fbSIlya Fursov and the initial condition 14*fa9584fbSIlya Fursov u(0,x) = 0 if 0 <= x <= 0.1, 15*fa9584fbSIlya Fursov = (x - 0.1)/0.15 if 0.1 < x < 0.25 16*fa9584fbSIlya Fursov = 1 if 0.25 <= x <= 1 17*fa9584fbSIlya Fursov We discretize the right-hand side using finite differences with 18*fa9584fbSIlya Fursov uniform grid spacing h: 19*fa9584fbSIlya Fursov u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) 20*fa9584fbSIlya Fursov 21*fa9584fbSIlya Fursov Reference: L. Shampine and S. Thompson, "Event Location for Ordinary Differential Equations", 22*fa9584fbSIlya Fursov http://www.radford.edu/~thompson/webddes/eventsweb.pdf 23*fa9584fbSIlya Fursov ------------------------------------------------------------------------- */ 24*fa9584fbSIlya Fursov 25*fa9584fbSIlya Fursov #include <petscdmda.h> 26*fa9584fbSIlya Fursov #include <petscts.h> 27*fa9584fbSIlya Fursov /* 28*fa9584fbSIlya Fursov User-defined application context - contains data needed by the 29*fa9584fbSIlya Fursov application-provided call-back routines. 30*fa9584fbSIlya Fursov */ 31*fa9584fbSIlya Fursov typedef struct { 32*fa9584fbSIlya Fursov PetscReal A; 33*fa9584fbSIlya Fursov PetscReal uc; 34*fa9584fbSIlya Fursov PetscInt *sw; 35*fa9584fbSIlya Fursov } AppCtx; 36*fa9584fbSIlya Fursov 37*fa9584fbSIlya Fursov PetscErrorCode InitialConditions(Vec U, DM da, AppCtx *app) 38*fa9584fbSIlya Fursov { 39*fa9584fbSIlya Fursov Vec xcoord; 40*fa9584fbSIlya Fursov PetscScalar *x, *u; 41*fa9584fbSIlya Fursov PetscInt lsize, M, xs, xm, i; 42*fa9584fbSIlya Fursov 43*fa9584fbSIlya Fursov PetscFunctionBeginUser; 44*fa9584fbSIlya Fursov PetscCall(DMGetCoordinates(da, &xcoord)); 45*fa9584fbSIlya Fursov PetscCall(DMDAVecGetArrayRead(da, xcoord, &x)); 46*fa9584fbSIlya Fursov 47*fa9584fbSIlya Fursov PetscCall(VecGetLocalSize(U, &lsize)); 48*fa9584fbSIlya Fursov PetscCall(PetscMalloc1(lsize, &app->sw)); 49*fa9584fbSIlya Fursov 50*fa9584fbSIlya Fursov PetscCall(DMDAVecGetArray(da, U, &u)); 51*fa9584fbSIlya Fursov 52*fa9584fbSIlya Fursov PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 53*fa9584fbSIlya Fursov PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 54*fa9584fbSIlya Fursov 55*fa9584fbSIlya Fursov for (i = xs; i < xs + xm; i++) { 56*fa9584fbSIlya Fursov if (x[i] <= 0.1) u[i] = 0.; 57*fa9584fbSIlya Fursov else if (x[i] > 0.1 && x[i] < 0.25) u[i] = (x[i] - 0.1) / 0.15; 58*fa9584fbSIlya Fursov else u[i] = 1.0; 59*fa9584fbSIlya Fursov 60*fa9584fbSIlya Fursov app->sw[i - xs] = 1; 61*fa9584fbSIlya Fursov } 62*fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArray(da, U, &u)); 63*fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArrayRead(da, xcoord, &x)); 64*fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 65*fa9584fbSIlya Fursov } 66*fa9584fbSIlya Fursov 67*fa9584fbSIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, void *ctx) 68*fa9584fbSIlya Fursov { 69*fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx; 70*fa9584fbSIlya Fursov const PetscScalar *u; 71*fa9584fbSIlya Fursov PetscInt i, lsize; 72*fa9584fbSIlya Fursov 73*fa9584fbSIlya Fursov PetscFunctionBeginUser; 74*fa9584fbSIlya Fursov PetscCall(VecGetLocalSize(U, &lsize)); 75*fa9584fbSIlya Fursov PetscCall(VecGetArrayRead(U, &u)); 76*fa9584fbSIlya Fursov for (i = 0; i < lsize; i++) fvalue[i] = PetscRealPart(u[i]) - app->uc; 77*fa9584fbSIlya Fursov PetscCall(VecRestoreArrayRead(U, &u)); 78*fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 79*fa9584fbSIlya Fursov } 80*fa9584fbSIlya Fursov 81*fa9584fbSIlya Fursov PetscErrorCode PostEventFunction(TS ts, PetscInt nevents_zero, PetscInt events_zero[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) 82*fa9584fbSIlya Fursov { 83*fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx; 84*fa9584fbSIlya Fursov PetscInt i, idx; 85*fa9584fbSIlya Fursov 86*fa9584fbSIlya Fursov PetscFunctionBeginUser; 87*fa9584fbSIlya Fursov for (i = 0; i < nevents_zero; i++) { 88*fa9584fbSIlya Fursov idx = events_zero[i]; 89*fa9584fbSIlya Fursov app->sw[idx] = 0; 90*fa9584fbSIlya Fursov } 91*fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 92*fa9584fbSIlya Fursov } 93*fa9584fbSIlya Fursov 94*fa9584fbSIlya Fursov /* 95*fa9584fbSIlya Fursov Defines the ODE passed to the ODE solver 96*fa9584fbSIlya Fursov */ 97*fa9584fbSIlya Fursov static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 98*fa9584fbSIlya Fursov { 99*fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx; 100*fa9584fbSIlya Fursov PetscScalar *f; 101*fa9584fbSIlya Fursov const PetscScalar *u, *udot; 102*fa9584fbSIlya Fursov DM da; 103*fa9584fbSIlya Fursov PetscInt M, xs, xm, i; 104*fa9584fbSIlya Fursov PetscReal h, h2; 105*fa9584fbSIlya Fursov Vec Ulocal; 106*fa9584fbSIlya Fursov 107*fa9584fbSIlya Fursov PetscFunctionBeginUser; 108*fa9584fbSIlya Fursov PetscCall(TSGetDM(ts, &da)); 109*fa9584fbSIlya Fursov 110*fa9584fbSIlya Fursov PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 111*fa9584fbSIlya Fursov PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 112*fa9584fbSIlya Fursov 113*fa9584fbSIlya Fursov PetscCall(DMGetLocalVector(da, &Ulocal)); 114*fa9584fbSIlya Fursov PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, Ulocal)); 115*fa9584fbSIlya Fursov PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, Ulocal)); 116*fa9584fbSIlya Fursov 117*fa9584fbSIlya Fursov h = 1.0 / (M - 1); 118*fa9584fbSIlya Fursov h2 = h * h; 119*fa9584fbSIlya Fursov PetscCall(DMDAVecGetArrayRead(da, Udot, &udot)); 120*fa9584fbSIlya Fursov PetscCall(DMDAVecGetArrayRead(da, Ulocal, &u)); 121*fa9584fbSIlya Fursov PetscCall(DMDAVecGetArray(da, F, &f)); 122*fa9584fbSIlya Fursov 123*fa9584fbSIlya Fursov for (i = xs; i < xs + xm; i++) { 124*fa9584fbSIlya Fursov if (i == 0) { 125*fa9584fbSIlya Fursov f[i] = u[i]; 126*fa9584fbSIlya Fursov } else if (i == M - 1) { 127*fa9584fbSIlya Fursov f[i] = (u[i] - u[i - 1]) / h; 128*fa9584fbSIlya Fursov } else { 129*fa9584fbSIlya Fursov f[i] = (u[i + 1] - 2 * u[i] + u[i - 1]) / h2 + app->sw[i - xs] * (-app->A * u[i]) - udot[i]; 130*fa9584fbSIlya Fursov } 131*fa9584fbSIlya Fursov } 132*fa9584fbSIlya Fursov 133*fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot)); 134*fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArrayRead(da, Ulocal, &u)); 135*fa9584fbSIlya Fursov PetscCall(DMDAVecRestoreArray(da, F, &f)); 136*fa9584fbSIlya Fursov PetscCall(DMRestoreLocalVector(da, &Ulocal)); 137*fa9584fbSIlya Fursov 138*fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 139*fa9584fbSIlya Fursov } 140*fa9584fbSIlya Fursov 141*fa9584fbSIlya Fursov /* 142*fa9584fbSIlya Fursov Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 143*fa9584fbSIlya Fursov */ 144*fa9584fbSIlya Fursov static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) 145*fa9584fbSIlya Fursov { 146*fa9584fbSIlya Fursov AppCtx *app = (AppCtx *)ctx; 147*fa9584fbSIlya Fursov DM da; 148*fa9584fbSIlya Fursov MatStencil row, col[3]; 149*fa9584fbSIlya Fursov PetscScalar v[3]; 150*fa9584fbSIlya Fursov PetscInt M, xs, xm, i; 151*fa9584fbSIlya Fursov PetscReal h, h2; 152*fa9584fbSIlya Fursov 153*fa9584fbSIlya Fursov PetscFunctionBeginUser; 154*fa9584fbSIlya Fursov PetscCall(TSGetDM(ts, &da)); 155*fa9584fbSIlya Fursov 156*fa9584fbSIlya Fursov PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 157*fa9584fbSIlya Fursov PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 158*fa9584fbSIlya Fursov 159*fa9584fbSIlya Fursov h = 1.0 / (M - 1); 160*fa9584fbSIlya Fursov h2 = h * h; 161*fa9584fbSIlya Fursov for (i = xs; i < xs + xm; i++) { 162*fa9584fbSIlya Fursov row.i = i; 163*fa9584fbSIlya Fursov if (i == 0) { 164*fa9584fbSIlya Fursov v[0] = 1.0; 165*fa9584fbSIlya Fursov PetscCall(MatSetValuesStencil(A, 1, &row, 1, &row, v, INSERT_VALUES)); 166*fa9584fbSIlya Fursov } else if (i == M - 1) { 167*fa9584fbSIlya Fursov col[0].i = i; 168*fa9584fbSIlya Fursov v[0] = 1 / h; 169*fa9584fbSIlya Fursov col[1].i = i - 1; 170*fa9584fbSIlya Fursov v[1] = -1 / h; 171*fa9584fbSIlya Fursov PetscCall(MatSetValuesStencil(A, 1, &row, 2, col, v, INSERT_VALUES)); 172*fa9584fbSIlya Fursov } else { 173*fa9584fbSIlya Fursov col[0].i = i + 1; 174*fa9584fbSIlya Fursov v[0] = 1 / h2; 175*fa9584fbSIlya Fursov col[1].i = i; 176*fa9584fbSIlya Fursov v[1] = -2 / h2 + app->sw[i - xs] * (-app->A) - a; 177*fa9584fbSIlya Fursov col[2].i = i - 1; 178*fa9584fbSIlya Fursov v[2] = 1 / h2; 179*fa9584fbSIlya Fursov PetscCall(MatSetValuesStencil(A, 1, &row, 3, col, v, INSERT_VALUES)); 180*fa9584fbSIlya Fursov } 181*fa9584fbSIlya Fursov } 182*fa9584fbSIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 183*fa9584fbSIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 184*fa9584fbSIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 185*fa9584fbSIlya Fursov } 186*fa9584fbSIlya Fursov 187*fa9584fbSIlya Fursov int main(int argc, char **argv) 188*fa9584fbSIlya Fursov { 189*fa9584fbSIlya Fursov TS ts; /* ODE integrator */ 190*fa9584fbSIlya Fursov Vec U; /* solution will be stored here */ 191*fa9584fbSIlya Fursov Mat J; /* Jacobian matrix */ 192*fa9584fbSIlya Fursov PetscInt n = 16; 193*fa9584fbSIlya Fursov AppCtx app; 194*fa9584fbSIlya Fursov DM da; 195*fa9584fbSIlya Fursov 196*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197*fa9584fbSIlya Fursov Initialize program 198*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 199*fa9584fbSIlya Fursov PetscFunctionBeginUser; 200*fa9584fbSIlya Fursov PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 201*fa9584fbSIlya Fursov 202*fa9584fbSIlya Fursov PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex22 options", ""); 203*fa9584fbSIlya Fursov { 204*fa9584fbSIlya Fursov app.A = 200000; 205*fa9584fbSIlya Fursov PetscCall(PetscOptionsReal("-A", "", "", app.A, &app.A, NULL)); 206*fa9584fbSIlya Fursov app.uc = 0.5; 207*fa9584fbSIlya Fursov PetscCall(PetscOptionsReal("-uc", "", "", app.uc, &app.uc, NULL)); 208*fa9584fbSIlya Fursov } 209*fa9584fbSIlya Fursov PetscOptionsEnd(); 210*fa9584fbSIlya Fursov 211*fa9584fbSIlya Fursov PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, -n, 1, 1, 0, &da)); 212*fa9584fbSIlya Fursov PetscCall(DMSetFromOptions(da)); 213*fa9584fbSIlya Fursov PetscCall(DMSetUp(da)); 214*fa9584fbSIlya Fursov PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0, 0, 0, 0)); 215*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 216*fa9584fbSIlya Fursov Create necessary matrix and vectors 217*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 218*fa9584fbSIlya Fursov PetscCall(DMCreateMatrix(da, &J)); 219*fa9584fbSIlya Fursov PetscCall(DMCreateGlobalVector(da, &U)); 220*fa9584fbSIlya Fursov 221*fa9584fbSIlya Fursov PetscCall(InitialConditions(U, da, &app)); 222*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 223*fa9584fbSIlya Fursov Create timestepping solver context 224*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 225*fa9584fbSIlya Fursov PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 226*fa9584fbSIlya Fursov PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 227*fa9584fbSIlya Fursov PetscCall(TSSetType(ts, TSROSW)); 228*fa9584fbSIlya Fursov PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, (void *)&app)); 229*fa9584fbSIlya Fursov PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobian)IJacobian, (void *)&app)); 230*fa9584fbSIlya Fursov 231*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 232*fa9584fbSIlya Fursov Set initial conditions 233*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 234*fa9584fbSIlya Fursov PetscCall(TSSetSolution(ts, U)); 235*fa9584fbSIlya Fursov 236*fa9584fbSIlya Fursov PetscCall(TSSetDM(ts, da)); 237*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 238*fa9584fbSIlya Fursov Set solver options 239*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 240*fa9584fbSIlya Fursov PetscCall(TSSetTimeStep(ts, 0.1)); 241*fa9584fbSIlya Fursov PetscCall(TSSetMaxTime(ts, 30.0)); 242*fa9584fbSIlya Fursov PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 243*fa9584fbSIlya Fursov PetscCall(TSSetFromOptions(ts)); 244*fa9584fbSIlya Fursov 245*fa9584fbSIlya Fursov PetscInt lsize; 246*fa9584fbSIlya Fursov PetscCall(VecGetLocalSize(U, &lsize)); 247*fa9584fbSIlya Fursov PetscInt *direction; 248*fa9584fbSIlya Fursov PetscBool *terminate; 249*fa9584fbSIlya Fursov PetscInt i; 250*fa9584fbSIlya Fursov PetscCall(PetscMalloc1(lsize, &direction)); 251*fa9584fbSIlya Fursov PetscCall(PetscMalloc1(lsize, &terminate)); 252*fa9584fbSIlya Fursov for (i = 0; i < lsize; i++) { 253*fa9584fbSIlya Fursov direction[i] = -1; 254*fa9584fbSIlya Fursov terminate[i] = PETSC_FALSE; 255*fa9584fbSIlya Fursov } 256*fa9584fbSIlya Fursov PetscCall(TSSetEventHandler(ts, lsize, direction, terminate, EventFunction, PostEventFunction, (void *)&app)); 257*fa9584fbSIlya Fursov 258*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 259*fa9584fbSIlya Fursov Run timestepping solver 260*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 261*fa9584fbSIlya Fursov PetscCall(TSSolve(ts, U)); 262*fa9584fbSIlya Fursov 263*fa9584fbSIlya Fursov /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264*fa9584fbSIlya Fursov Free work space. All PETSc objects should be destroyed when they are no longer needed. 265*fa9584fbSIlya Fursov - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 266*fa9584fbSIlya Fursov 267*fa9584fbSIlya Fursov PetscCall(MatDestroy(&J)); 268*fa9584fbSIlya Fursov PetscCall(VecDestroy(&U)); 269*fa9584fbSIlya Fursov PetscCall(DMDestroy(&da)); 270*fa9584fbSIlya Fursov PetscCall(TSDestroy(&ts)); 271*fa9584fbSIlya Fursov PetscCall(PetscFree(direction)); 272*fa9584fbSIlya Fursov PetscCall(PetscFree(terminate)); 273*fa9584fbSIlya Fursov 274*fa9584fbSIlya Fursov PetscCall(PetscFree(app.sw)); 275*fa9584fbSIlya Fursov PetscCall(PetscFinalize()); 276*fa9584fbSIlya Fursov return 0; 277*fa9584fbSIlya Fursov } 278