1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static const char help[] = "Solves PDE optimization problem using full-space method, interlaces state and adjoint variables.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscdm.h> 6*c4762a1bSJed Brown #include <petscdmda.h> 7*c4762a1bSJed Brown #include <petscdmredundant.h> 8*c4762a1bSJed Brown #include <petscdmcomposite.h> 9*c4762a1bSJed Brown #include <petscpf.h> 10*c4762a1bSJed Brown #include <petscsnes.h> 11*c4762a1bSJed Brown #include <petsc/private/dmimpl.h> 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown /* 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown w - design variables (what we change to get an optimal solution) 16*c4762a1bSJed Brown u - state variables (i.e. the PDE solution) 17*c4762a1bSJed Brown lambda - the Lagrange multipliers 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown U = (w [u_0 lambda_0 u_1 lambda_1 .....]) 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown fu, fw, flambda contain the gradient of L(w,u,lambda) 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown FU = (fw [fu_0 flambda_0 .....]) 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown In this example the PDE is 26*c4762a1bSJed Brown Uxx = 2, 27*c4762a1bSJed Brown u(0) = w(0), thus this is the free parameter 28*c4762a1bSJed Brown u(1) = 0 29*c4762a1bSJed Brown the function we wish to minimize is 30*c4762a1bSJed Brown \integral u^{2} 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown The exact solution for u is given by u(x) = x*x - 1.25*x + .25 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown Use the usual centered finite differences. 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown Note we treat the problem as non-linear though it happens to be linear 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown See ex21.c for the same code, but that does NOT interlaces the u and the lambda 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown The vectors u_lambda and fu_lambda contain the u and the lambda interlaced 41*c4762a1bSJed Brown */ 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown typedef struct { 44*c4762a1bSJed Brown PetscViewer u_lambda_viewer; 45*c4762a1bSJed Brown PetscViewer fu_lambda_viewer; 46*c4762a1bSJed Brown } UserCtx; 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown extern PetscErrorCode ComputeFunction(SNES,Vec,Vec,void*); 49*c4762a1bSJed Brown extern PetscErrorCode ComputeJacobian_MF(SNES,Vec,Mat,Mat,void*); 50*c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown /* 53*c4762a1bSJed Brown Uses full multigrid preconditioner with GMRES (with no preconditioner inside the GMRES) as the 54*c4762a1bSJed Brown smoother on all levels. This is because (1) in the matrix free case no matrix entries are 55*c4762a1bSJed Brown available for doing Jacobi or SOR preconditioning and (2) the explicit matrix case the diagonal 56*c4762a1bSJed Brown entry for the control variable is zero which means default SOR will not work. 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown */ 59*c4762a1bSJed Brown char common_options[] = "-ksp_type fgmres\ 60*c4762a1bSJed Brown -snes_grid_sequence 2 \ 61*c4762a1bSJed Brown -pc_type mg\ 62*c4762a1bSJed Brown -mg_levels_pc_type none \ 63*c4762a1bSJed Brown -mg_coarse_pc_type none \ 64*c4762a1bSJed Brown -pc_mg_type full \ 65*c4762a1bSJed Brown -mg_coarse_ksp_type gmres \ 66*c4762a1bSJed Brown -mg_levels_ksp_type gmres \ 67*c4762a1bSJed Brown -mg_coarse_ksp_max_it 6 \ 68*c4762a1bSJed Brown -mg_levels_ksp_max_it 3"; 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown char matrix_free_options[] = "-mat_mffd_compute_normu no \ 71*c4762a1bSJed Brown -mat_mffd_type wp"; 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown extern PetscErrorCode DMCreateMatrix_MF(DM,Mat*); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown int main(int argc,char **argv) 76*c4762a1bSJed Brown { 77*c4762a1bSJed Brown PetscErrorCode ierr; 78*c4762a1bSJed Brown UserCtx user; 79*c4762a1bSJed Brown DM red,da; 80*c4762a1bSJed Brown SNES snes; 81*c4762a1bSJed Brown DM packer; 82*c4762a1bSJed Brown PetscBool use_monitor = PETSC_FALSE; 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 85*c4762a1bSJed Brown ierr = PetscOptionsSetFromOptions(NULL);CHKERRQ(ierr); 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown /* Hardwire several options; can be changed at command line */ 88*c4762a1bSJed Brown ierr = PetscOptionsInsertString(NULL,common_options);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = PetscOptionsInsertString(NULL,matrix_free_options);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = PetscOptionsInsert(NULL,&argc,&argv,NULL);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-use_monitor",&use_monitor,PETSC_IGNORE);CHKERRQ(ierr); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown /* Create a global vector that includes a single redundant array and two da arrays */ 94*c4762a1bSJed Brown ierr = DMCompositeCreate(PETSC_COMM_WORLD,&packer);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = DMRedundantCreate(PETSC_COMM_WORLD,0,1,&red);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = DMSetOptionsPrefix(red,"red_");CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = DMCompositeAddDM(packer,red);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,5,2,1,NULL,&da);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = DMSetOptionsPrefix(red,"da_");CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"lambda");CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = DMCompositeAddDM(packer,(DM)da);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = DMSetApplicationContext(packer,&user);CHKERRQ(ierr); 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown packer->ops->creatematrix = DMCreateMatrix_MF; 108*c4762a1bSJed Brown 109*c4762a1bSJed Brown /* create nonlinear multi-level solver */ 110*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = SNESSetDM(snes,packer);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = SNESSetFunction(snes,NULL,ComputeFunction,NULL);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,NULL, NULL,ComputeJacobian_MF,NULL);CHKERRQ(ierr); 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown if (use_monitor) { 118*c4762a1bSJed Brown /* create graphics windows */ 119*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"u_lambda - state variables and Lagrange multipliers",-1,-1,-1,-1,&user.u_lambda_viewer);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"fu_lambda - derivate w.r.t. state variables and Lagrange multipliers",-1,-1,-1,-1,&user.fu_lambda_viewer);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = SNESMonitorSet(snes,Monitor,0,0);CHKERRQ(ierr); 122*c4762a1bSJed Brown } 123*c4762a1bSJed Brown 124*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,NULL);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown ierr = DMDestroy(&red);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = DMDestroy(&packer);CHKERRQ(ierr); 130*c4762a1bSJed Brown if (use_monitor) { 131*c4762a1bSJed Brown ierr = PetscViewerDestroy(&user.u_lambda_viewer);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = PetscViewerDestroy(&user.fu_lambda_viewer);CHKERRQ(ierr); 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown ierr = PetscFinalize(); 135*c4762a1bSJed Brown return ierr; 136*c4762a1bSJed Brown } 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown typedef struct { 139*c4762a1bSJed Brown PetscScalar u; 140*c4762a1bSJed Brown PetscScalar lambda; 141*c4762a1bSJed Brown } ULambda; 142*c4762a1bSJed Brown 143*c4762a1bSJed Brown /* 144*c4762a1bSJed Brown Evaluates FU = Gradiant(L(w,u,lambda)) 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown This local function acts on the ghosted version of U (accessed via DMCompositeGetLocalVectors() and 147*c4762a1bSJed Brown DMCompositeScatter()) BUT the global, nonghosted version of FU (via DMCompositeGetAccess()). 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown */ 150*c4762a1bSJed Brown PetscErrorCode ComputeFunction(SNES snes,Vec U,Vec FU,void *ctx) 151*c4762a1bSJed Brown { 152*c4762a1bSJed Brown PetscErrorCode ierr; 153*c4762a1bSJed Brown PetscInt xs,xm,i,N; 154*c4762a1bSJed Brown ULambda *u_lambda,*fu_lambda; 155*c4762a1bSJed Brown PetscScalar d,h,*w,*fw; 156*c4762a1bSJed Brown Vec vw,vfw,vu_lambda,vfu_lambda; 157*c4762a1bSJed Brown DM packer,red,da; 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown PetscFunctionBeginUser; 160*c4762a1bSJed Brown ierr = VecGetDM(U, &packer);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = DMCompositeGetEntries(packer,&red,&da);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(packer,&vw,&vu_lambda);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = DMCompositeScatter(packer,U,vw,vu_lambda);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = DMCompositeGetAccess(packer,FU,&vfw,&vfu_lambda);CHKERRQ(ierr); 165*c4762a1bSJed Brown 166*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = VecGetArray(vw,&w);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = VecGetArray(vfw,&fw);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,vu_lambda,&u_lambda);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,vfu_lambda,&fu_lambda);CHKERRQ(ierr); 172*c4762a1bSJed Brown d = N-1.0; 173*c4762a1bSJed Brown h = 1.0/d; 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown /* derivative of L() w.r.t. w */ 176*c4762a1bSJed Brown if (xs == 0) { /* only first processor computes this */ 177*c4762a1bSJed Brown fw[0] = -2.0*d*u_lambda[0].lambda; 178*c4762a1bSJed Brown } 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown /* derivative of L() w.r.t. u */ 181*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 182*c4762a1bSJed Brown if (i == 0) fu_lambda[0].lambda = h*u_lambda[0].u + 2.*d*u_lambda[0].lambda - d*u_lambda[1].lambda; 183*c4762a1bSJed 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; 184*c4762a1bSJed 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; 185*c4762a1bSJed 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; 186*c4762a1bSJed 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); 187*c4762a1bSJed Brown } 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* derivative of L() w.r.t. lambda */ 190*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 191*c4762a1bSJed Brown if (i == 0) fu_lambda[0].u = 2.0*d*(u_lambda[0].u - w[0]); 192*c4762a1bSJed Brown else if (i == N-1) fu_lambda[N-1].u = 2.0*d*u_lambda[N-1].u; 193*c4762a1bSJed 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); 194*c4762a1bSJed Brown } 195*c4762a1bSJed Brown 196*c4762a1bSJed Brown ierr = VecRestoreArray(vw,&w);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = VecRestoreArray(vfw,&fw);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,vu_lambda,&u_lambda);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,vfu_lambda,&fu_lambda);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(packer,&vw,&vu_lambda);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(packer,FU,&vfw,&vfu_lambda);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = PetscLogFlops(13.0*N);CHKERRQ(ierr); 203*c4762a1bSJed Brown PetscFunctionReturn(0); 204*c4762a1bSJed Brown } 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown /* 207*c4762a1bSJed Brown Computes the exact solution 208*c4762a1bSJed Brown */ 209*c4762a1bSJed Brown PetscErrorCode u_solution(void *dummy,PetscInt n,const PetscScalar *x,PetscScalar *u) 210*c4762a1bSJed Brown { 211*c4762a1bSJed Brown PetscInt i; 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown PetscFunctionBeginUser; 214*c4762a1bSJed Brown for (i=0; i<n; i++) u[2*i] = x[i]*x[i] - 1.25*x[i] + .25; 215*c4762a1bSJed Brown PetscFunctionReturn(0); 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown PetscErrorCode ExactSolution(DM packer,Vec U) 219*c4762a1bSJed Brown { 220*c4762a1bSJed Brown PF pf; 221*c4762a1bSJed Brown Vec x,u_global; 222*c4762a1bSJed Brown PetscScalar *w; 223*c4762a1bSJed Brown DM da; 224*c4762a1bSJed Brown PetscErrorCode ierr; 225*c4762a1bSJed Brown PetscInt m; 226*c4762a1bSJed Brown 227*c4762a1bSJed Brown PetscFunctionBeginUser; 228*c4762a1bSJed Brown ierr = DMCompositeGetEntries(packer,&m,&da);CHKERRQ(ierr); 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown ierr = PFCreate(PETSC_COMM_WORLD,1,2,&pf);CHKERRQ(ierr); 231*c4762a1bSJed Brown /* The cast through PETSC_UINTPTR_T is so that compilers will warn about casting to void * from void(*)(void) */ 232*c4762a1bSJed Brown ierr = PFSetType(pf,PFQUICK,(void*)(PETSC_UINTPTR_T)u_solution);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = DMGetCoordinates(da,&x);CHKERRQ(ierr); 234*c4762a1bSJed Brown if (!x) { 235*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = DMGetCoordinates(da,&x);CHKERRQ(ierr); 237*c4762a1bSJed Brown } 238*c4762a1bSJed Brown ierr = DMCompositeGetAccess(packer,U,&w,&u_global,0);CHKERRQ(ierr); 239*c4762a1bSJed Brown if (w) w[0] = .25; 240*c4762a1bSJed Brown ierr = PFApplyVec(pf,x,u_global);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = PFDestroy(&pf);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(packer,U,&w,&u_global,0);CHKERRQ(ierr); 243*c4762a1bSJed Brown PetscFunctionReturn(0); 244*c4762a1bSJed Brown } 245*c4762a1bSJed Brown 246*c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal rnorm,void *dummy) 247*c4762a1bSJed Brown { 248*c4762a1bSJed Brown UserCtx *user; 249*c4762a1bSJed Brown PetscErrorCode ierr; 250*c4762a1bSJed Brown PetscInt m,N; 251*c4762a1bSJed Brown PetscScalar *w,*dw; 252*c4762a1bSJed Brown Vec u_lambda,U,F,Uexact; 253*c4762a1bSJed Brown DM packer; 254*c4762a1bSJed Brown PetscReal norm; 255*c4762a1bSJed Brown DM da; 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown PetscFunctionBeginUser; 258*c4762a1bSJed Brown ierr = SNESGetDM(snes,&packer);CHKERRQ(ierr); 259*c4762a1bSJed Brown ierr = DMGetApplicationContext(packer,&user);CHKERRQ(ierr); 260*c4762a1bSJed Brown ierr = SNESGetSolution(snes,&U);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = DMCompositeGetAccess(packer,U,&w,&u_lambda);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = VecView(u_lambda,user->u_lambda_viewer);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(packer,U,&w,&u_lambda);CHKERRQ(ierr); 264*c4762a1bSJed Brown 265*c4762a1bSJed Brown ierr = SNESGetFunction(snes,&F,0,0);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = DMCompositeGetAccess(packer,F,&w,&u_lambda);CHKERRQ(ierr); 267*c4762a1bSJed Brown /* ierr = VecView(u_lambda,user->fu_lambda_viewer); */ 268*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(packer,U,&w,&u_lambda);CHKERRQ(ierr); 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown ierr = DMCompositeGetEntries(packer,&m,&da);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = VecDuplicate(U,&Uexact);CHKERRQ(ierr); 273*c4762a1bSJed Brown ierr = ExactSolution(packer,Uexact);CHKERRQ(ierr); 274*c4762a1bSJed Brown ierr = VecAXPY(Uexact,-1.0,U);CHKERRQ(ierr); 275*c4762a1bSJed Brown ierr = DMCompositeGetAccess(packer,Uexact,&dw,&u_lambda);CHKERRQ(ierr); 276*c4762a1bSJed Brown ierr = VecStrideNorm(u_lambda,0,NORM_2,&norm);CHKERRQ(ierr); 277*c4762a1bSJed Brown norm = norm/PetscSqrtReal((PetscReal)N-1.); 278*c4762a1bSJed Brown if (dw) ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Error at x = 0 %g\n",(double)norm,(double)PetscRealPart(dw[0]));CHKERRQ(ierr); 279*c4762a1bSJed Brown ierr = VecView(u_lambda,user->fu_lambda_viewer);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(packer,Uexact,&dw,&u_lambda);CHKERRQ(ierr); 281*c4762a1bSJed Brown ierr = VecDestroy(&Uexact);CHKERRQ(ierr); 282*c4762a1bSJed Brown PetscFunctionReturn(0); 283*c4762a1bSJed Brown } 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown PetscErrorCode DMCreateMatrix_MF(DM packer,Mat *A) 286*c4762a1bSJed Brown { 287*c4762a1bSJed Brown PetscErrorCode ierr; 288*c4762a1bSJed Brown Vec t; 289*c4762a1bSJed Brown PetscInt m; 290*c4762a1bSJed Brown 291*c4762a1bSJed Brown PetscFunctionBeginUser; 292*c4762a1bSJed Brown ierr = DMGetGlobalVector(packer,&t);CHKERRQ(ierr); 293*c4762a1bSJed Brown ierr = VecGetLocalSize(t,&m);CHKERRQ(ierr); 294*c4762a1bSJed Brown ierr = DMRestoreGlobalVector(packer,&t);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = MatCreateMFFD(PETSC_COMM_WORLD,m,m,PETSC_DETERMINE,PETSC_DETERMINE,A);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = MatSetUp(*A);CHKERRQ(ierr); 297*c4762a1bSJed Brown ierr = MatSetDM(*A,packer);CHKERRQ(ierr); 298*c4762a1bSJed Brown PetscFunctionReturn(0); 299*c4762a1bSJed Brown } 300*c4762a1bSJed Brown 301*c4762a1bSJed Brown PetscErrorCode ComputeJacobian_MF(SNES snes,Vec x,Mat A,Mat B,void *ctx) 302*c4762a1bSJed Brown { 303*c4762a1bSJed Brown PetscErrorCode ierr; 304*c4762a1bSJed Brown 305*c4762a1bSJed Brown PetscFunctionBeginUser; 306*c4762a1bSJed Brown ierr = MatMFFDSetFunction(A,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);CHKERRQ(ierr); 307*c4762a1bSJed Brown ierr = MatMFFDSetBase(A,x,NULL);CHKERRQ(ierr); 308*c4762a1bSJed Brown PetscFunctionReturn(0); 309*c4762a1bSJed Brown } 310*c4762a1bSJed Brown 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown /*TEST 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown test: 315*c4762a1bSJed Brown nsize: 2 316*c4762a1bSJed Brown args: -da_grid_x 10 -snes_converged_reason -ksp_converged_reason -snes_view 317*c4762a1bSJed Brown requires: !single 318*c4762a1bSJed Brown 319*c4762a1bSJed Brown TEST*/ 320