Actual source code: ex1f.F90

  1: !
  2: !
  3: !   Description: Demonstrates how users can augment the PETSc profiling by
  4: !                inserting their own event logging.
  5: !
  6: !/*T
  7: !   Concepts: PetscLog^user-defined event profiling (basic example);
  8: !   Concepts: PetscLog^activating/deactivating events for profiling (basic example);
  9: !   Processors: n
 10: !T*/
 11: ! -----------------------------------------------------------------------

 13:       program SchoolDay
 14: #include <petsc/finclude/petscsys.h>
 15: #include <petsc/finclude/petsclog.h>
 16:       use petscmpi  ! or mpi or mpi_f08
 17:       use petscsys
 18:       implicit none

 20: !====================================================================
 21: !     Local Variables

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

 28:       ! Codes
 29:       integer, parameter        :: BOY=1, GIRL=2, TEACHER=0
 30:       PetscMPIInt, parameter    :: tagMsg   = 1200;

 32:       ! Timers
 33:       PetscLogEvent :: Morning,  Afternoon
 34:       PetscLogEvent :: PlayBall, SkipRope
 35:       PetscLogEvent :: TidyClass
 36:       PetscLogEvent :: Lessons,  CorrectHomework
 37:       PetscClassId classid

 39:       ! Petsc-stuff
 40:       PetscErrorCode            :: ierr

 42:       ! MPI-stuff
 43:       PetscMPIInt              :: rank, size
 44:       PetscReal, allocatable    :: message(:,:)
 45:       integer                   :: item, maxItem
 46:       integer4                  :: status(MPI_STATUS_SIZE)
 47:       PetscMPIInt                  req

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

 63:       if (rank==0) then
 64:          role = TEACHER
 65:       else if (rank<0.4*size) then
 66:          role = GIRL
 67:       else
 68:          role = BOY
 69:       end if

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

 93: !====================================================================
 94: !     Go through the school day
 95:       call PetscLogEventBegin(Morning,ierr)

 97:          call PetscLogFlops(190000d0,ierr)
 98:          call PetscSleep(0.5*second,ierr)

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

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

134:          call PetscLogEventBegin(Lessons,ierr)
135:          call PetscLogFlops(120000d0,ierr)
136:          call PetscSleep(0.7*second, ierr)
137:          call PetscLogEventEnd(Lessons,ierr)

139:       call PetscLogEventEnd(Morning,ierr)

141:       call PetscLogEventBegin(Afternoon,ierr)

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

147:          item = rank*(10-rank)
148:          call MPI_Allreduce(item, maxItem, 1, MPI_INTEGER, MPI_MAX,                    &
149:      &                           PETSC_COMM_WORLD, ierr)

151:          call PetscLogFlops(58988d0,ierr)
152:          call PetscSleep(0.6*second,ierr)

154:          call PetscLogEventBegin(Lessons,ierr)
155:          call PetscLogFlops(123456d0,ierr)
156:          call PetscSleep(1*second, ierr)
157:          call PetscLogEventEnd(Lessons,ierr)

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

178:          call PetscLogEventBegin(Lessons,ierr)
179:          call PetscLogFlops(72344d0,ierr)
180:          call PetscSleep(0.5*second, ierr)
181:          call PetscLogEventEnd(Lessons,ierr)

183:       call PetscLogEventEnd(Afternoon,ierr)

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

205:       call PetscLogEventBegin(Lessons,ierr)
206:       call PetscLogFlops(120000d0,ierr)
207:       call PetscSleep(0.7*second, ierr)
208:       call PetscLogEventEnd(Lessons,ierr)

210:       call PetscSleep(0.25*second,ierr)

212:       call PetscLogEventBegin(Morning,ierr)

214:          call PetscLogFlops(190000d0,ierr)
215:          call PetscSleep(0.5*second,ierr)

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

232:          if (role==TEACHER) then
233:             call PetscLogEventBegin(TidyClass,ierr)
234:             call PetscLogFlops(600000d0,ierr)
235:             call PetscSleep(1.2*second, ierr)
236:             call PetscLogEventEnd(TidyClass,ierr)
237:          else if (role==BOY) then
238:             call PetscLogEventBegin(SkipRope,ierr)
239:             call PetscSleep(0.8*second, ierr)
240:             call PetscLogEventEnd(SkipRope,ierr)
241:          else
242:             call PetscLogEventBegin(PlayBall,ierr)
243:             call PetscSleep(0.9*second, ierr)
244:             call PetscLogEventEnd(PlayBall,ierr)
245:          end if

247:          call PetscLogEventBegin(Lessons,ierr)
248:          call PetscLogFlops(120000d0,ierr)
249:          call PetscSleep(0.7*second, ierr)
250:          call PetscLogEventEnd(Lessons,ierr)

252:       call PetscLogEventEnd(Morning,ierr)

254:       deallocate(message)

256:       call PetscFinalize(ierr)

258:       end program SchoolDay

260: !/*TEST
261: !
262: ! testset:
263: !   args: -log_view ascii:filename.txt
264: !   output_file: output/ex1f.out
265: !   test:
266: !     suffix: 1
267: !     nsize: 1
268: !   test:
269: !     suffix: 2
270: !     nsize: 2
271: !   test:
272: !     suffix: 3
273: !     nsize: 3
274: !
275: ! testset:
276: !   suffix: detail
277: !   args: -log_view ascii:filename.txt:ascii_info_detail
278: !   output_file: output/ex1f.out
279: !   test:
280: !     suffix: 1
281: !     nsize: 1
282: !   test:
283: !     suffix: 2
284: !     nsize: 2
285: !   test:
286: !     suffix: 3
287: !     nsize: 3
288: !
289: ! testset:
290: !   suffix: xml
291: !   args: -log_view ascii:filename.xml:ascii_xml
292: !   output_file: output/ex1f.out
293: !   test:
294: !     suffix: 1
295: !     nsize: 1
296: !   test:
297: !     suffix: 2
298: !     nsize: 2
299: !   test:
300: !     suffix: 3
301: !     nsize: 3
302: !
303: !TEST*/