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 74c4762a1bSJed Brown int main(int argc,char **argv) 75c4762a1bSJed Brown { 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 82*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Hardwire several options; can be changed at command line */ 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInsertString(NULL,common_options)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInsertString(NULL,matrix_free_options)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInsert(NULL,&argc,&argv,NULL)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_monitor",&use_monitor,PETSC_IGNORE)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Create a global vector that includes a single redundant array and two da arrays */ 915f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeCreate(PETSC_COMM_WORLD,&packer)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMRedundantCreate(PETSC_COMM_WORLD,0,1,&red)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(red,"red_")); 945f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(packer,red)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,2,1,NULL,&da)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(red,"da_")); 975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,0,"u")); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,1,"lambda")); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(packer,(DM)da)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(packer,&user)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown packer->ops->creatematrix = DMCreateMatrix_MF; 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* create nonlinear multi-level solver */ 1075f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes,packer)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,NULL,ComputeFunction,NULL)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,NULL, NULL,ComputeJacobian_MF,NULL)); 111c4762a1bSJed Brown 1125f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown if (use_monitor) { 115c4762a1bSJed Brown /* create graphics windows */ 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(SNESMonitorSet(snes,Monitor,0,0)); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 1215f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,NULL)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 123c4762a1bSJed Brown 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&red)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&packer)); 127c4762a1bSJed Brown if (use_monitor) { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&user.u_lambda_viewer)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&user.fu_lambda_viewer)); 130c4762a1bSJed Brown } 131*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 132*b122ec5aSJacob Faibussowitsch return 0; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown typedef struct { 136c4762a1bSJed Brown PetscScalar u; 137c4762a1bSJed Brown PetscScalar lambda; 138c4762a1bSJed Brown } ULambda; 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* 141c4762a1bSJed Brown Evaluates FU = Gradiant(L(w,u,lambda)) 142c4762a1bSJed Brown 143c4762a1bSJed Brown This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and 144c4762a1bSJed Brown DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeGetAccess()). 145c4762a1bSJed Brown 146c4762a1bSJed Brown */ 147c4762a1bSJed Brown PetscErrorCode ComputeFunction(SNES snes,Vec U,Vec FU,void *ctx) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown PetscInt xs,xm,i,N; 150c4762a1bSJed Brown ULambda *u_lambda,*fu_lambda; 151c4762a1bSJed Brown PetscScalar d,h,*w,*fw; 152c4762a1bSJed Brown Vec vw,vfw,vu_lambda,vfu_lambda; 153c4762a1bSJed Brown DM packer,red,da; 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscFunctionBeginUser; 1565f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(U, &packer)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetEntries(packer,&red,&da)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetLocalVectors(packer,&vw,&vu_lambda)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeScatter(packer,U,vw,vu_lambda)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(packer,FU,&vfw,&vfu_lambda)); 161c4762a1bSJed Brown 1625f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(vw,&w)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(vfw,&fw)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,vu_lambda,&u_lambda)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,vfu_lambda,&fu_lambda)); 168c4762a1bSJed Brown d = N-1.0; 169c4762a1bSJed Brown h = 1.0/d; 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* derivative of L() w.r.t. w */ 172c4762a1bSJed Brown if (xs == 0) { /* only first processor computes this */ 173c4762a1bSJed Brown fw[0] = -2.0*d*u_lambda[0].lambda; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* derivative of L() w.r.t. u */ 177c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 178c4762a1bSJed Brown if (i == 0) fu_lambda[0].lambda = h*u_lambda[0].u + 2.*d*u_lambda[0].lambda - d*u_lambda[1].lambda; 179c4762a1bSJed 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; 180c4762a1bSJed 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; 181c4762a1bSJed 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; 182c4762a1bSJed 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); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* derivative of L() w.r.t. lambda */ 186c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 187c4762a1bSJed Brown if (i == 0) fu_lambda[0].u = 2.0*d*(u_lambda[0].u - w[0]); 188c4762a1bSJed Brown else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u; 189c4762a1bSJed 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); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 1925f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(vw,&w)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(vfw,&fw)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,vu_lambda,&u_lambda)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,vfu_lambda,&fu_lambda)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreLocalVectors(packer,&vw,&vu_lambda)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccess(packer,FU,&vfw,&vfu_lambda)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(13.0*N)); 199c4762a1bSJed Brown PetscFunctionReturn(0); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* 203c4762a1bSJed Brown Computes the exact solution 204c4762a1bSJed Brown */ 205c4762a1bSJed Brown PetscErrorCode u_solution(void *dummy,PetscInt n,const PetscScalar *x,PetscScalar *u) 206c4762a1bSJed Brown { 207c4762a1bSJed Brown PetscInt i; 208c4762a1bSJed Brown 209c4762a1bSJed Brown PetscFunctionBeginUser; 210c4762a1bSJed Brown for (i=0; i<n; i++) u[2*i] = x[i]*x[i] - 1.25*x[i] + .25; 211c4762a1bSJed Brown PetscFunctionReturn(0); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscErrorCode ExactSolution(DM packer,Vec U) 215c4762a1bSJed Brown { 216c4762a1bSJed Brown PF pf; 217c4762a1bSJed Brown Vec x,u_global; 218c4762a1bSJed Brown PetscScalar *w; 219c4762a1bSJed Brown DM da; 220c4762a1bSJed Brown PetscInt m; 221c4762a1bSJed Brown 222c4762a1bSJed Brown PetscFunctionBeginUser; 2235f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetEntries(packer,&m,&da)); 224c4762a1bSJed Brown 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PFCreate(PETSC_COMM_WORLD,1,2,&pf)); 226c4762a1bSJed Brown /* The cast through PETSC_UINTPTR_T is so that compilers will warn about casting to void * from void(*)(void) */ 2275f80ce2aSJacob Faibussowitsch CHKERRQ(PFSetType(pf,PFQUICK,(void*)(PETSC_UINTPTR_T)u_solution)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da,&x)); 229c4762a1bSJed Brown if (!x) { 2305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da,&x)); 232c4762a1bSJed Brown } 2335f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(packer,U,&w,&u_global,0)); 234c4762a1bSJed Brown if (w) w[0] = .25; 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PFApplyVec(pf,x,u_global)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PFDestroy(&pf)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccess(packer,U,&w,&u_global,0)); 238c4762a1bSJed Brown PetscFunctionReturn(0); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy) 242c4762a1bSJed Brown { 243c4762a1bSJed Brown UserCtx *user; 244c4762a1bSJed Brown PetscInt m,N; 245c4762a1bSJed Brown PetscScalar *w,*dw; 246c4762a1bSJed Brown Vec u_lambda,U,F,Uexact; 247c4762a1bSJed Brown DM packer; 248c4762a1bSJed Brown PetscReal norm; 249c4762a1bSJed Brown DM da; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBeginUser; 2525f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&packer)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(packer,&user)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes,&U)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(packer,U,&w,&u_lambda)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u_lambda,user->u_lambda_viewer)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccess(packer,U,&w,&u_lambda)); 258c4762a1bSJed Brown 2595f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&F,0,0)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(packer,F,&w,&u_lambda)); 261c4762a1bSJed Brown /* ierr = VecView(u_lambda,user->fu_lambda_viewer); */ 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccess(packer,U,&w,&u_lambda)); 263c4762a1bSJed Brown 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetEntries(packer,&m,&da)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Uexact)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(packer,Uexact)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Uexact,-1.0,U)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccess(packer,Uexact,&dw,&u_lambda)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(VecStrideNorm(u_lambda,0,NORM_2,&norm)); 271c4762a1bSJed Brown norm = norm/PetscSqrtReal((PetscReal)N-1.); 2725f80ce2aSJacob Faibussowitsch if (dw) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Error at x = 0 %g\n",(double)norm,(double)PetscRealPart(dw[0]))); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u_lambda,user->fu_lambda_viewer)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccess(packer,Uexact,&dw,&u_lambda)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Uexact)); 276c4762a1bSJed Brown PetscFunctionReturn(0); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown PetscErrorCode DMCreateMatrix_MF(DM packer,Mat *A) 280c4762a1bSJed Brown { 281c4762a1bSJed Brown Vec t; 282c4762a1bSJed Brown PetscInt m; 283c4762a1bSJed Brown 284c4762a1bSJed Brown PetscFunctionBeginUser; 2855f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(packer,&t)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(t,&m)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(packer,&t)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMFFD(PETSC_COMM_WORLD,m,m,PETSC_DETERMINE,PETSC_DETERMINE,A)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(*A)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetDM(*A,packer)); 291c4762a1bSJed Brown PetscFunctionReturn(0); 292c4762a1bSJed Brown } 293c4762a1bSJed Brown 294c4762a1bSJed Brown PetscErrorCode ComputeJacobian_MF(SNES snes,Vec x,Mat A,Mat B,void *ctx) 295c4762a1bSJed Brown { 296c4762a1bSJed Brown PetscFunctionBeginUser; 2975f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetFunction(A,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetBase(A,x,NULL)); 299c4762a1bSJed Brown PetscFunctionReturn(0); 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown /*TEST 303c4762a1bSJed Brown 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown nsize: 2 306c4762a1bSJed Brown args: -da_grid_x 10 -snes_converged_reason -ksp_converged_reason -snes_view 307c4762a1bSJed Brown requires: !single 308c4762a1bSJed Brown 309c4762a1bSJed Brown TEST*/ 310