Actual source code: ex3f90.F90

  1: !
  2: !
  3: !   Description: Demonstrates how users can augment the PETSc profiling by
  4: !                inserting their own event logging.
  5: !
  6: ! -----------------------------------------------------------------------

  8:       program SchoolDay
  9: #include <petsc/finclude/petscsys.h>
 10: #include <petsc/finclude/petsclog.h>
 11:       use petscmpi  ! or mpi or mpi_f08
 12:       use petscsys
 13:       implicit none

 15: !====================================================================
 16: !     Local Variables

 18:       ! Settings:
 19:       integer, parameter        :: verbose=0               ! 0: silent, >=1 : increasing amount of debugging output
 20:       integer, parameter        :: msgLen = 30             ! number of reals which is sent with MPI_Isend
 21:       PetscReal, parameter      :: second=0.1;             ! time is sped up by a factor 10

 23:       ! Codes
 24:       integer, parameter        :: BOY=1, GIRL=2, TEACHER=0
 25:       PetscMPIInt, parameter    :: tagMsg   = 1200;

 27:       ! Timers
 28:       PetscLogEvent :: Morning,  Afternoon
 29:       PetscLogEvent :: PlayBall, SkipRope
 30:       PetscLogEvent :: TidyClass
 31:       PetscLogEvent :: Lessons,  CorrectHomework
 32:       PetscClassId classid

 34:       ! Petsc-stuff
 35:       PetscErrorCode            :: ierr

 37:       ! MPI-stuff
 38:       PetscMPIInt              :: rank, size
 39:       PetscReal, allocatable    :: message(:,:)
 40:       integer                   :: item, maxItem
 41:       integer4                  :: status(MPI_STATUS_SIZE)
 42:       PetscMPIInt                  req

 44:       ! Own stuff
 45:       integer4                  :: role                 ! is this process a BOY, a GIRL or a TEACHER?
 46:       integer4                  :: i, j
 47:       integer4,parameter        :: one=1
 48: !====================================================================
 49: !     Initializations
 50:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 51:       if (ierr .ne. 0) then
 52:         print*,'Unable to initialize PETSc'
 53:         stop
 54:       endif
 55:       call MPI_Comm_size(PETSC_COMM_WORLD, size,ierr)
 56:       call MPI_Comm_rank(PETSC_COMM_WORLD, rank,ierr)

 58:       if (rank==0) then
 59:          role = TEACHER
 60:       else if (rank<0.4*size) then
 61:          role = GIRL
 62:       else
 63:          role = BOY
 64:       end if

 66:       allocate(message(msgLen,msglen))
 67:       do i = 1,msgLen
 68:          do j  = 1,msgLen
 69:             message(i,j) = 10.0*j + i*1.0/(rank+one)
 70:          end do
 71:       end do
 72: !
 73: !====================================================================
 74: !     Create new user-defined events
 75:       classid = 0
 76:       call PetscLogEventRegister('Morning',         classid, Morning,   ierr)
 77:       call PetscLogEventRegister('Afternoon',       classid, Afternoon, ierr)
 78:       call PetscLogEventRegister('Play Ball',       classid, PlayBall,  ierr)
 79:       call PetscLogEventRegister('Skip Rope',       classid, SkipRope,  ierr)
 80:       call PetscLogEventRegister('Tidy Classroom',  classid, TidyClass, ierr)
 81:       call PetscLogEventRegister('Lessons',         classid, Lessons,   ierr)
 82:       call PetscLogEventRegister('Correct Homework',classid,CorrectHomework,          &
 83:      &                                                            ierr)
 84:       if (verbose>=1) then
 85:       print '(a,i0,a)','[',rank,'] SchoolDay events have been defined'
 86:       endif

 88: !====================================================================
 89: !     Go through the school day
 90:       call PetscLogEventBegin(Morning,ierr)

 92:          call PetscLogFlops(190000d0,ierr)
 93:          call PetscSleep(0.5*second,ierr)

 95:          call PetscLogEventBegin(Lessons,ierr)
 96:          call PetscLogFlops(23000d0,ierr)
 97:          call PetscSleep(1*second, ierr)
 98:          if (size>1) then
 99:          call MPI_Isend( message, msgLen, MPI_DOUBLE_PRECISION,                             &
100:      &                        mod(rank+1,size),                                             &
101:      &                        tagMsg+rank, PETSC_COMM_WORLD, req, ierr)
102:          call  MPI_Recv( message, msgLen, MPI_DOUBLE_PRECISION,                             &
103:      &                       mod(rank-1+size,size),                                         &
104:      &                  tagMsg+mod(rank-1+size,size), PETSC_COMM_WORLD,                     &
105:      &        status, ierr)
106:          call MPI_Wait(req,MPI_STATUS_IGNORE,ierr)
107:          end if
108:          call PetscLogEventEnd(Lessons,ierr)

