xref: /petsc/src/mat/tests/ex219f.F90 (revision d8606c274c09e255c003062beb17b1be973467bc)
1c4762a1bSJed Brownprogram newnonzero
2c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
3c4762a1bSJed Brown use petscmat
4c4762a1bSJed Brown implicit none
5c4762a1bSJed Brown
6c4762a1bSJed Brown Mat :: A
706946f3aSJose E. Roman PetscInt :: n,m,idxm(1),idxn(1),nl1,nl2,zero,one,i
8c4762a1bSJed Brown PetscScalar :: v(1)
9c4762a1bSJed Brown PetscErrorCode :: ierr
10c4762a1bSJed Brown
11*d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr))
1206946f3aSJose E. Roman zero = 0
13c4762a1bSJed Brown one = 1
14c4762a1bSJed Brown n=3
15c4762a1bSJed Brown m=n
16*d8606c27SBarry Smith PetscCallA(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,m,one,PETSC_NULL_INTEGER,zero,PETSC_NULL_INTEGER,A,ierr))
17c4762a1bSJed Brown
18*d8606c27SBarry Smith PetscCallA(MatGetOwnershipRange(A,nl1,nl2,ierr))
19c4762a1bSJed Brown do i=nl1,nl2-1
20c4762a1bSJed Brown   idxn(1)=i
21c4762a1bSJed Brown   idxm(1)=i
22c4762a1bSJed Brown   v(1)=1.0
23*d8606c27SBarry Smith   PetscCallA(MatSetValues(A,one,idxn,one,idxm, v,INSERT_VALUES,ierr))
24c4762a1bSJed Brown end do
25*d8606c27SBarry Smith PetscCallA(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr))
26*d8606c27SBarry Smith PetscCallA(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr))
27c4762a1bSJed Brown
28c4762a1bSJed Brown! Ignore any values set into new nonzero locations
29*d8606c27SBarry Smith PetscCallA(MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE,ierr))
30c4762a1bSJed Brown
31c4762a1bSJed Brown idxn(1)=0
32c4762a1bSJed Brown idxm(1)=n-1
33c4762a1bSJed Brown if ((idxn(1).ge.nl1).and.(idxn(1).le.nl2-1)) then
34c4762a1bSJed Brown   v(1)=2.0
35*d8606c27SBarry Smith   PetscCallA(MatSetValues(A,one,idxn,one,idxm, v,INSERT_VALUES,ierr))
36c4762a1bSJed Brown end if
37*d8606c27SBarry Smith PetscCallA(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr))
38*d8606c27SBarry Smith PetscCallA(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr))
39c4762a1bSJed Brown
40c4762a1bSJed Brown if ((idxn(1).ge.nl1).and.(idxn(1).le.nl2-1)) then
41*d8606c27SBarry Smith   PetscCallA(MatGetValues(A,one,idxn,one,idxm, v,ierr))
42c4762a1bSJed Brown   write(6,*) PetscRealPart(v)
43c4762a1bSJed Brown end if
44c4762a1bSJed Brown
45*d8606c27SBarry Smith PetscCallA(MatDestroy(A,ierr))
46*d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr))
47c4762a1bSJed Brown
48c4762a1bSJed Brown end program newnonzero
49c4762a1bSJed Brown
50c4762a1bSJed Brown!/*TEST
51c4762a1bSJed Brown!
52c4762a1bSJed Brown!     test:
53c4762a1bSJed Brown!       nsize: 2
54c4762a1bSJed Brown!       filter: Error:
55c4762a1bSJed Brown!
56c4762a1bSJed Brown!     test:
57dfd57a17SPierre Jolivet!       requires: defined(PETSC_USE_INFO)
58c4762a1bSJed Brown!       suffix: 2
59c4762a1bSJed Brown!       nsize: 2
60c4762a1bSJed Brown!       args: -info
61c4762a1bSJed Brown!       filter: grep "Skipping"
62c4762a1bSJed Brown!
63c4762a1bSJed Brown!TEST*/
64