Actual source code: ex30f.F90

  1: !
  2: !
  3: !  Tests parallel to parallel scatter where a to from index are
  4: !  duplicated
  5: #include <petsc/finclude/petscvec.h>
  6: program main
  7:   use petscvec
  8:   implicit none

 10:   PetscErrorCode ierr
 11:   PetscInt, parameter :: nlocal = 2, n = 8
 12:   PetscInt row
 13:   PetscMPIInt rank, size
 14:   PetscInt, parameter, dimension(8) :: &
 15:     from = [1, 5, 9, 13, 3, 7, 11, 15], &
 16:     to = [0, 0, 0, 0, 7, 7, 7, 7]
 17:   PetscScalar num
 18:   Vec v1, v2, v3
 19:   VecScatter scat1, scat2
 20:   IS fromis, tois
 21:   PetscCallA(PetscInitialize(ierr))
 22:   PetscCallMPIA(MPI_COMM_RANK(PETSC_COMM_WORLD, rank, ierr))
 23:   PetscCallMPIA(MPI_COMM_SIZE(PETSC_COMM_WORLD, size, ierr))
 24:   if (size /= 4) then
 25:     print *, 'Four processor test'
 26:     stop
 27:   end if

 29:   PetscCallA(VecCreateFromOptions(PETSC_COMM_WORLD, PETSC_NULL_CHARACTER, 1_PETSC_INT_KIND, nlocal*2, n*2, v1, ierr))
 30:   PetscCallA(VecCreateFromOptions(PETSC_COMM_WORLD, PETSC_NULL_CHARACTER, 1_PETSC_INT_KIND, nlocal, n, v2, ierr))
 31:   PetscCallA(VecCreateSeq(PETSC_COMM_SELF, n, v3, ierr))

 33:   num = 2.0
 34:   row = 1
 35:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 36:   row = 5
 37:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 38:   row = 9
 39:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 40:   row = 13
 41:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 42:   num = 1.0
 43:   row = 15
 44:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 45:   row = 3
 46:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 47:   row = 7
 48:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))
 49:   row = 11
 50:   PetscCallA(VecSetValue(v1, row, num, INSERT_VALUES, ierr))

 52:   PetscCallA(VecAssemblyBegin(v1, ierr))
 53:   PetscCallA(VecAssemblyEnd(v1, ierr))

 55:   num = 0.0
 56:   PetscCallA(VecScale(v2, num, ierr))
 57:   PetscCallA(VecScale(v3, num, ierr))

 59:   PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, 8_PETSC_INT_KIND, from, PETSC_COPY_VALUES, fromis, ierr))
 60:   PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, 8_PETSC_INT_KIND, to, PETSC_COPY_VALUES, tois, ierr))
 61:   PetscCallA(VecScatterCreate(v1, fromis, v2, tois, scat1, ierr))
 62:   PetscCallA(VecScatterCreate(v1, fromis, v3, tois, scat2, ierr))
 63:   PetscCallA(ISDestroy(fromis, ierr))
 64:   PetscCallA(ISDestroy(tois, ierr))

 66:   PetscCallA(VecScatterBegin(scat1, v1, v2, ADD_VALUES, SCATTER_FORWARD, ierr))
 67:   PetscCallA(VecScatterEnd(scat1, v1, v2, ADD_VALUES, SCATTER_FORWARD, ierr))

 69:   PetscCallA(VecScatterBegin(scat2, v1, v3, ADD_VALUES, SCATTER_FORWARD, ierr))
 70:   PetscCallA(VecScatterEnd(scat2, v1, v3, ADD_VALUES, SCATTER_FORWARD, ierr))

 72:   PetscCallA(PetscObjectSetName(v1, 'V1', ierr))
 73:   PetscCallA(VecView(v1, PETSC_VIEWER_STDOUT_WORLD, ierr))

 75:   PetscCallA(PetscObjectSetName(v2, 'V2', ierr))
 76:   PetscCallA(VecView(v2, PETSC_VIEWER_STDOUT_WORLD, ierr))

 78:   if (rank == 0) then
 79:     PetscCallA(PetscObjectSetName(v3, 'V3', ierr))
 80:     PetscCallA(VecView(v3, PETSC_VIEWER_STDOUT_SELF, ierr))
 81:   end if

 83:   PetscCallA(VecScatterDestroy(scat1, ierr))
 84:   PetscCallA(VecScatterDestroy(scat2, ierr))
 85:   PetscCallA(VecDestroy(v1, ierr))
 86:   PetscCallA(VecDestroy(v2, ierr))
 87:   PetscCallA(VecDestroy(v3, ierr))

 89:   PetscCallA(PetscFinalize(ierr))

 91: end

 93: !/*TEST
 94: !
 95: !     test:
 96: !       nsize: 4
 97: !
 98: !TEST*/