Actual source code: ex31f.F90

  1: !
  2: !   Tests TSRKGetTableau()
  3: !
  4: #include <petsc/finclude/petscts.h>

  6: program main
  7:   use petscts
  8:   implicit none
  9: !
 10:   TS ts
 11:   PetscInt s, p, sz
 12:   PetscBool FSAL
 13:   PetscReal, pointer :: A(:), b(:), c(:), bembed(:), binterp(:)
 14:   PetscErrorCode ierr
 15:   Vec x

 17:   PetscCallA(PetscInitialize(ierr))
 18:   PetscCallA(TSCreate(PETSC_COMM_WORLD, ts, ierr))
 19:   PetscCallA(TSSetType(ts, TSRK, ierr))
 20:   PetscCallA(TSSetFromOptions(ts, ierr))
 21:   PetscCallA(VecCreateSeq(PETSC_COMM_SELF, 1_PETSC_INT_KIND, x, ierr))
 22:   PetscCallA(TSSetSolution(ts, x, ierr))

 24:   PetscCallA(TSSetUp(ts, ierr))
 25:   PetscCallA(TSRKGetTableau(ts, s, A, b, c, bembed, p, binterp, FSAL, ierr))
 26:   sz = size(A)
 27:   PetscCallA(PetscRealView(sz, A, PETSC_VIEWER_STDOUT_SELF, ierr))
 28:   sz = size(b)
 29:   PetscCallA(PetscRealView(sz, b, PETSC_VIEWER_STDOUT_SELF, ierr))
 30:   sz = size(c)
 31:   PetscCallA(PetscRealView(sz, c, PETSC_VIEWER_STDOUT_SELF, ierr))
 32:   sz = size(bembed)
 33:   PetscCallA(PetscRealView(sz, bembed, PETSC_VIEWER_STDOUT_SELF, ierr))
 34:   sz = size(binterp)
 35:   PetscCallA(PetscRealView(sz, binterp, PETSC_VIEWER_STDOUT_SELF, ierr))
 36:   PetscCallA(TSRKRestoreTableau(ts, s, A, b, c, bembed, p, binterp, FSAL, ierr))

 38:   PetscCallA(TSRKSetType(ts, "5bs", ierr))
 39:   PetscCallA(TSSetUp(ts, ierr))
 40:   PetscCallA(TSRKGetTableau(ts, s, A, PETSC_NULL_REAL_POINTER, c, bembed, p, binterp, PETSC_NULL_BOOL, ierr))
 41:   sz = size(A)
 42:   PetscCallA(PetscRealView(sz, A, PETSC_VIEWER_STDOUT_SELF, ierr))
 43:   sz = size(c)
 44:   PetscCallA(PetscRealView(sz, c, PETSC_VIEWER_STDOUT_SELF, ierr))
 45:   sz = size(bembed)
 46:   PetscCallA(PetscRealView(sz, bembed, PETSC_VIEWER_STDOUT_SELF, ierr))
 47:   sz = size(binterp)
 48:   PetscCallA(PetscRealView(sz, binterp, PETSC_VIEWER_STDOUT_SELF, ierr))
 49:   PetscCallA(TSRKRestoreTableau(ts, s, A, PETSC_NULL_REAL_POINTER, c, bembed, p, binterp, PETSC_NULL_BOOL, ierr))

 51:   PetscCallA(VecDestroy(x, ierr))
 52:   PetscCallA(TSDestroy(ts, ierr))
 53:   PetscCallA(PetscFinalize(ierr))
 54: end
 55: !/*TEST
 56: !
 57: !    test:
 58: !
 59: !TEST*/