xref: /petsc/src/mat/tests/ex36f.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675) !
1d8606c27SBarry Smith!
2d8606c27SBarry Smith!
3d8606c27SBarry Smith!   This program demonstrates use of PETSc dense matrices.
4d8606c27SBarry Smith!
5c5e229c2SMartin Diehl#include <petsc/finclude/petscmat.h>
6*01fa2b5aSMartin Diehlmodule ex36fmodule
7*01fa2b5aSMartin Diehl  use petscmat
8d8606c27SBarry Smith  implicit none
9d8606c27SBarry Smith! -----------------------------------------------------------------
10d8606c27SBarry Smith!
11d8606c27SBarry Smith!  Demo1 -  This subroutine demonstrates the use of PETSc-allocated dense
12ce78bad3SBarry Smith!  matrix storage.  Here MatDenseGetArray() is used for direct access to the
1342ce371bSBarry Smith!  array that stores the dense matrix.
14d8606c27SBarry Smith!
152c5cda21SPierre Jolivet!  Note the use of PETSC_NULL_SCALAR_ARRAY in MatCreateSeqDense() to indicate that no
16d8606c27SBarry Smith!  storage is being provided by the user. (Do NOT pass a zero in that
17d8606c27SBarry Smith!  location.)
18d8606c27SBarry Smith!
19*01fa2b5aSMartin Diehlcontains
20d8606c27SBarry Smith  subroutine Demo1()
21d8606c27SBarry Smith
22d8606c27SBarry Smith    Mat A
23d8606c27SBarry Smith    PetscInt n, m
24d8606c27SBarry Smith    PetscErrorCode ierr
2542ce371bSBarry Smith    PetscScalar, pointer :: aa(:, :)
26d8606c27SBarry Smith
27d8606c27SBarry Smith    n = 4
28d8606c27SBarry Smith    m = 5
29d8606c27SBarry Smith
30d8606c27SBarry Smith!  Create matrix
31d8606c27SBarry Smith
32d8606c27SBarry Smith    PetscCall(MatCreate(PETSC_COMM_SELF, A, ierr))
33d8606c27SBarry Smith    PetscCall(MatSetSizes(A, m, n, m, n, ierr))
34d8606c27SBarry Smith    PetscCall(MatSetType(A, MATSEQDENSE, ierr))
35d8606c27SBarry Smith    PetscCall(MatSetUp(A, ierr))
36d8606c27SBarry Smith
37d8606c27SBarry Smith!  Access array storage
38ce78bad3SBarry Smith    PetscCall(MatDenseGetArray(A, aa, ierr))
39d8606c27SBarry Smith
40d8606c27SBarry Smith!  Set matrix values directly
4142ce371bSBarry Smith    PetscCall(FillUpMatrix(m, n, aa))
42d8606c27SBarry Smith
43ce78bad3SBarry Smith    PetscCall(MatDenseRestoreArray(A, aa, ierr))
44d8606c27SBarry Smith
45d8606c27SBarry Smith!  Finalize matrix assembly
46d8606c27SBarry Smith    PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
47d8606c27SBarry Smith    PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
48d8606c27SBarry Smith
49d8606c27SBarry Smith!  View matrix
50d8606c27SBarry Smith    PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF, ierr))
51d8606c27SBarry Smith
52d8606c27SBarry Smith!  Clean up
53d8606c27SBarry Smith    PetscCall(MatDestroy(A, ierr))
54d8606c27SBarry Smith  end
55d8606c27SBarry Smith
56d8606c27SBarry Smith! -----------------------------------------------------------------
57d8606c27SBarry Smith!
58d8606c27SBarry Smith!  Demo2 -  This subroutine demonstrates the use of user-allocated dense
59d8606c27SBarry Smith!  matrix storage.
60d8606c27SBarry Smith!
61d8606c27SBarry Smith  subroutine Demo2()
62d8606c27SBarry Smith
63d8606c27SBarry Smith    PetscInt n, m
64d8606c27SBarry Smith    PetscErrorCode ierr
65d8606c27SBarry Smith    parameter(m=5, n=4)
66d8606c27SBarry Smith    Mat A
67d8606c27SBarry Smith    PetscScalar aa(m, n)
68d8606c27SBarry Smith
69d8606c27SBarry Smith!  Create matrix
70d8606c27SBarry Smith    PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, n, aa, A, ierr))
71d8606c27SBarry Smith
72d8606c27SBarry Smith!  Set matrix values directly
7342ce371bSBarry Smith    PetscCall(FillUpMatrix(m, n, aa))
74d8606c27SBarry Smith
75d8606c27SBarry Smith!  Finalize matrix assembly
76d8606c27SBarry Smith    PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
77d8606c27SBarry Smith    PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
78d8606c27SBarry Smith
79d8606c27SBarry Smith!  View matrix
80d8606c27SBarry Smith    PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF, ierr))
81d8606c27SBarry Smith
82d8606c27SBarry Smith!  Clean up
83d8606c27SBarry Smith    PetscCall(MatDestroy(A, ierr))
84d8606c27SBarry Smith  end
85d8606c27SBarry Smith
86d8606c27SBarry Smith! -----------------------------------------------------------------
87d8606c27SBarry Smith
88d8606c27SBarry Smith  subroutine FillUpMatrix(m, n, X)
89d8606c27SBarry Smith    PetscInt m, n, i, j
90d8606c27SBarry Smith    PetscScalar X(m, n)
91d8606c27SBarry Smith
92ceeae6b5SMartin Diehl    do j = 1, n
93ceeae6b5SMartin Diehl      do i = 1, m
94d8606c27SBarry Smith        X(i, j) = 1.0/real(i + j - 1)
95ceeae6b5SMartin Diehl      end do
96ceeae6b5SMartin Diehl    end do
97*01fa2b5aSMartin Diehl    end module ex36fmodule
98*01fa2b5aSMartin Diehl
99*01fa2b5aSMartin Diehl    end program main
100*01fa2b5aSMartin Diehl    use ex36fmodule
101*01fa2b5aSMartin Diehl    implicit none
102*01fa2b5aSMartin Diehl
103*01fa2b5aSMartin Diehl    PetscErrorCode ierr
104*01fa2b5aSMartin Diehl
105*01fa2b5aSMartin Diehl    PetscCallA(PetscInitialize(ierr))
106*01fa2b5aSMartin Diehl
107*01fa2b5aSMartin Diehl!  Demo of PETSc-allocated dense matrix storage
108*01fa2b5aSMartin Diehl    call Demo1()
109*01fa2b5aSMartin Diehl
110*01fa2b5aSMartin Diehl!  Demo of user-allocated dense matrix storage
111*01fa2b5aSMartin Diehl    call Demo2()
112*01fa2b5aSMartin Diehl
113*01fa2b5aSMartin Diehl    PetscCallA(PetscFinalize(ierr))
114d8606c27SBarry Smith  end
115