xref: /petsc/src/mat/tutorials/ex20f.F90 (revision d52a580b706c59ca78066c1e38754e45b6b56e2b)
1!
2!     Demonstrates use of MatDuplicate() for a shell matrix with a context
3!
4#include "petsc/finclude/petscmat.h"
5module ex20fmodule
6  use petscmat
7  implicit none
8  type :: MatCtx
9    PetscReal :: lambda
10  end type MatCtx
11
12contains
13  subroutine MatDuplicate_F(F, opt, M, ierr)
14
15    Mat                  :: F, M
16    MatDuplicateOption   :: opt
17    PetscErrorCode       :: ierr
18    PetscInt             :: ml, nl
19    type(MatCtx), pointer :: ctxM, ctxF_pt
20
21    PetscCall(MatGetLocalSize(F, ml, nl, ierr))
22    PetscCall(MatShellGetContext(F, ctxF_pt, ierr))
23    allocate (ctxM)
24    ctxM%lambda = ctxF_pt%lambda
25    PetscCall(MatCreateShell(PETSC_COMM_WORLD, ml, nl, PETSC_DETERMINE, PETSC_DETERMINE, ctxM, M, ierr))
26    PetscCall(MatShellSetOperation(M, MATOP_DESTROY, MatDestroy_F, ierr))
27  end subroutine MatDuplicate_F
28
29  subroutine MatDestroy_F(F, ierr)
30
31    Mat                  :: F
32    PetscErrorCode       :: ierr
33    type(MatCtx), pointer :: ctxF_pt
34    PetscCall(MatShellGetContext(F, ctxF_pt, ierr))
35    deallocate (ctxF_pt)
36  end subroutine MatDestroy_F
37
38end module ex20fmodule
39
40! ----------------------------------------------------
41!                    main program
42! ----------------------------------------------------
43program main
44  use ex20fmodule
45  implicit none
46  Mat                  :: F, Fcopy
47  type(MatCtx)         :: ctxF
48  type(MatCtx), pointer :: ctxF_pt, ctxFcopy_pt
49  PetscErrorCode       :: ierr
50  PetscInt             :: n = 128
51
52  PetscCallA(PetscInitialize(ierr))
53  ctxF%lambda = 3.14d0
54  PetscCallA(MatCreateShell(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, n, ctxF, F, ierr))
55  PetscCallA(MatShellSetOperation(F, MATOP_DUPLICATE, MatDuplicate_F, ierr))
56  print *, 'ctxF%lambda = ', ctxF%lambda
57
58  PetscCallA(MatShellGetContext(F, ctxF_pt, ierr))
59  print *, 'ctxF_pt%lambda = ', ctxF_pt%lambda
60
61  PetscCallA(MatDuplicate(F, MAT_DO_NOT_COPY_VALUES, Fcopy, ierr))
62  PetscCallA(MatShellGetContext(Fcopy, ctxFcopy_pt, ierr))
63  print *, 'ctxFcopy_pt%lambda = ', ctxFcopy_pt%lambda
64
65  PetscCallA(MatDestroy(F, ierr))
66  PetscCallA(MatDestroy(Fcopy, ierr))
67  PetscCallA(PetscFinalize(ierr))
68end program main
69
70!/*TEST
71!
72!     build:
73!       requires: double
74!
75!     test:
76!
77!TEST*/
78