1c4762a1bSJed Brown static const char help[] = "p-Bratu nonlinear PDE in 2d.\n\ 2c4762a1bSJed Brown We solve the p-Laplacian (nonlinear diffusion) combined with\n\ 3c4762a1bSJed Brown the Bratu (solid fuel ignition) nonlinearity in a 2D rectangular\n\ 4c4762a1bSJed Brown domain, using distributed arrays (DAs) to partition the parallel grid.\n\ 5c4762a1bSJed Brown The command line options include:\n\ 6c4762a1bSJed Brown -p <2>: `p' in p-Laplacian term\n\ 7c4762a1bSJed Brown -epsilon <1e-05>: Strain-regularization in p-Laplacian\n\ 8c4762a1bSJed Brown -lambda <6>: Bratu parameter\n\ 9c4762a1bSJed Brown -blocks <bx,by>: number of coefficient interfaces in x and y direction\n\ 10c4762a1bSJed Brown -kappa <1e-3>: diffusivity in odd regions\n\ 11c4762a1bSJed Brown \n"; 12c4762a1bSJed Brown 13c4762a1bSJed Brown /*F 14c4762a1bSJed Brown The $p$-Bratu problem is a combination of the $p$-Laplacian (nonlinear diffusion) and the Brutu solid fuel ignition problem. 15c4762a1bSJed Brown This problem is modeled by the partial differential equation 16c4762a1bSJed Brown 17c4762a1bSJed Brown \begin{equation*} 18c4762a1bSJed Brown -\nabla\cdot (\eta \nabla u) - \lambda \exp(u) = 0 19c4762a1bSJed Brown \end{equation*} 20c4762a1bSJed Brown 21c4762a1bSJed Brown on $\Omega = (-1,1)^2$ with closure 22c4762a1bSJed Brown 23c4762a1bSJed Brown \begin{align*} 24c4762a1bSJed Brown \eta(\gamma) &= (\epsilon^2 + \gamma)^{(p-2)/2} & \gamma &= \frac 1 2 |\nabla u|^2 25c4762a1bSJed Brown \end{align*} 26c4762a1bSJed Brown 27c4762a1bSJed Brown and boundary conditions $u = 0$ for $(x,y) \in \partial \Omega$ 28c4762a1bSJed Brown 29c4762a1bSJed Brown A 9-point finite difference stencil is used to discretize 30c4762a1bSJed Brown the boundary value problem to obtain a nonlinear system of equations. 31c4762a1bSJed Brown This would be a 5-point stencil if not for the $p$-Laplacian's nonlinearity. 32c4762a1bSJed Brown F*/ 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown mpiexec -n 2 ./ex15 -snes_monitor -ksp_monitor log_summary 36c4762a1bSJed Brown */ 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* 39c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 40c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 41c4762a1bSJed Brown file automatically includes: 42c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors 43c4762a1bSJed Brown petscsys.h - system routines petscmat.h - matrices 44c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 45c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 46c4762a1bSJed Brown petscksp.h - linear solvers 47c4762a1bSJed Brown */ 48c4762a1bSJed Brown 49c4762a1bSJed Brown #include <petscdm.h> 50c4762a1bSJed Brown #include <petscdmda.h> 51c4762a1bSJed Brown #include <petscsnes.h> 52c4762a1bSJed Brown 53c4762a1bSJed Brown typedef enum {JAC_BRATU,JAC_PICARD,JAC_STAR,JAC_NEWTON} JacType; 54c4762a1bSJed Brown static const char *const JacTypes[] = {"BRATU","PICARD","STAR","NEWTON","JacType","JAC_",0}; 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* 57c4762a1bSJed Brown User-defined application context - contains data needed by the 58c4762a1bSJed Brown application-provided call-back routines, FormJacobianLocal() and 59c4762a1bSJed Brown FormFunctionLocal(). 60c4762a1bSJed Brown */ 61c4762a1bSJed Brown typedef struct { 62c4762a1bSJed Brown PetscReal lambda; /* Bratu parameter */ 63c4762a1bSJed Brown PetscReal p; /* Exponent in p-Laplacian */ 64c4762a1bSJed Brown PetscReal epsilon; /* Regularization */ 65c4762a1bSJed Brown PetscReal source; /* Source term */ 66c4762a1bSJed Brown JacType jtype; /* What type of Jacobian to assemble */ 67c4762a1bSJed Brown PetscBool picard; 68c4762a1bSJed Brown PetscInt blocks[2]; 69c4762a1bSJed Brown PetscReal kappa; 70c4762a1bSJed Brown PetscInt initial; /* initial conditions type */ 71c4762a1bSJed Brown } AppCtx; 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* 74c4762a1bSJed Brown User-defined routines 75c4762a1bSJed Brown */ 76c4762a1bSJed Brown static PetscErrorCode FormRHS(AppCtx*,DM,Vec); 77c4762a1bSJed Brown static PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec); 78c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*); 79c4762a1bSJed Brown static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*); 80c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,Mat,AppCtx*); 81c4762a1bSJed Brown static PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*); 82c4762a1bSJed Brown 83c4762a1bSJed Brown typedef struct _n_PreCheck *PreCheck; 84c4762a1bSJed Brown struct _n_PreCheck { 85c4762a1bSJed Brown MPI_Comm comm; 86c4762a1bSJed Brown PetscReal angle; 87c4762a1bSJed Brown Vec Ylast; 88c4762a1bSJed Brown PetscViewer monitor; 89c4762a1bSJed Brown }; 90c4762a1bSJed Brown PetscErrorCode PreCheckCreate(MPI_Comm,PreCheck*); 91c4762a1bSJed Brown PetscErrorCode PreCheckDestroy(PreCheck*); 92c4762a1bSJed Brown PetscErrorCode PreCheckFunction(SNESLineSearch,Vec,Vec,PetscBool*,void*); 93c4762a1bSJed Brown PetscErrorCode PreCheckSetFromOptions(PreCheck); 94c4762a1bSJed Brown 95c4762a1bSJed Brown int main(int argc,char **argv) 96c4762a1bSJed Brown { 97c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 98c4762a1bSJed Brown Vec x,r,b; /* solution, residual, rhs vectors */ 99c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 100c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 101c4762a1bSJed Brown SNESConvergedReason reason; /* Check convergence */ 102c4762a1bSJed Brown PetscBool alloc_star; /* Only allocate for the STAR stencil */ 103c4762a1bSJed Brown PetscBool write_output; 104c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "ex15.vts"; 105c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.; 10685b95db3SBarry Smith DM da; /* distributed array data structure */ 107c4762a1bSJed Brown PreCheck precheck = NULL; /* precheck context for version in this file */ 108c4762a1bSJed Brown PetscInt use_precheck; /* 0=none, 1=version in this file, 2=SNES-provided version */ 109c4762a1bSJed Brown PetscReal precheck_angle; /* When manually setting the SNES-provided precheck function */ 110c4762a1bSJed Brown PetscErrorCode ierr; 111c4762a1bSJed Brown SNESLineSearch linesearch; 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Initialize program 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 116c4762a1bSJed Brown 117c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 120c4762a1bSJed Brown Initialize problem parameters 121c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 122c4762a1bSJed Brown user.lambda = 0.0; user.p = 2.0; user.epsilon = 1e-5; user.source = 0.1; user.jtype = JAC_NEWTON;user.initial=-1; 123c4762a1bSJed Brown user.blocks[0] = 1; user.blocks[1] = 1; user.kappa = 1e-3; 124c4762a1bSJed Brown alloc_star = PETSC_FALSE; 125c4762a1bSJed Brown use_precheck = 0; precheck_angle = 10.; 126c4762a1bSJed Brown user.picard = PETSC_FALSE; 127c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"p-Bratu options",__FILE__);CHKERRQ(ierr); 128c4762a1bSJed Brown { 129c4762a1bSJed Brown PetscInt two=2; 130c4762a1bSJed Brown ierr = PetscOptionsReal("-lambda","Bratu parameter","",user.lambda,&user.lambda,NULL);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = PetscOptionsReal("-p","Exponent `p' in p-Laplacian","",user.p,&user.p,NULL);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = PetscOptionsReal("-epsilon","Strain-regularization in p-Laplacian","",user.epsilon,&user.epsilon,NULL);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = PetscOptionsReal("-source","Constant source term","",user.source,&user.source,NULL);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = PetscOptionsEnum("-jtype","Jacobian approximation to assemble","",JacTypes,(PetscEnum)user.jtype,(PetscEnum*)&user.jtype,NULL);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = PetscOptionsName("-picard","Solve with defect-correction Picard iteration","",&user.picard);CHKERRQ(ierr); 13685b95db3SBarry Smith if (user.picard) { 13785b95db3SBarry Smith user.jtype = JAC_PICARD; 1389ace16cdSJacob Faibussowitsch PetscAssertFalse(user.p != 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Picard iteration is only supported for p == 3"); 13985b95db3SBarry Smith /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */ 14085b95db3SBarry Smith /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */ 141c4762a1bSJed Brown ierr = PetscOptionsBool("-alloc_star","Allocate for STAR stencil (5-point)","",alloc_star,&alloc_star,NULL);CHKERRQ(ierr); 14285b95db3SBarry Smith } 143c4762a1bSJed Brown ierr = PetscOptionsInt("-precheck","Use a pre-check correction intended for use with Picard iteration 1=this version, 2=library","",use_precheck,&use_precheck,NULL);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = PetscOptionsInt("-initial","Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)","",user.initial,&user.initial,NULL);CHKERRQ(ierr); 145c4762a1bSJed Brown if (use_precheck == 2) { /* Using library version, get the angle */ 146c4762a1bSJed Brown ierr = PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck_angle,&precheck_angle,NULL);CHKERRQ(ierr); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown ierr = PetscOptionsIntArray("-blocks","number of coefficient interfaces in x and y direction","",user.blocks,&two,NULL);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = PetscOptionsReal("-kappa","diffusivity in odd regions","",user.kappa,&user.kappa,NULL);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = PetscOptionsString("-o","Output solution in vts format","",filename,filename,sizeof(filename),&write_output);CHKERRQ(ierr); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 153c4762a1bSJed Brown if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) { 154c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"WARNING: lambda %g out of range for p=2\n",(double)user.lambda);CHKERRQ(ierr); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 158c4762a1bSJed Brown Create nonlinear solver context 159c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 160c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 163c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 16585b95db3SBarry Smith ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,alloc_star ? DMDA_STENCIL_STAR : DMDA_STENCIL_BOX,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170c4762a1bSJed Brown Extract global vectors from DM; then duplicate for remaining 171c4762a1bSJed Brown vectors that are the same types 172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 173c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = VecDuplicate(x,&b);CHKERRQ(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 17885b95db3SBarry Smith User can override with: 179c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 180c4762a1bSJed Brown (unless user explicitly sets preconditioner) 181c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 182c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 183c4762a1bSJed Brown products within Newton-Krylov method 184c4762a1bSJed Brown 185c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188c4762a1bSJed Brown Set local function evaluation routine 189c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 190c4762a1bSJed Brown ierr = DMSetApplicationContext(da, &user);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = SNESSetDM(snes,da);CHKERRQ(ierr); 192c4762a1bSJed Brown if (user.picard) { 193c4762a1bSJed Brown /* 194c4762a1bSJed Brown This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access 195c4762a1bSJed Brown the SNES to set it 196c4762a1bSJed Brown */ 197c4762a1bSJed Brown ierr = DMDASNESSetPicardLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionPicardLocal, 198c4762a1bSJed Brown (PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = SNESSetJacobian(snes,NULL,NULL,SNESPicardComputeJacobian,&user);CHKERRQ(ierr); 201c4762a1bSJed Brown } else { 202c4762a1bSJed Brown ierr = DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = DMDASNESSetJacobianLocal(da,(PetscErrorCode (*)(DMDALocalInfo*,void*,Mat,Mat,void*))FormJacobianLocal,&user);CHKERRQ(ierr); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown Customize nonlinear solver; set runtime options 208c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 209c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = SNESSetNGS(snes,NonlinearGS,&user);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 212c4762a1bSJed Brown /* Set up the precheck context if requested */ 213c4762a1bSJed Brown if (use_precheck == 1) { /* Use the precheck routines in this file */ 214c4762a1bSJed Brown ierr = PreCheckCreate(PETSC_COMM_WORLD,&precheck);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = PreCheckSetFromOptions(precheck);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = SNESLineSearchSetPreCheck(linesearch,PreCheckFunction,precheck);CHKERRQ(ierr); 217c4762a1bSJed Brown } else if (use_precheck == 2) { /* Use the version provided by the library */ 218c4762a1bSJed Brown ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&precheck_angle);CHKERRQ(ierr); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222c4762a1bSJed Brown Evaluate initial guess 223c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 224c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 225c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 226c4762a1bSJed Brown this vector to zero by calling VecSet(). 227c4762a1bSJed Brown */ 228c4762a1bSJed Brown 229c4762a1bSJed Brown ierr = FormInitialGuess(&user,da,x);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = FormRHS(&user,da,b);CHKERRQ(ierr); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233c4762a1bSJed Brown Solve nonlinear system 234c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235c4762a1bSJed Brown ierr = SNESSolve(snes,b,x);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr); 238c4762a1bSJed Brown 239c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%s Number of nonlinear iterations = %D\n",SNESConvergedReasons[reason],its);CHKERRQ(ierr); 240c4762a1bSJed Brown 241c4762a1bSJed Brown if (write_output) { 242c4762a1bSJed Brown PetscViewer viewer; 243c4762a1bSJed Brown ierr = PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecView(x,viewer);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 249c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 250c4762a1bSJed Brown are no longer needed. 251c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 252c4762a1bSJed Brown 253c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = PreCheckDestroy(&precheck);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = PetscFinalize(); 260c4762a1bSJed Brown return ierr; 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 264c4762a1bSJed Brown /* 265c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 266c4762a1bSJed Brown 267c4762a1bSJed Brown Input Parameters: 268c4762a1bSJed Brown user - user-defined application context 269c4762a1bSJed Brown X - vector 270c4762a1bSJed Brown 271c4762a1bSJed Brown Output Parameter: 272c4762a1bSJed Brown X - vector 273c4762a1bSJed Brown */ 274c4762a1bSJed Brown static PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X) 275c4762a1bSJed Brown { 276c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 277c4762a1bSJed Brown PetscErrorCode ierr; 278c4762a1bSJed Brown PetscReal temp1,temp,hx,hy; 279c4762a1bSJed Brown PetscScalar **x; 280c4762a1bSJed Brown 281c4762a1bSJed Brown PetscFunctionBeginUser; 282c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 283c4762a1bSJed Brown 284c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 285c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 286c4762a1bSJed Brown temp1 = user->lambda / (user->lambda + 1.); 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* 289c4762a1bSJed Brown Get a pointer to vector data. 290c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 291c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 292c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 293c4762a1bSJed Brown the array. 294c4762a1bSJed Brown */ 295c4762a1bSJed Brown ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr); 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DA): 299c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points) 300c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points) 301c4762a1bSJed Brown 302c4762a1bSJed Brown */ 303c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* 306c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid 307c4762a1bSJed Brown */ 308c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 309c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j,My-j-1))*hy; 310c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 311c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 312c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 313c4762a1bSJed Brown x[j][i] = 0.0; 314c4762a1bSJed Brown } else { 315c4762a1bSJed Brown if (user->initial == -1) { 316c4762a1bSJed Brown if (user->lambda != 0) { 317c4762a1bSJed Brown x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp)); 318c4762a1bSJed Brown } else { 319c4762a1bSJed Brown /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting 320c4762a1bSJed Brown * with an exact solution. */ 321c4762a1bSJed Brown const PetscReal 322c4762a1bSJed Brown xx = 2*(PetscReal)i/(Mx-1) - 1, 323c4762a1bSJed Brown yy = 2*(PetscReal)j/(My-1) - 1; 324c4762a1bSJed Brown x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown } else if (user->initial == 0) { 327c4762a1bSJed Brown x[j][i] = 0.; 328c4762a1bSJed Brown } else if (user->initial == 1) { 329c4762a1bSJed Brown const PetscReal 330c4762a1bSJed Brown xx = 2*(PetscReal)i/(Mx-1) - 1, 331c4762a1bSJed Brown yy = 2*(PetscReal)j/(My-1) - 1; 332c4762a1bSJed Brown x[j][i] = (1 - xx*xx) * (1-yy*yy) * xx * yy; 333c4762a1bSJed Brown } else { 334c4762a1bSJed Brown if (user->lambda != 0) { 335c4762a1bSJed Brown x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp)); 336c4762a1bSJed Brown } else { 337c4762a1bSJed Brown x[j][i] = 0.5*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp)); 338c4762a1bSJed Brown } 339c4762a1bSJed Brown } 340c4762a1bSJed Brown } 341c4762a1bSJed Brown } 342c4762a1bSJed Brown } 343c4762a1bSJed Brown /* 344c4762a1bSJed Brown Restore vector 345c4762a1bSJed Brown */ 346c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr); 347c4762a1bSJed Brown PetscFunctionReturn(0); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown /* 351c4762a1bSJed Brown FormRHS - Forms constant RHS for the problem. 352c4762a1bSJed Brown 353c4762a1bSJed Brown Input Parameters: 354c4762a1bSJed Brown user - user-defined application context 355c4762a1bSJed Brown B - RHS vector 356c4762a1bSJed Brown 357c4762a1bSJed Brown Output Parameter: 358c4762a1bSJed Brown B - vector 359c4762a1bSJed Brown */ 360c4762a1bSJed Brown static PetscErrorCode FormRHS(AppCtx *user,DM da,Vec B) 361c4762a1bSJed Brown { 362c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 363c4762a1bSJed Brown PetscErrorCode ierr; 364c4762a1bSJed Brown PetscReal hx,hy; 365c4762a1bSJed Brown PetscScalar **b; 366c4762a1bSJed Brown 367c4762a1bSJed Brown PetscFunctionBeginUser; 368c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 369c4762a1bSJed Brown 370c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 371c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 372c4762a1bSJed Brown ierr = DMDAVecGetArray(da,B,&b);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 374c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 375c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 376c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 377c4762a1bSJed Brown b[j][i] = 0.0; 378c4762a1bSJed Brown } else { 379c4762a1bSJed Brown b[j][i] = hx*hy*user->source; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown } 382c4762a1bSJed Brown } 383c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,B,&b);CHKERRQ(ierr); 384c4762a1bSJed Brown PetscFunctionReturn(0); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387*9fbee547SJacob Faibussowitsch static inline PetscReal kappa(const AppCtx *ctx,PetscReal x,PetscReal y) 388c4762a1bSJed Brown { 389c4762a1bSJed Brown return (((PetscInt)(x*ctx->blocks[0])) + ((PetscInt)(y*ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0; 390c4762a1bSJed Brown } 391c4762a1bSJed Brown /* p-Laplacian diffusivity */ 392*9fbee547SJacob Faibussowitsch static inline PetscScalar eta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy) 393c4762a1bSJed Brown { 394c4762a1bSJed Brown return kappa(ctx,x,y) * PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-2.)); 395c4762a1bSJed Brown } 396*9fbee547SJacob Faibussowitsch static inline PetscScalar deta(const AppCtx *ctx,PetscReal x,PetscReal y,PetscScalar ux,PetscScalar uy) 397c4762a1bSJed Brown { 398c4762a1bSJed Brown return (ctx->p == 2) 399c4762a1bSJed Brown ? 0 400c4762a1bSJed Brown : kappa(ctx,x,y)*PetscPowScalar(PetscSqr(ctx->epsilon)+0.5*(ux*ux + uy*uy),0.5*(ctx->p-4)) * 0.5 * (ctx->p-2.); 401c4762a1bSJed Brown } 402c4762a1bSJed Brown 403c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 404c4762a1bSJed Brown /* 405c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x). 406c4762a1bSJed Brown */ 407c4762a1bSJed Brown static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user) 408c4762a1bSJed Brown { 409c4762a1bSJed Brown PetscReal hx,hy,dhx,dhy,sc; 410c4762a1bSJed Brown PetscInt i,j; 411c4762a1bSJed Brown PetscScalar eu; 412c4762a1bSJed Brown PetscErrorCode ierr; 413c4762a1bSJed Brown 414c4762a1bSJed Brown PetscFunctionBeginUser; 415c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 416c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 417c4762a1bSJed Brown sc = hx*hy*user->lambda; 418c4762a1bSJed Brown dhx = 1/hx; 419c4762a1bSJed Brown dhy = 1/hy; 420c4762a1bSJed Brown /* 421c4762a1bSJed Brown Compute function over the locally owned part of the grid 422c4762a1bSJed Brown */ 423c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 424c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 425c4762a1bSJed Brown PetscReal xx = i*hx,yy = j*hy; 426c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 427c4762a1bSJed Brown f[j][i] = x[j][i]; 428c4762a1bSJed Brown } else { 429c4762a1bSJed Brown const PetscScalar 430c4762a1bSJed Brown u = x[j][i], 431c4762a1bSJed Brown ux_E = dhx*(x[j][i+1]-x[j][i]), 432c4762a1bSJed Brown uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]), 433c4762a1bSJed Brown ux_W = dhx*(x[j][i]-x[j][i-1]), 434c4762a1bSJed Brown uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]), 435c4762a1bSJed Brown ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]), 436c4762a1bSJed Brown uy_N = dhy*(x[j+1][i]-x[j][i]), 437c4762a1bSJed Brown ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]), 438c4762a1bSJed Brown uy_S = dhy*(x[j][i]-x[j-1][i]), 439c4762a1bSJed Brown e_E = eta(user,xx,yy,ux_E,uy_E), 440c4762a1bSJed Brown e_W = eta(user,xx,yy,ux_W,uy_W), 441c4762a1bSJed Brown e_N = eta(user,xx,yy,ux_N,uy_N), 442c4762a1bSJed Brown e_S = eta(user,xx,yy,ux_S,uy_S), 443c4762a1bSJed Brown uxx = -hy * (e_E*ux_E - e_W*ux_W), 444c4762a1bSJed Brown uyy = -hx * (e_N*uy_N - e_S*uy_S); 445c4762a1bSJed Brown if (sc) eu = PetscExpScalar(u); 446c4762a1bSJed Brown else eu = 0.; 4470e3d61c9SBarry Smith /* For p=2, these terms decay to: 4480e3d61c9SBarry Smith uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx 4490e3d61c9SBarry Smith uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy 4500e3d61c9SBarry Smith */ 451c4762a1bSJed Brown f[j][i] = uxx + uyy - sc*eu; 452c4762a1bSJed Brown } 453c4762a1bSJed Brown } 454c4762a1bSJed Brown } 455c4762a1bSJed Brown ierr = PetscLogFlops(info->xm*info->ym*(72.0));CHKERRQ(ierr); 456c4762a1bSJed Brown PetscFunctionReturn(0); 457c4762a1bSJed Brown } 458c4762a1bSJed Brown 459c4762a1bSJed Brown /* 460c4762a1bSJed Brown This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation, 461c4762a1bSJed Brown that is FormFunction applies A(x) x - b(x) while this applies b(x) because for Picard we think of it as solving A(x) x = b(x) 462c4762a1bSJed Brown 463c4762a1bSJed Brown */ 464c4762a1bSJed Brown static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user) 465c4762a1bSJed Brown { 466c4762a1bSJed Brown PetscReal hx,hy,sc; 467c4762a1bSJed Brown PetscInt i,j; 468c4762a1bSJed Brown PetscErrorCode ierr; 469c4762a1bSJed Brown 470c4762a1bSJed Brown PetscFunctionBeginUser; 471c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 472c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 473c4762a1bSJed Brown sc = hx*hy*user->lambda; 474c4762a1bSJed Brown /* 475c4762a1bSJed Brown Compute function over the locally owned part of the grid 476c4762a1bSJed Brown */ 477c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 478c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 479c4762a1bSJed Brown if (!(i == 0 || j == 0 || i == info->mx-1 || j == info->my-1)) { 480c4762a1bSJed Brown const PetscScalar u = x[j][i]; 481c4762a1bSJed Brown f[j][i] = sc*PetscExpScalar(u); 4820df40c35SBarry Smith } else { 4830df40c35SBarry Smith f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */ 484c4762a1bSJed Brown } 485c4762a1bSJed Brown } 486c4762a1bSJed Brown } 487c4762a1bSJed Brown ierr = PetscLogFlops(info->xm*info->ym*2.0);CHKERRQ(ierr); 488c4762a1bSJed Brown PetscFunctionReturn(0); 489c4762a1bSJed Brown } 490c4762a1bSJed Brown 491c4762a1bSJed Brown /* 492c4762a1bSJed Brown FormJacobianLocal - Evaluates Jacobian matrix. 493c4762a1bSJed Brown */ 494c4762a1bSJed Brown static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat J,Mat B,AppCtx *user) 495c4762a1bSJed Brown { 496c4762a1bSJed Brown PetscErrorCode ierr; 497c4762a1bSJed Brown PetscInt i,j; 498c4762a1bSJed Brown MatStencil col[9],row; 499c4762a1bSJed Brown PetscScalar v[9]; 500c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx,dhx,dhy,sc; 501c4762a1bSJed Brown 502c4762a1bSJed Brown PetscFunctionBeginUser; 5030df40c35SBarry Smith ierr = MatZeroEntries(B);CHKERRQ(ierr); 504c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 505c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 506c4762a1bSJed Brown sc = hx*hy*user->lambda; 507c4762a1bSJed Brown dhx = 1/hx; 508c4762a1bSJed Brown dhy = 1/hy; 509c4762a1bSJed Brown hxdhy = hx/hy; 510c4762a1bSJed Brown hydhx = hy/hx; 511c4762a1bSJed Brown 512c4762a1bSJed Brown /* 513c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 514c4762a1bSJed Brown - PETSc parallel matrix formats are partitioned by 515c4762a1bSJed Brown contiguous chunks of rows across the processors. 516c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 517c4762a1bSJed Brown locally (but any non-local elements will be sent to the 518c4762a1bSJed Brown appropriate processor during matrix assembly). 519c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 520c4762a1bSJed Brown */ 521c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 522c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 523c4762a1bSJed Brown PetscReal xx = i*hx,yy = j*hy; 524c4762a1bSJed Brown row.j = j; row.i = i; 525c4762a1bSJed Brown /* boundary points */ 526c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 527c4762a1bSJed Brown v[0] = 1.0; 528c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr); 529c4762a1bSJed Brown } else { 530c4762a1bSJed Brown /* interior grid points */ 531c4762a1bSJed Brown const PetscScalar 532c4762a1bSJed Brown ux_E = dhx*(x[j][i+1]-x[j][i]), 533c4762a1bSJed Brown uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]), 534c4762a1bSJed Brown ux_W = dhx*(x[j][i]-x[j][i-1]), 535c4762a1bSJed Brown uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]), 536c4762a1bSJed Brown ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]), 537c4762a1bSJed Brown uy_N = dhy*(x[j+1][i]-x[j][i]), 538c4762a1bSJed Brown ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]), 539c4762a1bSJed Brown uy_S = dhy*(x[j][i]-x[j-1][i]), 540c4762a1bSJed Brown u = x[j][i], 541c4762a1bSJed Brown e_E = eta(user,xx,yy,ux_E,uy_E), 542c4762a1bSJed Brown e_W = eta(user,xx,yy,ux_W,uy_W), 543c4762a1bSJed Brown e_N = eta(user,xx,yy,ux_N,uy_N), 544c4762a1bSJed Brown e_S = eta(user,xx,yy,ux_S,uy_S), 545c4762a1bSJed Brown de_E = deta(user,xx,yy,ux_E,uy_E), 546c4762a1bSJed Brown de_W = deta(user,xx,yy,ux_W,uy_W), 547c4762a1bSJed Brown de_N = deta(user,xx,yy,ux_N,uy_N), 548c4762a1bSJed Brown de_S = deta(user,xx,yy,ux_S,uy_S), 549c4762a1bSJed Brown skew_E = de_E*ux_E*uy_E, 550c4762a1bSJed Brown skew_W = de_W*ux_W*uy_W, 551c4762a1bSJed Brown skew_N = de_N*ux_N*uy_N, 552c4762a1bSJed Brown skew_S = de_S*ux_S*uy_S, 553c4762a1bSJed Brown cross_EW = 0.25*(skew_E - skew_W), 554c4762a1bSJed Brown cross_NS = 0.25*(skew_N - skew_S), 555c4762a1bSJed Brown newt_E = e_E+de_E*PetscSqr(ux_E), 556c4762a1bSJed Brown newt_W = e_W+de_W*PetscSqr(ux_W), 557c4762a1bSJed Brown newt_N = e_N+de_N*PetscSqr(uy_N), 558c4762a1bSJed Brown newt_S = e_S+de_S*PetscSqr(uy_S); 559c4762a1bSJed Brown /* interior grid points */ 560c4762a1bSJed Brown switch (user->jtype) { 561c4762a1bSJed Brown case JAC_BRATU: 562c4762a1bSJed Brown /* Jacobian from p=2 */ 563c4762a1bSJed Brown v[0] = -hxdhy; col[0].j = j-1; col[0].i = i; 564c4762a1bSJed Brown v[1] = -hydhx; col[1].j = j; col[1].i = i-1; 565c4762a1bSJed Brown v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(u); col[2].j = row.j; col[2].i = row.i; 566c4762a1bSJed Brown v[3] = -hydhx; col[3].j = j; col[3].i = i+1; 567c4762a1bSJed Brown v[4] = -hxdhy; col[4].j = j+1; col[4].i = i; 568c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); 569c4762a1bSJed Brown break; 570c4762a1bSJed Brown case JAC_PICARD: 571c4762a1bSJed Brown /* Jacobian arising from Picard linearization */ 572c4762a1bSJed Brown v[0] = -hxdhy*e_S; col[0].j = j-1; col[0].i = i; 573c4762a1bSJed Brown v[1] = -hydhx*e_W; col[1].j = j; col[1].i = i-1; 574c4762a1bSJed Brown v[2] = (e_W+e_E)*hydhx + (e_S+e_N)*hxdhy; col[2].j = row.j; col[2].i = row.i; 575c4762a1bSJed Brown v[3] = -hydhx*e_E; col[3].j = j; col[3].i = i+1; 576c4762a1bSJed Brown v[4] = -hxdhy*e_N; col[4].j = j+1; col[4].i = i; 577c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); 578c4762a1bSJed Brown break; 579c4762a1bSJed Brown case JAC_STAR: 580c4762a1bSJed Brown /* Full Jacobian, but only a star stencil */ 581c4762a1bSJed Brown col[0].j = j-1; col[0].i = i; 582c4762a1bSJed Brown col[1].j = j; col[1].i = i-1; 583c4762a1bSJed Brown col[2].j = j; col[2].i = i; 584c4762a1bSJed Brown col[3].j = j; col[3].i = i+1; 585c4762a1bSJed Brown col[4].j = j+1; col[4].i = i; 586c4762a1bSJed Brown v[0] = -hxdhy*newt_S + cross_EW; 587c4762a1bSJed Brown v[1] = -hydhx*newt_W + cross_NS; 588c4762a1bSJed Brown v[2] = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u); 589c4762a1bSJed Brown v[3] = -hydhx*newt_E - cross_NS; 590c4762a1bSJed Brown v[4] = -hxdhy*newt_N - cross_EW; 591c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); 592c4762a1bSJed Brown break; 593c4762a1bSJed Brown case JAC_NEWTON: 5940e3d61c9SBarry Smith /* The Jacobian is 5950e3d61c9SBarry Smith 5960e3d61c9SBarry Smith -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u 5970e3d61c9SBarry Smith 5980e3d61c9SBarry Smith */ 599c4762a1bSJed Brown col[0].j = j-1; col[0].i = i-1; 600c4762a1bSJed Brown col[1].j = j-1; col[1].i = i; 601c4762a1bSJed Brown col[2].j = j-1; col[2].i = i+1; 602c4762a1bSJed Brown col[3].j = j; col[3].i = i-1; 603c4762a1bSJed Brown col[4].j = j; col[4].i = i; 604c4762a1bSJed Brown col[5].j = j; col[5].i = i+1; 605c4762a1bSJed Brown col[6].j = j+1; col[6].i = i-1; 606c4762a1bSJed Brown col[7].j = j+1; col[7].i = i; 607c4762a1bSJed Brown col[8].j = j+1; col[8].i = i+1; 608c4762a1bSJed Brown v[0] = -0.25*(skew_S + skew_W); 609c4762a1bSJed Brown v[1] = -hxdhy*newt_S + cross_EW; 610c4762a1bSJed Brown v[2] = 0.25*(skew_S + skew_E); 611c4762a1bSJed Brown v[3] = -hydhx*newt_W + cross_NS; 612c4762a1bSJed Brown v[4] = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*PetscExpScalar(u); 613c4762a1bSJed Brown v[5] = -hydhx*newt_E - cross_NS; 614c4762a1bSJed Brown v[6] = 0.25*(skew_N + skew_W); 615c4762a1bSJed Brown v[7] = -hxdhy*newt_N - cross_EW; 616c4762a1bSJed Brown v[8] = -0.25*(skew_N + skew_E); 617c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&row,9,col,v,INSERT_VALUES);CHKERRQ(ierr); 618c4762a1bSJed Brown break; 619c4762a1bSJed Brown default: 62098921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)info->da),PETSC_ERR_SUP,"Jacobian type %d not implemented",user->jtype); 621c4762a1bSJed Brown } 622c4762a1bSJed Brown } 623c4762a1bSJed Brown } 624c4762a1bSJed Brown } 625c4762a1bSJed Brown 626c4762a1bSJed Brown /* 627c4762a1bSJed Brown Assemble matrix, using the 2-step process: 628c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 629c4762a1bSJed Brown */ 630c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 631c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 632c4762a1bSJed Brown 633c4762a1bSJed Brown if (J != B) { 634c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 635c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown 638c4762a1bSJed Brown /* 639c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 640c4762a1bSJed Brown matrix. If we do, it will generate an error. 641c4762a1bSJed Brown */ 642c4762a1bSJed Brown if (user->jtype == JAC_NEWTON) { 643c4762a1bSJed Brown ierr = PetscLogFlops(info->xm*info->ym*(131.0));CHKERRQ(ierr); 644c4762a1bSJed Brown } 645c4762a1bSJed Brown PetscFunctionReturn(0); 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648c4762a1bSJed Brown /*********************************************************** 649c4762a1bSJed Brown * PreCheck implementation 650c4762a1bSJed Brown ***********************************************************/ 651c4762a1bSJed Brown PetscErrorCode PreCheckSetFromOptions(PreCheck precheck) 652c4762a1bSJed Brown { 653c4762a1bSJed Brown PetscErrorCode ierr; 654c4762a1bSJed Brown PetscBool flg; 655c4762a1bSJed Brown 656c4762a1bSJed Brown PetscFunctionBeginUser; 657c4762a1bSJed Brown ierr = PetscOptionsBegin(precheck->comm,NULL,"PreCheck Options","none");CHKERRQ(ierr); 658c4762a1bSJed Brown ierr = PetscOptionsReal("-precheck_angle","Angle in degrees between successive search directions necessary to activate step correction","",precheck->angle,&precheck->angle,NULL);CHKERRQ(ierr); 659c4762a1bSJed Brown flg = PETSC_FALSE; 660c4762a1bSJed Brown ierr = PetscOptionsBool("-precheck_monitor","Monitor choices made by precheck routine","",flg,&flg,NULL);CHKERRQ(ierr); 661c4762a1bSJed Brown if (flg) { 662c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(precheck->comm,"stdout",&precheck->monitor);CHKERRQ(ierr); 663c4762a1bSJed Brown } 664c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 665c4762a1bSJed Brown PetscFunctionReturn(0); 666c4762a1bSJed Brown } 667c4762a1bSJed Brown 668c4762a1bSJed Brown /* 669c4762a1bSJed Brown Compare the direction of the current and previous step, modify the current step accordingly 670c4762a1bSJed Brown */ 671c4762a1bSJed Brown PetscErrorCode PreCheckFunction(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed, void *ctx) 672c4762a1bSJed Brown { 673c4762a1bSJed Brown PetscErrorCode ierr; 674c4762a1bSJed Brown PreCheck precheck; 675c4762a1bSJed Brown Vec Ylast; 676c4762a1bSJed Brown PetscScalar dot; 677c4762a1bSJed Brown PetscInt iter; 678c4762a1bSJed Brown PetscReal ynorm,ylastnorm,theta,angle_radians; 679c4762a1bSJed Brown SNES snes; 680c4762a1bSJed Brown 681c4762a1bSJed Brown PetscFunctionBeginUser; 682c4762a1bSJed Brown ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 683c4762a1bSJed Brown precheck = (PreCheck)ctx; 684c4762a1bSJed Brown if (!precheck->Ylast) {ierr = VecDuplicate(Y,&precheck->Ylast);CHKERRQ(ierr);} 685c4762a1bSJed Brown Ylast = precheck->Ylast; 686c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 687c4762a1bSJed Brown if (iter < 1) { 688c4762a1bSJed Brown ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 689c4762a1bSJed Brown *changed = PETSC_FALSE; 690c4762a1bSJed Brown PetscFunctionReturn(0); 691c4762a1bSJed Brown } 692c4762a1bSJed Brown 693c4762a1bSJed Brown ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 694c4762a1bSJed Brown ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 695c4762a1bSJed Brown ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 696c4762a1bSJed Brown /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 697c4762a1bSJed Brown theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 698c4762a1bSJed Brown angle_radians = precheck->angle * PETSC_PI / 180.; 699c4762a1bSJed Brown if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 700c4762a1bSJed Brown /* Modify the step Y */ 701c4762a1bSJed Brown PetscReal alpha,ydiffnorm; 702c4762a1bSJed Brown ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 703c4762a1bSJed Brown ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 704c4762a1bSJed Brown alpha = ylastnorm / ydiffnorm; 705c4762a1bSJed Brown ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 706c4762a1bSJed Brown ierr = VecScale(Y,alpha);CHKERRQ(ierr); 707c4762a1bSJed Brown if (precheck->monitor) { 708c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(precheck->monitor,"Angle %E degrees less than threshold %g, corrected step by alpha=%g\n",(double)(theta*180./PETSC_PI),(double)precheck->angle,(double)alpha);CHKERRQ(ierr); 709c4762a1bSJed Brown } 710c4762a1bSJed Brown } else { 711c4762a1bSJed Brown ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 712c4762a1bSJed Brown *changed = PETSC_FALSE; 713c4762a1bSJed Brown if (precheck->monitor) { 714c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(precheck->monitor,"Angle %E degrees exceeds threshold %g, no correction applied\n",(double)(theta*180./PETSC_PI),(double)precheck->angle);CHKERRQ(ierr); 715c4762a1bSJed Brown } 716c4762a1bSJed Brown } 717c4762a1bSJed Brown PetscFunctionReturn(0); 718c4762a1bSJed Brown } 719c4762a1bSJed Brown 720c4762a1bSJed Brown PetscErrorCode PreCheckDestroy(PreCheck *precheck) 721c4762a1bSJed Brown { 722c4762a1bSJed Brown PetscErrorCode ierr; 723c4762a1bSJed Brown 724c4762a1bSJed Brown PetscFunctionBeginUser; 725c4762a1bSJed Brown if (!*precheck) PetscFunctionReturn(0); 726c4762a1bSJed Brown ierr = VecDestroy(&(*precheck)->Ylast);CHKERRQ(ierr); 727c4762a1bSJed Brown ierr = PetscViewerDestroy(&(*precheck)->monitor);CHKERRQ(ierr); 728c4762a1bSJed Brown ierr = PetscFree(*precheck);CHKERRQ(ierr); 729c4762a1bSJed Brown PetscFunctionReturn(0); 730c4762a1bSJed Brown } 731c4762a1bSJed Brown 732c4762a1bSJed Brown PetscErrorCode PreCheckCreate(MPI_Comm comm,PreCheck *precheck) 733c4762a1bSJed Brown { 734c4762a1bSJed Brown PetscErrorCode ierr; 735c4762a1bSJed Brown 736c4762a1bSJed Brown PetscFunctionBeginUser; 737c4762a1bSJed Brown ierr = PetscNew(precheck);CHKERRQ(ierr); 738c4762a1bSJed Brown 739c4762a1bSJed Brown (*precheck)->comm = comm; 740c4762a1bSJed Brown (*precheck)->angle = 10.; /* only active if angle is less than 10 degrees */ 741c4762a1bSJed Brown PetscFunctionReturn(0); 742c4762a1bSJed Brown } 743c4762a1bSJed Brown 744c4762a1bSJed Brown /* 745c4762a1bSJed Brown Applies some sweeps on nonlinear Gauss-Seidel on each process 746c4762a1bSJed Brown 747c4762a1bSJed Brown */ 748c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx) 749c4762a1bSJed Brown { 750c4762a1bSJed Brown PetscInt i,j,k,xs,ys,xm,ym,its,tot_its,sweeps,l,m; 751c4762a1bSJed Brown PetscErrorCode ierr; 752c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx,dhx,dhy,sc; 753c4762a1bSJed Brown PetscScalar **x,**b,bij,F,F0=0,J,y,u,eu; 754c4762a1bSJed Brown PetscReal atol,rtol,stol; 755c4762a1bSJed Brown DM da; 756c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 757c4762a1bSJed Brown Vec localX,localB; 758c4762a1bSJed Brown DMDALocalInfo info; 759c4762a1bSJed Brown 760c4762a1bSJed Brown PetscFunctionBeginUser; 761c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 762c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 763c4762a1bSJed Brown 764c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 765c4762a1bSJed Brown hy = 1.0/(PetscReal)(info.my-1); 766c4762a1bSJed Brown sc = hx*hy*user->lambda; 767c4762a1bSJed Brown dhx = 1/hx; 768c4762a1bSJed Brown dhy = 1/hy; 769c4762a1bSJed Brown hxdhy = hx/hy; 770c4762a1bSJed Brown hydhx = hy/hx; 771c4762a1bSJed Brown 772c4762a1bSJed Brown tot_its = 0; 773c4762a1bSJed Brown ierr = SNESNGSGetSweeps(snes,&sweeps);CHKERRQ(ierr); 774c4762a1bSJed Brown ierr = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 775c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 776c4762a1bSJed Brown if (B) { 777c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localB);CHKERRQ(ierr); 778c4762a1bSJed Brown } 779c4762a1bSJed Brown if (B) { 780c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB);CHKERRQ(ierr); 781c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB);CHKERRQ(ierr); 782c4762a1bSJed Brown } 78360d4fc61SSatish Balay if (B) {ierr = DMDAVecGetArrayRead(da,localB,&b);CHKERRQ(ierr);} 784c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 785c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 786c4762a1bSJed Brown ierr = DMDAVecGetArray(da,localX,&x);CHKERRQ(ierr); 787c4762a1bSJed Brown for (l=0; l<sweeps; l++) { 788c4762a1bSJed Brown /* 789c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DMDA): 790c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points) 791c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points) 792c4762a1bSJed Brown */ 793c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 794c4762a1bSJed Brown for (m=0; m<2; m++) { 795c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 796c4762a1bSJed Brown for (i=xs+(m+j)%2; i<xs+xm; i+=2) { 797c4762a1bSJed Brown PetscReal xx = i*hx,yy = j*hy; 798c4762a1bSJed Brown if (B) bij = b[j][i]; 799c4762a1bSJed Brown else bij = 0.; 800c4762a1bSJed Brown 801c4762a1bSJed Brown if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) { 802c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 803c4762a1bSJed Brown x[j][i] = 0.0 + bij; 804c4762a1bSJed Brown } else { 805c4762a1bSJed Brown const PetscScalar 806c4762a1bSJed Brown u_E = x[j][i+1], 807c4762a1bSJed Brown u_W = x[j][i-1], 808c4762a1bSJed Brown u_N = x[j+1][i], 809c4762a1bSJed Brown u_S = x[j-1][i]; 810c4762a1bSJed Brown const PetscScalar 811c4762a1bSJed Brown uy_E = 0.25*dhy*(x[j+1][i]+x[j+1][i+1]-x[j-1][i]-x[j-1][i+1]), 812c4762a1bSJed Brown uy_W = 0.25*dhy*(x[j+1][i-1]+x[j+1][i]-x[j-1][i-1]-x[j-1][i]), 813c4762a1bSJed Brown ux_N = 0.25*dhx*(x[j][i+1]+x[j+1][i+1]-x[j][i-1]-x[j+1][i-1]), 814c4762a1bSJed Brown ux_S = 0.25*dhx*(x[j-1][i+1]+x[j][i+1]-x[j-1][i-1]-x[j][i-1]); 815c4762a1bSJed Brown u = x[j][i]; 816c4762a1bSJed Brown for (k=0; k<its; k++) { 817c4762a1bSJed Brown const PetscScalar 818c4762a1bSJed Brown ux_E = dhx*(u_E-u), 819c4762a1bSJed Brown ux_W = dhx*(u-u_W), 820c4762a1bSJed Brown uy_N = dhy*(u_N-u), 821c4762a1bSJed Brown uy_S = dhy*(u-u_S), 822c4762a1bSJed Brown e_E = eta(user,xx,yy,ux_E,uy_E), 823c4762a1bSJed Brown e_W = eta(user,xx,yy,ux_W,uy_W), 824c4762a1bSJed Brown e_N = eta(user,xx,yy,ux_N,uy_N), 825c4762a1bSJed Brown e_S = eta(user,xx,yy,ux_S,uy_S), 826c4762a1bSJed Brown de_E = deta(user,xx,yy,ux_E,uy_E), 827c4762a1bSJed Brown de_W = deta(user,xx,yy,ux_W,uy_W), 828c4762a1bSJed Brown de_N = deta(user,xx,yy,ux_N,uy_N), 829c4762a1bSJed Brown de_S = deta(user,xx,yy,ux_S,uy_S), 830c4762a1bSJed Brown newt_E = e_E+de_E*PetscSqr(ux_E), 831c4762a1bSJed Brown newt_W = e_W+de_W*PetscSqr(ux_W), 832c4762a1bSJed Brown newt_N = e_N+de_N*PetscSqr(uy_N), 833c4762a1bSJed Brown newt_S = e_S+de_S*PetscSqr(uy_S), 834c4762a1bSJed Brown uxx = -hy * (e_E*ux_E - e_W*ux_W), 835c4762a1bSJed Brown uyy = -hx * (e_N*uy_N - e_S*uy_S); 836c4762a1bSJed Brown 837c4762a1bSJed Brown if (sc) eu = PetscExpScalar(u); 838c4762a1bSJed Brown else eu = 0; 839c4762a1bSJed Brown 840c4762a1bSJed Brown F = uxx + uyy - sc*eu - bij; 841c4762a1bSJed Brown if (k == 0) F0 = F; 842c4762a1bSJed Brown J = hxdhy*(newt_N + newt_S) + hydhx*(newt_E + newt_W) - sc*eu; 843c4762a1bSJed Brown y = F/J; 844c4762a1bSJed Brown u -= y; 845c4762a1bSJed Brown tot_its++; 846c4762a1bSJed Brown if (atol > PetscAbsReal(PetscRealPart(F)) || 847c4762a1bSJed Brown rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || 848c4762a1bSJed Brown stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) { 849c4762a1bSJed Brown break; 850c4762a1bSJed Brown } 851c4762a1bSJed Brown } 852c4762a1bSJed Brown x[j][i] = u; 853c4762a1bSJed Brown } 854c4762a1bSJed Brown } 855c4762a1bSJed Brown } 856c4762a1bSJed Brown } 857c4762a1bSJed Brown /* 858c4762a1bSJed Brown x Restore vector 859c4762a1bSJed Brown */ 860c4762a1bSJed Brown } 861c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,localX,&x);CHKERRQ(ierr); 862c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);CHKERRQ(ierr); 863c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);CHKERRQ(ierr); 864c4762a1bSJed Brown ierr = PetscLogFlops(tot_its*(118.0));CHKERRQ(ierr); 865c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 866c4762a1bSJed Brown if (B) { 867c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localB,&b);CHKERRQ(ierr); 868c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localB);CHKERRQ(ierr); 869c4762a1bSJed Brown } 870c4762a1bSJed Brown PetscFunctionReturn(0); 871c4762a1bSJed Brown } 872c4762a1bSJed Brown 873c4762a1bSJed Brown /*TEST 874c4762a1bSJed Brown 875c4762a1bSJed Brown test: 876c4762a1bSJed Brown nsize: 2 877c4762a1bSJed Brown args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON 878c4762a1bSJed Brown requires: !single 879c4762a1bSJed Brown 880c4762a1bSJed Brown test: 881c4762a1bSJed Brown suffix: 2 882c4762a1bSJed Brown nsize: 2 883c4762a1bSJed Brown args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1 884c4762a1bSJed Brown requires: !single 885c4762a1bSJed Brown 886c4762a1bSJed Brown test: 887c4762a1bSJed Brown suffix: 3 888c4762a1bSJed Brown nsize: 2 88985b95db3SBarry Smith args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3 89085b95db3SBarry Smith requires: !single 89185b95db3SBarry Smith 89285b95db3SBarry Smith test: 89385b95db3SBarry Smith suffix: 3_star 89485b95db3SBarry Smith nsize: 2 89585b95db3SBarry Smith args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3 -alloc_star 89685b95db3SBarry Smith output_file: output/ex15_3.out 897c4762a1bSJed Brown requires: !single 898c4762a1bSJed Brown 899c4762a1bSJed Brown test: 900c4762a1bSJed Brown suffix: 4 901c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtonls -npc_snes_type ngs -snes_npc_side left -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -ksp_monitor_short -pc_type none 902c4762a1bSJed Brown requires: !single 903c4762a1bSJed Brown 904c4762a1bSJed Brown test: 905c4762a1bSJed Brown suffix: lag_jac 906c4762a1bSJed Brown nsize: 4 907c4762a1bSJed Brown args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 6.0 -lambda 0 -jtype NEWTON -snes_type ngmres -npc_snes_type newtonls -npc_snes_lag_jacobian 5 -npc_pc_type asm -npc_ksp_converged_reason -npc_snes_lag_jacobian_persists 908c4762a1bSJed Brown requires: !single 909c4762a1bSJed Brown 910c4762a1bSJed Brown test: 911c4762a1bSJed Brown suffix: lag_pc 912c4762a1bSJed Brown nsize: 4 913c4762a1bSJed Brown args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 6.0 -lambda 0 -jtype NEWTON -snes_type ngmres -npc_snes_type newtonls -npc_snes_lag_preconditioner 5 -npc_pc_type asm -npc_ksp_converged_reason -npc_snes_lag_preconditioner_persists 914c4762a1bSJed Brown requires: !single 915c4762a1bSJed Brown 916c4762a1bSJed Brown test: 917c4762a1bSJed Brown suffix: nleqerr 918c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtonls -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -snes_linesearch_monitor -pc_type lu -snes_linesearch_type nleqerr 919c4762a1bSJed Brown requires: !single 920c4762a1bSJed Brown 921bbc1464cSBarry Smith test: 922bbc1464cSBarry Smith suffix: mf 923bbc1464cSBarry Smith args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_mf_operator 924bbc1464cSBarry Smith requires: !single 925bbc1464cSBarry Smith 926bbc1464cSBarry Smith test: 927bbc1464cSBarry Smith suffix: mf_picard 928bbc1464cSBarry Smith args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_mf_operator -picard 929bbc1464cSBarry Smith requires: !single 930bbc1464cSBarry Smith output_file: output/ex15_mf.out 931bbc1464cSBarry Smith 93285b95db3SBarry Smith test: 93385b95db3SBarry Smith suffix: fd_picard 93485b95db3SBarry Smith args: -snes_monitor_short -pc_type lu -da_refine 2 -p 3 -ksp_rtol 1.e-12 -snes_fd -picard 93585b95db3SBarry Smith requires: !single 93685b95db3SBarry Smith 93785b95db3SBarry Smith test: 93885b95db3SBarry Smith suffix: fd_color_picard 93985b95db3SBarry Smith args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_fd_color -picard 94085b95db3SBarry Smith requires: !single 94185b95db3SBarry Smith output_file: output/ex15_mf.out 94285b95db3SBarry Smith 943c4762a1bSJed Brown TEST*/ 944