32 public :: prof_papi_rapstart
33 public :: prof_papi_rapstop
34 public :: prof_papi_rapreport
37 public :: prof_valcheck
39 interface prof_valcheck
40 module procedure prof_valcheck_sp_1d
41 module procedure prof_valcheck_sp_2d
42 module procedure prof_valcheck_sp_3d
43 module procedure prof_valcheck_sp_4d
44 module procedure prof_valcheck_sp_5d
45 module procedure prof_valcheck_sp_6d
46 module procedure prof_valcheck_dp_1d
47 module procedure prof_valcheck_dp_2d
48 module procedure prof_valcheck_dp_3d
49 module procedure prof_valcheck_dp_4d
50 module procedure prof_valcheck_dp_5d
51 module procedure prof_valcheck_dp_6d
52 end interface prof_valcheck
68 integer,
private,
parameter :: prof_rapnlimit = 300
69 character(len=H_SHORT),
private :: prof_prefix =
'' 70 integer,
private :: prof_rapnmax = 0
71 character(len=H_SHORT),
private :: prof_rapname(prof_rapnlimit)
72 integer,
private :: prof_grpnmax = 0
73 character(len=H_SHORT),
private :: prof_grpname(prof_rapnlimit)
74 integer,
private :: prof_grpid (prof_rapnlimit)
75 real(DP),
private :: prof_raptstr(prof_rapnlimit)
76 real(DP),
private :: prof_rapttot(prof_rapnlimit)
77 integer,
private :: prof_rapnstr(prof_rapnlimit)
78 integer,
private :: prof_rapnend(prof_rapnlimit)
79 integer,
private :: prof_raplevel(prof_rapnlimit)
81 integer,
private,
parameter :: prof_default_rap_level = 2
82 integer,
private :: prof_rap_level = 2
83 logical,
private :: prof_mpi_barrier = .false.
86 integer(DP),
private :: prof_papi_flops = 0
87 real(SP),
private :: prof_papi_real_time = 0.0
88 real(SP),
private :: prof_papi_proc_time = 0.0
89 real(SP),
private :: prof_papi_mflops = 0.0
90 integer,
private :: prof_papi_check
93 character(len=7),
private :: prof_header
94 character(len=16),
private :: prof_item
95 real(DP),
private :: prof_max
96 real(DP),
private :: prof_min
97 real(DP),
private :: prof_sum
107 namelist / param_prof / &
114 log_info(
"PRF_setup",*)
'Setup' 120 log_info(
"PROF_setup",*)
'Not found namelist. Default used.' 121 elseif( ierr > 0 )
then 122 log_error(
"PROF_setup",*)
'Not appropriate names in namelist PARAM_PROF. Check!' 128 log_info(
"PROF_setup",*)
'Rap output level = ', prof_rap_level
129 log_info(
"PROF_setup",*)
'Add MPI_barrier in every rap? = ', prof_mpi_barrier
141 character(len=*),
intent(in) :: prefxname
145 if ( prefxname ==
'' )
then 148 prof_prefix = trim(prefxname)//
'_' 162 character(len=*),
intent(in) :: rapname_base
163 integer,
intent(in),
optional :: level
165 character(len=H_SHORT) :: rapname
177 if (
present(level) )
then 180 level_ = prof_default_rap_level
183 if( level_ > prof_rap_level )
return 185 rapname = trim(prof_prefix)//trim(rapname_base)
187 id = get_rapid( rapname, level_ )
192 prof_rapnstr(id) = prof_rapnstr(id) + 1
198 call fapp_start( trim(prof_grpname(get_grpid(rapname))), id, level_ )
201 call start_collection( rapname )
215 character(len=*),
intent(in) :: rapname_base
216 integer,
intent(in),
optional :: level
218 character(len=H_SHORT) :: rapname
230 if (
present(level) )
then 231 if( level > prof_rap_level )
return 234 rapname = trim(prof_prefix)//trim(rapname_base)
236 id = get_rapid( rapname, level_ )
238 if( level_ > prof_rap_level )
return 242 prof_rapttot(id) = prof_rapttot(id) + (
prc_mpitime()-prof_raptstr(id) )
243 prof_rapnend(id) = prof_rapnend(id) + 1
246 call stop_collection( rapname )
249 call fapp_stop( trim(prof_grpname(prof_grpid(id))), id, level_ )
263 real(DP) :: avgvar(prof_rapnlimit)
264 real(DP) :: maxvar(prof_rapnlimit)
265 real(DP) :: minvar(prof_rapnlimit)
266 integer :: maxidx(prof_rapnlimit)
267 integer :: minidx(prof_rapnlimit)
273 do id = 1, prof_rapnmax
274 if ( prof_rapnstr(id) /= prof_rapnend(id) )
then 275 log_warn(
"PROF_rapreport",*)
'Mismatch Report',id,prof_rapname(id),prof_rapnstr(id),prof_rapnend(id)
280 log_info(
"PROF_rapreport",
'(1x,A,I2,A)')
'Computational Time Report (Rap level = ', prof_rap_level,
')' 284 do gid = 1, prof_rapnmax
285 do id = 1, prof_rapnmax
286 if ( prof_raplevel(id) <= prof_rap_level &
287 .AND. prof_grpid(id) == gid )
then 288 log_info_cont(
'(1x,A,I3.3,A,A,A,F10.3,A,I9)') &
289 'ID=',id,
' : ',prof_rapname(id),
' T=',prof_rapttot(id),
' N=',prof_rapnstr(id)
297 maxvar(1:prof_rapnmax), &
298 minvar(1:prof_rapnmax), &
299 maxidx(1:prof_rapnmax), &
300 minidx(1:prof_rapnmax), &
301 prof_rapttot(1:prof_rapnmax) )
306 write(*,*)
'INFO [PROF_rapreport] Computational Time Report' 313 do gid = 1, prof_rapnmax
314 do id = 1, prof_rapnmax
315 if ( prof_raplevel(id) <= prof_rap_level &
316 .AND. prof_grpid(id) == gid &
318 write(fid,
'(6x,A,I3.3,3A,F10.3,A,F10.3,A,I5,2A,F10.3,A,I5,2A,I9)') &
319 'ID=',id,
' : ',prof_rapname(id), &
320 ' T(avg)=',avgvar(id), &
321 ', T(max)=',maxvar(id),
'[',maxidx(id),
']', &
322 ', T(min)=',minvar(id),
'[',minidx(id),
']', &
323 ', N=',prof_rapnstr(id)
336 subroutine prof_papi_rapstart
340 call papif_flops( prof_papi_real_time, prof_papi_proc_time, prof_papi_flops, prof_papi_mflops, prof_papi_check )
343 end subroutine prof_papi_rapstart
347 subroutine prof_papi_rapstop
351 call papif_flops( prof_papi_real_time, prof_papi_proc_time, prof_papi_flops, prof_papi_mflops, prof_papi_check )
354 end subroutine prof_papi_rapstop
358 subroutine prof_papi_rapreport
365 real(DP) :: avgvar(3)
366 real(DP) :: maxvar(3)
367 real(DP) :: minvar(3)
372 real(DP) :: PROF_PAPI_gflop
373 real(DP) :: statistics(3)
376 prof_papi_gflop =
real(PROF_PAPI_flops,kind=8) / 1024.0_DP**3
381 log_info(
"PROF_PAPI_rapreport",*)
'PAPI Report [Local PE information]' 382 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)')
'Real time [sec] : ', prof_papi_real_time
383 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)')
'CPU time [sec] : ', prof_papi_proc_time
384 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)')
'FLOP [GFLOP] : ', prof_papi_gflop
385 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)')
'FLOPS by PAPI [GFLOPS] : ', prof_papi_mflops/1024.0_dp
386 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)')
'FLOP / CPU Time [GFLOPS] : ', prof_papi_gflop/prof_papi_proc_time
389 statistics(1) =
real(prof_papi_real_time,kind=8)
390 statistics(2) =
real(prof_papi_proc_time,kind=8)
391 statistics(3) = prof_papi_gflop
400 zerosw = 0.5_dp - sign(0.5_dp,maxvar(2)-1.d-12)
403 log_info(
"PROF_PAPI_rapreport",*)
'PAPI Report' 404 log_info(
"PROF_PAPI_rapreport",
'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
405 'Real time [sec]',
' T(avg)=',avgvar(1), &
406 ', T(max)=',maxvar(1),
'[',maxidx(1),
']',
', T(min)=',minvar(1),
'[',minidx(1),
']' 407 log_info(
"PROF_PAPI_rapreport",
'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
408 'CPU time [sec]',
' T(avg)=',avgvar(2), &
409 ', T(max)=',maxvar(2),
'[',maxidx(2),
']',
', T(min)=',minvar(2),
'[',minidx(2),
']' 410 log_info(
"PROF_PAPI_rapreport",
'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
411 'FLOP [GFLOP]',
' N(avg)=',avgvar(3), &
412 ', N(max)=',maxvar(3),
'[',maxidx(3),
']',
', N(min)=',minvar(3),
'[',minidx(3),
']' 414 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3,A,I6,A)') &
416 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)') &
417 'FLOPS [GFLOPS] : ', avgvar(3)*
prc_nprocs * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
418 log_info(
"PROF_PAPI_rapreport",
'(1x,A,F15.3)') &
419 'FLOPS per PE [GFLOPS] : ', avgvar(3) * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
425 write(*,*)
'*** PAPI Report' 426 write(*,
'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
427 '*** Real time [sec]',
' T(avg)=',avgvar(1), &
428 ', T(max)=',maxvar(1),
'[',maxidx(1),
']',
', T(min)=',minvar(1),
'[',minidx(1),
']' 429 write(*,
'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
430 '*** CPU time [sec]',
' T(avg)=',avgvar(2), &
431 ', T(max)=',maxvar(2),
'[',maxidx(2),
']',
', T(min)=',minvar(2),
'[',minidx(2),
']' 432 write(*,
'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
433 '*** FLOP [GFLOP]',
' N(avg)=',avgvar(3), &
434 ', N(max)=',maxvar(3),
'[',maxidx(3),
']',
', N(min)=',minvar(3),
'[',minidx(3),
']' 436 write(*,
'(1x,A,F15.3,A,I6,A)') &
438 write(*,
'(1x,A,F15.3)') &
439 '*** FLOPS [GFLOPS] : ', avgvar(3)*
prc_nprocs * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
440 write(*,
'(1x,A,F15.3)') &
441 '*** FLOPS per PE [GFLOPS] : ', avgvar(3) * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
447 end subroutine prof_papi_rapreport
452 function get_rapid( rapname, level )
result(id)
455 character(len=*),
intent(in) :: rapname
456 integer,
intent(inout) :: level
458 character (len=H_SHORT) :: trapname
463 trapname = trim(rapname)
465 do id = 1, prof_rapnmax
466 if ( trapname == prof_rapname(id) )
then 467 level = prof_raplevel(id)
472 prof_rapnmax = prof_rapnmax + 1
474 prof_rapname(id) = trapname
478 prof_rapttot(id) = 0.0_dp
480 prof_grpid(id) = get_grpid(trapname)
481 prof_raplevel(id) = level
484 end function get_rapid
488 function get_grpid( rapname )
result(gid)
491 character(len=*),
intent(in) :: rapname
493 character(len=H_SHORT) :: grpname
499 idx = index(rapname,
" ")
501 grpname = rapname(1:idx-1)
506 do gid = 1, prof_grpnmax
507 if( grpname == prof_grpname(gid) )
return 510 prof_grpnmax = prof_grpnmax + 1
512 prof_grpname(gid) = grpname
515 end function get_grpid
518 subroutine prof_valcheck_sp_1d( &
524 character(len=*),
intent(in) :: header
525 character(len=*),
intent(in) :: varname
526 real(SP),
intent(in) :: var(:)
529 prof_header = trim(header)
530 prof_item = trim(varname)
531 prof_max =
real(maxval(var),kind=
dp)
532 prof_min =
real(minval(var),kind=
dp)
533 prof_sum =
real(sum (var),kind=
dp)
534 log_info(
"PROF_valcheck_SP_1D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
535 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
538 end subroutine prof_valcheck_sp_1d
541 subroutine prof_valcheck_sp_2d( &
547 character(len=*),
intent(in) :: header
548 character(len=*),
intent(in) :: varname
549 real(SP),
intent(in) :: var(:,:)
552 prof_header = trim(header)
553 prof_item = trim(varname)
554 prof_max =
real(maxval(var),kind=
dp)
555 prof_min =
real(minval(var),kind=
dp)
556 prof_sum =
real(sum (var),kind=
dp)
557 log_info(
"PROF_valcheck_SP_2D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
558 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
561 end subroutine prof_valcheck_sp_2d
564 subroutine prof_valcheck_sp_3d( &
570 character(len=*),
intent(in) :: header
571 character(len=*),
intent(in) :: varname
572 real(SP),
intent(in) :: var(:,:,:)
575 prof_header = trim(header)
576 prof_item = trim(varname)
577 prof_max =
real(maxval(var),kind=
dp)
578 prof_min =
real(minval(var),kind=
dp)
579 prof_sum =
real(sum (var),kind=
dp)
580 log_info(
"PROF_valcheck_SP_3D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
581 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
584 end subroutine prof_valcheck_sp_3d
587 subroutine prof_valcheck_sp_4d( &
593 character(len=*),
intent(in) :: header
594 character(len=*),
intent(in) :: varname
595 real(SP),
intent(in) :: var(:,:,:,:)
598 prof_header = trim(header)
599 prof_item = trim(varname)
600 prof_max =
real(maxval(var),kind=
dp)
601 prof_min =
real(minval(var),kind=
dp)
602 prof_sum =
real(sum (var),kind=
dp)
603 log_info(
"PROF_valcheck_SP_4D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
604 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
607 end subroutine prof_valcheck_sp_4d
610 subroutine prof_valcheck_sp_5d( &
616 character(len=*),
intent(in) :: header
617 character(len=*),
intent(in) :: varname
618 real(SP),
intent(in) :: var(:,:,:,:,:)
621 prof_header = trim(header)
622 prof_item = trim(varname)
623 prof_max =
real(maxval(var),kind=
dp)
624 prof_min =
real(minval(var),kind=
dp)
625 prof_sum =
real(sum (var),kind=
dp)
626 log_info(
"PROF_valcheck_SP_5D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
627 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
630 end subroutine prof_valcheck_sp_5d
633 subroutine prof_valcheck_sp_6d( &
639 character(len=*),
intent(in) :: header
640 character(len=*),
intent(in) :: varname
641 real(SP),
intent(in) :: var(:,:,:,:,:,:)
644 prof_header = trim(header)
645 prof_item = trim(varname)
646 prof_max =
real(maxval(var),kind=
dp)
647 prof_min =
real(minval(var),kind=
dp)
648 prof_sum =
real(sum (var),kind=
dp)
649 log_info(
"PROF_valcheck_SP_6D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
650 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
653 end subroutine prof_valcheck_sp_6d
656 subroutine prof_valcheck_dp_1d( &
662 character(len=*),
intent(in) :: header
663 character(len=*),
intent(in) :: varname
664 real(DP),
intent(in) :: var(:)
667 prof_header = trim(header)
668 prof_item = trim(varname)
669 prof_max =
real(maxval(var),kind=
dp)
670 prof_min =
real(minval(var),kind=
dp)
671 prof_sum =
real(sum (var),kind=
dp)
672 log_info(
"PROF_valcheck_DP_1D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
673 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
676 end subroutine prof_valcheck_dp_1d
679 subroutine prof_valcheck_dp_2d( &
685 character(len=*),
intent(in) :: header
686 character(len=*),
intent(in) :: varname
687 real(DP),
intent(in) :: var(:,:)
690 prof_header = trim(header)
691 prof_item = trim(varname)
692 prof_max =
real(maxval(var),kind=
dp)
693 prof_min =
real(minval(var),kind=
dp)
694 prof_sum =
real(sum (var),kind=
dp)
695 log_info(
"PROF_valcheck_DP_2D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
696 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
699 end subroutine prof_valcheck_dp_2d
702 subroutine prof_valcheck_dp_3d( &
708 character(len=*),
intent(in) :: header
709 character(len=*),
intent(in) :: varname
710 real(DP),
intent(in) :: var(:,:,:)
713 prof_header = trim(header)
714 prof_item = trim(varname)
715 prof_max =
real(maxval(var),kind=
dp)
716 prof_min =
real(minval(var),kind=
dp)
717 prof_sum =
real(sum (var),kind=
dp)
718 log_info(
"PROF_valcheck_DP_3D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
719 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
722 end subroutine prof_valcheck_dp_3d
725 subroutine prof_valcheck_dp_4d( &
731 character(len=*),
intent(in) :: header
732 character(len=*),
intent(in) :: varname
733 real(DP),
intent(in) :: var(:,:,:,:)
736 prof_header = trim(header)
737 prof_item = trim(varname)
738 prof_max =
real(maxval(var),kind=
dp)
739 prof_min =
real(minval(var),kind=
dp)
740 prof_sum =
real(sum (var),kind=
dp)
741 log_info(
"PROF_valcheck_DP_4D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
742 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
745 end subroutine prof_valcheck_dp_4d
748 subroutine prof_valcheck_dp_5d( &
754 character(len=*),
intent(in) :: header
755 character(len=*),
intent(in) :: varname
756 real(DP),
intent(in) :: var(:,:,:,:,:)
759 prof_header = trim(header)
760 prof_item = trim(varname)
761 prof_max =
real(maxval(var),kind=
dp)
762 prof_min =
real(minval(var),kind=
dp)
763 prof_sum =
real(sum (var),kind=
dp)
764 log_info(
"PROF_valcheck_DP_5D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
765 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
768 end subroutine prof_valcheck_dp_5d
771 subroutine prof_valcheck_dp_6d( &
777 character(len=*),
intent(in) :: header
778 character(len=*),
intent(in) :: varname
779 real(DP),
intent(in) :: var(:,:,:,:,:,:)
782 prof_header = trim(header)
783 prof_item = trim(varname)
784 prof_max =
real(maxval(var),kind=
dp)
785 prof_min =
real(minval(var),kind=
dp)
786 prof_sum =
real(sum (var),kind=
dp)
787 log_info(
"PROF_valcheck_DP_6D",
'(1x,A,A7,A,A16,3(A,ES24.16))') &
788 '+',prof_header,
'[',prof_item,
'] max=',prof_max,
',min=',prof_min,
',sum=',prof_sum
791 end subroutine prof_valcheck_dp_6d
subroutine, public prof_setup
subroutine, public prof_setprefx(prefxname)
logical, public io_log_suppress
suppress all of log output?
integer, public io_fid_conf
Config file ID.
integer, public prc_nprocs
myrank in local communicator
logical, public io_l
output log or not? (this process)
integer, parameter, public io_fid_stdout
logical, public io_log_allnode
output log for each node?
real(dp) function, public prc_mpitime()
Get MPI time.
subroutine, public prc_abort
Abort Process.
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
subroutine, public prc_mpibarrier
Barrier MPI.
logical, public prc_ismaster
master process in local communicator?
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
integer, public io_fid_log
Log file ID.
subroutine, public prc_mpitimestat(avgvar, maxvar, minvar, maxidx, minidx, var)
Calc global statistics for timer.
integer, parameter, public dp
subroutine, public prof_rapreport
Report raptime.