Actual source code: ex1f.F90

  1: !
  2: !  Description: This example solves a nonlinear system on 1 processor with SNES.
  3: !  We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
  4: !  domain.  The command line options include:
  5: !    -par <parameter>, where <parameter> indicates the nonlinearity of the problem
  6: !       problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)
  7: !    -mx <xg>, where <xg> = number of grid points in the x-direction
  8: !    -my <yg>, where <yg> = number of grid points in the y-direction
  9: !

 11: !
 12: !  --------------------------------------------------------------------------
 13: !
 14: !  Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 15: !  the partial differential equation
 16: !
 17: !          -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 18: !
 19: !  with boundary conditions
 20: !
 21: !           u = 0  for  x = 0, x = 1, y = 0, y = 1.
 22: !
 23: !  A finite difference approximation with the usual 5-point stencil
 24: !  is used to discretize the boundary value problem to obtain a nonlinear
 25: !  system of equations.
 26: !
 27: !  The parallel version of this code is snes/tutorials/ex5f.F
 28: !
 29: !  --------------------------------------------------------------------------
 30: #include <petsc/finclude/petscsnes.h>
 31: #include <petsc/finclude/petscdraw.h>
 32: module ex1fmodule
 33:   use petscsnes
 34:   implicit none
 35:   PetscReal lambda
 36:   PetscInt mx, my
 37:   PetscBool fd_coloring
 38: contains
 39:   subroutine postcheck(snes, x, y, w, changed_y, changed_w, ctx, ierr)
 40:     SNES snes
 41:     PetscReal norm
 42:     Vec tmp, x, y, w
 43:     PetscBool changed_w, changed_y
 44:     PetscErrorCode ierr
 45:     PetscInt ctx
 46:     PetscScalar, parameter :: mone = -1.0
 47:     MPIU_Comm comm

 49:     character(len=PETSC_MAX_PATH_LEN) :: outputString

 51:     PetscCallA(PetscObjectGetComm(snes, comm, ierr))
 52:     PetscCallA(VecDuplicate(x, tmp, ierr))
 53:     PetscCallA(VecWAXPY(tmp, mone, x, w, ierr))
 54:     PetscCallA(VecNorm(tmp, NORM_2, norm, ierr))
 55:     PetscCallA(VecDestroy(tmp, ierr))
 56:     write (outputString, *) norm
 57:     PetscCallA(PetscPrintf(comm, 'Norm of search step '//trim(outputString)//'\n', ierr))
 58:   end

 60: ! ---------------------------------------------------------------------
 61: !
 62: !  FormInitialGuess - Forms initial approximation.
 63: !
 64: !  Input Parameter:
 65: !  X - vector
 66: !
 67: !  Output Parameters:
 68: !  X - vector
 69: !  ierr - error code
 70: !
 71: !  Notes:
 72: !  This routine serves as a wrapper for the lower-level routine
 73: !  "ApplicationInitialGuess", where the actual computations are
 74: !  done using the standard Fortran style of treating the local
 75: !  vector data as a multidimensional array over the local mesh.
 76: !  This routine merely accesses the local vector data via
 77: !  VecGetArray() and VecRestoreArray().
 78: !
 79:   subroutine FormInitialGuess(X, ierr)
 80:     Vec X
 81:     PetscErrorCode, intent(out) :: ierr
 82: !   Declarations for use with local arrays:
 83:     PetscScalar, pointer :: lx_v(:)

 85: !   Get a pointer to vector data.
 86: !     - VecGetArray() returns a pointer to the data array.
 87: !     - You MUST call VecRestoreArray() when you no longer need access to
 88: !       the array.
 89:     PetscCallA(VecGetArray(X, lx_v, ierr))

 91: !   Compute initial guess
 92:     PetscCallA(ApplicationInitialGuess(lx_v, ierr))

 94: !   Restore vector
 95:     PetscCallA(VecRestoreArray(X, lx_v, ierr))

 97:   end

 99: !  ApplicationInitialGuess - Computes initial approximation, called by
100: !  the higher level routine FormInitialGuess().
101: !
102: !  Input Parameter:
103: !  x - local vector data
104: !
105: !  Output Parameters:
106: !  f - local vector data, f(x)
107: !  ierr - error code
108: !
109: !  Notes:
110: !  This routine uses standard Fortran-style computations over a 2-dim array.
111: !
112:   subroutine ApplicationInitialGuess(x, ierr)
113:     PetscScalar, intent(out) :: x(mx, my)
114:     PetscErrorCode, intent(out) :: ierr
115: !   Local variables:
116:     PetscInt i, j
117:     PetscReal temp1, temp, hx, hy

119: !   Set parameters
120:     hx = 1.0_PETSC_REAL_KIND/(mx - 1)
121:     hy = 1.0_PETSC_REAL_KIND/(my - 1)
122:     temp1 = lambda/(lambda + 1.0_PETSC_REAL_KIND)

124:     do j = 1, my
125:       temp = min(j - 1, my - j)*hy
126:       do i = 1, mx
127:         if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then
128:           x(i, j) = 0.0
129:         else
130:           x(i, j) = temp1*sqrt(min(min(i - 1, mx - i)*hx, temp))
131:         end if
132:       end do
133:     end do
134:     ierr = 0
135:   end

137: ! ---------------------------------------------------------------------
138: !
139: !  FormFunction - Evaluates nonlinear function, F(x).
140: !
141: !  Input Parameters:
142: !  snes  - the SNES context
143: !  X     - input vector
144: !  dummy - optional user-defined context, as set by SNESSetFunction()
145: !          (not used here)
146: !
147: !  Output Parameter:
148: !  F     - vector with newly computed function
149: !
150: !  Notes:
151: !  This routine serves as a wrapper for the lower-level routine
152: !  "ApplicationFunction", where the actual computations are
153: !  done using the standard Fortran style of treating the local
154: !  vector data as a multidimensional array over the local mesh.
155: !  This routine merely accesses the local vector data via
156: !  VecGetArray() and VecRestoreArray().
157: !
158:   subroutine FormFunction(snes, X, F, fdcoloring, ierr)
159: !   Input/output variables:
160:     SNES snes
161:     Vec X, F
162:     PetscErrorCode, intent(out) :: ierr
163:     MatFDColoring fdcoloring
164: !   Declarations for use with local arrays:
165:     PetscScalar, pointer :: lx_v(:), lf_v(:)
166:     PetscInt, pointer :: indices(:)

168: !  Get pointers to vector data.
169: !    - VecGetArray() returns a pointer to the data array.
170: !    - You MUST call VecRestoreArray() when you no longer need access to
171: !      the array.
172:     PetscCallA(VecGetArrayRead(X, lx_v, ierr))
173:     PetscCallA(VecGetArray(F, lf_v, ierr))

175: !   Compute function
176:     PetscCallA(ApplicationFunction(lx_v, lf_v, ierr))

178: !   Restore vectors
179:     PetscCallA(VecRestoreArrayRead(X, lx_v, ierr))
180:     PetscCallA(VecRestoreArray(F, lf_v, ierr))

182:     PetscCallA(PetscLogFlops(11.0d0*mx*my, ierr))
183: !
184: !   fd_coloring is a global variable used here ONLY to test the
185: !   calls to MatFDColoringGetPerturbedColumns() and  MatFDColoringRestorePerturbedColumns()
186:     if (fd_coloring) then
187:       PetscCallA(MatFDColoringGetPerturbedColumns(fdcoloring, PETSC_NULL_INTEGER, indices, ierr))
188:       print *, 'Indices from GetPerturbedColumns'
189:       write (*, 1000) indices
190: 1000  format(50i4)
191:       PetscCallA(MatFDColoringRestorePerturbedColumns(fdcoloring, PETSC_NULL_INTEGER, indices, ierr))
192:     end if
193:   end

195: ! ---------------------------------------------------------------------
196: !
197: !  ApplicationFunction - Computes nonlinear function, called by
198: !  the higher level routine FormFunction().
199: !
200: !  Input Parameter:
201: !  x    - local vector data
202: !
203: !  Output Parameters:
204: !  f    - local vector data, f(x)
205: !  ierr - error code
206: !
207: !  Notes:
208: !  This routine uses standard Fortran-style computations over a 2-dim array.
209: !
210:   subroutine ApplicationFunction(x, f, ierr)
211:     PetscScalar, intent(in) :: x(mx, my)
212:     PetscScalar, intent(out) :: f(mx, my)
213:     PetscErrorCode, intent(out) :: ierr
214: !   Local variables:
215:     PetscScalar, parameter :: one = 1.0, two = 2.0
216:     PetscScalar hx, hy, hxdhy, hydhx, sc
217:     PetscScalar u, uxx, uyy
218:     PetscInt i, j

220:     hx = one/(mx - 1)
221:     hy = one/(my - 1)
222:     sc = hx*hy*lambda
223:     hxdhy = hx/hy
224:     hydhx = hy/hx

226: !  Compute function

228:     do j = 1, my
229:       do i = 1, mx
230:         if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then
231:           f(i, j) = x(i, j)
232:         else
233:           u = x(i, j)
234:           uxx = hydhx*(two*u - x(i - 1, j) - x(i + 1, j))
235:           uyy = hxdhy*(two*u - x(i, j - 1) - x(i, j + 1))
236:           f(i, j) = uxx + uyy - sc*exp(u)
237:         end if
238:       end do
239:     end do
240:     ierr = 0
241:   end

243: ! ---------------------------------------------------------------------
244: !
245: !  FormJacobian - Evaluates Jacobian matrix.
246: !
247: !  Input Parameters:
248: !  snes    - the SNES context
249: !  x       - input vector
250: !  dummy   - optional user-defined context, as set by SNESSetJacobian()
251: !            (not used here)
252: !
253: !  Output Parameters:
254: !  jac      - Jacobian matrix
255: !  jac_prec - optionally different matrix used to construct the preconditioner (not used here)
256: !
257: !  Notes:
258: !  This routine serves as a wrapper for the lower-level routine
259: !  "ApplicationJacobian", where the actual computations are
260: !  done using the standard Fortran style of treating the local
261: !  vector data as a multidimensional array over the local mesh.
262: !  This routine merely accesses the local vector data via
263: !  VecGetArray() and VecRestoreArray().
264: !
265:   subroutine FormJacobian(snes, X, jac, jac_prec, dummy, ierr)
266: !   Input/output variables:
267:     SNES snes
268:     Vec X
269:     Mat jac, jac_prec
270:     PetscErrorCode, intent(out) :: ierr
271:     integer dummy
272: !   Declarations for use with local array:
273:     PetscScalar, pointer :: lx_v(:)

275: !   Get a pointer to vector data
276:     PetscCallA(VecGetArrayRead(X, lx_v, ierr))

278: !   Compute Jacobian entries
279:     PetscCallA(ApplicationJacobian(lx_v, jac, jac_prec, ierr))

281: !   Restore vector
282:     PetscCallA(VecRestoreArrayRead(X, lx_v, ierr))

284: !   Assemble matrix
285:     PetscCallA(MatAssemblyBegin(jac_prec, MAT_FINAL_ASSEMBLY, ierr))
286:     PetscCallA(MatAssemblyEnd(jac_prec, MAT_FINAL_ASSEMBLY, ierr))
287:   end

289: ! ---------------------------------------------------------------------
290: !
291: !  ApplicationJacobian - Computes Jacobian matrix, called by
292: !  the higher level routine FormJacobian().
293: !
294: !  Input Parameters:
295: !  x        - local vector data
296: !
297: !  Output Parameters:
298: !  jac      - Jacobian matrix
299: !  jac_prec - optionally different matrix used to construct the preconditioner (not used here)
300: !  ierr     - error code
301: !
302: !  Notes:
303: !  This routine uses standard Fortran-style computations over a 2-dim array.
304: !
305:   subroutine ApplicationJacobian(x, jac, jac_prec, ierr)
306: !   Input/output variables:
307:     PetscScalar, intent(in) :: x(mx, my)
308:     Mat jac, jac_prec
309:     PetscErrorCode, intent(out) :: ierr
310: !   Local variables:
311:     PetscInt i, j, row(1), col(5)
312:     PetscScalar, parameter :: one = 1.0, two = 2.0
313:     PetscScalar hx, hy, hxdhy, hydhx, sc, v(5)

315: !  Set parameters
316:     hx = one/(mx - 1)
317:     hy = one/(my - 1)
318:     sc = hx*hy
319:     hxdhy = hx/hy
320:     hydhx = hy/hx

322: !  Compute entries of the Jacobian matrix
323: !   - Here, we set all entries for a particular row at once.
324: !   - Note that MatSetValues() uses 0-based row and column numbers
325: !     in Fortran as well as in C.

327:     do j = 1, my
328:       row(1) = (j - 1)*mx - 1
329:       do i = 1, mx
330:         row(1) = row(1) + 1
331: !           boundary points
332:         if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then
333:           PetscCallA(MatSetValues(jac_prec, 1_PETSC_INT_KIND, row, 1_PETSC_INT_KIND, row, [one], INSERT_VALUES, ierr))
334: !           interior grid points
335:         else
336:           v(1) = -hxdhy
337:           v(2) = -hydhx
338:           v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i, j))
339:           v(4) = -hydhx
340:           v(5) = -hxdhy
341:           col(1) = row(1) - mx
342:           col(2) = row(1) - 1
343:           col(3) = row(1)
344:           col(4) = row(1) + 1
345:           col(5) = row(1) + mx
346:           PetscCallA(MatSetValues(jac_prec, 1_PETSC_INT_KIND, row, 5_PETSC_INT_KIND, col, v, INSERT_VALUES, ierr))
347:         end if
348:       end do
349:     end do

