xref: /petsc/src/mat/tests/ex196f90.F90 (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown!
2*c4762a1bSJed Brown!
3*c4762a1bSJed Brown!   This program demonstrates use of MatSeqAIJGetArrayF90()
4*c4762a1bSJed Brown!
5*c4762a1bSJed Brown      program main
6*c4762a1bSJed Brown
7*c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
8*c4762a1bSJed Brown      use petscmat
9*c4762a1bSJed Brown      implicit none
10*c4762a1bSJed Brown
11*c4762a1bSJed Brown      Mat      A
12*c4762a1bSJed Brown      PetscErrorCode ierr
13*c4762a1bSJed Brown      PetscViewer   v
14*c4762a1bSJed Brown      PetscScalar, pointer :: aa(:)
15*c4762a1bSJed Brown      character*(256)  f
16*c4762a1bSJed Brown      PetscBool flg
17*c4762a1bSJed Brown
18*c4762a1bSJed Brown      call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
19*c4762a1bSJed Brown      if (ierr .ne. 0) then
20*c4762a1bSJed Brown        print*,'Unable to initialize PETSc'
21*c4762a1bSJed Brown        stop
22*c4762a1bSJed Brown      endif
23*c4762a1bSJed Brown
24*c4762a1bSJed Brown      call PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-f',f,flg,ierr);CHKERRA(ierr)
25*c4762a1bSJed Brown      call PetscViewerBinaryOpen(PETSC_COMM_WORLD,f,FILE_MODE_READ,v,ierr);CHKERRA(ierr)
26*c4762a1bSJed Brown
27*c4762a1bSJed Brown      call MatCreate(PETSC_COMM_WORLD,A,ierr);CHKERRA(ierr)
28*c4762a1bSJed Brown      call MatSetType(A, MATSEQAIJ,ierr);CHKERRA(ierr)
29*c4762a1bSJed Brown      call MatLoad(A,v,ierr);CHKERRA(ierr)
30*c4762a1bSJed Brown
31*c4762a1bSJed Brown      call MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr)
32*c4762a1bSJed Brown
33*c4762a1bSJed Brown      call MatSeqAIJGetArrayF90(A,aa,ierr);CHKERRA(ierr)
34*c4762a1bSJed Brown      print*,aa(3)
35*c4762a1bSJed Brown
36*c4762a1bSJed Brown      call MatDestroy(A,ierr);CHKERRA(ierr)
37*c4762a1bSJed Brown      call PetscViewerDestroy(v,ierr);CHKERRA(ierr)
38*c4762a1bSJed Brown
39*c4762a1bSJed Brown      call PetscFinalize(ierr)
40*c4762a1bSJed Brown      end
41*c4762a1bSJed Brown
42*c4762a1bSJed Brown!/*TEST
43*c4762a1bSJed Brown!
44*c4762a1bSJed Brown!   test:
45*c4762a1bSJed Brown!      args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -malloc_dump
46*c4762a1bSJed Brown!      requires: !complex double !define(PETSC_USE_64BIT_INDICES)
47*c4762a1bSJed Brown!
48*c4762a1bSJed Brown!TEST*/
49