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