351:   end

353: end module ex1fmodule
354: program main
355:   use petscdraw
356:   use petscsnes
357:   use ex1fmodule
358:   implicit none
359: !
360: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
361: !                   Variable declarations
362: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
363: !
364: !  Variables:
365: !     snes        - nonlinear solver
366: !     x, r        - solution, residual vectors
367: !     J           - Jacobian matrix
368: !     its         - iterations for convergence
369: !     matrix_free - flag - 1 indicates matrix-free version
370: !     lambda      - nonlinearity parameter
371: !     draw        - drawing context
372: !
373:   SNES snes
374:   MatColoring mc
375:   Vec x, r
376:   PetscDraw draw
377:   Mat J
378:   PetscBool matrix_free, flg
379:   PetscErrorCode ierr
380:   PetscInt its, N
381:   MatFDColoring fdcoloring
382:   PetscMPIInt size, rank
383:   PetscReal, parameter :: lambda_min = 0.0_PETSC_REAL_KIND, lambda_max = 6.81_PETSC_REAL_KIND
384:   ISColoring iscoloring
385:   PetscBool pc
386:   integer4 imx, imy
387:   character(len=PETSC_MAX_PATH_LEN) :: outputString
388:   PetscScalar, pointer :: lx_v(:)
389:   integer4, parameter :: xl = 300, yl = 0, width = 300, height = 300

