1c4421ceaSFande Kong 2c4421ceaSFande Kong static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\ 3c4421ceaSFande Kong This example employs a user-defined reasonview routine.\n\n"; 4c4421ceaSFande Kong 5c4421ceaSFande Kong /*T 6c4421ceaSFande Kong Concepts: SNES^basic uniprocessor example 7c4421ceaSFande Kong Concepts: SNES^setting a user-defined reasonview routine 8c4421ceaSFande Kong Processors: 1 9c4421ceaSFande Kong T*/ 10c4421ceaSFande Kong 11c4421ceaSFande Kong #include <petscsnes.h> 12c4421ceaSFande Kong 13c4421ceaSFande Kong /* 14c4421ceaSFande Kong User-defined routines 15c4421ceaSFande Kong */ 16c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 17c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 18c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode FormInitialGuess(Vec); 19c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MySNESConvergedReasonView(SNES,void*); 20c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode MyKSPConvergedReasonView(KSP,void*); 21c4421ceaSFande Kong 22c4421ceaSFande Kong /* 23c4421ceaSFande Kong User-defined context for monitoring 24c4421ceaSFande Kong */ 25c4421ceaSFande Kong typedef struct { 26c4421ceaSFande Kong PetscViewer viewer; 27c4421ceaSFande Kong } ReasonViewCtx; 28c4421ceaSFande Kong 29c4421ceaSFande Kong int main(int argc,char **argv) 30c4421ceaSFande Kong { 31c4421ceaSFande Kong SNES snes; /* SNES context */ 32c4421ceaSFande Kong KSP ksp; /* KSP context */ 33c4421ceaSFande Kong Vec x,r,F,U; /* vectors */ 34c4421ceaSFande Kong Mat J; /* Jacobian matrix */ 35c4421ceaSFande Kong ReasonViewCtx monP; /* monitoring context */ 36c4421ceaSFande Kong PetscErrorCode ierr; 37c4421ceaSFande Kong PetscInt its,n = 5,i; 38c4421ceaSFande Kong PetscMPIInt size; 39c4421ceaSFande Kong PetscScalar h,xp,v; 40c4421ceaSFande Kong MPI_Comm comm; 41c4421ceaSFande Kong 42c4421ceaSFande Kong ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 43c4421ceaSFande Kong ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 44*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 45c4421ceaSFande Kong ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 46c4421ceaSFande Kong h = 1.0/(n-1); 47c4421ceaSFande Kong comm = PETSC_COMM_WORLD; 48c4421ceaSFande Kong /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49c4421ceaSFande Kong Create nonlinear solver context 50c4421ceaSFande Kong - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 51c4421ceaSFande Kong 52c4421ceaSFande Kong ierr = SNESCreate(comm,&snes);CHKERRQ(ierr); 53c4421ceaSFande Kong 54c4421ceaSFande Kong /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 55c4421ceaSFande Kong Create vector data structures; set function evaluation routine 56c4421ceaSFande Kong - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 57c4421ceaSFande Kong /* 58c4421ceaSFande Kong Note that we form 1 vector from scratch and then duplicate as needed. 59c4421ceaSFande Kong */ 60c4421ceaSFande Kong ierr = VecCreate(comm,&x);CHKERRQ(ierr); 61c4421ceaSFande Kong ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr); 62c4421ceaSFande Kong ierr = VecSetFromOptions(x);CHKERRQ(ierr); 63c4421ceaSFande Kong ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 64c4421ceaSFande Kong ierr = VecDuplicate(x,&F);CHKERRQ(ierr); 65c4421ceaSFande Kong ierr = VecDuplicate(x,&U);CHKERRQ(ierr); 66c4421ceaSFande Kong 67c4421ceaSFande Kong /* 68c4421ceaSFande Kong Set function evaluation routine and vector 69c4421ceaSFande Kong */ 70c4421ceaSFande Kong ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr); 71c4421ceaSFande Kong 72c4421ceaSFande Kong /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73c4421ceaSFande Kong Create matrix data structure; set Jacobian evaluation routine 74c4421ceaSFande Kong - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 75c4421ceaSFande Kong 76c4421ceaSFande Kong ierr = MatCreate(comm,&J);CHKERRQ(ierr); 77c4421ceaSFande Kong ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); 78c4421ceaSFande Kong ierr = MatSetFromOptions(J);CHKERRQ(ierr); 79c4421ceaSFande Kong ierr = MatSeqAIJSetPreallocation(J,3,NULL);CHKERRQ(ierr); 80c4421ceaSFande Kong 81c4421ceaSFande Kong /* 82c4421ceaSFande Kong Set Jacobian matrix data structure and default Jacobian evaluation 83c4421ceaSFande Kong routine. User can override with: 84c4421ceaSFande Kong -snes_fd : default finite differencing approximation of Jacobian 85c4421ceaSFande Kong -snes_mf : matrix-free Newton-Krylov method with no preconditioning 86c4421ceaSFande Kong (unless user explicitly sets preconditioner) 87c4421ceaSFande Kong -snes_mf_operator : form preconditioning matrix as set by the user, 88c4421ceaSFande Kong but use matrix-free approx for Jacobian-vector 89c4421ceaSFande Kong products within Newton-Krylov method 90c4421ceaSFande Kong */ 91c4421ceaSFande Kong 92c4421ceaSFande Kong ierr = SNESSetJacobian(snes,J,J,FormJacobian,NULL);CHKERRQ(ierr); 93c4421ceaSFande Kong 94c4421ceaSFande Kong /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4421ceaSFande Kong Customize nonlinear solver; set runtime options 96c4421ceaSFande Kong - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 97c4421ceaSFande Kong 98c4421ceaSFande Kong /* 99c4421ceaSFande Kong Set an optional user-defined reasonview routine 100c4421ceaSFande Kong */ 101c4421ceaSFande Kong ierr = PetscViewerASCIIGetStdout(comm,&monP.viewer);CHKERRQ(ierr); 102c4421ceaSFande Kong /* Just make sure we can not repeat addding the same function 103c4421ceaSFande Kong * PETSc will be able to igore the repeated function 104c4421ceaSFande Kong */ 105c4421ceaSFande Kong for (i=0; i<4; i++) { 106c4421ceaSFande Kong ierr = SNESConvergedReasonViewSet(snes,MySNESConvergedReasonView,&monP,0);CHKERRQ(ierr); 107c4421ceaSFande Kong } 108c4421ceaSFande Kong ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 109c4421ceaSFande Kong /* Just make sure we can not repeat addding the same function 110c4421ceaSFande Kong * PETSc will be able to igore the repeated function 111c4421ceaSFande Kong */ 112c4421ceaSFande Kong for (i=0; i<4; i++) { 113c4421ceaSFande Kong ierr = KSPConvergedReasonViewSet(ksp,MyKSPConvergedReasonView,&monP,0);CHKERRQ(ierr); 114c4421ceaSFande Kong } 115c4421ceaSFande Kong /* 116c4421ceaSFande Kong Set SNES/KSP/KSP/PC runtime options, e.g., 117c4421ceaSFande Kong -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 118c4421ceaSFande Kong */ 119c4421ceaSFande Kong ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 120c4421ceaSFande Kong 121c4421ceaSFande Kong /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4421ceaSFande Kong Initialize application: 123c4421ceaSFande Kong Store right-hand-side of PDE and exact solution 124c4421ceaSFande Kong - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125c4421ceaSFande Kong 126c4421ceaSFande Kong xp = 0.0; 127c4421ceaSFande Kong for (i=0; i<n; i++) { 128c4421ceaSFande Kong v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 129c4421ceaSFande Kong ierr = VecSetValues(F,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 130c4421ceaSFande Kong v = xp*xp*xp; 131c4421ceaSFande Kong ierr = VecSetValues(U,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 132c4421ceaSFande Kong xp += h; 133c4421ceaSFande Kong } 134c4421ceaSFande Kong 135c4421ceaSFande Kong /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4421ceaSFande Kong Evaluate initial guess; then solve nonlinear system 137c4421ceaSFande Kong - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4421ceaSFande Kong /* 139c4421ceaSFande Kong Note: The user should initialize the vector, x, with the initial guess 140c4421ceaSFande Kong for the nonlinear solver prior to calling SNESSolve(). In particular, 141c4421ceaSFande Kong to employ an initial guess of zero, the user should explicitly set 142c4421ceaSFande Kong this vector to zero by calling VecSet(). 143c4421ceaSFande Kong */ 144c4421ceaSFande Kong ierr = FormInitialGuess(x);CHKERRQ(ierr); 145c4421ceaSFande Kong ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 146c4421ceaSFande Kong ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 147c4421ceaSFande Kong 148c4421ceaSFande Kong /* 149c4421ceaSFande Kong Free work space. All PETSc objects should be destroyed when they 150c4421ceaSFande Kong are no longer needed. 151c4421ceaSFande Kong */ 152c4421ceaSFande Kong ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 153c4421ceaSFande Kong ierr = VecDestroy(&U);CHKERRQ(ierr); ierr = VecDestroy(&F);CHKERRQ(ierr); 154c4421ceaSFande Kong ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 155c4421ceaSFande Kong /*ierr = PetscViewerDestroy(&monP.viewer);CHKERRQ(ierr);*/ 156c4421ceaSFande Kong ierr = PetscFinalize(); 157c4421ceaSFande Kong return ierr; 158c4421ceaSFande Kong } 159c4421ceaSFande Kong /* ------------------------------------------------------------------- */ 160c4421ceaSFande Kong /* 161c4421ceaSFande Kong FormInitialGuess - Computes initial guess. 162c4421ceaSFande Kong 163c4421ceaSFande Kong Input/Output Parameter: 164c4421ceaSFande Kong . x - the solution vector 165c4421ceaSFande Kong */ 166c4421ceaSFande Kong PetscErrorCode FormInitialGuess(Vec x) 167c4421ceaSFande Kong { 168c4421ceaSFande Kong PetscErrorCode ierr; 169c4421ceaSFande Kong PetscScalar pfive = .50; 170c4421ceaSFande Kong ierr = VecSet(x,pfive);CHKERRQ(ierr); 171c4421ceaSFande Kong return 0; 172c4421ceaSFande Kong } 173c4421ceaSFande Kong /* ------------------------------------------------------------------- */ 174c4421ceaSFande Kong /* 175c4421ceaSFande Kong FormFunction - Evaluates nonlinear function, F(x). 176c4421ceaSFande Kong 177c4421ceaSFande Kong Input Parameters: 178c4421ceaSFande Kong . snes - the SNES context 179c4421ceaSFande Kong . x - input vector 180c4421ceaSFande Kong . ctx - optional user-defined context, as set by SNESSetFunction() 181c4421ceaSFande Kong 182c4421ceaSFande Kong Output Parameter: 183c4421ceaSFande Kong . f - function vector 184c4421ceaSFande Kong 185c4421ceaSFande Kong Note: 186c4421ceaSFande Kong The user-defined context can contain any application-specific data 187c4421ceaSFande Kong needed for the function evaluation (such as various parameters, work 188c4421ceaSFande Kong vectors, and grid information). In this program the context is just 189c4421ceaSFande Kong a vector containing the right-hand-side of the discretized PDE. 190c4421ceaSFande Kong */ 191c4421ceaSFande Kong 192c4421ceaSFande Kong PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx) 193c4421ceaSFande Kong { 194c4421ceaSFande Kong Vec g = (Vec)ctx; 195c4421ceaSFande Kong const PetscScalar *xx,*gg; 196c4421ceaSFande Kong PetscScalar *ff,d; 197c4421ceaSFande Kong PetscErrorCode ierr; 198c4421ceaSFande Kong PetscInt i,n; 199c4421ceaSFande Kong 200c4421ceaSFande Kong /* 201c4421ceaSFande Kong Get pointers to vector data. 202c4421ceaSFande Kong - For default PETSc vectors, VecGetArray() returns a pointer to 203c4421ceaSFande Kong the data array. Otherwise, the routine is implementation dependent. 204c4421ceaSFande Kong - You MUST call VecRestoreArray() when you no longer need access to 205c4421ceaSFande Kong the array. 206c4421ceaSFande Kong */ 207c4421ceaSFande Kong ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 208c4421ceaSFande Kong ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 209c4421ceaSFande Kong ierr = VecGetArrayRead(g,&gg);CHKERRQ(ierr); 210c4421ceaSFande Kong 211c4421ceaSFande Kong /* 212c4421ceaSFande Kong Compute function 213c4421ceaSFande Kong */ 214c4421ceaSFande Kong ierr = VecGetSize(x,&n);CHKERRQ(ierr); 215c4421ceaSFande Kong d = (PetscReal)(n - 1); d = d*d; 216c4421ceaSFande Kong ff[0] = xx[0]; 217c4421ceaSFande Kong 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]; 218c4421ceaSFande Kong ff[n-1] = xx[n-1] - 1.0; 219c4421ceaSFande Kong 220c4421ceaSFande Kong /* 221c4421ceaSFande Kong Restore vectors 222c4421ceaSFande Kong */ 223c4421ceaSFande Kong ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 224c4421ceaSFande Kong ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 225c4421ceaSFande Kong ierr = VecRestoreArrayRead(g,&gg);CHKERRQ(ierr); 226c4421ceaSFande Kong return 0; 227c4421ceaSFande Kong } 228c4421ceaSFande Kong /* ------------------------------------------------------------------- */ 229c4421ceaSFande Kong /* 230c4421ceaSFande Kong FormJacobian - Evaluates Jacobian matrix. 231c4421ceaSFande Kong 232c4421ceaSFande Kong Input Parameters: 233c4421ceaSFande Kong . snes - the SNES context 234c4421ceaSFande Kong . x - input vector 235c4421ceaSFande Kong . dummy - optional user-defined context (not used here) 236c4421ceaSFande Kong 237c4421ceaSFande Kong Output Parameters: 238c4421ceaSFande Kong . jac - Jacobian matrix 239c4421ceaSFande Kong . B - optionally different preconditioning matrix 240c4421ceaSFande Kong 241c4421ceaSFande Kong */ 242c4421ceaSFande Kong 243c4421ceaSFande Kong PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 244c4421ceaSFande Kong { 245c4421ceaSFande Kong const PetscScalar *xx; 246c4421ceaSFande Kong PetscScalar A[3],d; 247c4421ceaSFande Kong PetscErrorCode ierr; 248c4421ceaSFande Kong PetscInt i,n,j[3]; 249c4421ceaSFande Kong 250c4421ceaSFande Kong /* 251c4421ceaSFande Kong Get pointer to vector data 252c4421ceaSFande Kong */ 253c4421ceaSFande Kong ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 254c4421ceaSFande Kong 255c4421ceaSFande Kong /* 256c4421ceaSFande Kong Compute Jacobian entries and insert into matrix. 257c4421ceaSFande Kong - Note that in this case we set all elements for a particular 258c4421ceaSFande Kong row at once. 259c4421ceaSFande Kong */ 260c4421ceaSFande Kong ierr = VecGetSize(x,&n);CHKERRQ(ierr); 261c4421ceaSFande Kong d = (PetscReal)(n - 1); d = d*d; 262c4421ceaSFande Kong 263c4421ceaSFande Kong /* 264c4421ceaSFande Kong Interior grid points 265c4421ceaSFande Kong */ 266c4421ceaSFande Kong for (i=1; i<n-1; i++) { 267c4421ceaSFande Kong j[0] = i - 1; j[1] = i; j[2] = i + 1; 268c4421ceaSFande Kong A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 269c4421ceaSFande Kong ierr = MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr); 270c4421ceaSFande Kong } 271c4421ceaSFande Kong 272c4421ceaSFande Kong /* 273c4421ceaSFande Kong Boundary points 274c4421ceaSFande Kong */ 275c4421ceaSFande Kong i = 0; A[0] = 1.0; 276c4421ceaSFande Kong 277c4421ceaSFande Kong ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr); 278c4421ceaSFande Kong 279c4421ceaSFande Kong i = n-1; A[0] = 1.0; 280c4421ceaSFande Kong 281c4421ceaSFande Kong ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr); 282c4421ceaSFande Kong 283c4421ceaSFande Kong /* 284c4421ceaSFande Kong Restore vector 285c4421ceaSFande Kong */ 286c4421ceaSFande Kong ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 287c4421ceaSFande Kong 288c4421ceaSFande Kong /* 289c4421ceaSFande Kong Assemble matrix 290c4421ceaSFande Kong */ 291c4421ceaSFande Kong ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 292c4421ceaSFande Kong ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 293c4421ceaSFande Kong if (jac != B) { 294c4421ceaSFande Kong ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 295c4421ceaSFande Kong ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 296c4421ceaSFande Kong } 297c4421ceaSFande Kong return 0; 298c4421ceaSFande Kong } 299c4421ceaSFande Kong 300c4421ceaSFande Kong PetscErrorCode MySNESConvergedReasonView(SNES snes,void *ctx) 301c4421ceaSFande Kong { 302c4421ceaSFande Kong PetscErrorCode ierr; 303c4421ceaSFande Kong ReasonViewCtx *monP = (ReasonViewCtx*) ctx; 304c4421ceaSFande Kong PetscViewer viewer = monP->viewer; 305c4421ceaSFande Kong SNESConvergedReason reason; 306c4421ceaSFande Kong const char *strreason; 307c4421ceaSFande Kong 308c4421ceaSFande Kong ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr); 309c4421ceaSFande Kong ierr = SNESGetConvergedReasonString(snes,&strreason);CHKERRQ(ierr); 310c4421ceaSFande Kong ierr = PetscViewerASCIIPrintf(viewer,"Customized SNES converged reason view\n");CHKERRQ(ierr); 311c4421ceaSFande Kong ierr = PetscViewerASCIIAddTab(viewer,1);CHKERRQ(ierr); 312c4421ceaSFande Kong if (reason > 0) { 313c4421ceaSFande Kong ierr = PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",strreason);CHKERRQ(ierr); 314c4421ceaSFande Kong } else if (reason <= 0) { 315c4421ceaSFande Kong ierr = PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",strreason);CHKERRQ(ierr); 316c4421ceaSFande Kong } 317c4421ceaSFande Kong ierr = PetscViewerASCIISubtractTab(viewer,1);CHKERRQ(ierr); 318c4421ceaSFande Kong return 0; 319c4421ceaSFande Kong } 320c4421ceaSFande Kong 321c4421ceaSFande Kong PetscErrorCode MyKSPConvergedReasonView(KSP ksp,void *ctx) 322c4421ceaSFande Kong { 323c4421ceaSFande Kong PetscErrorCode ierr; 324c4421ceaSFande Kong ReasonViewCtx *monP = (ReasonViewCtx*) ctx; 325c4421ceaSFande Kong PetscViewer viewer = monP->viewer; 326c4421ceaSFande Kong KSPConvergedReason reason; 327c4421ceaSFande Kong const char *reasonstr; 328c4421ceaSFande Kong 329c4421ceaSFande Kong ierr = KSPGetConvergedReason(ksp,&reason);CHKERRQ(ierr); 330c4421ceaSFande Kong ierr = KSPGetConvergedReasonString(ksp,&reasonstr);CHKERRQ(ierr); 331c4421ceaSFande Kong ierr = PetscViewerASCIIAddTab(viewer,2);CHKERRQ(ierr); 332c4421ceaSFande Kong ierr = PetscViewerASCIIPrintf(viewer,"Customized KSP converged reason view\n");CHKERRQ(ierr); 333c4421ceaSFande Kong ierr = PetscViewerASCIIAddTab(viewer,1);CHKERRQ(ierr); 334c4421ceaSFande Kong if (reason > 0) { 335c4421ceaSFande Kong ierr = PetscViewerASCIIPrintf(viewer,"Converged due to %s\n",reasonstr);CHKERRQ(ierr); 336c4421ceaSFande Kong } else if (reason <= 0) { 337c4421ceaSFande Kong ierr = PetscViewerASCIIPrintf(viewer,"Did not converge due to %s\n",reasonstr);CHKERRQ(ierr); 338c4421ceaSFande Kong } 339c4421ceaSFande Kong ierr = PetscViewerASCIISubtractTab(viewer,3);CHKERRQ(ierr); 340c4421ceaSFande Kong return 0; 341c4421ceaSFande Kong } 342c4421ceaSFande Kong 343c4421ceaSFande Kong /*TEST 344c4421ceaSFande Kong 345c4421ceaSFande Kong test: 346c4421ceaSFande Kong suffix: 1 347c4421ceaSFande Kong nsize: 1 348bc44df21SStefano Zampini filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" 349c4421ceaSFande Kong 350c4421ceaSFande Kong test: 351c4421ceaSFande Kong suffix: 2 352c4421ceaSFande Kong nsize: 1 353c4421ceaSFande Kong args: -ksp_converged_reason_view_cancel 354bc44df21SStefano Zampini filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" 355c4421ceaSFande Kong 356c4421ceaSFande Kong test: 357c4421ceaSFande Kong suffix: 3 358c4421ceaSFande Kong nsize: 1 359c4421ceaSFande Kong args: -ksp_converged_reason_view_cancel -ksp_converged_reason 360bc44df21SStefano Zampini filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" 361c4421ceaSFande Kong 362c4421ceaSFande Kong test: 363c4421ceaSFande Kong suffix: 4 364c4421ceaSFande Kong nsize: 1 365c4421ceaSFande Kong args: -snes_converged_reason_view_cancel 366bc44df21SStefano Zampini filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" 367c4421ceaSFande Kong 368c4421ceaSFande Kong test: 369c4421ceaSFande Kong suffix: 5 370c4421ceaSFande Kong nsize: 1 371c4421ceaSFande Kong args: -snes_converged_reason_view_cancel -snes_converged_reason 372bc44df21SStefano Zampini filter: sed -e "s/CONVERGED_ATOL/CONVERGED_RTOL/g" 373c4421ceaSFande Kong 374c4421ceaSFande Kong TEST*/ 375