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