391: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
392: !  Initialize program
393: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

395:   PetscCallA(PetscInitialize(ierr))
396:   PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
397:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))

399:   PetscCheckA(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, 'This is a uniprocessor example only')

401: !  Initialize problem parameters
402:   mx = 4
403:   PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-mx', mx, flg, ierr))
404:   my = 4
405:   PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-my', my, flg, ierr))
406:   lambda = 6.0
407:   PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-par', lambda, flg, ierr))
408:   PetscCheckA(lambda < lambda_max .and. lambda > lambda_min, PETSC_COMM_SELF, PETSC_ERR_USER, 'Lambda out of range ')
409:   N = mx*my
410:   pc = PETSC_FALSE
411:   PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-pc', pc, PETSC_NULL_BOOL, ierr))

413: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
414: !  Create nonlinear solver context
415: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

417:   PetscCallA(SNESCreate(PETSC_COMM_WORLD, snes, ierr))

419:   if (pc .eqv. PETSC_TRUE) then
420:     PetscCallA(SNESSetType(snes, SNESNEWTONTR, ierr))
421:     PetscCallA(SNESNewtonTRSetPostCheck(snes, postcheck, snes, ierr))
422:   end if

424: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
425: !  Create vector data structures; set function evaluation routine
426: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

428:   PetscCallA(VecCreate(PETSC_COMM_WORLD, x, ierr))
429:   PetscCallA(VecSetSizes(x, PETSC_DECIDE, N, ierr))
430:   PetscCallA(VecSetFromOptions(x, ierr))
431:   PetscCallA(VecDuplicate(x, r, ierr))

