xref: /petsc/src/mat/ftn-mod/petscmat.h90 (revision 0337bfe0b9dcc77abc5d44df0b7f57cdcdf2ff74)
1
2
3      Interface
4        Subroutine MatFDColoringRestorePerturbedColumns(i,len,array,ierr)
5      use, intrinsic :: ISO_C_binding
6          import tMatFDColoring
7          PetscInt, pointer :: array(:)
8          PetscInt len
9          PetscErrorCode ierr
10         MatFDColoring      i
11        End Subroutine
12      End Interface
13
14        interface MatDenseGetArray
15        Subroutine MatDenseGetArray1d(v,array,ierr)
16      use, intrinsic :: ISO_C_binding
17        import tMat
18          PetscScalar, pointer :: array(:)
19          PetscErrorCode ierr
20          Mat     v
21        End Subroutine
22        Subroutine MatDenseGetArray2d(v,array,ierr)
23      use, intrinsic :: ISO_C_binding
24         import tMat
25          PetscScalar, pointer :: array(:,:)
26          PetscErrorCode ierr
27          Mat     v
28        End Subroutine
29        end interface
30
31        interface MatDenseRestoreArray
32        Subroutine MatDenseRestoreArray1d(v,array,ierr)
33      use, intrinsic :: ISO_C_binding
34         import tMat
35         PetscScalar, pointer :: array(:)
36          PetscErrorCode ierr
37          Mat     v
38        End Subroutine
39        Subroutine MatDenseRestoreArray2d(v,array,ierr)
40      use, intrinsic :: ISO_C_binding
41         import tMat
42         PetscScalar, pointer :: array(:,:)
43          PetscErrorCode ierr
44          Mat     v
45        End Subroutine
46        end interface
47
48        interface MatDenseGetArrayRead
49        Subroutine MatDenseGetArrayRead1d(v,array,ierr)
50      use, intrinsic :: ISO_C_binding
51         import tMat
52         PetscScalar, pointer :: array(:)
53          PetscErrorCode ierr
54          Mat     v
55        End Subroutine
56        Subroutine MatDenseGetArrayRead2d(v,array,ierr)
57      use, intrinsic :: ISO_C_binding
58         import tMat
59         PetscScalar, pointer :: array(:,:)
60          PetscErrorCode ierr
61          Mat     v
62        End Subroutine
63        end interface
64
65        interface MatDenseRestoreArrayRead
66        Subroutine MatDenseRestoreArrayRead1d(v,array,ierr)
67      use, intrinsic :: ISO_C_binding
68         import tMat
69         PetscScalar, pointer :: array(:)
70          PetscErrorCode ierr
71          Mat     v
72        End Subroutine
73        Subroutine MatDenseRestoreArrayRead2d(v,array,ierr)
74      use, intrinsic :: ISO_C_binding
75         import tMat
76         PetscScalar, pointer :: array(:,:)
77          PetscErrorCode ierr
78          Mat     v
79        End Subroutine
80        end interface
81
82        interface MatDenseGetArrayWrite
83        Subroutine MatDenseGetArrayWrite1d(v,array,ierr)
84      use, intrinsic :: ISO_C_binding
85         import tMat
86         PetscScalar, pointer :: array(:)
87          PetscErrorCode ierr
88          Mat     v
89        End Subroutine
90        Subroutine MatDenseGetArrayWrite2d(v,array,ierr)
91      use, intrinsic :: ISO_C_binding
92         import tMat
93         PetscScalar, pointer :: array(:,:)
94          PetscErrorCode ierr
95          Mat     v
96        End Subroutine
97        end interface
98
99        interface MatDenseRestoreArrayWrite
100        Subroutine MatDenseRestoreArrayWrite1d(v,array,ierr)
101      use, intrinsic :: ISO_C_binding
102         import tMat
103         PetscScalar, pointer :: array(:)
104          PetscErrorCode ierr
105          Mat     v
106        End Subroutine
107        Subroutine MatDenseRestoreArrayWrite2d(v,array,ierr)
108      use, intrinsic :: ISO_C_binding
109         import tMat
110         PetscScalar, pointer :: array(:,:)
111          PetscErrorCode ierr
112          Mat     v
113        End Subroutine
114        end interface
115