xref: /petsc/src/mat/tests/ex63f.F90 (revision d8606c274c09e255c003062beb17b1be973467bc)
1*d8606c27SBarry Smith!
2*d8606c27SBarry Smith!
3*d8606c27SBarry Smith!   This program tests storage of PETSc Dense matrix.
4*d8606c27SBarry Smith!   It Creates a Dense matrix, and stores it into a file,
5*d8606c27SBarry Smith!   and then reads it back in as a SeqDense and MPIDense
6*d8606c27SBarry Smith!   matrix, and prints out the contents.
7*d8606c27SBarry Smith!
8*d8606c27SBarry Smith      program main
9*d8606c27SBarry Smith#include <petsc/finclude/petscmat.h>
10*d8606c27SBarry Smith      use petscmat
11*d8606c27SBarry Smith      implicit none
12*d8606c27SBarry Smith
13*d8606c27SBarry Smith      PetscErrorCode ierr
14*d8606c27SBarry Smith      PetscInt row,col,ten
15*d8606c27SBarry Smith      PetscMPIInt rank
16*d8606c27SBarry Smith      PetscScalar  v
17*d8606c27SBarry Smith      Mat     A
18*d8606c27SBarry Smith      PetscViewer  view
19*d8606c27SBarry Smith
20*d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
21*d8606c27SBarry Smith
22*d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
23*d8606c27SBarry Smith!
24*d8606c27SBarry Smith!     Proc-0 Create a seq-dense matrix and write it to a file
25*d8606c27SBarry Smith!
26*d8606c27SBarry Smith      if (rank .eq. 0) then
27*d8606c27SBarry Smith         ten = 10
28*d8606c27SBarry Smith         PetscCallA(MatCreateSeqDense(PETSC_COMM_SELF,ten,ten,PETSC_NULL_SCALAR,A,ierr))
29*d8606c27SBarry Smith         v = 1.0
30*d8606c27SBarry Smith         do row=0,9
31*d8606c27SBarry Smith            do col=0,9
32*d8606c27SBarry Smith               PetscCallA(MatSetValue(A,row,col,v,INSERT_VALUES,ierr))
33*d8606c27SBarry Smith               v = v + 1.0
34*d8606c27SBarry Smith            enddo
35*d8606c27SBarry Smith         enddo
36*d8606c27SBarry Smith
37*d8606c27SBarry Smith         PetscCallA(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr))
38*d8606c27SBarry Smith         PetscCallA(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr))
39*d8606c27SBarry Smith
40*d8606c27SBarry Smith         PetscCallA(PetscObjectSetName(A,'Original Matrix',ierr))
41*d8606c27SBarry Smith         PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_SELF,ierr))
42*d8606c27SBarry Smith!
43*d8606c27SBarry Smith!        Now Write this matrix to a binary file
44*d8606c27SBarry Smith!
45*d8606c27SBarry Smith         PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF,'dense.mat',FILE_MODE_WRITE,view,ierr))
46*d8606c27SBarry Smith         PetscCallA(MatView(A,view,ierr))
47*d8606c27SBarry Smith         PetscCallA(PetscViewerDestroy(view,ierr))
48*d8606c27SBarry Smith         PetscCallA(MatDestroy(A,ierr))
49*d8606c27SBarry Smith!
50*d8606c27SBarry Smith!        Read this matrix into a SeqDense matrix
51*d8606c27SBarry Smith
52*d8606c27SBarry Smith         PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_SELF,'dense.mat',FILE_MODE_READ,view,ierr))
53*d8606c27SBarry Smith         PetscCallA(MatCreate(PETSC_COMM_SELF,A,ierr))
54*d8606c27SBarry Smith         PetscCallA(MatSetType(A, MATSEQDENSE,ierr))
55*d8606c27SBarry Smith         PetscCallA(MatLoad(A,view,ierr))
56*d8606c27SBarry Smith
57*d8606c27SBarry Smith         PetscCallA(PetscObjectSetName(A,'SeqDense Matrix read in from file',ierr))
58*d8606c27SBarry Smith         PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_SELF,ierr))
59*d8606c27SBarry Smith         PetscCallA(MatDestroy(A,ierr))
60*d8606c27SBarry Smith         PetscCallA(PetscViewerDestroy(view,ierr))
61*d8606c27SBarry Smith      endif
62*d8606c27SBarry Smith
63*d8606c27SBarry Smith!
64*d8606c27SBarry Smith!     Use a barrier, so that the procs do not try opening the file before
65*d8606c27SBarry Smith!     it is created.
66*d8606c27SBarry Smith!
67*d8606c27SBarry Smith      PetscCallMPIA(MPI_Barrier(PETSC_COMM_WORLD,ierr))
68*d8606c27SBarry Smith
69*d8606c27SBarry Smith      PetscCallA(PetscViewerBinaryOpen(PETSC_COMM_WORLD,'dense.mat',FILE_MODE_READ,view,ierr))
70*d8606c27SBarry Smith      PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
71*d8606c27SBarry Smith      PetscCallA(MatSetType(A, MATMPIDENSE,ierr))
72*d8606c27SBarry Smith      PetscCallA(MatLoad(A,view,ierr))
73*d8606c27SBarry Smith
74*d8606c27SBarry Smith      PetscCallA(PetscObjectSetName(A, 'MPIDense Matrix read in from file',ierr))
75*d8606c27SBarry Smith      PetscCallA(MatView(A,PETSC_VIEWER_STDOUT_WORLD,ierr))
76*d8606c27SBarry Smith      PetscCallA(MatDestroy(A,ierr))
77*d8606c27SBarry Smith      PetscCallA(PetscViewerDestroy(view,ierr))
78*d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
79*d8606c27SBarry Smith      end
80*d8606c27SBarry Smith
81*d8606c27SBarry Smith!/*TEST
82*d8606c27SBarry Smith!
83*d8606c27SBarry Smith!   test:
84*d8606c27SBarry Smith!      nsize: 2
85*d8606c27SBarry Smith!      output_file: output/ex63_1.out
86*d8606c27SBarry Smith!
87*d8606c27SBarry Smith!TEST*/
88