Actual source code: ex20f.F90
1: !
2: ! Demonstrates use of MatDuplicate() for a shell matrix with a context
3: !
4: #include "petsc/finclude/petscmat.h"
5: module ex20fmodule
6: use petscmat
7: implicit none
8: type :: MatCtx
9: PetscReal :: lambda
10: end type MatCtx
12: contains
13: subroutine MatDuplicate_F(F, opt, M, ierr)
15: Mat :: F, M
16: MatDuplicateOption :: opt
17: PetscErrorCode :: ierr
18: PetscInt :: ml, nl
19: type(MatCtx), pointer :: ctxM, ctxF_pt
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
29: subroutine MatDestroy_F(F, ierr)
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
38: end module ex20fmodule
40: ! ----------------------------------------------------
41: ! main program
42: ! ----------------------------------------------------
43: program 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
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
58: PetscCallA(MatShellGetContext(F, ctxF_pt, ierr))
59: print *, 'ctxF_pt%lambda = ', ctxF_pt%lambda
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
65: PetscCallA(MatDestroy(F, ierr))
66: PetscCallA(MatDestroy(Fcopy, ierr))
67: PetscCallA(PetscFinalize(ierr))
68: end program main
70: !/*TEST
71: !
72: ! build:
73: ! requires: double
74: !
75: ! test:
76: !
77: !TEST*/