433: !  Set function evaluation routine and vector.  Whenever the nonlinear
434: !  solver needs to evaluate the nonlinear function, it will call this
435: !  routine.
436: !   - Note that the final routine argument is the user-defined
437: !     context that provides application-specific data for the
438: !     function evaluation routine.

440:   PetscCallA(SNESSetFunction(snes, r, FormFunction, fdcoloring, ierr))

442: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
443: !  Create matrix data structure; set Jacobian evaluation routine
444: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

446: !  Create matrix. Here we only approximately preallocate storage space
447: !  for the Jacobian.  See the users manual for a discussion of better
448: !  techniques for preallocating matrix memory.

450:   PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-snes_mf', matrix_free, ierr))
451:   if (.not. matrix_free) then
452:     PetscCallA(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, J, ierr))
453:   end if

455: !
456: !     This option will cause the Jacobian to be computed via finite differences
457: !    efficiently using a coloring of the columns of the matrix.
458: !
459:   fd_coloring = .false.
460:   PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-snes_fd_coloring', fd_coloring, ierr))
461:   if (fd_coloring) then

463: !
464: !      This initializes the nonzero structure of the Jacobian. This is artificial
465: !      because clearly if we had a routine to compute the Jacobian we won't need
466: !      to use finite differences.
467: !
468:     PetscCallA(FormJacobian(snes, x, J, J, 0, ierr))
469: !
470: !       Color the matrix, i.e. determine groups of columns that share no common
471: !      rows. These columns in the Jacobian can all be computed simultaneously.
472: !
473:     PetscCallA(MatColoringCreate(J, mc, ierr))
474:     PetscCallA(MatColoringSetType(mc, MATCOLORINGNATURAL, ierr))
475:     PetscCallA(MatColoringSetFromOptions(mc, ierr))
476:     PetscCallA(MatColoringApply(mc, iscoloring, ierr))
477:     PetscCallA(MatColoringDestroy(mc, ierr))
478: !
479: !       Create the data structure that SNESComputeJacobianDefaultColor() uses
480: !       to compute the actual Jacobians via finite differences.
481: !
482:     PetscCallA(MatFDColoringCreate(J, iscoloring, fdcoloring, ierr))
483:     PetscCallA(MatFDColoringSetFunction(fdcoloring, FormFunction, fdcoloring, ierr))
484:     PetscCallA(MatFDColoringSetFromOptions(fdcoloring, ierr))
485:     PetscCallA(MatFDColoringSetUp(J, iscoloring, fdcoloring, ierr))
486: !
487: !        Tell SNES to use the routine SNESComputeJacobianDefaultColor()
488: !      to compute Jacobians.
489: !
490:     PetscCallA(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring, ierr))
491:     PetscCallA(ISColoringDestroy(iscoloring, ierr))