110:          if (role==TEACHER) then
111:             call PetscLogEventBegin(TidyClass,ierr)
112:             call PetscLogFlops(600000d0,ierr)
113:             call PetscSleep(0.6*second, ierr)
114:                call PetscLogEventBegin(CorrectHomework,ierr)
115:                call PetscLogFlops(234700d0,ierr)
116:                call PetscSleep(0.4*second, ierr)
117:                call PetscLogEventEnd(CorrectHomework,ierr)
118:             call PetscLogEventEnd(TidyClass,ierr)
119:          else if (role==BOY) then
120:             call PetscLogEventBegin(SkipRope,ierr)
121:             call PetscSleep(0.8*second, ierr)
122:             call PetscLogEventEnd(SkipRope,ierr)
123:          else
124:             call PetscLogEventBegin(PlayBall,ierr)
125:             call PetscSleep(0.9*second, ierr)
126:             call PetscLogEventEnd(PlayBall,ierr)
127:          end if

129:          call PetscLogEventBegin(Lessons,ierr)
130:          call PetscLogFlops(120000d0,ierr)
131:          call PetscSleep(0.7*second, ierr)
132:          call PetscLogEventEnd(Lessons,ierr)

134:       call PetscLogEventEnd(Morning,ierr)

136:       call PetscLogEventBegin(Afternoon,ierr)

138:          item = rank*(3-rank)
139:          call MPI_Allreduce(item, maxItem, 1, MPI_INTEGER, MPI_MAX,                    &
140:      &                           PETSC_COMM_WORLD, ierr)

142:          item = rank*(10-rank)
143:          call MPI_Allreduce(item, maxItem, 1, MPI_INTEGER, MPI_MAX,                    &
144:      &                           PETSC_COMM_WORLD, ierr)

146:          call PetscLogFlops(58988d0,ierr)
147:          call PetscSleep(0.6*second,ierr)

149:          call PetscLogEventBegin(Lessons,ierr)
150:          call PetscLogFlops(123456d0,ierr)
151:          call PetscSleep(1*second, ierr)
152:          call PetscLogEventEnd(Lessons,ierr)

154:          if (role==TEACHER) then
155:             call PetscLogEventBegin(TidyClass,ierr)
156:             call PetscLogFlops(17800d0,ierr)
157:             call PetscSleep(1.1*second, ierr)
158:             call PetscLogEventBegin(Lessons,ierr)
159:             call PetscLogFlops(72344d0,ierr)
160:             call PetscSleep(0.5*second, ierr)
161:             call PetscLogEventEnd(Lessons,ierr)
162:             call PetscLogEventEnd(TidyClass,ierr)
163:          else if (role==GIRL) then
164:             call PetscLogEventBegin(SkipRope,ierr)
165:             call PetscSleep(0.7*second, ierr)
166:             call PetscLogEventEnd(SkipRope,ierr)
167:          else
168:             call PetscLogEventBegin(PlayBall,ierr)
169:             call PetscSleep(0.8*second, ierr)
170:             call PetscLogEventEnd(PlayBall,ierr)
171:          end if

173:          call PetscLogEventBegin(Lessons,ierr)
174:          call PetscLogFlops(72344d0,ierr)
175:          call PetscSleep(0.5*second, ierr)
176:          call PetscLogEventEnd(Lessons,ierr)

178:       call PetscLogEventEnd(Afternoon,ierr)

180:       if (.false.) then
181:          continue
182:       else if (role==TEACHER) then
183:          call PetscLogEventBegin(TidyClass,ierr)
184:          call PetscLogFlops(612300d0,ierr)
185:          call PetscSleep(1.1*second, ierr)
186:          call PetscLogEventEnd(TidyClass,ierr)
187:          call PetscLogEventBegin(CorrectHomework,ierr)
188:          call PetscLogFlops(234700d0,ierr)
189:          call PetscSleep(1.1*second, ierr)
190:          call PetscLogEventEnd(CorrectHomework,ierr)
191:       else
192:          call PetscLogEventBegin(SkipRope,ierr)
193:          call PetscSleep(0.7*second, ierr)
194:          call PetscLogEventEnd(SkipRope,ierr)
195:          call PetscLogEventBegin(PlayBall,ierr)
196:          call PetscSleep(0.8*second, ierr)
197:          call PetscLogEventEnd(PlayBall,ierr)
198:       end if

