1c4762a1bSJed Brown static const char help[] = "An elastic wave equation driven by Dieterich-Ruina friction\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown This whole derivation comes from Erickson, Birnir, and Lavallee [2010]. The model comes from the continuum limit in Carlson and Langer [1989], 4c4762a1bSJed Brown 5c4762a1bSJed Brown u_{tt} = c^2 u_{xx} - \tilde\gamma^2 u - (\gamma^2 / \xi) (\theta + \ln(u_t + 1)) 6c4762a1bSJed Brown \theta_t = -(u_t + 1) (\theta + (1 + \epsilon) \ln(u_t +1)) 7c4762a1bSJed Brown 8c4762a1bSJed Brown which can be reduced to a first order system, 9c4762a1bSJed Brown 10c4762a1bSJed Brown u_t = v 11c4762a1bSJed Brown v_t = c^2 u_{xx} - \tilde\gamma^2 u - (\gamma^2 / \xi)(\theta + ln(v + 1))) 12c4762a1bSJed Brown \theta_t = -(v + 1) (\theta + (1 + \epsilon) \ln(v+1)) 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown #include <petscts.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown PetscScalar u, v, th; 21c4762a1bSJed Brown } Field; 22c4762a1bSJed Brown 23c4762a1bSJed Brown typedef struct _User *User; 24c4762a1bSJed Brown struct _User { 25c4762a1bSJed Brown PetscReal epsilon; /* inverse of seismic ratio, B-A / A */ 26c4762a1bSJed Brown PetscReal gamma; /* wave frequency for interblock coupling */ 27c4762a1bSJed Brown PetscReal gammaTilde; /* wave frequency for coupling to plate */ 28c4762a1bSJed Brown PetscReal xi; /* interblock spring constant */ 29c4762a1bSJed Brown PetscReal c; /* wavespeed */ 30c4762a1bSJed Brown }; 31c4762a1bSJed Brown 32*2a8381b2SBarry Smith static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec U, Vec F, PetscCtx ctx) 33d71ae5a4SJacob Faibussowitsch { 34c4762a1bSJed Brown User user = (User)ctx; 35c4762a1bSJed Brown DM dm, cdm; 36c4762a1bSJed Brown DMDALocalInfo info; 37c4762a1bSJed Brown Vec C; 38c4762a1bSJed Brown Field *f; 39c4762a1bSJed Brown const Field *u; 40c4762a1bSJed Brown const PetscScalar *x; 41c4762a1bSJed Brown PetscInt i; 42c4762a1bSJed Brown 43c4762a1bSJed Brown PetscFunctionBeginUser; 449566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 459566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &C)); 479566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 489566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(dm, U, (void *)&u)); 499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, F, &f)); 509566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cdm, C, (void *)&x)); 51c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i) { 52c4762a1bSJed Brown const PetscScalar hx = i + 1 == info.xs + info.xm ? x[i] - x[i - 1] : x[i + 1] - x[i]; 53c4762a1bSJed Brown 54c4762a1bSJed Brown f[i].u = hx * (u[i].v); 55c4762a1bSJed Brown f[i].v = -hx * (PetscSqr(user->gammaTilde) * u[i].u + (PetscSqr(user->gamma) / user->xi) * (u[i].th + PetscLogScalar(u[i].v + 1))); 56c4762a1bSJed Brown f[i].th = -hx * (u[i].v + 1) * (u[i].th + (1 + user->epsilon) * PetscLogScalar(u[i].v + 1)); 57c4762a1bSJed Brown } 589566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(dm, U, (void *)&u)); 599566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, F, &f)); 609566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cdm, C, (void *)&x)); 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64*2a8381b2SBarry Smith static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx) 65d71ae5a4SJacob Faibussowitsch { 66c4762a1bSJed Brown User user = (User)ctx; 67c4762a1bSJed Brown DM dm, cdm; 68c4762a1bSJed Brown DMDALocalInfo info; 69c4762a1bSJed Brown Vec Uloc, C; 70c4762a1bSJed Brown Field *u, *udot, *f; 71c4762a1bSJed Brown PetscScalar *x; 72c4762a1bSJed Brown PetscInt i; 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscFunctionBeginUser; 759566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 769566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 789566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &C)); 799566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Uloc)); 809566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, U, INSERT_VALUES, Uloc)); 819566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, U, INSERT_VALUES, Uloc)); 829566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(dm, Uloc, &u)); 839566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(dm, Udot, &udot)); 849566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, F, &f)); 859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cdm, C, &x)); 86c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i) { 87c4762a1bSJed Brown if (i == 0) { 88c4762a1bSJed Brown const PetscScalar hx = x[i + 1] - x[i]; 89c4762a1bSJed Brown f[i].u = hx * udot[i].u; 90c4762a1bSJed Brown f[i].v = hx * udot[i].v - PetscSqr(user->c) * (u[i + 1].u - u[i].u) / hx; 91c4762a1bSJed Brown f[i].th = hx * udot[i].th; 92c4762a1bSJed Brown } else if (i == info.mx - 1) { 93c4762a1bSJed Brown const PetscScalar hx = x[i] - x[i - 1]; 94c4762a1bSJed Brown f[i].u = hx * udot[i].u; 95c4762a1bSJed Brown f[i].v = hx * udot[i].v - PetscSqr(user->c) * (u[i - 1].u - u[i].u) / hx; 96c4762a1bSJed Brown f[i].th = hx * udot[i].th; 97c4762a1bSJed Brown } else { 98c4762a1bSJed Brown const PetscScalar hx = x[i + 1] - x[i]; 99c4762a1bSJed Brown f[i].u = hx * udot[i].u; 100c4762a1bSJed Brown f[i].v = hx * udot[i].v - PetscSqr(user->c) * (u[i - 1].u - 2. * u[i].u + u[i + 1].u) / hx; 101c4762a1bSJed Brown f[i].th = hx * udot[i].th; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown } 1049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(dm, Uloc, &u)); 1059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(dm, Udot, &udot)); 1069566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, F, &f)); 1079566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cdm, C, &x)); 1089566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Uloc)); 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* IJacobian - Compute IJacobian = dF/dU + a dF/dUdot */ 113*2a8381b2SBarry Smith PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, PetscCtx ctx) 114d71ae5a4SJacob Faibussowitsch { 115c4762a1bSJed Brown User user = (User)ctx; 116c4762a1bSJed Brown DM dm, cdm; 117c4762a1bSJed Brown DMDALocalInfo info; 118c4762a1bSJed Brown Vec C; 119c4762a1bSJed Brown Field *u, *udot; 120c4762a1bSJed Brown PetscScalar *x; 121c4762a1bSJed Brown PetscInt i; 122c4762a1bSJed Brown 123c4762a1bSJed Brown PetscFunctionBeginUser; 1249566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1259566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 1269566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1279566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &C)); 1289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(dm, U, &u)); 1299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(dm, Udot, &udot)); 1309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cdm, C, &x)); 131c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i) { 132c4762a1bSJed Brown if (i == 0) { 133c4762a1bSJed Brown const PetscScalar hx = x[i + 1] - x[i]; 134c4762a1bSJed Brown const PetscInt row = i, col[] = {i, i + 1}; 135c4762a1bSJed Brown const PetscScalar dxx0 = PetscSqr(user->c) / hx, dxxR = -PetscSqr(user->c) / hx; 1369371c9d4SSatish Balay const PetscScalar vals[3][2][3] = { 1379371c9d4SSatish Balay {{a * hx, 0, 0}, {0, 0, 0} }, 138c4762a1bSJed Brown {{0, a * hx + dxx0, 0}, {0, dxxR, 0}}, 1399371c9d4SSatish Balay {{0, 0, a * hx}, {0, 0, 0} } 1409371c9d4SSatish Balay }; 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES)); 143c4762a1bSJed Brown } else if (i == info.mx - 1) { 144c4762a1bSJed Brown const PetscScalar hx = x[i + 1] - x[i]; 145c4762a1bSJed Brown const PetscInt row = i, col[] = {i - 1, i}; 146c4762a1bSJed Brown const PetscScalar dxxL = -PetscSqr(user->c) / hx, dxx0 = PetscSqr(user->c) / hx; 1479371c9d4SSatish Balay const PetscScalar vals[3][2][3] = { 1489371c9d4SSatish Balay {{0, 0, 0}, {a * hx, 0, 0} }, 149c4762a1bSJed Brown {{0, dxxL, 0}, {0, a * hx + dxx0, 0}}, 1509371c9d4SSatish Balay {{0, 0, 0}, {0, 0, a * hx} } 1519371c9d4SSatish Balay }; 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES)); 154c4762a1bSJed Brown } else { 155c4762a1bSJed Brown const PetscScalar hx = x[i + 1] - x[i]; 156c4762a1bSJed Brown const PetscInt row = i, col[] = {i - 1, i, i + 1}; 157c4762a1bSJed Brown const PetscScalar dxxL = -PetscSqr(user->c) / hx, dxx0 = 2. * PetscSqr(user->c) / hx, dxxR = -PetscSqr(user->c) / hx; 1589371c9d4SSatish Balay const PetscScalar vals[3][3][3] = { 1599371c9d4SSatish Balay {{0, 0, 0}, {a * hx, 0, 0}, {0, 0, 0} }, 160c4762a1bSJed Brown {{0, dxxL, 0}, {0, a * hx + dxx0, 0}, {0, dxxR, 0}}, 1619371c9d4SSatish Balay {{0, 0, 0}, {0, 0, a * hx}, {0, 0, 0} } 1629371c9d4SSatish Balay }; 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES)); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 1679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(dm, U, &u)); 1689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(dm, Udot, &udot)); 1699566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cdm, C, &x)); 1709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 1719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 172c4762a1bSJed Brown if (J != Jpre) { 1739566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1749566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 175c4762a1bSJed Brown } 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179*2a8381b2SBarry Smith PetscErrorCode FormInitialSolution(TS ts, Vec U, PetscCtx ctx) 180d71ae5a4SJacob Faibussowitsch { 181c4762a1bSJed Brown /* User user = (User) ctx; */ 182c4762a1bSJed Brown DM dm, cdm; 183c4762a1bSJed Brown DMDALocalInfo info; 184c4762a1bSJed Brown Vec C; 185c4762a1bSJed Brown Field *u; 186c4762a1bSJed Brown PetscScalar *x; 187c4762a1bSJed Brown const PetscReal sigma = 1.0; 188c4762a1bSJed Brown PetscInt i; 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscFunctionBeginUser; 1919566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1929566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1939566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &C)); 1949566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 1959566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, U, &u)); 1969566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cdm, C, &x)); 197c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i) { 198c4762a1bSJed Brown u[i].u = 1.5 * PetscExpScalar(-PetscSqr(x[i] - 10) / PetscSqr(sigma)); 199c4762a1bSJed Brown u[i].v = 0.0; 200c4762a1bSJed Brown u[i].th = 0.0; 201c4762a1bSJed Brown } 2029566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, U, &u)); 2039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cdm, C, &x)); 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 208d71ae5a4SJacob Faibussowitsch { 209c4762a1bSJed Brown DM dm; 210c4762a1bSJed Brown TS ts; 211c4762a1bSJed Brown Vec X; 212c4762a1bSJed Brown Mat J; 213c4762a1bSJed Brown PetscInt steps, mx; 214c4762a1bSJed Brown PetscReal ftime, hx, dt; 215c4762a1bSJed Brown TSConvergedReason reason; 216c4762a1bSJed Brown struct _User user; 217c4762a1bSJed Brown 218327415f7SBarry Smith PetscFunctionBeginUser; 2199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2209566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 3, 1, NULL, &dm)); 2219566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 2229566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 2239566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dm, 0.0, 20.0, 0.0, 0.0, 0.0, 0.0)); 2249566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &X)); 225c4762a1bSJed Brown 226d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Dynamic Friction Options", ""); 227c4762a1bSJed Brown { 228c4762a1bSJed Brown user.epsilon = 0.1; 229c4762a1bSJed Brown user.gamma = 0.5; 230c4762a1bSJed Brown user.gammaTilde = 0.5; 231c4762a1bSJed Brown user.xi = 0.5; 232c4762a1bSJed Brown user.c = 0.5; 2339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-epsilon", "Inverse of seismic ratio", "", user.epsilon, &user.epsilon, NULL)); 2349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-gamma", "Wave frequency for interblock coupling", "", user.gamma, &user.gamma, NULL)); 2359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-gamma_tilde", "Wave frequency for plate coupling", "", user.gammaTilde, &user.gammaTilde, NULL)); 2369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xi", "Interblock spring constant", "", user.xi, &user.xi, NULL)); 2379566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-c", "Wavespeed", "", user.c, &user.c, NULL)); 238c4762a1bSJed Brown } 239d0609cedSBarry Smith PetscOptionsEnd(); 240c4762a1bSJed Brown 2419566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2429566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 2439566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user)); 2449566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user)); 2459566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm, MATAIJ)); 2469566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J)); 2479566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown ftime = 800.0; 2509566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 2519566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 2529566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts, X, &user)); 2539566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X)); 2549566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &mx)); 255c4762a1bSJed Brown hx = 20.0 / (PetscReal)(mx - 1); 256c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / PetscSqr(user.c); /* Diffusive stability limit */ 2579566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 2589566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 2619566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2629566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 2639566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 26463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 2689566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2709566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 271b122ec5aSJacob Faibussowitsch return 0; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown /*TEST 275c4762a1bSJed Brown 276c4762a1bSJed Brown build: 277f56ea12dSJed Brown requires: !single !complex 278c4762a1bSJed Brown 279c4762a1bSJed Brown test: 280c4762a1bSJed Brown TODO: broken, was not nightly tested, SNES solve eventually fails, -snes_test_jacobian indicates Jacobian is wrong, but even -snes_mf_operator fails 281c4762a1bSJed Brown 282c4762a1bSJed Brown TEST*/ 283