xref: /petsc/src/ts/tutorials/ex34.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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*9371c9d4SSatish Balay static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) {
33c4762a1bSJed Brown   User               user = (User)ctx;
34c4762a1bSJed Brown   DM                 dm, cdm;
35c4762a1bSJed Brown   DMDALocalInfo      info;
36c4762a1bSJed Brown   Vec                C;
37c4762a1bSJed Brown   Field             *f;
38c4762a1bSJed Brown   const Field       *u;
39c4762a1bSJed Brown   const PetscScalar *x;
40c4762a1bSJed Brown   PetscInt           i;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBeginUser;
439566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
449566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
459566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &C));
469566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
479566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(dm, U, (void *)&u));
489566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, F, &f));
499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cdm, C, (void *)&x));
50c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; ++i) {
51c4762a1bSJed Brown     const PetscScalar hx = i + 1 == info.xs + info.xm ? x[i] - x[i - 1] : x[i + 1] - x[i];
52c4762a1bSJed Brown 
53c4762a1bSJed Brown     f[i].u  = hx * (u[i].v);
54c4762a1bSJed Brown     f[i].v  = -hx * (PetscSqr(user->gammaTilde) * u[i].u + (PetscSqr(user->gamma) / user->xi) * (u[i].th + PetscLogScalar(u[i].v + 1)));
55c4762a1bSJed Brown     f[i].th = -hx * (u[i].v + 1) * (u[i].th + (1 + user->epsilon) * PetscLogScalar(u[i].v + 1));
56c4762a1bSJed Brown   }
579566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(dm, U, (void *)&u));
589566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, F, &f));
599566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cdm, C, (void *)&x));
60c4762a1bSJed Brown   PetscFunctionReturn(0);
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
63*9371c9d4SSatish Balay static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
64c4762a1bSJed Brown   User          user = (User)ctx;
65c4762a1bSJed Brown   DM            dm, cdm;
66c4762a1bSJed Brown   DMDALocalInfo info;
67c4762a1bSJed Brown   Vec           Uloc, C;
68c4762a1bSJed Brown   Field        *u, *udot, *f;
69c4762a1bSJed Brown   PetscScalar  *x;
70c4762a1bSJed Brown   PetscInt      i;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   PetscFunctionBeginUser;
739566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
749566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
759566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
769566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &C));
779566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Uloc));
789566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, U, INSERT_VALUES, Uloc));
799566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, U, INSERT_VALUES, Uloc));
809566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(dm, Uloc, &u));
819566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(dm, Udot, &udot));
829566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, F, &f));
839566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cdm, C, &x));
84c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; ++i) {
85c4762a1bSJed Brown     if (i == 0) {
86c4762a1bSJed Brown       const PetscScalar hx = x[i + 1] - x[i];
87c4762a1bSJed Brown       f[i].u               = hx * udot[i].u;
88c4762a1bSJed Brown       f[i].v               = hx * udot[i].v - PetscSqr(user->c) * (u[i + 1].u - u[i].u) / hx;
89c4762a1bSJed Brown       f[i].th              = hx * udot[i].th;
90c4762a1bSJed Brown     } else if (i == info.mx - 1) {
91c4762a1bSJed Brown       const PetscScalar hx = x[i] - x[i - 1];
92c4762a1bSJed Brown       f[i].u               = hx * udot[i].u;
93c4762a1bSJed Brown       f[i].v               = hx * udot[i].v - PetscSqr(user->c) * (u[i - 1].u - u[i].u) / hx;
94c4762a1bSJed Brown       f[i].th              = hx * udot[i].th;
95c4762a1bSJed Brown     } else {
96c4762a1bSJed Brown       const PetscScalar hx = x[i + 1] - x[i];
97c4762a1bSJed Brown       f[i].u               = hx * udot[i].u;
98c4762a1bSJed Brown       f[i].v               = hx * udot[i].v - PetscSqr(user->c) * (u[i - 1].u - 2. * u[i].u + u[i + 1].u) / hx;
99c4762a1bSJed Brown       f[i].th              = hx * udot[i].th;
100c4762a1bSJed Brown     }
101c4762a1bSJed Brown   }
1029566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(dm, Uloc, &u));
1039566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(dm, Udot, &udot));
1049566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, F, &f));
1059566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cdm, C, &x));
1069566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Uloc));
107c4762a1bSJed Brown   PetscFunctionReturn(0);
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown /* IJacobian - Compute IJacobian = dF/dU + a dF/dUdot */
111*9371c9d4SSatish Balay PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, void *ctx) {
112c4762a1bSJed Brown   User          user = (User)ctx;
113c4762a1bSJed Brown   DM            dm, cdm;
114c4762a1bSJed Brown   DMDALocalInfo info;
115c4762a1bSJed Brown   Vec           C;
116c4762a1bSJed Brown   Field        *u, *udot;
117c4762a1bSJed Brown   PetscScalar  *x;
118c4762a1bSJed Brown   PetscInt      i;
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   PetscFunctionBeginUser;
1219566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
1229566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
1239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
1249566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &C));
1259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(dm, U, &u));
1269566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(dm, Udot, &udot));
1279566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cdm, C, &x));
128c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; ++i) {
129c4762a1bSJed Brown     if (i == 0) {
130c4762a1bSJed Brown       const PetscScalar hx  = x[i + 1] - x[i];
131c4762a1bSJed Brown       const PetscInt    row = i, col[] = {i, i + 1};
132c4762a1bSJed Brown       const PetscScalar dxx0 = PetscSqr(user->c) / hx, dxxR = -PetscSqr(user->c) / hx;
133*9371c9d4SSatish Balay       const PetscScalar vals[3][2][3] = {
134*9371c9d4SSatish Balay         {{a * hx, 0, 0},        {0, 0, 0}   },
135c4762a1bSJed Brown         {{0, a * hx + dxx0, 0}, {0, dxxR, 0}},
136*9371c9d4SSatish Balay         {{0, 0, a * hx},        {0, 0, 0}   }
137*9371c9d4SSatish Balay       };
138c4762a1bSJed Brown 
1399566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES));
140c4762a1bSJed Brown     } else if (i == info.mx - 1) {
141c4762a1bSJed Brown       const PetscScalar hx  = x[i + 1] - x[i];
142c4762a1bSJed Brown       const PetscInt    row = i, col[] = {i - 1, i};
143c4762a1bSJed Brown       const PetscScalar dxxL = -PetscSqr(user->c) / hx, dxx0 = PetscSqr(user->c) / hx;
144*9371c9d4SSatish Balay       const PetscScalar vals[3][2][3] = {
145*9371c9d4SSatish Balay         {{0, 0, 0},    {a * hx, 0, 0}       },
146c4762a1bSJed Brown         {{0, dxxL, 0}, {0, a * hx + dxx0, 0}},
147*9371c9d4SSatish Balay         {{0, 0, 0},    {0, 0, a * hx}       }
148*9371c9d4SSatish Balay       };
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES));
151c4762a1bSJed Brown     } else {
152c4762a1bSJed Brown       const PetscScalar hx  = x[i + 1] - x[i];
153c4762a1bSJed Brown       const PetscInt    row = i, col[] = {i - 1, i, i + 1};
154c4762a1bSJed Brown       const PetscScalar dxxL = -PetscSqr(user->c) / hx, dxx0 = 2. * PetscSqr(user->c) / hx, dxxR = -PetscSqr(user->c) / hx;
155*9371c9d4SSatish Balay       const PetscScalar vals[3][3][3] = {
156*9371c9d4SSatish Balay         {{0, 0, 0},    {a * hx, 0, 0},        {0, 0, 0}   },
157c4762a1bSJed Brown         {{0, dxxL, 0}, {0, a * hx + dxx0, 0}, {0, dxxR, 0}},
158*9371c9d4SSatish Balay         {{0, 0, 0},    {0, 0, a * hx},        {0, 0, 0}   }
159*9371c9d4SSatish Balay       };
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES));
162c4762a1bSJed Brown     }
163c4762a1bSJed Brown   }
1649566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(dm, U, &u));
1659566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(dm, Udot, &udot));
1669566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cdm, C, &x));
1679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
1689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
169c4762a1bSJed Brown   if (J != Jpre) {
1709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1719566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176*9371c9d4SSatish Balay PetscErrorCode FormInitialSolution(TS ts, Vec U, void *ctx) {
177c4762a1bSJed Brown   /* User            user = (User) ctx; */
178c4762a1bSJed Brown   DM              dm, cdm;
179c4762a1bSJed Brown   DMDALocalInfo   info;
180c4762a1bSJed Brown   Vec             C;
181c4762a1bSJed Brown   Field          *u;
182c4762a1bSJed Brown   PetscScalar    *x;
183c4762a1bSJed Brown   const PetscReal sigma = 1.0;
184c4762a1bSJed Brown   PetscInt        i;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
1889566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
1899566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &C));
1909566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
1919566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, U, &u));
1929566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cdm, C, &x));
193c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; ++i) {
194c4762a1bSJed Brown     u[i].u  = 1.5 * PetscExpScalar(-PetscSqr(x[i] - 10) / PetscSqr(sigma));
195c4762a1bSJed Brown     u[i].v  = 0.0;
196c4762a1bSJed Brown     u[i].th = 0.0;
197c4762a1bSJed Brown   }
1989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, U, &u));
1999566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cdm, C, &x));
200c4762a1bSJed Brown   PetscFunctionReturn(0);
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203*9371c9d4SSatish Balay int main(int argc, char **argv) {
204c4762a1bSJed Brown   DM                dm;
205c4762a1bSJed Brown   TS                ts;
206c4762a1bSJed Brown   Vec               X;
207c4762a1bSJed Brown   Mat               J;
208c4762a1bSJed Brown   PetscInt          steps, mx;
209c4762a1bSJed Brown   PetscReal         ftime, hx, dt;
210c4762a1bSJed Brown   TSConvergedReason reason;
211c4762a1bSJed Brown   struct _User      user;
212c4762a1bSJed Brown 
213327415f7SBarry Smith   PetscFunctionBeginUser;
2149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2159566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 3, 1, NULL, &dm));
2169566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
2179566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
2189566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dm, 0.0, 20.0, 0.0, 0.0, 0.0, 0.0));
2199566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &X));
220c4762a1bSJed Brown 
221d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Dynamic Friction Options", "");
222c4762a1bSJed Brown   {
223c4762a1bSJed Brown     user.epsilon    = 0.1;
224c4762a1bSJed Brown     user.gamma      = 0.5;
225c4762a1bSJed Brown     user.gammaTilde = 0.5;
226c4762a1bSJed Brown     user.xi         = 0.5;
227c4762a1bSJed Brown     user.c          = 0.5;
2289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-epsilon", "Inverse of seismic ratio", "", user.epsilon, &user.epsilon, NULL));
2299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-gamma", "Wave frequency for interblock coupling", "", user.gamma, &user.gamma, NULL));
2309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-gamma_tilde", "Wave frequency for plate coupling", "", user.gammaTilde, &user.gammaTilde, NULL));
2319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-xi", "Interblock spring constant", "", user.xi, &user.xi, NULL));
2329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-c", "Wavespeed", "", user.c, &user.c, NULL));
233c4762a1bSJed Brown   }
234d0609cedSBarry Smith   PetscOptionsEnd();
235c4762a1bSJed Brown 
2369566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2379566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
2389566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user));
2399566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user));
2409566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(dm, MATAIJ));
2419566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
2429566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   ftime = 800.0;
2459566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
2469566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2479566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts, X, &user));
2489566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
2499566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &mx));
250c4762a1bSJed Brown   hx = 20.0 / (PetscReal)(mx - 1);
251c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / PetscSqr(user.c); /* Diffusive stability limit */
2529566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
2539566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
2569566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2579566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
2589566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
25963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps));
260c4762a1bSJed Brown 
2619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
2629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
2639566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2649566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
266b122ec5aSJacob Faibussowitsch   return 0;
267c4762a1bSJed Brown }
268c4762a1bSJed Brown 
269c4762a1bSJed Brown /*TEST
270c4762a1bSJed Brown 
271c4762a1bSJed Brown     build:
272f56ea12dSJed Brown       requires: !single  !complex
273c4762a1bSJed Brown 
274c4762a1bSJed Brown     test:
275c4762a1bSJed Brown       TODO: broken, was not nightly tested, SNES solve eventually fails, -snes_test_jacobian indicates Jacobian is wrong, but even -snes_mf_operator fails
276c4762a1bSJed Brown 
277c4762a1bSJed Brown TEST*/
278