Actual source code: ex15f.F90

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: program main

  3:  #include <petsc/finclude/petscvec.h>
  4:  #include <petsc/finclude/petscmat.h>

  6:   use petscvec
  7:   use petscmat

  9:   implicit none

 11:   Mat             :: A
 12:   MatPartitioning :: part
 13:   IS              :: is
 14:   PetscInt        :: r,myStart,myEnd
 15:   PetscInt,parameter :: N = 10
 16:   PetscErrorCode  :: ierr
 17:   PetscScalar,pointer,dimension(:) :: vals
 18:   PetscInt,pointer,dimension(:) :: cols
 19:   PetscBool :: flg

 21:   call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 22:   if (ierr /= 0) then
 23:    print*,'PetscInitialize failed'
 24:    stop
 25:   endif

 27:   call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-N",N,flg,ierr);CHKERRA(ierr)
 28:   call MatCreate(PETSC_COMM_WORLD, A,ierr);CHKERRA(ierr)
 29:   call MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N,ierr);CHKERRA(ierr)
 30:   call MatSetFromOptions(A,ierr);CHKERRA(ierr)
 31:   call MatSeqAIJSetPreallocation(A, 3, PETSC_NULL_INTEGER,ierr);CHKERRA(ierr)
 32:   call MatMPIAIJSetPreallocation(A, 3, PETSC_NULL_INTEGER, 2, PETSC_NULL_INTEGER,ierr);CHKERRA(ierr)

 34:   !/* Create a linear mesh */
 35:   call MatGetOwnershipRange(A, myStart, myEnd,ierr);CHKERRA(ierr)
 36: 
 37:   do r=myStart,myEnd-1
 38:     if (r == 0) then
 39:      allocate(vals(0:1))
 40:      vals = 1.0
 41:      allocate(cols(0:1),source=[r,r+1])
 42:      call MatSetValues(A, 1, r, 2, cols, vals, INSERT_VALUES,ierr);CHKERRA(ierr)
 43:      deallocate(cols)
 44:      deallocate(vals)
 45:     else if (r == N-1) then
 46:      allocate(vals(0:1))
 47:      vals = 1.0
 48:      allocate(cols(0:1),source=[r-1,r])
 49:      call MatSetValues(A, 1, r, 2, cols, vals, INSERT_VALUES,ierr);CHKERRA(ierr)
 50:      deallocate(cols)
 51:      deallocate(vals)
 52:     else
 53:      allocate(vals(0:2))
 54:      vals = 1.0
 55:      allocate(cols(0:2),source=[r-1,r,r+1])
 56:      call MatSetValues(A, 1, r, 3, cols, vals, INSERT_VALUES,ierr);CHKERRA(ierr)
 57:      deallocate(cols)
 58:      deallocate(vals)
 59:     end if
 60:   end do
 61:   call MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY,ierr);CHKERRA(ierr)
 62:   call MatAssemblyend(A, MAT_FINAL_ASSEMBLY,ierr);CHKERRA(ierr)
 63:   call MatPartitioningCreate(PETSC_COMM_WORLD, part,ierr);CHKERRA(ierr)
 64:   call MatPartitioningSetAdjacency(part, A,ierr);CHKERRA(ierr)
 65:   call MatPartitioningSetFromOptions(part,ierr);CHKERRA(ierr)
 66:   call MatPartitioningApply(part, is,ierr);CHKERRA(ierr)
 67:   call ISView(is, PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr)
 68:   call ISDestroy(is,ierr);CHKERRA(ierr)
 69:   call MatPartitioningDestroy(part,ierr);CHKERRA(ierr)
 70:   call MatDestroy(A,ierr);CHKERRA(ierr)
 71:   call PetscFinalize(ierr);CHKERRA(ierr)
 72: 
 73: end program