493:   else if (.not. matrix_free) then

495: !  Set Jacobian matrix data structure and default Jacobian evaluation
496: !  routine.  Whenever the nonlinear solver needs to compute the
497: !  Jacobian matrix, it will call this routine.
498: !   - Note that the final routine argument is the user-defined
499: !     context that provides application-specific data for the
500: !     Jacobian evaluation routine.
501: !   - The user can override with:
502: !      -snes_fd : default finite differencing approximation of Jacobian
503: !      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
504: !                 (unless user explicitly sets preconditioner)
505: !      -snes_mf_operator : form matrix from which to construct the preconditioner as set by the user,
506: !                          but use matrix-free approx for Jacobian-vector
507: !                          products within Newton-Krylov method
508: !
509:     PetscCallA(SNESSetJacobian(snes, J, J, FormJacobian, 0, ierr))
510:   end if

512: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
513: !  Customize nonlinear solver; set runtime options
514: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

516: !  Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)

518:   PetscCallA(SNESSetFromOptions(snes, ierr))

520: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
521: !  Evaluate initial guess; then solve nonlinear system.
522: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

524: !  Note: The user should initialize the vector, x, with the initial guess
525: !  for the nonlinear solver prior to calling SNESSolve().  In particular,
526: !  to employ an initial guess of zero, the user should explicitly set
527: !  this vector to zero by calling VecSet().