200:       call PetscLogEventBegin(Lessons,ierr)
201:       call PetscLogFlops(120000d0,ierr)
202:       call PetscSleep(0.7*second, ierr)
203:       call PetscLogEventEnd(Lessons,ierr)

205:       call PetscSleep(0.25*second,ierr)

207:       call PetscLogEventBegin(Morning,ierr)

209:          call PetscLogFlops(190000d0,ierr)
210:          call PetscSleep(0.5*second,ierr)

212:          call PetscLogEventBegin(Lessons,ierr)
213:          call PetscLogFlops(23000d0,ierr)
214:          call PetscSleep(1*second, ierr)
215:          if (size>1) then
216:          call MPI_Isend( message, msgLen, MPI_DOUBLE_PRECISION,                             &
217:      &                        mod(rank+1,size),                                             &
218:      &                   tagMsg+rank, PETSC_COMM_WORLD, req, ierr)
219:          call MPI_Recv( message, msgLen, MPI_DOUBLE_PRECISION,                              &
220:      &                  mod(rank-1+size,size),                                              &
221:      &                  tagMsg+mod(rank-1+size,size), PETSC_COMM_WORLD,                     &
222:      &                   status, ierr)
223:          call MPI_Wait(req,MPI_STATUS_IGNORE,ierr)
224:          end if
225:          call PetscLogEventEnd(Lessons,ierr)

227:          if (role==TEACHER) then
228:             call PetscLogEventBegin(TidyClass,ierr)
229:             call PetscLogFlops(600000d0,ierr)
230:             call PetscSleep(1.2*second, ierr)
231:             call PetscLogEventEnd(TidyClass,ierr)
232:          else if (role==BOY) then
233:             call PetscLogEventBegin(SkipRope,ierr)
234:             call PetscSleep(0.8*second, ierr)
235:             call PetscLogEventEnd(SkipRope,ierr)
236:          else
237:             call PetscLogEventBegin(PlayBall,ierr)
238:             call PetscSleep(0.9*second, ierr)
239:             call PetscLogEventEnd(PlayBall,ierr)
240:          end if

242:          call PetscLogEventBegin(Lessons,ierr)
243:          call PetscLogFlops(120000d0,ierr)
244:          call PetscSleep(0.7*second, ierr)
245:          call PetscLogEventEnd(Lessons,ierr)

247:       call PetscLogEventEnd(Morning,ierr)

249:       deallocate(message)

251:       call PetscFinalize(ierr)

253:       end program SchoolDay

255: !/*TEST
256: !
257: ! testset:
258: !   suffix: no_log
259: !   requires: !defined(PETSC_USE_LOG)
260: !   test:
261: !     suffix: ascii
262: !     args: -log_view ascii:filename.txt
263: !   test:
264: !     suffix: detail
265: !     args: -log_view ascii:filename.txt:ascii_info_detail
266: !   test:
267: !     suffix: xml
268: !     args: -log_view ascii:filename.xml:ascii_xml
269: !
270: ! testset:
271: !   args: -log_view ascii:filename.txt
272: !   output_file: output/ex3f90.out
273: !   requires: defined(PETSC_USE_LOG)
274: !   test:
275: !     suffix: 1
276: !     nsize: 1
277: !   test:
278: !     suffix: 2
279: !     nsize: 2
280: !   test:
281: !     suffix: 3
282: !     nsize: 3
283: !
284: ! testset:
285: !   suffix: detail
286: !   args: -log_view ascii:filename.txt:ascii_info_detail
287: !   output_file: output/ex3f90.out
288: !   requires: defined(PETSC_USE_LOG)
289: !   test:
290: !     suffix: 1
291: !     nsize: 1
292: !   test:
293: !     suffix: 2
294: !     nsize: 2
295: !   test:
296: !     suffix: 3
297: !     nsize: 3
298: !
299: ! testset:
300: !   suffix: xml
301: !   args: -log_view ascii:filename.xml:ascii_xml
302: !   output_file: output/ex3f90.out
303: !   requires: defined(PETSC_USE_LOG)
304: !   test:
305: !     suffix: 1
306: !     nsize: 1
307: !   test:
308: !     suffix: 2
309: !     nsize: 2
310: !   test:
311: !     suffix: 3
312: !     nsize: 3
313: !
314: !TEST*/