1c4762a1bSJed Brown 2c4762a1bSJed Brown static const char help[] = "Solves PDE optimization problem using full-space method, interlaces state and adjoint variables.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdmredundant.h> 7c4762a1bSJed Brown #include <petscdmcomposite.h> 8c4762a1bSJed Brown #include <petscpf.h> 9c4762a1bSJed Brown #include <petscsnes.h> 10c4762a1bSJed Brown #include <petsc/private/dmimpl.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown 14c4762a1bSJed Brown w - design variables (what we change to get an optimal solution) 15c4762a1bSJed Brown u - state variables (i.e. the PDE solution) 16c4762a1bSJed Brown lambda - the Lagrange multipliers 17c4762a1bSJed Brown 18c4762a1bSJed Brown U = (w [u_0 lambda_0 u_1 lambda_1 .....]) 19c4762a1bSJed Brown 20c4762a1bSJed Brown fu, fw, flambda contain the gradient of L(w,u,lambda) 21c4762a1bSJed Brown 22c4762a1bSJed Brown FU = (fw [fu_0 flambda_0 .....]) 23c4762a1bSJed Brown 24c4762a1bSJed Brown In this example the PDE is 25c4762a1bSJed Brown Uxx = 2, 26c4762a1bSJed Brown u(0) = w(0), thus this is the free parameter 27c4762a1bSJed Brown u(1) = 0 28c4762a1bSJed Brown the function we wish to minimize is 29c4762a1bSJed Brown \integral u^{2} 30c4762a1bSJed Brown 31c4762a1bSJed Brown The exact solution for u is given by u(x) = x*x - 1.25*x + .25 32c4762a1bSJed Brown 33c4762a1bSJed Brown Use the usual centered finite differences. 34c4762a1bSJed Brown 35c4762a1bSJed Brown Note we treat the problem as non-linear though it happens to be linear 36c4762a1bSJed Brown 37c4762a1bSJed Brown See ex21.c for the same code, but that does NOT interlaces the u and the lambda 38c4762a1bSJed Brown 39c4762a1bSJed Brown The vectors u_lambda and fu_lambda contain the u and the lambda interlaced 40c4762a1bSJed Brown */ 41c4762a1bSJed Brown 42c4762a1bSJed Brown typedef struct { 43c4762a1bSJed Brown PetscViewer u_lambda_viewer; 44c4762a1bSJed Brown PetscViewer fu_lambda_viewer; 45c4762a1bSJed Brown } UserCtx; 46c4762a1bSJed Brown 47c4762a1bSJed Brown extern PetscErrorCode ComputeFunction(SNES, Vec, Vec, void *); 48c4762a1bSJed Brown extern PetscErrorCode ComputeJacobian_MF(SNES, Vec, Mat, Mat, void *); 49c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown Uses full multigrid preconditioner with GMRES (with no preconditioner inside the GMRES) as the 53c4762a1bSJed Brown smoother on all levels. This is because (1) in the matrix free case no matrix entries are 54c4762a1bSJed Brown available for doing Jacobi or SOR preconditioning and (2) the explicit matrix case the diagonal 55c4762a1bSJed Brown entry for the control variable is zero which means default SOR will not work. 56c4762a1bSJed Brown 57c4762a1bSJed Brown */ 58c4762a1bSJed Brown char common_options[] = "-ksp_type fgmres\ 59c4762a1bSJed Brown -snes_grid_sequence 2 \ 60c4762a1bSJed Brown -pc_type mg\ 61c4762a1bSJed Brown -mg_levels_pc_type none \ 62c4762a1bSJed Brown -mg_coarse_pc_type none \ 63c4762a1bSJed Brown -pc_mg_type full \ 64c4762a1bSJed Brown -mg_coarse_ksp_type gmres \ 65c4762a1bSJed Brown -mg_levels_ksp_type gmres \ 66c4762a1bSJed Brown -mg_coarse_ksp_max_it 6 \ 67c4762a1bSJed Brown -mg_levels_ksp_max_it 3"; 68c4762a1bSJed Brown 69c4762a1bSJed Brown char matrix_free_options[] = "-mat_mffd_compute_normu no \ 70c4762a1bSJed Brown -mat_mffd_type wp"; 71c4762a1bSJed Brown 72c4762a1bSJed Brown extern PetscErrorCode DMCreateMatrix_MF(DM, Mat *); 73c4762a1bSJed Brown 74d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 75d71ae5a4SJacob Faibussowitsch { 76c4762a1bSJed Brown UserCtx user; 77c4762a1bSJed Brown DM red, da; 78c4762a1bSJed Brown SNES snes; 79c4762a1bSJed Brown DM packer; 80c4762a1bSJed Brown PetscBool use_monitor = PETSC_FALSE; 81c4762a1bSJed Brown 82327415f7SBarry Smith PetscFunctionBeginUser; 839566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Hardwire several options; can be changed at command line */ 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsertString(NULL, common_options)); 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsertString(NULL, matrix_free_options)); 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &argc, &argv, NULL)); 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_monitor", &use_monitor, PETSC_IGNORE)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Create a global vector that includes a single redundant array and two da arrays */ 929566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &packer)); 939566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(PETSC_COMM_WORLD, 0, 1, &red)); 949566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(red, "red_")); 959566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(packer, red)); 969566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 5, 2, 1, NULL, &da)); 979566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(red, "da_")); 989566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 999566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1009566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u")); 1019566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "lambda")); 1029566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(packer, (DM)da)); 1039566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(packer, &user)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown packer->ops->creatematrix = DMCreateMatrix_MF; 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* create nonlinear multi-level solver */ 1089566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 1099566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, packer)); 1109566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, ComputeFunction, NULL)); 1119566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, ComputeJacobian_MF, NULL)); 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown if (use_monitor) { 116c4762a1bSJed Brown /* create graphics windows */ 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "u_lambda - state variables and Lagrange multipliers", -1, -1, -1, -1, &user.u_lambda_viewer)); 118*9d2ffcd0SPierre Jolivet PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "fu_lambda - derivative w.r.t. state variables and Lagrange multipliers", -1, -1, -1, -1, &user.fu_lambda_viewer)); 1199566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(snes, Monitor, 0, 0)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, NULL)); 1239566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(DMDestroy(&red)); 1269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&packer)); 128c4762a1bSJed Brown if (use_monitor) { 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&user.u_lambda_viewer)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&user.fu_lambda_viewer)); 131c4762a1bSJed Brown } 1329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 133b122ec5aSJacob Faibussowitsch return 0; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown typedef struct { 137c4762a1bSJed Brown PetscScalar u; 138c4762a1bSJed Brown PetscScalar lambda; 139c4762a1bSJed Brown } ULambda; 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* 142d5b43468SJose E. Roman Evaluates FU = Gradient(L(w,u,lambda)) 143c4762a1bSJed Brown 144c4762a1bSJed Brown This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and 145c4762a1bSJed Brown DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeGetAccess()). 146c4762a1bSJed Brown 147c4762a1bSJed Brown */ 148d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeFunction(SNES snes, Vec U, Vec FU, void *ctx) 149d71ae5a4SJacob Faibussowitsch { 150c4762a1bSJed Brown PetscInt xs, xm, i, N; 151c4762a1bSJed Brown ULambda *u_lambda, *fu_lambda; 152c4762a1bSJed Brown PetscScalar d, h, *w, *fw; 153c4762a1bSJed Brown Vec vw, vfw, vu_lambda, vfu_lambda; 154c4762a1bSJed Brown DM packer, red, da; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBeginUser; 1579566063dSJacob Faibussowitsch PetscCall(VecGetDM(U, &packer)); 1589566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(packer, &red, &da)); 1599566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(packer, &vw, &vu_lambda)); 1609566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(packer, U, vw, vu_lambda)); 1619566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(packer, FU, &vfw, &vfu_lambda)); 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 1649566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &N, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 1659566063dSJacob Faibussowitsch PetscCall(VecGetArray(vw, &w)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetArray(vfw, &fw)); 1679566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, vu_lambda, &u_lambda)); 1689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, vfu_lambda, &fu_lambda)); 169c4762a1bSJed Brown d = N - 1.0; 170c4762a1bSJed Brown h = 1.0 / d; 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* derivative of L() w.r.t. w */ 173c4762a1bSJed Brown if (xs == 0) { /* only first processor computes this */ 174c4762a1bSJed Brown fw[0] = -2.0 * d * u_lambda[0].lambda; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* derivative of L() w.r.t. u */ 178c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 179c4762a1bSJed Brown if (i == 0) fu_lambda[0].lambda = h * u_lambda[0].u + 2. * d * u_lambda[0].lambda - d * u_lambda[1].lambda; 180c4762a1bSJed Brown else if (i == 1) fu_lambda[1].lambda = 2. * h * u_lambda[1].u + 2. * d * u_lambda[1].lambda - d * u_lambda[2].lambda; 181c4762a1bSJed Brown else if (i == N - 1) fu_lambda[N - 1].lambda = h * u_lambda[N - 1].u + 2. * d * u_lambda[N - 1].lambda - d * u_lambda[N - 2].lambda; 182c4762a1bSJed Brown else if (i == N - 2) fu_lambda[N - 2].lambda = 2. * h * u_lambda[N - 2].u + 2. * d * u_lambda[N - 2].lambda - d * u_lambda[N - 3].lambda; 183c4762a1bSJed Brown else fu_lambda[i].lambda = 2. * h * u_lambda[i].u - d * (u_lambda[i + 1].lambda - 2.0 * u_lambda[i].lambda + u_lambda[i - 1].lambda); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* derivative of L() w.r.t. lambda */ 187c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 188c4762a1bSJed Brown if (i == 0) fu_lambda[0].u = 2.0 * d * (u_lambda[0].u - w[0]); 189c4762a1bSJed Brown else if (i == N - 1) fu_lambda[N - 1].u = 2.0 * d * u_lambda[N - 1].u; 190c4762a1bSJed Brown else fu_lambda[i].u = -(d * (u_lambda[i + 1].u - 2.0 * u_lambda[i].u + u_lambda[i - 1].u) - 2.0 * h); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vw, &w)); 1949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vfw, &fw)); 1959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, vu_lambda, &u_lambda)); 1969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, vfu_lambda, &fu_lambda)); 1979566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(packer, &vw, &vu_lambda)); 1989566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(packer, FU, &vfw, &vfu_lambda)); 1999566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(13.0 * N)); 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* 204c4762a1bSJed Brown Computes the exact solution 205c4762a1bSJed Brown */ 206d71ae5a4SJacob Faibussowitsch PetscErrorCode u_solution(void *dummy, PetscInt n, const PetscScalar *x, PetscScalar *u) 207d71ae5a4SJacob Faibussowitsch { 208c4762a1bSJed Brown PetscInt i; 209c4762a1bSJed Brown 210c4762a1bSJed Brown PetscFunctionBeginUser; 211c4762a1bSJed Brown for (i = 0; i < n; i++) u[2 * i] = x[i] * x[i] - 1.25 * x[i] + .25; 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(DM packer, Vec U) 216d71ae5a4SJacob Faibussowitsch { 217c4762a1bSJed Brown PF pf; 218c4762a1bSJed Brown Vec x, u_global; 219c4762a1bSJed Brown PetscScalar *w; 220c4762a1bSJed Brown DM da; 221c4762a1bSJed Brown PetscInt m; 222c4762a1bSJed Brown 223c4762a1bSJed Brown PetscFunctionBeginUser; 2249566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(packer, &m, &da)); 225c4762a1bSJed Brown 2269566063dSJacob Faibussowitsch PetscCall(PFCreate(PETSC_COMM_WORLD, 1, 2, &pf)); 227c4762a1bSJed Brown /* The cast through PETSC_UINTPTR_T is so that compilers will warn about casting to void * from void(*)(void) */ 2289566063dSJacob Faibussowitsch PetscCall(PFSetType(pf, PFQUICK, (void *)(PETSC_UINTPTR_T)u_solution)); 2299566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &x)); 230c4762a1bSJed Brown if (!x) { 2319566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 2329566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &x)); 233c4762a1bSJed Brown } 2349566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(packer, U, &w, &u_global, 0)); 235c4762a1bSJed Brown if (w) w[0] = .25; 2369566063dSJacob Faibussowitsch PetscCall(PFApplyVec(pf, x, u_global)); 2379566063dSJacob Faibussowitsch PetscCall(PFDestroy(&pf)); 2389566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(packer, U, &w, &u_global, 0)); 2393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal rnorm, void *dummy) 243d71ae5a4SJacob Faibussowitsch { 244c4762a1bSJed Brown UserCtx *user; 245c4762a1bSJed Brown PetscInt m, N; 246c4762a1bSJed Brown PetscScalar *w, *dw; 247c4762a1bSJed Brown Vec u_lambda, U, F, Uexact; 248c4762a1bSJed Brown DM packer; 249c4762a1bSJed Brown PetscReal norm; 250c4762a1bSJed Brown DM da; 251c4762a1bSJed Brown 252c4762a1bSJed Brown PetscFunctionBeginUser; 2539566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &packer)); 2549566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(packer, &user)); 2559566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &U)); 2569566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(packer, U, &w, &u_lambda)); 2579566063dSJacob Faibussowitsch PetscCall(VecView(u_lambda, user->u_lambda_viewer)); 2589566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(packer, U, &w, &u_lambda)); 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &F, 0, 0)); 2619566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(packer, F, &w, &u_lambda)); 262c4762a1bSJed Brown /* ierr = VecView(u_lambda,user->fu_lambda_viewer); */ 2639566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(packer, U, &w, &u_lambda)); 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(packer, &m, &da)); 2669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &N, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 2679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &Uexact)); 2689566063dSJacob Faibussowitsch PetscCall(ExactSolution(packer, Uexact)); 2699566063dSJacob Faibussowitsch PetscCall(VecAXPY(Uexact, -1.0, U)); 2709566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(packer, Uexact, &dw, &u_lambda)); 2719566063dSJacob Faibussowitsch PetscCall(VecStrideNorm(u_lambda, 0, NORM_2, &norm)); 272c4762a1bSJed Brown norm = norm / PetscSqrtReal((PetscReal)N - 1.); 2739566063dSJacob Faibussowitsch if (dw) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Error at x = 0 %g\n", (double)norm, (double)PetscRealPart(dw[0]))); 2749566063dSJacob Faibussowitsch PetscCall(VecView(u_lambda, user->fu_lambda_viewer)); 2759566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(packer, Uexact, &dw, &u_lambda)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Uexact)); 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateMatrix_MF(DM packer, Mat *A) 281d71ae5a4SJacob Faibussowitsch { 282c4762a1bSJed Brown Vec t; 283c4762a1bSJed Brown PetscInt m; 284c4762a1bSJed Brown 285c4762a1bSJed Brown PetscFunctionBeginUser; 2869566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(packer, &t)); 2879566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(t, &m)); 2889566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(packer, &t)); 2899566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PETSC_COMM_WORLD, m, m, PETSC_DETERMINE, PETSC_DETERMINE, A)); 2909566063dSJacob Faibussowitsch PetscCall(MatSetUp(*A)); 2919566063dSJacob Faibussowitsch PetscCall(MatSetDM(*A, packer)); 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeJacobian_MF(SNES snes, Vec x, Mat A, Mat B, void *ctx) 296d71ae5a4SJacob Faibussowitsch { 297c4762a1bSJed Brown PetscFunctionBeginUser; 2989566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(A, (PetscErrorCode(*)(void *, Vec, Vec))SNESComputeFunction, snes)); 2999566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(A, x, NULL)); 3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown /*TEST 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown nsize: 2 307c4762a1bSJed Brown args: -da_grid_x 10 -snes_converged_reason -ksp_converged_reason -snes_view 308c4762a1bSJed Brown requires: !single 309c4762a1bSJed Brown 310c4762a1bSJed Brown TEST*/ 311