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 PetscErrorCode ierr; 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscFunctionBeginUser; 45c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(dm, U, (void*)&u);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = DMDAVecGetArray(dm, F, &f);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cdm, C, (void*)&x);CHKERRQ(ierr); 52c4762a1bSJed Brown for (i = info.xs; i < info.xs+info.xm; ++i) { 53c4762a1bSJed Brown const PetscScalar hx = i+1 == info.xs+info.xm ? x[i] - x[i-1] : x[i+1] - x[i]; 54c4762a1bSJed Brown 55c4762a1bSJed Brown f[i].u = hx*(u[i].v); 56c4762a1bSJed Brown f[i].v = -hx*(PetscSqr(user->gammaTilde)*u[i].u + (PetscSqr(user->gamma) / user->xi)*(u[i].th + PetscLogScalar(u[i].v + 1))); 57c4762a1bSJed Brown f[i].th = -hx*(u[i].v + 1)*(u[i].th + (1 + user->epsilon)*PetscLogScalar(u[i].v + 1)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(dm, U, (void*)&u);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = DMDAVecRestoreArray(dm, F, &f);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cdm, C, (void*)&x);CHKERRQ(ierr); 62c4762a1bSJed Brown PetscFunctionReturn(0); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 66c4762a1bSJed Brown { 67c4762a1bSJed Brown User user = (User) ctx; 68c4762a1bSJed Brown DM dm, cdm; 69c4762a1bSJed Brown DMDALocalInfo info; 70c4762a1bSJed Brown Vec Uloc, C; 71c4762a1bSJed Brown Field *u, *udot, *f; 72c4762a1bSJed Brown PetscScalar *x; 73c4762a1bSJed Brown PetscInt i; 74c4762a1bSJed Brown PetscErrorCode ierr; 75c4762a1bSJed Brown 76c4762a1bSJed Brown PetscFunctionBeginUser; 77c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &Uloc);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, Uloc);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, Uloc);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(dm, Uloc, &u);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(dm, Udot, &udot);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = DMDAVecGetArray(dm, F, &f);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cdm, C, &x);CHKERRQ(ierr); 88c4762a1bSJed Brown for (i = info.xs; i < info.xs+info.xm; ++i) { 89c4762a1bSJed Brown if (i == 0) { 90c4762a1bSJed Brown const PetscScalar hx = x[i+1] - x[i]; 91c4762a1bSJed Brown f[i].u = hx * udot[i].u; 92c4762a1bSJed Brown f[i].v = hx * udot[i].v - PetscSqr(user->c) * (u[i+1].u - u[i].u) / hx; 93c4762a1bSJed Brown f[i].th = hx * udot[i].th; 94c4762a1bSJed Brown } else if (i == info.mx-1) { 95c4762a1bSJed Brown const PetscScalar hx = x[i] - x[i-1]; 96c4762a1bSJed Brown f[i].u = hx * udot[i].u; 97c4762a1bSJed Brown f[i].v = hx * udot[i].v - PetscSqr(user->c) * (u[i-1].u - u[i].u) / hx; 98c4762a1bSJed Brown f[i].th = hx * udot[i].th; 99c4762a1bSJed Brown } else { 100c4762a1bSJed Brown const PetscScalar hx = x[i+1] - x[i]; 101c4762a1bSJed Brown f[i].u = hx * udot[i].u; 102c4762a1bSJed Brown f[i].v = hx * udot[i].v - PetscSqr(user->c) * (u[i-1].u - 2.*u[i].u + u[i+1].u) / hx; 103c4762a1bSJed Brown f[i].th = hx * udot[i].th; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown } 106c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(dm, Uloc, &u);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(dm, Udot, &udot);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = DMDAVecRestoreArray(dm, F, &f);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cdm, C, &x);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &Uloc);CHKERRQ(ierr); 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* IJacobian - Compute IJacobian = dF/dU + a dF/dUdot */ 115c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, void *ctx) 116c4762a1bSJed Brown { 117c4762a1bSJed Brown User user = (User) ctx; 118c4762a1bSJed Brown DM dm, cdm; 119c4762a1bSJed Brown DMDALocalInfo info; 120c4762a1bSJed Brown Vec C; 121c4762a1bSJed Brown Field *u, *udot; 122c4762a1bSJed Brown PetscScalar *x; 123c4762a1bSJed Brown PetscInt i; 124c4762a1bSJed Brown PetscErrorCode ierr; 125c4762a1bSJed Brown 126c4762a1bSJed Brown PetscFunctionBeginUser; 127c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(dm, U, &u);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(dm, Udot, &udot);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cdm, C, &x);CHKERRQ(ierr); 134c4762a1bSJed Brown for (i = info.xs; i < info.xs+info.xm; ++i) { 135c4762a1bSJed Brown if (i == 0) { 136c4762a1bSJed Brown const PetscScalar hx = x[i+1] - x[i]; 137c4762a1bSJed Brown const PetscInt row = i, col[] = {i,i+1}; 138c4762a1bSJed Brown const PetscScalar dxx0 = PetscSqr(user->c)/hx,dxxR = -PetscSqr(user->c)/hx; 139c4762a1bSJed Brown const PetscScalar vals[3][2][3] = {{{a*hx, 0,0},{0,0, 0}}, 140c4762a1bSJed Brown {{0,a*hx+dxx0,0},{0,dxxR,0}}, 141c4762a1bSJed Brown {{0,0, a*hx},{0,0, 0}}}; 142c4762a1bSJed Brown 143c4762a1bSJed Brown ierr = MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES);CHKERRQ(ierr); 144c4762a1bSJed Brown } else if (i == info.mx-1) { 145c4762a1bSJed Brown const PetscScalar hx = x[i+1] - x[i]; 146c4762a1bSJed Brown const PetscInt row = i, col[] = {i-1,i}; 147c4762a1bSJed Brown const PetscScalar dxxL = -PetscSqr(user->c)/hx, dxx0 = PetscSqr(user->c)/hx; 148c4762a1bSJed Brown const PetscScalar vals[3][2][3] = {{{0,0, 0},{a*hx, 0,0}}, 149c4762a1bSJed Brown {{0,dxxL,0},{0,a*hx+dxx0,0}}, 150c4762a1bSJed Brown {{0,0, 0},{0,0, a*hx}}}; 151c4762a1bSJed Brown 152c4762a1bSJed Brown ierr = MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES);CHKERRQ(ierr); 153c4762a1bSJed Brown } else { 154c4762a1bSJed Brown const PetscScalar hx = x[i+1] - x[i]; 155c4762a1bSJed Brown const PetscInt row = i, col[] = {i-1,i,i+1}; 156c4762a1bSJed Brown const PetscScalar dxxL = -PetscSqr(user->c)/hx, dxx0 = 2.*PetscSqr(user->c)/hx,dxxR = -PetscSqr(user->c)/hx; 157c4762a1bSJed Brown const PetscScalar vals[3][3][3] = {{{0,0, 0},{a*hx, 0,0},{0,0, 0}}, 158c4762a1bSJed Brown {{0,dxxL,0},{0,a*hx+dxx0,0},{0,dxxR,0}}, 159c4762a1bSJed Brown {{0,0, 0},{0,0, a*hx},{0,0, 0}}}; 160c4762a1bSJed Brown 161c4762a1bSJed Brown ierr = MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES);CHKERRQ(ierr); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown } 164c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(dm, U, &u);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(dm, Udot, &udot);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cdm, C, &x);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169c4762a1bSJed Brown if (J != Jpre) { 170c4762a1bSJed Brown ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown PetscFunctionReturn(0); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts, Vec U, void *ctx) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown /* User user = (User) ctx; */ 179c4762a1bSJed Brown DM dm, cdm; 180c4762a1bSJed Brown DMDALocalInfo info; 181c4762a1bSJed Brown Vec C; 182c4762a1bSJed Brown Field *u; 183c4762a1bSJed Brown PetscScalar *x; 184c4762a1bSJed Brown const PetscReal sigma = 1.0; 185c4762a1bSJed Brown PetscInt i; 186c4762a1bSJed Brown PetscErrorCode ierr; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBeginUser; 189c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = DMDAVecGetArray(dm, U, &u);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cdm, C, &x);CHKERRQ(ierr); 195c4762a1bSJed Brown for (i = info.xs; i < info.xs+info.xm; ++i) { 196c4762a1bSJed Brown u[i].u = 1.5 * PetscExpScalar(-PetscSqr(x[i] - 10)/PetscSqr(sigma)); 197c4762a1bSJed Brown u[i].v = 0.0; 198c4762a1bSJed Brown u[i].th = 0.0; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown ierr = DMDAVecRestoreArray(dm, U, &u);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cdm, C, &x);CHKERRQ(ierr); 202c4762a1bSJed Brown PetscFunctionReturn(0); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown int main(int argc, char **argv) 206c4762a1bSJed Brown { 207c4762a1bSJed Brown DM dm; 208c4762a1bSJed Brown TS ts; 209c4762a1bSJed Brown Vec X; 210c4762a1bSJed Brown Mat J; 211c4762a1bSJed Brown PetscInt steps, mx; 212c4762a1bSJed Brown PetscReal ftime, hx, dt; 213c4762a1bSJed Brown TSConvergedReason reason; 214c4762a1bSJed Brown struct _User user; 215c4762a1bSJed Brown PetscErrorCode ierr; 216c4762a1bSJed Brown 217c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 218c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 3, 1, NULL, &dm);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(dm, 0.0, 20.0, 0.0, 0.0, 0.0, 0.0);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &X);CHKERRQ(ierr); 223c4762a1bSJed Brown 224c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Dynamic Friction Options", ""); 225c4762a1bSJed Brown { 226c4762a1bSJed Brown user.epsilon = 0.1; 227c4762a1bSJed Brown user.gamma = 0.5; 228c4762a1bSJed Brown user.gammaTilde = 0.5; 229c4762a1bSJed Brown user.xi = 0.5; 230c4762a1bSJed Brown user.c = 0.5; 231c4762a1bSJed Brown ierr = PetscOptionsReal("-epsilon", "Inverse of seismic ratio", "", user.epsilon, &user.epsilon, NULL);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = PetscOptionsReal("-gamma", "Wave frequency for interblock coupling", "", user.gamma, &user.gamma, NULL);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = PetscOptionsReal("-gamma_tilde", "Wave frequency for plate coupling", "", user.gammaTilde, &user.gammaTilde, NULL);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = PetscOptionsReal("-xi", "Interblock spring constant", "", user.xi, &user.xi, NULL);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = PetscOptionsReal("-c", "Wavespeed", "", user.c, &user.c, NULL);CHKERRQ(ierr); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 238c4762a1bSJed Brown 239c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = TSSetRHSFunction(ts, NULL, FormRHSFunction, &user);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = TSSetIFunction(ts, NULL, FormIFunction, &user);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = DMSetMatType(dm, MATAIJ);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = TSSetIJacobian(ts, J, J, FormIJacobian, &user);CHKERRQ(ierr); 246c4762a1bSJed Brown 247c4762a1bSJed Brown ftime = 800.0; 248c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = FormInitialSolution(ts, X, &user);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = TSSetSolution(ts, X);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = VecGetSize(X, &mx);CHKERRQ(ierr); 253c4762a1bSJed Brown hx = 20.0/(PetscReal)(mx-1); 254c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / PetscSqr(user.c); /* Diffusive stability limit */ 255c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 257c4762a1bSJed Brown 258c4762a1bSJed Brown ierr = TSSolve(ts, X);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = TSGetStepNumber(ts, &steps);CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = TSGetConvergedReason(ts, &reason);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double)ftime, steps);CHKERRQ(ierr); 263c4762a1bSJed Brown 264c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 265c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 268c4762a1bSJed Brown ierr = PetscFinalize(); 269c4762a1bSJed Brown return ierr; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown /*TEST 273c4762a1bSJed Brown 274c4762a1bSJed Brown build: 275*f56ea12dSJed Brown requires: !single !complex 276c4762a1bSJed Brown 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown TODO: broken, was not nightly tested, SNES solve eventually fails, -snes_test_jacobian indicates Jacobian is wrong, but even -snes_mf_operator fails 279c4762a1bSJed Brown 280c4762a1bSJed Brown TEST*/ 281