Actual source code: ex58f.F90

  1: !
  2: !
  3: !   This program demonstrates use of MatGetRow() and MatGetRowMaxAbs() from Fortran
  4: !
  5:       program main
  6: #include <petsc/finclude/petscmat.h>
  7:       use petscmat
  8:       implicit none

 10:       Mat      A
 11:       PetscErrorCode ierr
 12:       PetscInt M,N
 13:       PetscViewer   v
 14:       Vec           rowmax
 15:       PetscBool flg
 16:       character*(256)  f

 18:       PetscCallA(PetscInitialize(ierr))

 20:       PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-f',f,flg,ierr))
 21:       PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,f,FILE_MODE_READ,v,ierr))

 23:       PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
 24:       PetscCallA(MatSetType(A, MATSEQAIJ,ierr))
 25:       PetscCallA(MatLoad(A,v,ierr))

 27:       PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr))

 29: !
 30: !     Test MatGetRowMaxAbs()
 31:       PetscCallA(MatGetSize(A,M,N,ierr))
 32:       PetscCallA(VecCreate(PETSC_COMM_WORLD,rowmax,ierr))
 33:       PetscCallA(VecSetSizes(rowmax,M,M,ierr))
 34:       PetscCallA(VecSetFromOptions(rowmax,ierr))

 36:       PetscCallA(MatGetRowMaxAbs(A,rowmax,PETSC_NULL_INTEGER,ierr))
 37:       PetscCallA(VecView(rowmax,PETSC_VIEWER_STDOUT_WORLD,ierr))

 39:       PetscCallA(MatGetRowMax(A,rowmax,PETSC_NULL_INTEGER,ierr))
 40:       PetscCallA(VecView(rowmax,PETSC_VIEWER_STDOUT_WORLD,ierr))

 42:       PetscCallA(MatGetRowMinAbs(A,rowmax,PETSC_NULL_INTEGER,ierr))
 43:       PetscCallA(VecView(rowmax,PETSC_VIEWER_STDOUT_WORLD,ierr))

 45:       PetscCallA(MatGetRowMin(A,rowmax,PETSC_NULL_INTEGER,ierr))
 46:       PetscCallA(VecView(rowmax,PETSC_VIEWER_STDOUT_WORLD,ierr))

 48:       PetscCallA(MatDestroy(A,ierr))
 49:       PetscCallA(PetscViewerDestroy(v,ierr))
 50:       PetscCallA(VecDestroy(rowmax,ierr))

 52:       PetscCallA(PetscFinalize(ierr))
 53:       end

 55: !/*TEST
 56: !
 57: !     test:
 58: !       args: -f ${DATAFILESPATH}/matrices/tiny
 59: !       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
 60: !
 61: !TEST*/