1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\ 3c4762a1bSJed Brown This example employs a user-defined monitoring routine.\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /*T 6c4762a1bSJed Brown Concepts: SNES^basic uniprocessor example 7c4762a1bSJed Brown Concepts: SNES^setting a user-defined monitoring routine 8c4762a1bSJed Brown Processors: 1 9c4762a1bSJed Brown T*/ 10c4762a1bSJed Brown 11c4762a1bSJed Brown /* 12c4762a1bSJed Brown Include "petscdraw.h" so that we can use PETSc drawing routines. 13c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 14c4762a1bSJed Brown file automatically includes: 15c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 16c4762a1bSJed Brown petscmat.h - matrices 17c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 18c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 19c4762a1bSJed Brown petscksp.h - linear solvers 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown 22c4762a1bSJed Brown #include <petscsnes.h> 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* 25c4762a1bSJed Brown User-defined routines 26c4762a1bSJed Brown */ 27c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 28c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 29c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec); 30c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* 33c4762a1bSJed Brown User-defined context for monitoring 34c4762a1bSJed Brown */ 35c4762a1bSJed Brown typedef struct { 36c4762a1bSJed Brown PetscViewer viewer; 37c4762a1bSJed Brown } MonitorCtx; 38c4762a1bSJed Brown 39c4762a1bSJed Brown int main(int argc,char **argv) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown SNES snes; /* SNES context */ 42c4762a1bSJed Brown Vec x,r,F,U; /* vectors */ 43c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 44c4762a1bSJed Brown MonitorCtx monP; /* monitoring context */ 45c4762a1bSJed Brown PetscErrorCode ierr; 46c4762a1bSJed Brown PetscInt its,n = 5,i,maxit,maxf; 47c4762a1bSJed Brown PetscMPIInt size; 48c4762a1bSJed Brown PetscScalar h,xp,v,none = -1.0; 49c4762a1bSJed Brown PetscReal abstol,rtol,stol,norm; 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 52ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 53*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 54c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 55c4762a1bSJed Brown h = 1.0/(n-1); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 58c4762a1bSJed Brown Create nonlinear solver context 59c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 60c4762a1bSJed Brown 61c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 64c4762a1bSJed Brown Create vector data structures; set function evaluation routine 65c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown Note that we form 1 vector from scratch and then duplicate as needed. 68c4762a1bSJed Brown */ 69c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = VecDuplicate(x,&F);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = VecDuplicate(x,&U);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Set function evaluation routine and vector 78c4762a1bSJed Brown */ 79c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 83c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(J,3,NULL);CHKERRQ(ierr); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 92c4762a1bSJed Brown routine. User can override with: 93c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 94c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 95c4762a1bSJed Brown (unless user explicitly sets preconditioner) 96c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 97c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 98c4762a1bSJed Brown products within Newton-Krylov method 99c4762a1bSJed Brown */ 100c4762a1bSJed Brown 101c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian,NULL);CHKERRQ(ierr); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104c4762a1bSJed Brown Customize nonlinear solver; set runtime options 105c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* 108c4762a1bSJed Brown Set an optional user-defined monitoring routine 109c4762a1bSJed Brown */ 110c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = SNESMonitorSet(snes,Monitor,&monP,0);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* 114c4762a1bSJed Brown Set names for some vectors to facilitate monitoring (optional) 115c4762a1bSJed Brown */ 116c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* 120c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 121c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 122c4762a1bSJed Brown */ 123c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Print parameters used for convergence testing (optional) ... just 127c4762a1bSJed Brown to demonstrate this routine; this information is also printed with 128c4762a1bSJed Brown the option -snes_view 129c4762a1bSJed Brown */ 130c4762a1bSJed Brown ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf);CHKERRQ(ierr); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 134c4762a1bSJed Brown Initialize application: 135c4762a1bSJed Brown Store right-hand-side of PDE and exact solution 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137c4762a1bSJed Brown 138c4762a1bSJed Brown xp = 0.0; 139c4762a1bSJed Brown for (i=0; i<n; i++) { 140c4762a1bSJed Brown v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 141c4762a1bSJed Brown ierr = VecSetValues(F,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 142c4762a1bSJed Brown v = xp*xp*xp; 143c4762a1bSJed Brown ierr = VecSetValues(U,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 144c4762a1bSJed Brown xp += h; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 149c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 150c4762a1bSJed Brown /* 151c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 152c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 153c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 154c4762a1bSJed Brown this vector to zero by calling VecSet(). 155c4762a1bSJed Brown */ 156c4762a1bSJed Brown ierr = FormInitialGuess(x);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162c4762a1bSJed Brown Check solution and clean up 163c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* 166c4762a1bSJed Brown Check the error 167c4762a1bSJed Brown */ 168c4762a1bSJed Brown ierr = VecAXPY(x,none,U);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);CHKERRQ(ierr); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 174c4762a1bSJed Brown are no longer needed. 175c4762a1bSJed Brown */ 176c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); ierr = VecDestroy(&F);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = PetscViewerDestroy(&monP.viewer);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = PetscFinalize(); 181c4762a1bSJed Brown return ierr; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 184c4762a1bSJed Brown /* 185c4762a1bSJed Brown FormInitialGuess - Computes initial guess. 186c4762a1bSJed Brown 187c4762a1bSJed Brown Input/Output Parameter: 188c4762a1bSJed Brown . x - the solution vector 189c4762a1bSJed Brown */ 190c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x) 191c4762a1bSJed Brown { 192c4762a1bSJed Brown PetscErrorCode ierr; 193c4762a1bSJed Brown PetscScalar pfive = .50; 194c4762a1bSJed Brown ierr = VecSet(x,pfive);CHKERRQ(ierr); 195c4762a1bSJed Brown return 0; 196c4762a1bSJed Brown } 197c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 198c4762a1bSJed Brown /* 199c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 200c4762a1bSJed Brown 201c4762a1bSJed Brown Input Parameters: 202c4762a1bSJed Brown . snes - the SNES context 203c4762a1bSJed Brown . x - input vector 204c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction() 205c4762a1bSJed Brown 206c4762a1bSJed Brown Output Parameter: 207c4762a1bSJed Brown . f - function vector 208c4762a1bSJed Brown 209c4762a1bSJed Brown Note: 210c4762a1bSJed Brown The user-defined context can contain any application-specific data 211c4762a1bSJed Brown needed for the function evaluation (such as various parameters, work 212c4762a1bSJed Brown vectors, and grid information). In this program the context is just 213c4762a1bSJed Brown a vector containing the right-hand-side of the discretized PDE. 214c4762a1bSJed Brown */ 215c4762a1bSJed Brown 216c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx) 217c4762a1bSJed Brown { 218c4762a1bSJed Brown Vec g = (Vec)ctx; 219c4762a1bSJed Brown const PetscScalar *xx,*gg; 220c4762a1bSJed Brown PetscScalar *ff,d; 221c4762a1bSJed Brown PetscErrorCode ierr; 222c4762a1bSJed Brown PetscInt i,n; 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* 225c4762a1bSJed Brown Get pointers to vector data. 226c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 227c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 228c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 229c4762a1bSJed Brown the array. 230c4762a1bSJed Brown */ 231c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = VecGetArrayRead(g,&gg);CHKERRQ(ierr); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* 236c4762a1bSJed Brown Compute function 237c4762a1bSJed Brown */ 238c4762a1bSJed Brown ierr = VecGetSize(x,&n);CHKERRQ(ierr); 239c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 240c4762a1bSJed Brown ff[0] = xx[0]; 241c4762a1bSJed Brown for (i=1; i<n-1; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - gg[i]; 242c4762a1bSJed Brown ff[n-1] = xx[n-1] - 1.0; 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* 245c4762a1bSJed Brown Restore vectors 246c4762a1bSJed Brown */ 247c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = VecRestoreArrayRead(g,&gg);CHKERRQ(ierr); 250c4762a1bSJed Brown return 0; 251c4762a1bSJed Brown } 252c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 255c4762a1bSJed Brown 256c4762a1bSJed Brown Input Parameters: 257c4762a1bSJed Brown . snes - the SNES context 258c4762a1bSJed Brown . x - input vector 259c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 260c4762a1bSJed Brown 261c4762a1bSJed Brown Output Parameters: 262c4762a1bSJed Brown . jac - Jacobian matrix 263c4762a1bSJed Brown . B - optionally different preconditioning matrix 264c4762a1bSJed Brown 265c4762a1bSJed Brown */ 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 268c4762a1bSJed Brown { 269c4762a1bSJed Brown const PetscScalar *xx; 270c4762a1bSJed Brown PetscScalar A[3],d; 271c4762a1bSJed Brown PetscErrorCode ierr; 272c4762a1bSJed Brown PetscInt i,n,j[3]; 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* 275c4762a1bSJed Brown Get pointer to vector data 276c4762a1bSJed Brown */ 277c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 281c4762a1bSJed Brown - Note that in this case we set all elements for a particular 282c4762a1bSJed Brown row at once. 283c4762a1bSJed Brown */ 284c4762a1bSJed Brown ierr = VecGetSize(x,&n);CHKERRQ(ierr); 285c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 286c4762a1bSJed Brown 287c4762a1bSJed Brown /* 288c4762a1bSJed Brown Interior grid points 289c4762a1bSJed Brown */ 290c4762a1bSJed Brown for (i=1; i<n-1; i++) { 291c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 292c4762a1bSJed Brown A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 293c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* 297c4762a1bSJed Brown Boundary points 298c4762a1bSJed Brown */ 299c4762a1bSJed Brown i = 0; A[0] = 1.0; 300c4762a1bSJed Brown 301c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr); 302c4762a1bSJed Brown 303c4762a1bSJed Brown i = n-1; A[0] = 1.0; 304c4762a1bSJed Brown 305c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr); 306c4762a1bSJed Brown 307c4762a1bSJed Brown /* 308c4762a1bSJed Brown Restore vector 309c4762a1bSJed Brown */ 310c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* 313c4762a1bSJed Brown Assemble matrix 314c4762a1bSJed Brown */ 315c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 317c4762a1bSJed Brown if (jac != B) { 318c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 319c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown return 0; 322c4762a1bSJed Brown } 323c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 324c4762a1bSJed Brown /* 325c4762a1bSJed Brown Monitor - User-defined monitoring routine that views the 326c4762a1bSJed Brown current iterate with an x-window plot. 327c4762a1bSJed Brown 328c4762a1bSJed Brown Input Parameters: 329c4762a1bSJed Brown snes - the SNES context 330c4762a1bSJed Brown its - iteration number 331c4762a1bSJed Brown norm - 2-norm function value (may be estimated) 332c4762a1bSJed Brown ctx - optional user-defined context for private data for the 333c4762a1bSJed Brown monitor routine, as set by SNESMonitorSet() 334c4762a1bSJed Brown 335c4762a1bSJed Brown Note: 336c4762a1bSJed Brown See the manpage for PetscViewerDrawOpen() for useful runtime options, 337c4762a1bSJed Brown such as -nox to deactivate all x-window output. 338c4762a1bSJed Brown */ 339c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx) 340c4762a1bSJed Brown { 341c4762a1bSJed Brown PetscErrorCode ierr; 342c4762a1bSJed Brown MonitorCtx *monP = (MonitorCtx*) ctx; 343c4762a1bSJed Brown Vec x; 344c4762a1bSJed Brown 345c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"iter = %D, SNES Function norm %g\n",its,(double)fnorm);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr); 347c4762a1bSJed Brown ierr = VecView(x,monP->viewer);CHKERRQ(ierr); 348c4762a1bSJed Brown return 0; 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown /*TEST 352c4762a1bSJed Brown 353c4762a1bSJed Brown test: 354c4762a1bSJed Brown args: -nox -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always 355c4762a1bSJed Brown 356c4762a1bSJed Brown test: 357c4762a1bSJed Brown suffix: 2 358c4762a1bSJed Brown args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontr -snes_view 359c4762a1bSJed Brown requires: !single 360c4762a1bSJed Brown 361c4762a1bSJed Brown test: 362c4762a1bSJed Brown suffix: 3 363c4762a1bSJed Brown args: -nox -malloc no -options_left no -snes_monitor_cancel -snes_monitor_short -snes_view -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always 364c4762a1bSJed Brown 36541ba4c6cSHeeho Park test: 36641ba4c6cSHeeho Park suffix: 4 36741ba4c6cSHeeho Park args: -nox -snes_monitor_cancel -snes_monitor_short -snes_type newtontrdc -snes_view 36841ba4c6cSHeeho Park requires: !single 36941ba4c6cSHeeho Park 370c4762a1bSJed Brown TEST*/ 371