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