Actual source code: ex219f.F90
1: #include <petsc/finclude/petscmat.h>
2: program newnonzero
3: use petscmat
4: implicit none
6: Mat :: A
7: PetscInt, parameter :: n = 3, m = n
8: PetscInt :: idxm(1), idxn(1), nl1, nl2, i
9: PetscScalar :: v(1), value(1), values(2)
10: PetscErrorCode :: ierr
11: IS :: is
12: ISLocalToGlobalMapping :: ismap
14: PetscCallA(PetscInitialize(ierr))
15: PetscCallA(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, m, 1_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, 0_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, A, ierr))
17: PetscCallA(MatGetOwnershipRange(A, nl1, nl2, ierr))
18: do i = nl1, nl2 - 1
19: idxn(1) = i
20: idxm(1) = i
21: v(1) = 1.0
22: PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, idxn, 1_PETSC_INT_KIND, idxm, v, INSERT_VALUES, ierr))
23: end do
24: PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
25: PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
27: ! Ignore any values set into new nonzero locations
28: PetscCallA(MatSetOption(A, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE, ierr))
30: idxn(1) = 0
31: idxm(1) = n - 1
32: if ((idxn(1) >= nl1) .and. (idxn(1) <= nl2 - 1)) then
33: v(1) = 2.0
34: PetscCallA(MatSetValues(A, 1_PETSC_INT_KIND, idxn, 1_PETSC_INT_KIND, idxm, v, INSERT_VALUES, ierr))
35: end if
36: PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
37: PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
39: if ((idxn(1) >= nl1) .and. (idxn(1) <= nl2 - 1)) then
40: PetscCallA(MatGetValues(A, 1_PETSC_INT_KIND, idxn, 1_PETSC_INT_KIND, idxm, v, ierr))
41: write (6, *) PetscRealPart(v)
42: end if
44: PetscCallA(ISCreateStride(PETSC_COMM_WORLD, nl2 - nl1, nl1, 1_PETSC_INT_KIND, is, ierr))
45: PetscCallA(ISLocalToGlobalMappingCreateIS(is, ismap, ierr))
46: PetscCallA(MatSetLocalToGlobalMapping(A, ismap, ismap, ierr))
47: PetscCallA(ISLocalToGlobalMappingDestroy(ismap, ierr))
48: PetscCallA(ISDestroy(is, ierr))
49: PetscCallA(MatGetValuesLocal(A, 1_PETSC_INT_KIND, [0_PETSC_INT_KIND], 1_PETSC_INT_KIND, [0_PETSC_INT_KIND], value, ierr))
50: PetscCallA(MatGetValuesLocal(A, 1_PETSC_INT_KIND, [0_PETSC_INT_KIND], 1_PETSC_INT_KIND, [0_PETSC_INT_KIND], values, ierr))
51: idxn(1) = 0
52: PetscCallA(MatGetValuesLocal(A, 1_PETSC_INT_KIND, idxn, 1_PETSC_INT_KIND, [0_PETSC_INT_KIND], values, ierr))
53: PetscCallA(MatGetValuesLocal(A, 1_PETSC_INT_KIND, idxn, 1_PETSC_INT_KIND, idxn, values, ierr))
55: PetscCallA(MatDestroy(A, ierr))
56: PetscCallA(PetscFinalize(ierr))
58: end program newnonzero
60: !/*TEST
61: !
62: ! test:
63: ! nsize: 2
64: ! filter: Error:
65: !
66: ! test:
67: ! requires: defined(PETSC_USE_INFO)
68: ! suffix: 2
69: ! nsize: 2
70: ! args: -info
71: ! filter: grep "Skipping"
72: !
73: !TEST*/