xref: /petsc/src/ts/tutorials/ex34.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
32c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx)
33c4762a1bSJed Brown {
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;
445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(dm, &cdm));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &C));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm, &info));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(dm,  U, (void*)&u));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,  F, &f));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(dm,  U, (void*)&u));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,  F, &f));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cdm, C, (void*)&x));
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
65c4762a1bSJed Brown {
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;
755f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm, &info));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(dm, &cdm));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &C));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm, &Uloc));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(dm, U, INSERT_VALUES, Uloc));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(dm, U, INSERT_VALUES, Uloc));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(dm,  Uloc, &u));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(dm,  Udot, &udot));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,  F,    &f));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(dm,  Uloc, &u));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(dm,  Udot, &udot));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,  F,    &f));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cdm, C,    &x));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm, &Uloc));
109c4762a1bSJed Brown   PetscFunctionReturn(0);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown /* IJacobian - Compute IJacobian = dF/dU + a dF/dUdot */
113c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, void *ctx)
114c4762a1bSJed Brown {
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;
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm, &info));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(dm, &cdm));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &C));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(dm,  U,    &u));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(dm,  Udot, &udot));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
136c4762a1bSJed Brown       const PetscScalar vals[3][2][3] = {{{a*hx,     0,0},{0,0,   0}},
137c4762a1bSJed Brown                                          {{0,a*hx+dxx0,0},{0,dxxR,0}},
138c4762a1bSJed Brown                                          {{0,0,     a*hx},{0,0,   0}}};
139c4762a1bSJed Brown 
1405f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES));
141c4762a1bSJed Brown     } else if (i == info.mx-1) {
142c4762a1bSJed Brown       const PetscScalar hx            = x[i+1] - x[i];
143c4762a1bSJed Brown       const PetscInt    row           = i, col[] = {i-1,i};
144c4762a1bSJed Brown       const PetscScalar dxxL          = -PetscSqr(user->c)/hx, dxx0 = PetscSqr(user->c)/hx;
145c4762a1bSJed Brown       const PetscScalar vals[3][2][3] = {{{0,0,   0},{a*hx,     0,0}},
146c4762a1bSJed Brown                                          {{0,dxxL,0},{0,a*hx+dxx0,0}},
147c4762a1bSJed Brown                                          {{0,0,   0},{0,0,     a*hx}}};
148c4762a1bSJed Brown 
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES));
150c4762a1bSJed Brown     } else {
151c4762a1bSJed Brown       const PetscScalar hx            = x[i+1] - x[i];
152c4762a1bSJed Brown       const PetscInt    row           = i, col[] = {i-1,i,i+1};
153c4762a1bSJed Brown       const PetscScalar dxxL          = -PetscSqr(user->c)/hx, dxx0 = 2.*PetscSqr(user->c)/hx,dxxR = -PetscSqr(user->c)/hx;
154c4762a1bSJed Brown       const PetscScalar vals[3][3][3] = {{{0,0,   0},{a*hx,     0,0},{0,0,   0}},
155c4762a1bSJed Brown                                          {{0,dxxL,0},{0,a*hx+dxx0,0},{0,dxxR,0}},
156c4762a1bSJed Brown                                          {{0,0,   0},{0,0,     a*hx},{0,0,   0}}};
157c4762a1bSJed Brown 
1585f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES));
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   }
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(dm,  U,    &u));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(dm,  Udot, &udot));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cdm, C,    &x));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
166c4762a1bSJed Brown   if (J != Jpre) {
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
169c4762a1bSJed Brown   }
170c4762a1bSJed Brown   PetscFunctionReturn(0);
171c4762a1bSJed Brown }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts, Vec U, void *ctx)
174c4762a1bSJed Brown {
175c4762a1bSJed Brown   /* User            user = (User) ctx; */
176c4762a1bSJed Brown   DM              dm, cdm;
177c4762a1bSJed Brown   DMDALocalInfo   info;
178c4762a1bSJed Brown   Vec             C;
179c4762a1bSJed Brown   Field          *u;
180c4762a1bSJed Brown   PetscScalar    *x;
181c4762a1bSJed Brown   const PetscReal sigma = 1.0;
182c4762a1bSJed Brown   PetscInt        i;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   PetscFunctionBeginUser;
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(dm, &cdm));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &C));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm, &info));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,  U, &u));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cdm, C, &x));
191c4762a1bSJed Brown   for (i = info.xs; i < info.xs+info.xm; ++i) {
192c4762a1bSJed Brown     u[i].u  = 1.5 * PetscExpScalar(-PetscSqr(x[i] - 10)/PetscSqr(sigma));
193c4762a1bSJed Brown     u[i].v  = 0.0;
194c4762a1bSJed Brown     u[i].th = 0.0;
195c4762a1bSJed Brown   }
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,  U, &u));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cdm, C, &x));
198c4762a1bSJed Brown   PetscFunctionReturn(0);
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown int main(int argc, char **argv)
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   DM                dm;
204c4762a1bSJed Brown   TS                ts;
205c4762a1bSJed Brown   Vec               X;
206c4762a1bSJed Brown   Mat               J;
207c4762a1bSJed Brown   PetscInt          steps, mx;
208c4762a1bSJed Brown   PetscReal         ftime, hx, dt;
209c4762a1bSJed Brown   TSConvergedReason reason;
210c4762a1bSJed Brown   struct _User      user;
211c4762a1bSJed Brown   PetscErrorCode    ierr;
212c4762a1bSJed Brown 
213*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 3, 1, NULL, &dm));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(dm));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(dm, 0.0, 20.0, 0.0, 0.0, 0.0, 0.0));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &X));
219c4762a1bSJed Brown 
2201e1ea65dSPierre Jolivet   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Dynamic Friction Options", "");CHKERRQ(ierr);
221c4762a1bSJed Brown   {
222c4762a1bSJed Brown     user.epsilon    = 0.1;
223c4762a1bSJed Brown     user.gamma      = 0.5;
224c4762a1bSJed Brown     user.gammaTilde = 0.5;
225c4762a1bSJed Brown     user.xi         = 0.5;
226c4762a1bSJed Brown     user.c          = 0.5;
2275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-epsilon", "Inverse of seismic ratio", "", user.epsilon, &user.epsilon, NULL));
2285f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-gamma", "Wave frequency for interblock coupling", "", user.gamma, &user.gamma, NULL));
2295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-gamma_tilde", "Wave frequency for plate coupling", "", user.gammaTilde, &user.gammaTilde, NULL));
2305f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-xi", "Interblock spring constant", "", user.xi, &user.xi, NULL));
2315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-c", "Wavespeed", "", user.c, &user.c, NULL));
232c4762a1bSJed Brown   }
233c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
234c4762a1bSJed Brown 
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts));
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts, dm));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts, NULL, FormIFunction, &user));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetMatType(dm, MATAIJ));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(dm, &J));
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts, J, J, FormIJacobian, &user));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   ftime = 800.0;
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,ftime));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialSolution(ts, X, &user));
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts, X));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(X, &mx));
249c4762a1bSJed Brown   hx   = 20.0/(PetscReal)(mx-1);
250c4762a1bSJed Brown   dt   = 0.4 * PetscSqr(hx) / PetscSqr(user.c); /* Diffusive stability limit */
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
253c4762a1bSJed Brown 
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts, X));
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts, &ftime));
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts, &steps));
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetConvergedReason(ts, &reason));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double)ftime, steps));
259c4762a1bSJed Brown 
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
2615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
2625f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
264*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
265*b122ec5aSJacob Faibussowitsch   return 0;
266c4762a1bSJed Brown }
267c4762a1bSJed Brown 
268c4762a1bSJed Brown /*TEST
269c4762a1bSJed Brown 
270c4762a1bSJed Brown     build:
271f56ea12dSJed Brown       requires: !single  !complex
272c4762a1bSJed Brown 
273c4762a1bSJed Brown     test:
274c4762a1bSJed Brown       TODO: broken, was not nightly tested, SNES solve eventually fails, -snes_test_jacobian indicates Jacobian is wrong, but even -snes_mf_operator fails
275c4762a1bSJed Brown 
276c4762a1bSJed Brown TEST*/
277