xref: /petsc/src/vec/is/sf/tutorials/ex1f.F90 (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown
2*c4762a1bSJed Brown!    Description: A star forest is a simple tree with one root and zero or more leaves.
3*c4762a1bSJed Brown!    Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa.
4*c4762a1bSJed Brown!     This example creates a star forest, communicates values using the graph  views the graph, then destroys it.
5*c4762a1bSJed Brown!
6*c4762a1bSJed Brown!     This is a copy of ex1.c but currently only tests the broadcast operation
7*c4762a1bSJed Brown
8*c4762a1bSJed Brown      program main
9*c4762a1bSJed Brown#include <petsc/finclude/petscvec.h>
10*c4762a1bSJed Brown      use petscvec
11*c4762a1bSJed Brown      implicit none
12*c4762a1bSJed Brown
13*c4762a1bSJed Brown      PetscErrorCode                ierr
14*c4762a1bSJed Brown      PetscInt                      i,nroots,nrootsalloc,nleaves,nleavesalloc,mine(6),stride
15*c4762a1bSJed Brown      type(PetscSFNode)             remote(6)
16*c4762a1bSJed Brown      PetscMPIInt                   rank,size
17*c4762a1bSJed Brown      PetscSF                       sf
18*c4762a1bSJed Brown      PetscInt                      rootdata(6),leafdata(6)
19*c4762a1bSJed Brown
20*c4762a1bSJed Brown! used with PetscSFGetGraph()
21*c4762a1bSJed Brown      type(PetscSFNode), pointer :: gremote(:)
22*c4762a1bSJed Brown      PetscInt, pointer ::          gmine(:)
23*c4762a1bSJed Brown      PetscInt                      gnroots,gnleaves;
24*c4762a1bSJed Brown
25*c4762a1bSJed Brown      call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
26*c4762a1bSJed Brown      if (ierr .ne. 0) then
27*c4762a1bSJed Brown        print*,'Unable to initialize PETSc'
28*c4762a1bSJed Brown        stop
29*c4762a1bSJed Brown      endif
30*c4762a1bSJed Brown      stride = 2
31*c4762a1bSJed Brown      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr);CHKERRA(ierr);
32*c4762a1bSJed Brown      call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr);CHKERRA(ierr);
33*c4762a1bSJed Brown
34*c4762a1bSJed Brown      if (rank == 0) then
35*c4762a1bSJed Brown         nroots = 3
36*c4762a1bSJed Brown      else
37*c4762a1bSJed Brown         nroots = 2
38*c4762a1bSJed Brown      endif
39*c4762a1bSJed Brown      nrootsalloc  = nroots * stride;
40*c4762a1bSJed Brown      if (rank > 0) then
41*c4762a1bSJed Brown         nleaves = 3
42*c4762a1bSJed Brown      else
43*c4762a1bSJed Brown         nleaves = 2
44*c4762a1bSJed Brown      endif
45*c4762a1bSJed Brown      nleavesalloc = nleaves * stride
46*c4762a1bSJed Brown      if (stride > 1) then
47*c4762a1bSJed Brown         do i=1,nleaves
48*c4762a1bSJed Brown            mine(i) = stride * (i-1)
49*c4762a1bSJed Brown         enddo
50*c4762a1bSJed Brown      endif
51*c4762a1bSJed Brown
52*c4762a1bSJed Brown! Left periodic neighbor
53*c4762a1bSJed Brown      remote(1)%rank  = modulo(rank+size-1,size)
54*c4762a1bSJed Brown      remote(1)%index = 1 * stride
55*c4762a1bSJed Brown! Right periodic neighbor
56*c4762a1bSJed Brown      remote(2)%rank  = modulo(rank+1,size)
57*c4762a1bSJed Brown      remote(2)%index = 0 * stride
58*c4762a1bSJed Brown      if (rank > 0) then !               All processes reference rank 0, index
59*c4762a1bSJed Brown         remote(3)%rank  = 0
60*c4762a1bSJed Brown         remote(3)%index = 2 * stride
61*c4762a1bSJed Brown      endif
62*c4762a1bSJed Brown
63*c4762a1bSJed Brown!  Create a star forest for communication
64*c4762a1bSJed Brown      call PetscSFCreate(PETSC_COMM_WORLD,sf,ierr);CHKERRA(ierr);
65*c4762a1bSJed Brown      call PetscSFSetFromOptions(sf,ierr);CHKERRA(ierr);
66*c4762a1bSJed Brown      call PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr);CHKERRA(ierr);
67*c4762a1bSJed Brown      call PetscSFSetUp(sf,ierr);CHKERRA(ierr);
68*c4762a1bSJed Brown
69*c4762a1bSJed Brown!   View graph, mostly useful for debugging purposes.
70*c4762a1bSJed Brown      call PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr);CHKERRA(ierr);
71*c4762a1bSJed Brown      call PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr);
72*c4762a1bSJed Brown      call PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr);
73*c4762a1bSJed Brown
74*c4762a1bSJed Brown!   Allocate space for send and recieve buffers. This example communicates PetscInt, but other types, including
75*c4762a1bSJed Brown!     * user-defined structures, could also be used.
76*c4762a1bSJed Brown!     Set rootdata buffer to be broadcast
77*c4762a1bSJed Brown      do i=1,nrootsalloc
78*c4762a1bSJed Brown         rootdata(i) = -1
79*c4762a1bSJed Brown      enddo
80*c4762a1bSJed Brown      do i=1,nroots
81*c4762a1bSJed Brown         rootdata(1 + (i-1)*stride) = 100*(rank+1) + i - 1
82*c4762a1bSJed Brown      enddo
83*c4762a1bSJed Brown
84*c4762a1bSJed Brown!     Initialize local buffer, these values are never used.
85*c4762a1bSJed Brown      do i=1,nleavesalloc
86*c4762a1bSJed Brown         leafdata(i) = -1
87*c4762a1bSJed Brown      enddo
88*c4762a1bSJed Brown
89*c4762a1bSJed Brown!     Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls.
90*c4762a1bSJed Brown      call PetscSFBcastBegin(sf,MPIU_INTEGER,rootdata,leafdata,ierr);CHKERRA(ierr);
91*c4762a1bSJed Brown      call PetscSFBcastEnd(sf,MPIU_INTEGER,rootdata,leafdata,ierr);CHKERRA(ierr);
92*c4762a1bSJed Brown      call PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n",ierr);CHKERRA(ierr);
93*c4762a1bSJed Brown      call PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr);
94*c4762a1bSJed Brown      call PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n",ierr);CHKERRA(ierr);
95*c4762a1bSJed Brown      call PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr);
96*c4762a1bSJed Brown
97*c4762a1bSJed Brown      call PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr);CHKERRA(ierr);
98*c4762a1bSJed Brown      if (gnleaves .ne. nleaves) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()'); endif
99*c4762a1bSJed Brown      do i=1,nleaves
100*c4762a1bSJed Brown         if (gmine(i) .ne. mine(i)) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()'); endif
101*c4762a1bSJed Brown      enddo
102*c4762a1bSJed Brown      do i=1,nleaves
103*c4762a1bSJed Brown         if (gremote(i)%index .ne. remote(i)%index) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Leaf from PetscSFGetGraph() does not match that set with PetscSFSetGraph()'); endif
104*c4762a1bSJed Brown      enddo
105*c4762a1bSJed Brown
106*c4762a1bSJed Brown      deallocate(gremote)
107*c4762a1bSJed Brown!    Clean storage for star forest.
108*c4762a1bSJed Brown      call PetscSFDestroy(sf,ierr);CHKERRA(ierr);
109*c4762a1bSJed Brown      call PetscFinalize(ierr);
110*c4762a1bSJed Brown  end
111*c4762a1bSJed Brown
112*c4762a1bSJed Brown!/*TEST
113*c4762a1bSJed Brown!  build:
114*c4762a1bSJed Brown!    requires: define(PETSC_HAVE_FORTRAN_TYPE_STAR)
115*c4762a1bSJed Brown!
116*c4762a1bSJed Brown!  test:
117*c4762a1bSJed Brown!    nsize: 3
118*c4762a1bSJed Brown!
119*c4762a1bSJed Brown!TEST*/
120