529:   PetscCallA(FormInitialGuess(x, ierr))
530:   PetscCallA(SNESSolve(snes, PETSC_NULL_VEC, x, ierr))
531:   PetscCallA(SNESGetIterationNumber(snes, its, ierr))
532:   write (outputString, *) its
533:   PetscCallA(PetscPrintf(PETSC_COMM_WORLD, 'Number of SNES iterations = '//trim(outputString)//'\n', ierr))

535: !  PetscDraw contour plot of solution

537:   PetscCallA(PetscDrawCreate(PETSC_COMM_WORLD, PETSC_NULL_CHARACTER, 'Solution', xl, yl, width, height, draw, ierr))
538:   PetscCallA(PetscDrawSetFromOptions(draw, ierr))

540:   PetscCallA(VecGetArrayRead(x, lx_v, ierr))
541:   imx = int(mx, kind=kind(imx))
542:   imy = int(my, kind=kind(imy))
543:   PetscCallA(PetscDrawTensorContour(draw, imx, imy, PETSC_NULL_REAL_ARRAY, PETSC_NULL_REAL_ARRAY, lx_v, ierr))
544:   PetscCallA(VecRestoreArrayRead(x, lx_v, ierr))

546: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
547: !  Free work space.  All PETSc objects should be destroyed when they
548: !  are no longer needed.
549: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

551:   if (.not. matrix_free) PetscCallA(MatDestroy(J, ierr))
552:   if (fd_coloring) PetscCallA(MatFDColoringDestroy(fdcoloring, ierr))

554:   PetscCallA(VecDestroy(x, ierr))
555:   PetscCallA(VecDestroy(r, ierr))
556:   PetscCallA(SNESDestroy(snes, ierr))
557:   PetscCallA(PetscDrawDestroy(draw, ierr))
558:   PetscCallA(PetscFinalize(ierr))
559: end
560: !
561: !/*TEST
562: !
563: !   build:
564: !      requires: !single !complex
565: !
566: !   test:
567: !      args: -snes_monitor_short -nox -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
568: !
569: !   test:
570: !      suffix: 2
571: !      args: -snes_monitor_short -nox -snes_fd -ksp_gmres_cgs_refinement_type refine_always
572: !
573: !   test:
574: !      suffix: 3
575: !      args: -snes_monitor_short -nox -snes_fd_coloring -mat_coloring_type sl -ksp_gmres_cgs_refinement_type refine_always
576: !      filter: sort -b
577: !      filter_output: sort -b
578: !
579: !   test:
580: !     suffix: 4
581: !     args: -pc -par 6.807 -nox
582: !
583: !TEST*/