1c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 2d.\n\ 2c4762a1bSJed Brown We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\ 3c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\ 4c4762a1bSJed Brown The command line options include:\n\ 5c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 6c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n\ 7c4762a1bSJed Brown -m_par/n_par <parameter>, where <parameter> indicates an integer\n \ 8c4762a1bSJed Brown that MMS3 will be evaluated with 2^m_par, 2^n_par"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* ------------------------------------------------------------------------ 11c4762a1bSJed Brown 12c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by 13c4762a1bSJed Brown the partial differential equation 14c4762a1bSJed Brown 15c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 16c4762a1bSJed Brown 17c4762a1bSJed Brown with boundary conditions 18c4762a1bSJed Brown 19c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1. 20c4762a1bSJed Brown 21c4762a1bSJed Brown A finite difference approximation with the usual 5-point stencil 22c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear 23c4762a1bSJed Brown system of equations. 24c4762a1bSJed Brown 25c4762a1bSJed Brown This example shows how geometric multigrid can be run transparently with a nonlinear solver so long 26c4762a1bSJed Brown as SNESSetDM() is provided. Example usage 27c4762a1bSJed Brown 28c4762a1bSJed Brown ./ex5 -pc_type mg -ksp_monitor -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 29c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full 30c4762a1bSJed Brown 31c4762a1bSJed Brown or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of 32c4762a1bSJed Brown multigrid levels, it will be determined automatically based on the number of refinements done) 33c4762a1bSJed Brown 34c4762a1bSJed Brown ./ex5 -pc_type mg -ksp_monitor -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 35c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full 36c4762a1bSJed Brown 37c4762a1bSJed Brown ------------------------------------------------------------------------- */ 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* 40c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 41c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 42c4762a1bSJed Brown */ 43c4762a1bSJed Brown #include <petscdm.h> 44c4762a1bSJed Brown #include <petscdmda.h> 45c4762a1bSJed Brown #include <petscsnes.h> 46c4762a1bSJed Brown #include <petscmatlab.h> 47c4762a1bSJed Brown #include <petsc/private/snesimpl.h> /* For SNES_Solve event */ 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown User-defined application context - contains data needed by the 51c4762a1bSJed Brown application-provided call-back routines, FormJacobianLocal() and 52c4762a1bSJed Brown FormFunctionLocal(). 53c4762a1bSJed Brown */ 54c4762a1bSJed Brown typedef struct AppCtx AppCtx; 55c4762a1bSJed Brown struct AppCtx { 56c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 57c4762a1bSJed Brown PetscInt m,n; /* MMS3 parameters */ 58c4762a1bSJed Brown PetscErrorCode (*mms_solution)(AppCtx*,const DMDACoor2d*,PetscScalar*); 59c4762a1bSJed Brown PetscErrorCode (*mms_forcing)(AppCtx*,const DMDACoor2d*,PetscScalar*); 60c4762a1bSJed Brown }; 61c4762a1bSJed Brown 62*227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */ 63c4762a1bSJed Brown /* 64*227f9e03SJunchao Zhang FormInitialGuess - Forms initial approximation. 65*227f9e03SJunchao Zhang 66*227f9e03SJunchao Zhang Input Parameters: 67*227f9e03SJunchao Zhang da - The DM 68*227f9e03SJunchao Zhang user - user-defined application context 69*227f9e03SJunchao Zhang 70*227f9e03SJunchao Zhang Output Parameter: 71*227f9e03SJunchao Zhang X - vector 72c4762a1bSJed Brown */ 73*227f9e03SJunchao Zhang static PetscErrorCode FormInitialGuess(DM da,AppCtx *user,Vec X) 74*227f9e03SJunchao Zhang { 75*227f9e03SJunchao Zhang PetscInt i,j,Mx,My,xs,ys,xm,ym; 76*227f9e03SJunchao Zhang PetscReal lambda,temp1,temp,hx,hy; 77*227f9e03SJunchao Zhang PetscScalar **x; 78*227f9e03SJunchao Zhang 79*227f9e03SJunchao Zhang PetscFunctionBeginUser; 80*227f9e03SJunchao Zhang PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 81*227f9e03SJunchao Zhang 82*227f9e03SJunchao Zhang lambda = user->param; 83*227f9e03SJunchao Zhang hx = 1.0/(PetscReal)(Mx-1); 84*227f9e03SJunchao Zhang hy = 1.0/(PetscReal)(My-1); 85*227f9e03SJunchao Zhang temp1 = lambda/(lambda + 1.0); 86*227f9e03SJunchao Zhang 87*227f9e03SJunchao Zhang /* 88*227f9e03SJunchao Zhang Get a pointer to vector data. 89*227f9e03SJunchao Zhang - For default PETSc vectors, VecGetArray() returns a pointer to 90*227f9e03SJunchao Zhang the data array. Otherwise, the routine is implementation dependent. 91*227f9e03SJunchao Zhang - You MUST call VecRestoreArray() when you no longer need access to 92*227f9e03SJunchao Zhang the array. 93*227f9e03SJunchao Zhang */ 94*227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(da,X,&x)); 95*227f9e03SJunchao Zhang 96*227f9e03SJunchao Zhang /* 97*227f9e03SJunchao Zhang Get local grid boundaries (for 2-dimensional DMDA): 98*227f9e03SJunchao Zhang xs, ys - starting grid indices (no ghost points) 99*227f9e03SJunchao Zhang xm, ym - widths of local grid (no ghost points) 100*227f9e03SJunchao Zhang 101*227f9e03SJunchao Zhang */ 102*227f9e03SJunchao Zhang PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 103*227f9e03SJunchao Zhang 104*227f9e03SJunchao Zhang /* 105*227f9e03SJunchao Zhang Compute initial guess over the locally owned part of the grid 106*227f9e03SJunchao Zhang */ 107*227f9e03SJunchao Zhang for (j=ys; j<ys+ym; j++) { 108*227f9e03SJunchao Zhang temp = (PetscReal)(PetscMin(j,My-j-1))*hy; 109*227f9e03SJunchao Zhang for (i=xs; i<xs+xm; i++) { 110*227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 111*227f9e03SJunchao Zhang /* boundary conditions are all zero Dirichlet */ 112*227f9e03SJunchao Zhang x[j][i] = 0.0; 113*227f9e03SJunchao Zhang } else { 114*227f9e03SJunchao Zhang x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp)); 115*227f9e03SJunchao Zhang } 116*227f9e03SJunchao Zhang } 117*227f9e03SJunchao Zhang } 118*227f9e03SJunchao Zhang 119*227f9e03SJunchao Zhang /* 120*227f9e03SJunchao Zhang Restore vector 121*227f9e03SJunchao Zhang */ 122*227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da,X,&x)); 123*227f9e03SJunchao Zhang PetscFunctionReturn(0); 124*227f9e03SJunchao Zhang } 125*227f9e03SJunchao Zhang 126*227f9e03SJunchao Zhang /* 127*227f9e03SJunchao Zhang FormExactSolution - Forms MMS solution 128*227f9e03SJunchao Zhang 129*227f9e03SJunchao Zhang Input Parameters: 130*227f9e03SJunchao Zhang da - The DM 131*227f9e03SJunchao Zhang user - user-defined application context 132*227f9e03SJunchao Zhang 133*227f9e03SJunchao Zhang Output Parameter: 134*227f9e03SJunchao Zhang X - vector 135*227f9e03SJunchao Zhang */ 136*227f9e03SJunchao Zhang static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U) 137*227f9e03SJunchao Zhang { 138*227f9e03SJunchao Zhang DM coordDA; 139*227f9e03SJunchao Zhang Vec coordinates; 140*227f9e03SJunchao Zhang DMDACoor2d **coords; 141*227f9e03SJunchao Zhang PetscScalar **u; 142*227f9e03SJunchao Zhang PetscInt xs, ys, xm, ym, i, j; 143*227f9e03SJunchao Zhang 144*227f9e03SJunchao Zhang PetscFunctionBeginUser; 145*227f9e03SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 146*227f9e03SJunchao Zhang PetscCall(DMGetCoordinateDM(da, &coordDA)); 147*227f9e03SJunchao Zhang PetscCall(DMGetCoordinates(da, &coordinates)); 148*227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 149*227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(da, U, &u)); 150*227f9e03SJunchao Zhang for (j = ys; j < ys+ym; ++j) { 151*227f9e03SJunchao Zhang for (i = xs; i < xs+xm; ++i) { 152*227f9e03SJunchao Zhang user->mms_solution(user,&coords[j][i],&u[j][i]); 153*227f9e03SJunchao Zhang } 154*227f9e03SJunchao Zhang } 155*227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, U, &u)); 156*227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 157*227f9e03SJunchao Zhang PetscFunctionReturn(0); 158*227f9e03SJunchao Zhang } 159*227f9e03SJunchao Zhang 160*227f9e03SJunchao Zhang static PetscErrorCode ZeroBCSolution(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 161*227f9e03SJunchao Zhang { 162*227f9e03SJunchao Zhang u[0] = 0.; 163*227f9e03SJunchao Zhang return 0; 164*227f9e03SJunchao Zhang } 165*227f9e03SJunchao Zhang 166*227f9e03SJunchao Zhang /* The functions below evaluate the MMS solution u(x,y) and associated forcing 167*227f9e03SJunchao Zhang 168*227f9e03SJunchao Zhang f(x,y) = -u_xx - u_yy - lambda exp(u) 169*227f9e03SJunchao Zhang 170*227f9e03SJunchao Zhang such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term. 171*227f9e03SJunchao Zhang */ 172*227f9e03SJunchao Zhang static PetscErrorCode MMSSolution1(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 173*227f9e03SJunchao Zhang { 174*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 175*227f9e03SJunchao Zhang u[0] = x*(1 - x)*y*(1 - y); 176*227f9e03SJunchao Zhang PetscLogFlops(5); 177*227f9e03SJunchao Zhang return 0; 178*227f9e03SJunchao Zhang } 179*227f9e03SJunchao Zhang static PetscErrorCode MMSForcing1(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 180*227f9e03SJunchao Zhang { 181*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 182*227f9e03SJunchao Zhang f[0] = 2*x*(1 - x) + 2*y*(1 - y) - user->param*PetscExpReal(x*(1 - x)*y*(1 - y)); 183*227f9e03SJunchao Zhang return 0; 184*227f9e03SJunchao Zhang } 185*227f9e03SJunchao Zhang 186*227f9e03SJunchao Zhang static PetscErrorCode MMSSolution2(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 187*227f9e03SJunchao Zhang { 188*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 189*227f9e03SJunchao Zhang u[0] = PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y); 190*227f9e03SJunchao Zhang PetscLogFlops(5); 191*227f9e03SJunchao Zhang return 0; 192*227f9e03SJunchao Zhang } 193*227f9e03SJunchao Zhang static PetscErrorCode MMSForcing2(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 194*227f9e03SJunchao Zhang { 195*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 196*227f9e03SJunchao Zhang f[0] = 2*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y) - user->param*PetscExpReal(PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y)); 197*227f9e03SJunchao Zhang return 0; 198*227f9e03SJunchao Zhang } 199*227f9e03SJunchao Zhang 200*227f9e03SJunchao Zhang static PetscErrorCode MMSSolution3(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 201*227f9e03SJunchao Zhang { 202*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 203*227f9e03SJunchao Zhang u[0] = PetscSinReal(user->m*PETSC_PI*x*(1-y))*PetscSinReal(user->n*PETSC_PI*y*(1-x)); 204*227f9e03SJunchao Zhang PetscLogFlops(5); 205*227f9e03SJunchao Zhang return 0; 206*227f9e03SJunchao Zhang } 207*227f9e03SJunchao Zhang static PetscErrorCode MMSForcing3(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 208*227f9e03SJunchao Zhang { 209*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 210*227f9e03SJunchao Zhang PetscReal m = user->m, n = user->n, lambda = user->param; 211*227f9e03SJunchao Zhang f[0] = (-(PetscExpReal(PetscSinReal(m*PETSC_PI*x*(1 - y))*PetscSinReal(n*PETSC_PI*(1 - x)*y))*lambda) 212*227f9e03SJunchao Zhang + PetscSqr(PETSC_PI)*(-2*m*n*((-1 + x)*x + (-1 + y)*y)*PetscCosReal(m*PETSC_PI*x*(-1 + y))*PetscCosReal(n*PETSC_PI*(-1 + x)*y) 213*227f9e03SJunchao Zhang + (PetscSqr(m)*(PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n)*(PetscSqr(-1 + x) + PetscSqr(y))) 214*227f9e03SJunchao Zhang *PetscSinReal(m*PETSC_PI*x*(-1 + y))*PetscSinReal(n*PETSC_PI*(-1 + x)*y))); 215*227f9e03SJunchao Zhang return 0; 216*227f9e03SJunchao Zhang } 217*227f9e03SJunchao Zhang 218*227f9e03SJunchao Zhang static PetscErrorCode MMSSolution4(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 219*227f9e03SJunchao Zhang { 220*227f9e03SJunchao Zhang const PetscReal Lx = 1.,Ly = 1.; 221*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 222*227f9e03SJunchao Zhang u[0] = (PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y)); 223*227f9e03SJunchao Zhang PetscLogFlops(9); 224*227f9e03SJunchao Zhang return 0; 225*227f9e03SJunchao Zhang } 226*227f9e03SJunchao Zhang static PetscErrorCode MMSForcing4(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 227*227f9e03SJunchao Zhang { 228*227f9e03SJunchao Zhang const PetscReal Lx = 1.,Ly = 1.; 229*227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 230*227f9e03SJunchao Zhang f[0] = (2*PetscSqr(x)*(PetscSqr(x)-PetscSqr(Lx))*(PetscSqr(Ly)-6*PetscSqr(y)) 231*227f9e03SJunchao Zhang + 2*PetscSqr(y)*(PetscSqr(Lx)-6*PetscSqr(x))*(PetscSqr(y)-PetscSqr(Ly)) 232*227f9e03SJunchao Zhang - user->param*PetscExpReal((PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y)))); 233*227f9e03SJunchao Zhang return 0; 234*227f9e03SJunchao Zhang } 235*227f9e03SJunchao Zhang 236*227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */ 237*227f9e03SJunchao Zhang /* 238*227f9e03SJunchao Zhang FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch 239*227f9e03SJunchao Zhang 240*227f9e03SJunchao Zhang */ 241*227f9e03SJunchao Zhang static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user) 242*227f9e03SJunchao Zhang { 243*227f9e03SJunchao Zhang PetscInt i,j; 244*227f9e03SJunchao Zhang PetscReal lambda,hx,hy,hxdhy,hydhx; 245*227f9e03SJunchao Zhang PetscScalar u,ue,uw,un,us,uxx,uyy,mms_solution,mms_forcing; 246*227f9e03SJunchao Zhang DMDACoor2d c; 247*227f9e03SJunchao Zhang 248*227f9e03SJunchao Zhang PetscFunctionBeginUser; 249*227f9e03SJunchao Zhang lambda = user->param; 250*227f9e03SJunchao Zhang hx = 1.0/(PetscReal)(info->mx-1); 251*227f9e03SJunchao Zhang hy = 1.0/(PetscReal)(info->my-1); 252*227f9e03SJunchao Zhang hxdhy = hx/hy; 253*227f9e03SJunchao Zhang hydhx = hy/hx; 254*227f9e03SJunchao Zhang /* 255*227f9e03SJunchao Zhang Compute function over the locally owned part of the grid 256*227f9e03SJunchao Zhang */ 257*227f9e03SJunchao Zhang for (j=info->ys; j<info->ys+info->ym; j++) { 258*227f9e03SJunchao Zhang for (i=info->xs; i<info->xs+info->xm; i++) { 259*227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 260*227f9e03SJunchao Zhang c.x = i*hx; c.y = j*hy; 261*227f9e03SJunchao Zhang PetscCall(user->mms_solution(user,&c,&mms_solution)); 262*227f9e03SJunchao Zhang f[j][i] = 2.0*(hydhx+hxdhy)*(x[j][i] - mms_solution); 263*227f9e03SJunchao Zhang } else { 264*227f9e03SJunchao Zhang u = x[j][i]; 265*227f9e03SJunchao Zhang uw = x[j][i-1]; 266*227f9e03SJunchao Zhang ue = x[j][i+1]; 267*227f9e03SJunchao Zhang un = x[j-1][i]; 268*227f9e03SJunchao Zhang us = x[j+1][i]; 269*227f9e03SJunchao Zhang 270*227f9e03SJunchao Zhang /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */ 271*227f9e03SJunchao Zhang if (i-1 == 0) {c.x = (i-1)*hx; c.y = j*hy; PetscCall(user->mms_solution(user,&c,&uw));} 272*227f9e03SJunchao Zhang if (i+1 == info->mx-1) {c.x = (i+1)*hx; c.y = j*hy; PetscCall(user->mms_solution(user,&c,&ue));} 273*227f9e03SJunchao Zhang if (j-1 == 0) {c.x = i*hx; c.y = (j-1)*hy; PetscCall(user->mms_solution(user,&c,&un));} 274*227f9e03SJunchao Zhang if (j+1 == info->my-1) {c.x = i*hx; c.y = (j+1)*hy; PetscCall(user->mms_solution(user,&c,&us));} 275*227f9e03SJunchao Zhang 276*227f9e03SJunchao Zhang uxx = (2.0*u - uw - ue)*hydhx; 277*227f9e03SJunchao Zhang uyy = (2.0*u - un - us)*hxdhy; 278*227f9e03SJunchao Zhang mms_forcing = 0; 279*227f9e03SJunchao Zhang c.x = i*hx; c.y = j*hy; 280*227f9e03SJunchao Zhang if (user->mms_forcing) PetscCall(user->mms_forcing(user,&c,&mms_forcing)); 281*227f9e03SJunchao Zhang f[j][i] = uxx + uyy - hx*hy*(lambda*PetscExpScalar(u) + mms_forcing); 282*227f9e03SJunchao Zhang } 283*227f9e03SJunchao Zhang } 284*227f9e03SJunchao Zhang } 285*227f9e03SJunchao Zhang PetscCall(PetscLogFlops(11.0*info->ym*info->xm)); 286*227f9e03SJunchao Zhang PetscFunctionReturn(0); 287*227f9e03SJunchao Zhang } 288*227f9e03SJunchao Zhang 289*227f9e03SJunchao Zhang /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */ 290*227f9e03SJunchao Zhang static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info,PetscScalar **x,PetscReal *obj,AppCtx *user) 291*227f9e03SJunchao Zhang { 292*227f9e03SJunchao Zhang PetscInt i,j; 293*227f9e03SJunchao Zhang PetscReal lambda,hx,hy,hxdhy,hydhx,sc,lobj=0; 294*227f9e03SJunchao Zhang PetscScalar u,ue,uw,un,us,uxux,uyuy; 295*227f9e03SJunchao Zhang MPI_Comm comm; 296*227f9e03SJunchao Zhang 297*227f9e03SJunchao Zhang PetscFunctionBeginUser; 298*227f9e03SJunchao Zhang *obj = 0; 299*227f9e03SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)info->da,&comm)); 300*227f9e03SJunchao Zhang lambda = user->param; 301*227f9e03SJunchao Zhang hx = 1.0/(PetscReal)(info->mx-1); 302*227f9e03SJunchao Zhang hy = 1.0/(PetscReal)(info->my-1); 303*227f9e03SJunchao Zhang sc = hx*hy*lambda; 304*227f9e03SJunchao Zhang hxdhy = hx/hy; 305*227f9e03SJunchao Zhang hydhx = hy/hx; 306*227f9e03SJunchao Zhang /* 307*227f9e03SJunchao Zhang Compute function over the locally owned part of the grid 308*227f9e03SJunchao Zhang */ 309*227f9e03SJunchao Zhang for (j=info->ys; j<info->ys+info->ym; j++) { 310*227f9e03SJunchao Zhang for (i=info->xs; i<info->xs+info->xm; i++) { 311*227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 312*227f9e03SJunchao Zhang lobj += PetscRealPart((hydhx + hxdhy)*x[j][i]*x[j][i]); 313*227f9e03SJunchao Zhang } else { 314*227f9e03SJunchao Zhang u = x[j][i]; 315*227f9e03SJunchao Zhang uw = x[j][i-1]; 316*227f9e03SJunchao Zhang ue = x[j][i+1]; 317*227f9e03SJunchao Zhang un = x[j-1][i]; 318*227f9e03SJunchao Zhang us = x[j+1][i]; 319*227f9e03SJunchao Zhang 320*227f9e03SJunchao Zhang if (i-1 == 0) uw = 0.; 321*227f9e03SJunchao Zhang if (i+1 == info->mx-1) ue = 0.; 322*227f9e03SJunchao Zhang if (j-1 == 0) un = 0.; 323*227f9e03SJunchao Zhang if (j+1 == info->my-1) us = 0.; 324*227f9e03SJunchao Zhang 325*227f9e03SJunchao Zhang /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */ 326*227f9e03SJunchao Zhang 327*227f9e03SJunchao Zhang uxux = u*(2.*u - ue - uw)*hydhx; 328*227f9e03SJunchao Zhang uyuy = u*(2.*u - un - us)*hxdhy; 329*227f9e03SJunchao Zhang 330*227f9e03SJunchao Zhang lobj += PetscRealPart(0.5*(uxux + uyuy) - sc*PetscExpScalar(u)); 331*227f9e03SJunchao Zhang } 332*227f9e03SJunchao Zhang } 333*227f9e03SJunchao Zhang } 334*227f9e03SJunchao Zhang PetscCall(PetscLogFlops(12.0*info->ym*info->xm)); 335*227f9e03SJunchao Zhang PetscCallMPI(MPI_Allreduce(&lobj,obj,1,MPIU_REAL,MPIU_SUM,comm)); 336*227f9e03SJunchao Zhang PetscFunctionReturn(0); 337*227f9e03SJunchao Zhang } 338*227f9e03SJunchao Zhang 339*227f9e03SJunchao Zhang /* 340*227f9e03SJunchao Zhang FormJacobianLocal - Evaluates Jacobian matrix on local process patch 341*227f9e03SJunchao Zhang */ 342*227f9e03SJunchao Zhang static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,Mat jacpre,AppCtx *user) 343*227f9e03SJunchao Zhang { 344*227f9e03SJunchao Zhang PetscInt i,j,k; 345*227f9e03SJunchao Zhang MatStencil col[5],row; 346*227f9e03SJunchao Zhang PetscScalar lambda,v[5],hx,hy,hxdhy,hydhx,sc; 347*227f9e03SJunchao Zhang DM coordDA; 348*227f9e03SJunchao Zhang Vec coordinates; 349*227f9e03SJunchao Zhang DMDACoor2d **coords; 350*227f9e03SJunchao Zhang 351*227f9e03SJunchao Zhang PetscFunctionBeginUser; 352*227f9e03SJunchao Zhang lambda = user->param; 353*227f9e03SJunchao Zhang /* Extract coordinates */ 354*227f9e03SJunchao Zhang PetscCall(DMGetCoordinateDM(info->da, &coordDA)); 355*227f9e03SJunchao Zhang PetscCall(DMGetCoordinates(info->da, &coordinates)); 356*227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 357*227f9e03SJunchao Zhang hx = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs+1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0; 358*227f9e03SJunchao Zhang hy = info->ym > 1 ? PetscRealPart(coords[info->ys+1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0; 359*227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 360*227f9e03SJunchao Zhang hxdhy = hx/hy; 361*227f9e03SJunchao Zhang hydhx = hy/hx; 362*227f9e03SJunchao Zhang sc = hx*hy*lambda; 363*227f9e03SJunchao Zhang 364*227f9e03SJunchao Zhang /* 365*227f9e03SJunchao Zhang Compute entries for the locally owned part of the Jacobian. 366*227f9e03SJunchao Zhang - Currently, all PETSc parallel matrix formats are partitioned by 367*227f9e03SJunchao Zhang contiguous chunks of rows across the processors. 368*227f9e03SJunchao Zhang - Each processor needs to insert only elements that it owns 369*227f9e03SJunchao Zhang locally (but any non-local elements will be sent to the 370*227f9e03SJunchao Zhang appropriate processor during matrix assembly). 371*227f9e03SJunchao Zhang - Here, we set all entries for a particular row at once. 372*227f9e03SJunchao Zhang - We can set matrix entries either using either 373*227f9e03SJunchao Zhang MatSetValuesLocal() or MatSetValues(), as discussed above. 374*227f9e03SJunchao Zhang */ 375*227f9e03SJunchao Zhang for (j=info->ys; j<info->ys+info->ym; j++) { 376*227f9e03SJunchao Zhang for (i=info->xs; i<info->xs+info->xm; i++) { 377*227f9e03SJunchao Zhang row.j = j; row.i = i; 378*227f9e03SJunchao Zhang /* boundary points */ 379*227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 380*227f9e03SJunchao Zhang v[0] = 2.0*(hydhx + hxdhy); 381*227f9e03SJunchao Zhang PetscCall(MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES)); 382*227f9e03SJunchao Zhang } else { 383*227f9e03SJunchao Zhang k = 0; 384*227f9e03SJunchao Zhang /* interior grid points */ 385*227f9e03SJunchao Zhang if (j-1 != 0) { 386*227f9e03SJunchao Zhang v[k] = -hxdhy; 387*227f9e03SJunchao Zhang col[k].j = j - 1; col[k].i = i; 388*227f9e03SJunchao Zhang k++; 389*227f9e03SJunchao Zhang } 390*227f9e03SJunchao Zhang if (i-1 != 0) { 391*227f9e03SJunchao Zhang v[k] = -hydhx; 392*227f9e03SJunchao Zhang col[k].j = j; col[k].i = i-1; 393*227f9e03SJunchao Zhang k++; 394*227f9e03SJunchao Zhang } 395*227f9e03SJunchao Zhang 396*227f9e03SJunchao Zhang v[k] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[k].j = row.j; col[k].i = row.i; k++; 397*227f9e03SJunchao Zhang 398*227f9e03SJunchao Zhang if (i+1 != info->mx-1) { 399*227f9e03SJunchao Zhang v[k] = -hydhx; 400*227f9e03SJunchao Zhang col[k].j = j; col[k].i = i+1; 401*227f9e03SJunchao Zhang k++; 402*227f9e03SJunchao Zhang } 403*227f9e03SJunchao Zhang if (j+1 != info->mx-1) { 404*227f9e03SJunchao Zhang v[k] = -hxdhy; 405*227f9e03SJunchao Zhang col[k].j = j + 1; col[k].i = i; 406*227f9e03SJunchao Zhang k++; 407*227f9e03SJunchao Zhang } 408*227f9e03SJunchao Zhang PetscCall(MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES)); 409*227f9e03SJunchao Zhang } 410*227f9e03SJunchao Zhang } 411*227f9e03SJunchao Zhang } 412*227f9e03SJunchao Zhang 413*227f9e03SJunchao Zhang /* 414*227f9e03SJunchao Zhang Assemble matrix, using the 2-step process: 415*227f9e03SJunchao Zhang MatAssemblyBegin(), MatAssemblyEnd(). 416*227f9e03SJunchao Zhang */ 417*227f9e03SJunchao Zhang PetscCall(MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY)); 418*227f9e03SJunchao Zhang PetscCall(MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY)); 419*227f9e03SJunchao Zhang /* 420*227f9e03SJunchao Zhang Tell the matrix we will never add a new nonzero location to the 421*227f9e03SJunchao Zhang matrix. If we do, it will generate an error. 422*227f9e03SJunchao Zhang */ 423*227f9e03SJunchao Zhang PetscCall(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 424*227f9e03SJunchao Zhang PetscFunctionReturn(0); 425*227f9e03SJunchao Zhang } 426*227f9e03SJunchao Zhang 427*227f9e03SJunchao Zhang static PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr) 428*227f9e03SJunchao Zhang { 429*227f9e03SJunchao Zhang #if PetscDefined(HAVE_MATLAB_ENGINE) 430*227f9e03SJunchao Zhang AppCtx *user = (AppCtx*)ptr; 431*227f9e03SJunchao Zhang PetscInt Mx,My; 432*227f9e03SJunchao Zhang PetscReal lambda,hx,hy; 433*227f9e03SJunchao Zhang Vec localX,localF; 434*227f9e03SJunchao Zhang MPI_Comm comm; 435*227f9e03SJunchao Zhang DM da; 436*227f9e03SJunchao Zhang 437*227f9e03SJunchao Zhang PetscFunctionBeginUser; 438*227f9e03SJunchao Zhang PetscCall(SNESGetDM(snes,&da)); 439*227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da,&localX)); 440*227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da,&localF)); 441*227f9e03SJunchao Zhang PetscCall(PetscObjectSetName((PetscObject)localX,"localX")); 442*227f9e03SJunchao Zhang PetscCall(PetscObjectSetName((PetscObject)localF,"localF")); 443*227f9e03SJunchao Zhang PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 444*227f9e03SJunchao Zhang 445*227f9e03SJunchao Zhang lambda = user->param; 446*227f9e03SJunchao Zhang hx = 1.0/(PetscReal)(Mx-1); 447*227f9e03SJunchao Zhang hy = 1.0/(PetscReal)(My-1); 448*227f9e03SJunchao Zhang 449*227f9e03SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 450*227f9e03SJunchao Zhang /* 451*227f9e03SJunchao Zhang Scatter ghost points to local vector,using the 2-step process 452*227f9e03SJunchao Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 453*227f9e03SJunchao Zhang By placing code between these two statements, computations can be 454*227f9e03SJunchao Zhang done while messages are in transition. 455*227f9e03SJunchao Zhang */ 456*227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 457*227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 458*227f9e03SJunchao Zhang PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX)); 459*227f9e03SJunchao Zhang PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",(double)hx,(double)hy,(double)lambda)); 460*227f9e03SJunchao Zhang PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF)); 461*227f9e03SJunchao Zhang 462*227f9e03SJunchao Zhang /* 463*227f9e03SJunchao Zhang Insert values into global vector 464*227f9e03SJunchao Zhang */ 465*227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalBegin(da,localF,INSERT_VALUES,F)); 466*227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalEnd(da,localF,INSERT_VALUES,F)); 467*227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da,&localX)); 468*227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da,&localF)); 469*227f9e03SJunchao Zhang PetscFunctionReturn(0); 470*227f9e03SJunchao Zhang #else 471*227f9e03SJunchao Zhang return 0; /* Never called */ 472*227f9e03SJunchao Zhang #endif 473*227f9e03SJunchao Zhang } 474*227f9e03SJunchao Zhang 475*227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */ 476*227f9e03SJunchao Zhang /* 477*227f9e03SJunchao Zhang Applies some sweeps on nonlinear Gauss-Seidel on each process 478*227f9e03SJunchao Zhang 479*227f9e03SJunchao Zhang */ 480*227f9e03SJunchao Zhang static PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx) 481*227f9e03SJunchao Zhang { 482*227f9e03SJunchao Zhang PetscInt i,j,k,Mx,My,xs,ys,xm,ym,its,tot_its,sweeps,l; 483*227f9e03SJunchao Zhang PetscReal lambda,hx,hy,hxdhy,hydhx,sc; 484*227f9e03SJunchao Zhang PetscScalar **x,**b,bij,F,F0=0,J,u,un,us,ue,eu,uw,uxx,uyy,y; 485*227f9e03SJunchao Zhang PetscReal atol,rtol,stol; 486*227f9e03SJunchao Zhang DM da; 487*227f9e03SJunchao Zhang AppCtx *user; 488*227f9e03SJunchao Zhang Vec localX,localB; 489*227f9e03SJunchao Zhang 490*227f9e03SJunchao Zhang PetscFunctionBeginUser; 491*227f9e03SJunchao Zhang tot_its = 0; 492*227f9e03SJunchao Zhang PetscCall(SNESNGSGetSweeps(snes,&sweeps)); 493*227f9e03SJunchao Zhang PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its)); 494*227f9e03SJunchao Zhang PetscCall(SNESGetDM(snes,&da)); 495*227f9e03SJunchao Zhang PetscCall(DMGetApplicationContext(da,&user)); 496*227f9e03SJunchao Zhang 497*227f9e03SJunchao Zhang PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 498*227f9e03SJunchao Zhang 499*227f9e03SJunchao Zhang lambda = user->param; 500*227f9e03SJunchao Zhang hx = 1.0/(PetscReal)(Mx-1); 501*227f9e03SJunchao Zhang hy = 1.0/(PetscReal)(My-1); 502*227f9e03SJunchao Zhang sc = hx*hy*lambda; 503*227f9e03SJunchao Zhang hxdhy = hx/hy; 504*227f9e03SJunchao Zhang hydhx = hy/hx; 505*227f9e03SJunchao Zhang 506*227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da,&localX)); 507*227f9e03SJunchao Zhang if (B) { 508*227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da,&localB)); 509*227f9e03SJunchao Zhang } 510*227f9e03SJunchao Zhang for (l=0; l<sweeps; l++) { 511*227f9e03SJunchao Zhang 512*227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 513*227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 514*227f9e03SJunchao Zhang if (B) { 515*227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB)); 516*227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB)); 517*227f9e03SJunchao Zhang } 518*227f9e03SJunchao Zhang /* 519*227f9e03SJunchao Zhang Get a pointer to vector data. 520*227f9e03SJunchao Zhang - For default PETSc vectors, VecGetArray() returns a pointer to 521*227f9e03SJunchao Zhang the data array. Otherwise, the routine is implementation dependent. 522*227f9e03SJunchao Zhang - You MUST call VecRestoreArray() when you no longer need access to 523*227f9e03SJunchao Zhang the array. 524*227f9e03SJunchao Zhang */ 525*227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(da,localX,&x)); 526*227f9e03SJunchao Zhang if (B) PetscCall(DMDAVecGetArray(da,localB,&b)); 527*227f9e03SJunchao Zhang /* 528*227f9e03SJunchao Zhang Get local grid boundaries (for 2-dimensional DMDA): 529*227f9e03SJunchao Zhang xs, ys - starting grid indices (no ghost points) 530*227f9e03SJunchao Zhang xm, ym - widths of local grid (no ghost points) 531*227f9e03SJunchao Zhang */ 532*227f9e03SJunchao Zhang PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 533*227f9e03SJunchao Zhang 534*227f9e03SJunchao Zhang for (j=ys; j<ys+ym; j++) { 535*227f9e03SJunchao Zhang for (i=xs; i<xs+xm; i++) { 536*227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 537*227f9e03SJunchao Zhang /* boundary conditions are all zero Dirichlet */ 538*227f9e03SJunchao Zhang x[j][i] = 0.0; 539*227f9e03SJunchao Zhang } else { 540*227f9e03SJunchao Zhang if (B) bij = b[j][i]; 541*227f9e03SJunchao Zhang else bij = 0.; 542*227f9e03SJunchao Zhang 543*227f9e03SJunchao Zhang u = x[j][i]; 544*227f9e03SJunchao Zhang un = x[j-1][i]; 545*227f9e03SJunchao Zhang us = x[j+1][i]; 546*227f9e03SJunchao Zhang ue = x[j][i-1]; 547*227f9e03SJunchao Zhang uw = x[j][i+1]; 548*227f9e03SJunchao Zhang 549*227f9e03SJunchao Zhang for (k=0; k<its; k++) { 550*227f9e03SJunchao Zhang eu = PetscExpScalar(u); 551*227f9e03SJunchao Zhang uxx = (2.0*u - ue - uw)*hydhx; 552*227f9e03SJunchao Zhang uyy = (2.0*u - un - us)*hxdhy; 553*227f9e03SJunchao Zhang F = uxx + uyy - sc*eu - bij; 554*227f9e03SJunchao Zhang if (k == 0) F0 = F; 555*227f9e03SJunchao Zhang J = 2.0*(hydhx + hxdhy) - sc*eu; 556*227f9e03SJunchao Zhang y = F/J; 557*227f9e03SJunchao Zhang u -= y; 558*227f9e03SJunchao Zhang tot_its++; 559*227f9e03SJunchao Zhang 560*227f9e03SJunchao Zhang if (atol > PetscAbsReal(PetscRealPart(F)) || 561*227f9e03SJunchao Zhang rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || 562*227f9e03SJunchao Zhang stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) { 563*227f9e03SJunchao Zhang break; 564*227f9e03SJunchao Zhang } 565*227f9e03SJunchao Zhang } 566*227f9e03SJunchao Zhang x[j][i] = u; 567*227f9e03SJunchao Zhang } 568*227f9e03SJunchao Zhang } 569*227f9e03SJunchao Zhang } 570*227f9e03SJunchao Zhang /* 571*227f9e03SJunchao Zhang Restore vector 572*227f9e03SJunchao Zhang */ 573*227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da,localX,&x)); 574*227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X)); 575*227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X)); 576*227f9e03SJunchao Zhang } 577*227f9e03SJunchao Zhang PetscCall(PetscLogFlops(tot_its*(21.0))); 578*227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da,&localX)); 579*227f9e03SJunchao Zhang if (B) { 580*227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da,localB,&b)); 581*227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da,&localB)); 582*227f9e03SJunchao Zhang } 583*227f9e03SJunchao Zhang PetscFunctionReturn(0); 584*227f9e03SJunchao Zhang } 585c4762a1bSJed Brown 586c4762a1bSJed Brown int main(int argc,char **argv) 587c4762a1bSJed Brown { 588c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 589c4762a1bSJed Brown Vec x; /* solution vector */ 590c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 591c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 592c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81; 593c4762a1bSJed Brown PetscReal bratu_lambda_min = 0.; 594c4762a1bSJed Brown PetscInt MMS = 0; 595c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 596c4762a1bSJed Brown DM da; 597c4762a1bSJed Brown Vec r = NULL; 598c4762a1bSJed Brown KSP ksp; 599c4762a1bSJed Brown PetscInt lits,slits; 600c4762a1bSJed Brown 601c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 602c4762a1bSJed Brown Initialize program 603c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 604c4762a1bSJed Brown 6059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 606c4762a1bSJed Brown 607c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 608c4762a1bSJed Brown Initialize problem parameters 609c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 610c4762a1bSJed Brown user.param = 6.0; 6119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL)); 612e00437b9SBarry Smith PetscCheck(user.param <= bratu_lambda_max && user.param >= bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda, %g, is out of range, [%g, %g]", (double)user.param, (double)bratu_lambda_min, (double)bratu_lambda_max); 6139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mms",&MMS,NULL)); 614c4762a1bSJed Brown if (MMS == 3) { 615c4762a1bSJed Brown PetscInt mPar = 2, nPar = 1; 6169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_par",&mPar,NULL)); 6179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n_par",&nPar,NULL)); 618c4762a1bSJed Brown user.m = PetscPowInt(2,mPar); 619c4762a1bSJed Brown user.n = PetscPowInt(2,nPar); 620c4762a1bSJed Brown } 621c4762a1bSJed Brown 622c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 623c4762a1bSJed Brown Create nonlinear solver context 624c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6259566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 6269566063dSJacob Faibussowitsch PetscCall(SNESSetCountersReset(snes,PETSC_FALSE)); 6279566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes, NonlinearGS, NULL)); 628c4762a1bSJed Brown 629c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 630c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 631c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6329566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); 6339566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 6349566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 6359566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 6369566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da,&user)); 6379566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes,da)); 638c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 639c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 640c4762a1bSJed Brown vectors that are the same types 641c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6429566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 643c4762a1bSJed Brown 644c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 645c4762a1bSJed Brown Set local function evaluation routine 646c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 647c4762a1bSJed Brown user.mms_solution = ZeroBCSolution; 648c4762a1bSJed Brown switch (MMS) { 6492f613bf5SBarry Smith case 0: user.mms_solution = NULL; user.mms_forcing = NULL; 650c4762a1bSJed Brown case 1: user.mms_solution = MMSSolution1; user.mms_forcing = MMSForcing1; break; 651c4762a1bSJed Brown case 2: user.mms_solution = MMSSolution2; user.mms_forcing = MMSForcing2; break; 652c4762a1bSJed Brown case 3: user.mms_solution = MMSSolution3; user.mms_forcing = MMSForcing3; break; 653c4762a1bSJed Brown case 4: user.mms_solution = MMSSolution4; user.mms_forcing = MMSForcing4; break; 65463a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Unknown MMS type %" PetscInt_FMT,MMS); 655c4762a1bSJed Brown } 6569566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,&user)); 6579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-fd",&flg,NULL)); 658c4762a1bSJed Brown if (!flg) { 6599566063dSJacob Faibussowitsch PetscCall(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user)); 660c4762a1bSJed Brown } 661c4762a1bSJed Brown 6629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-obj",&flg,NULL)); 663c4762a1bSJed Brown if (flg) { 6649566063dSJacob Faibussowitsch PetscCall(DMDASNESSetObjectiveLocal(da,(DMDASNESObjective)FormObjectiveLocal,&user)); 665c4762a1bSJed Brown } 666c4762a1bSJed Brown 6674e1ad211SJed Brown if (PetscDefined(HAVE_MATLAB_ENGINE)) { 6684e1ad211SJed Brown PetscBool matlab_function = PETSC_FALSE; 6699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-matlab_function",&matlab_function,0)); 670c4762a1bSJed Brown if (matlab_function) { 6719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 6729566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunctionMatlab,&user)); 673c4762a1bSJed Brown } 6744e1ad211SJed Brown } 675c4762a1bSJed Brown 676c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 677c4762a1bSJed Brown Customize nonlinear solver; set runtime options 678c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6799566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 680c4762a1bSJed Brown 6819566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(da,&user,x)); 682c4762a1bSJed Brown 683c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 684c4762a1bSJed Brown Solve nonlinear system 685c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6869566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 6879566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 688c4762a1bSJed Brown 6899566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes,&slits)); 6909566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 6919566063dSJacob Faibussowitsch PetscCall(KSPGetTotalIterations(ksp,&lits)); 69263a3b9bcSJacob Faibussowitsch PetscCheck(lits == slits,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Number of total linear iterations reported by SNES %" PetscInt_FMT " does not match reported by KSP %" PetscInt_FMT,slits,lits); 693c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 694c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 695c4762a1bSJed Brown 696c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 697c4762a1bSJed Brown If using MMS, check the l_2 error 698c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 699c4762a1bSJed Brown if (MMS) { 700c4762a1bSJed Brown Vec e; 701c4762a1bSJed Brown PetscReal errorl2, errorinf; 702c4762a1bSJed Brown PetscInt N; 703c4762a1bSJed Brown 7049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &e)); 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) x, NULL, "-sol_view")); 7069566063dSJacob Faibussowitsch PetscCall(FormExactSolution(da, &user, e)); 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) e, NULL, "-exact_view")); 7089566063dSJacob Faibussowitsch PetscCall(VecAXPY(e, -1.0, x)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) e, NULL, "-error_view")); 7109566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_2, &errorl2)); 7119566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &errorinf)); 7129566063dSJacob Faibussowitsch PetscCall(VecGetSize(e, &N)); 71363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2/PetscSqrtReal((PetscReal)N)), (double) errorinf)); 7149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 7159566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N)); 7169566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2/PetscSqrtReal(N))); 717c4762a1bSJed Brown } 718c4762a1bSJed Brown 719c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 720c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 721c4762a1bSJed Brown are no longer needed. 722c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 7239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 7249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 7259566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 7269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 7279566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 728b122ec5aSJacob Faibussowitsch return 0; 729c4762a1bSJed Brown } 730c4762a1bSJed Brown 731c4762a1bSJed Brown /*TEST 732c4762a1bSJed Brown 733c4762a1bSJed Brown test: 734c4762a1bSJed Brown suffix: asm_0 735c4762a1bSJed Brown requires: !single 736c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu 737c4762a1bSJed Brown 738c4762a1bSJed Brown test: 739c4762a1bSJed Brown suffix: msm_0 740c4762a1bSJed Brown requires: !single 741c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu 742c4762a1bSJed Brown 743c4762a1bSJed Brown test: 744c4762a1bSJed Brown suffix: asm_1 745c4762a1bSJed Brown requires: !single 746c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 747c4762a1bSJed Brown 748c4762a1bSJed Brown test: 749c4762a1bSJed Brown suffix: msm_1 750c4762a1bSJed Brown requires: !single 751c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 752c4762a1bSJed Brown 753c4762a1bSJed Brown test: 754c4762a1bSJed Brown suffix: asm_2 755c4762a1bSJed Brown requires: !single 756c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 757c4762a1bSJed Brown 758c4762a1bSJed Brown test: 759c4762a1bSJed Brown suffix: msm_2 760c4762a1bSJed Brown requires: !single 761c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 762c4762a1bSJed Brown 763c4762a1bSJed Brown test: 764c4762a1bSJed Brown suffix: asm_3 765c4762a1bSJed Brown requires: !single 766c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 767c4762a1bSJed Brown 768c4762a1bSJed Brown test: 769c4762a1bSJed Brown suffix: msm_3 770c4762a1bSJed Brown requires: !single 771c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 772c4762a1bSJed Brown 773c4762a1bSJed Brown test: 774c4762a1bSJed Brown suffix: asm_4 775c4762a1bSJed Brown requires: !single 776c4762a1bSJed Brown nsize: 2 777c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 778c4762a1bSJed Brown 779c4762a1bSJed Brown test: 780c4762a1bSJed Brown suffix: msm_4 781c4762a1bSJed Brown requires: !single 782c4762a1bSJed Brown nsize: 2 783c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 784c4762a1bSJed Brown 785c4762a1bSJed Brown test: 786c4762a1bSJed Brown suffix: asm_5 787c4762a1bSJed Brown requires: !single 788c4762a1bSJed Brown nsize: 2 789c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 790c4762a1bSJed Brown 791c4762a1bSJed Brown test: 792c4762a1bSJed Brown suffix: msm_5 793c4762a1bSJed Brown requires: !single 794c4762a1bSJed Brown nsize: 2 795c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 796c4762a1bSJed Brown 797c4762a1bSJed Brown test: 798c4762a1bSJed Brown args: -snes_rtol 1.e-5 -pc_type mg -ksp_monitor_short -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 -mg_levels_ksp_monitor_short -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full 799c4762a1bSJed Brown requires: !single 800c4762a1bSJed Brown 801c4762a1bSJed Brown test: 802c4762a1bSJed Brown suffix: 2 803c4762a1bSJed Brown args: -pc_type mg -ksp_converged_reason -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full -ksp_atol -1. 804c4762a1bSJed Brown 805c4762a1bSJed Brown test: 806c4762a1bSJed Brown suffix: 3 807c4762a1bSJed Brown nsize: 2 808c4762a1bSJed Brown args: -snes_grid_sequence 2 -snes_mf_operator -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -snes_rtol 1.e-2 809c4762a1bSJed Brown filter: grep -v "otal number of function evaluations" 810c4762a1bSJed Brown 811c4762a1bSJed Brown test: 812c4762a1bSJed Brown suffix: 4 813c4762a1bSJed Brown nsize: 2 814c4762a1bSJed Brown args: -snes_grid_sequence 2 -snes_monitor_short -ksp_converged_reason -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -ksp_atol -1 815c4762a1bSJed Brown 816c4762a1bSJed Brown test: 817c4762a1bSJed Brown suffix: 5_anderson 818c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type anderson 819c4762a1bSJed Brown 820c4762a1bSJed Brown test: 821c4762a1bSJed Brown suffix: 5_aspin 822c4762a1bSJed Brown nsize: 4 823c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type aspin -snes_view 824c4762a1bSJed Brown 825c4762a1bSJed Brown test: 826c4762a1bSJed Brown suffix: 5_broyden 827c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_qn_type broyden -snes_qn_m 10 828c4762a1bSJed Brown 829c4762a1bSJed Brown test: 830c4762a1bSJed Brown suffix: 5_fas 831c4762a1bSJed Brown args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 832c4762a1bSJed Brown requires: !single 833c4762a1bSJed Brown 834c4762a1bSJed Brown test: 835c4762a1bSJed Brown suffix: 5_fas_additive 836c4762a1bSJed Brown args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 -snes_fas_type additive -snes_max_it 50 837c4762a1bSJed Brown 838c4762a1bSJed Brown test: 839c4762a1bSJed Brown suffix: 5_fas_monitor 840c4762a1bSJed Brown args: -da_refine 1 -snes_type fas -snes_fas_monitor 841c4762a1bSJed Brown requires: !single 842c4762a1bSJed Brown 843c4762a1bSJed Brown test: 844c4762a1bSJed Brown suffix: 5_ls 845c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls 846c4762a1bSJed Brown 847c4762a1bSJed Brown test: 848c4762a1bSJed Brown suffix: 5_ls_sell_sor 849c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls -dm_mat_type sell -pc_type sor 850c4762a1bSJed Brown output_file: output/ex5_5_ls.out 851c4762a1bSJed Brown 852c4762a1bSJed Brown test: 853c4762a1bSJed Brown suffix: 5_nasm 854c4762a1bSJed Brown nsize: 4 855c4762a1bSJed Brown args: -snes_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type nasm -snes_nasm_type restrict -snes_max_it 10 856c4762a1bSJed Brown 857c4762a1bSJed Brown test: 858c4762a1bSJed Brown suffix: 5_ncg 859c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ncg -snes_ncg_type fr 860c4762a1bSJed Brown 861c4762a1bSJed Brown test: 862c4762a1bSJed Brown suffix: 5_newton_asm_dmda 863c4762a1bSJed Brown nsize: 4 864c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type asm -pc_asm_dm_subdomains -malloc_dump 865c4762a1bSJed Brown requires: !single 866c4762a1bSJed Brown 867c4762a1bSJed Brown test: 868c4762a1bSJed Brown suffix: 5_newton_gasm_dmda 869c4762a1bSJed Brown nsize: 4 870c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type gasm -malloc_dump 871c4762a1bSJed Brown requires: !single 872c4762a1bSJed Brown 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown suffix: 5_ngmres 875c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 876c4762a1bSJed Brown 877c4762a1bSJed Brown test: 878c4762a1bSJed Brown suffix: 5_ngmres_fas 879c4762a1bSJed Brown args: -snes_rtol 1.e-4 -snes_type ngmres -npc_fas_coarse_snes_max_it 1 -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_pc_type lu -npc_fas_coarse_ksp_type preonly -snes_ngmres_m 10 -snes_monitor_short -npc_snes_max_it 1 -npc_snes_type fas -npc_fas_coarse_ksp_type richardson -da_refine 6 880c4762a1bSJed Brown 881c4762a1bSJed Brown test: 882c4762a1bSJed Brown suffix: 5_ngmres_ngs 883c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -npc_snes_type ngs -npc_snes_max_it 1 884c4762a1bSJed Brown 885c4762a1bSJed Brown test: 886c4762a1bSJed Brown suffix: 5_ngmres_nrichardson 887c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 -npc_snes_type nrichardson -npc_snes_max_it 3 888c4762a1bSJed Brown output_file: output/ex5_5_ngmres_richardson.out 889c4762a1bSJed Brown 890c4762a1bSJed Brown test: 891c4762a1bSJed Brown suffix: 5_nrichardson 892c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type nrichardson 893c4762a1bSJed Brown 894c4762a1bSJed Brown test: 895c4762a1bSJed Brown suffix: 5_qn 896c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_linesearch_type cp -snes_qn_m 10 897c4762a1bSJed Brown 898c4762a1bSJed Brown test: 899c4762a1bSJed Brown suffix: 6 900c4762a1bSJed Brown nsize: 4 901c4762a1bSJed Brown args: -snes_converged_reason -ksp_converged_reason -da_grid_x 129 -da_grid_y 129 -pc_type mg -pc_mg_levels 8 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.5,0,1.1 -mg_levels_ksp_max_it 2 902c4762a1bSJed Brown 903c4762a1bSJed Brown test: 904c4762a1bSJed Brown requires: complex !single 905c4762a1bSJed Brown suffix: complex 906c4762a1bSJed Brown args: -snes_mf_operator -mat_mffd_complex -snes_monitor 907c4762a1bSJed Brown 908c4762a1bSJed Brown TEST*/ 909