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