source: LMDZ6/branches/Optimisation_LMDZ/libf/phy_common/profiling_physic.F90 @ 3709

Last change on this file since 3709 was 3706, checked in by adurocher, 4 years ago

Added timers for physiq and display physic profiling

File size: 9.2 KB
Line 
1! Conditional compilation : use scorep user task api
2#ifdef CPP_SCOREP
3
4#include "scorep/SCOREP_User.inc"
5#define SCOREP_DEFINE_HANDLE_ARRAY( handle ) SCOREP_USER_REGION_HANDLE, allocatable :: handle(:)
6#define SCOREP_ALLOCATE_HANDLE_ARRAY( handle, size ) allocate( scorep_handle(size) )
7
8#else
9
10#define SCOREP_DEFINE_HANDLE_ARRAY(handle)
11#define SCOREP_ALLOCATE_HANDLE_ARRAY(handle,size)
12#define SCOREP_USER_REGION_INIT(handle,name,region_type)
13#define SCOREP_USER_REGION_ENTER(handle)
14#define SCOREP_USER_REGION_END(handle)
15
16#endif
17
18MODULE profiling_physic_mod
19  use mod_phys_lmdz_omp_data, only : omp_size, omp_rank
20  use dict_int_mod
21  IMPLICIT NONE
22  SAVE
23  PRIVATE
24 
25  ! Shared variables
26  CHARACTER(*), parameter :: modname = "profiling_physic_mod"
27  LOGICAL, PARAMETER :: check_accuracy = .true.
28  INTEGER, PARAMETER :: max_id=20, max_depth=10
29  INTEGER :: nb_id ! Number of timers
30  CHARACTER(15), DIMENSION(max_id) :: name ! Timer names
31  type(dict) :: ids
32 
33  ! Threadlocal variables
34  INTEGER, ALLOCATABLE :: depth(:)!(omp_size) ! Number of current nested timers
35  INTEGER, ALLOCATABLE :: current_id(:,:)!(omp_size, max_depth) ! Current nested timer ids
36  REAL, ALLOCATABLE :: chrono(:,:)!(omp_size, max_depth) ! Last timer start time
37  REAL, ALLOCATABLE :: elapsed(:,:)!(omp_size, max_id) ! Cumulative time for each timer
38  SCOREP_DEFINE_HANDLE_ARRAY(scorep_handle)
39 
40  PUBLIC :: enter_profile, exit_profile, print_profile
41
42CONTAINS
43  ! Be sure to init mod_phys_lmdz_omp_data before profiling
44  SUBROUTINE init_profiling
45    nb_id=0
46    allocate(depth(0:omp_size-1))
47    allocate(current_id(0:omp_size-1, max_depth))
48    allocate(chrono(0:omp_size-1, max_depth))
49    allocate(elapsed(0:omp_size-1, max_id))
50    depth(:)=0
51    SCOREP_ALLOCATE_HANDLE_ARRAY(scorep_handle, max_id)
52    call reset_profiling()
53  END SUBROUTINE init_profiling
54 
55  SUBROUTINE reset_profiling
56    if( any( depth /= 0 ) ) call abort_physic(modname, "Impossible to reset profiling : a profiling region is still open",1)   
57    elapsed(:,:)=0   
58  END SUBROUTINE
59
60  SUBROUTINE register_id(thename, id)
61    CHARACTER(*), INTENT(IN) :: thename
62    INTEGER, INTENT(OUT) :: id
63    nb_id = nb_id+1
64    if(nb_id > max_id) call abort_physic(modname, "Too many timers",1)
65    id = nb_id
66    name(id)=thename
67    SCOREP_USER_REGION_INIT( scorep_handle(id), thename, SCOREP_USER_REGION_TYPE_COMMON )
68  END SUBROUTINE register_id
69
70  FUNCTION get_elapsed(start)
71    INTEGER(kind=8) :: count, count_rate
72    REAL :: start,get_elapsed
73    CALL SYSTEM_CLOCK(count,count_rate)
74    if(check_accuracy .and. count < 10) print *, "Warning : Profiling elapsed time too short to be accurately measured"
75    get_elapsed = (1.*count)/(1.*count_rate) - start
76    IF(get_elapsed<0.) get_elapsed=0.
77  END FUNCTION get_elapsed
78
79  ! Enter profiling region with name 'name'
80  ! If 'name' doesn't exist, new id is created
81  ! Profiling region associated to name must be inactive (dosen't exist yet or 'exit_profile'd)
82  ! First call for 'name' must be called for every thread in the parallel region
83  SUBROUTINE enter_profile(name)
84    CHARACTER(*) :: name
85    integer :: id
86    logical, save :: first = .true.
87
88    ! Initialize profiling at first call
89    if(first) then
90      call checked_omp_barrier("Timers init")
91      !$omp master
92      first = .false.
93      call init_profiling
94      !$omp end master
95      !$omp barrier
96    end if   
97    ! Create prof_id at first call with 'name'
98    if( .not. exists(ids, name) ) then
99       call checked_omp_barrier(name)
100       !$omp master
101       call register_id(name,id)
102       call insert_or_assign(ids, name, id)
103       !$omp end master
104       !$omp barrier
105    end if
106    call enter_profile_aux(get_val(ids, name))
107  END SUBROUTINE
108 
109  SUBROUTINE enter_profile_aux(id)
110    INTEGER, INTENT(IN) :: id
111    depth(omp_rank) = depth(omp_rank)+1
112    chrono(omp_rank, depth(omp_rank)) = get_elapsed(0.)
113    current_id(omp_rank,depth(omp_rank)) = id
114    SCOREP_USER_REGION_ENTER(scorep_handle(id))
115  END SUBROUTINE
116 
117  ! Exit profiling region with name 'name'
118  ! Region must exist and be active on current thread
119  SUBROUTINE exit_profile(name)
120    CHARACTER(*) :: name   
121    if(.not. exists(ids, name)) call abort_physic(modname, "Exit unknown timer")
122    call exit_profile_aux(get_val(ids,name))
123  END SUBROUTINE
124   
125  SUBROUTINE exit_profile_aux(id)
126    INTEGER, INTENT(IN) :: id
127    INTEGER :: parent_id
128    REAL :: my_elapsed
129   
130    IF(depth(omp_rank)<=0) call abort_physic(modname, "exit_profile called without a matching enter_profile (depth=0)",1)
131    IF(id /= current_id(omp_rank, depth(omp_rank))) call abort_physic(modname, "wrong timer id : exit_profile id doesn't match enter_profile",1)
132   
133    SCOREP_USER_REGION_END(scorep_handle(id))
134   
135    my_elapsed = get_elapsed(chrono(omp_rank, depth(omp_rank)))
136    ! add elapsed to current profile
137    elapsed(omp_rank,id) = elapsed(omp_rank,id) + my_elapsed
138    depth(omp_rank) = depth(omp_rank)-1
139    ! and substract from parent profile
140    IF(depth(omp_rank)>0) THEN
141       parent_id = current_id(omp_rank,depth(omp_rank))
142       elapsed(omp_rank,parent_id) = elapsed(omp_rank,parent_id) - my_elapsed
143    END IF
144  END SUBROUTINE
145
146  ! Collective operation : must be called on every process/thread
147  SUBROUTINE print_profile
148    !use mpi_mod
149    !use mpipara
150    use mod_phys_lmdz_mpi_data, only : is_mpi_root, mpi_size
151    use mod_phys_lmdz_mpi_transfert , only : reduce_sum_mpi, reduce_min_mpi, reduce_max_mpi
152    INTEGER :: i
153    REAL :: mean_total, min_total, max_total, mean_total_local, min_total_local, max_total_local
154    REAL :: omp_mean_time_local, omp_min_time_local, omp_max_time_local
155    REAL :: mean_time, max_process, min_process, percent
156
157    call enter_profile("print_profile")
158
159    !$OMP MASTER
160    ! mean_total : mean of sum of all timers on all threads and all procs
161    mean_total_local = SUM(elapsed(:,1:nb_id))/omp_size/mpi_size
162    !call MPI_Reduce(mean_total_local, mean_total, 1, MPI_DOUBLE, MPI_SUM, 0, comm_icosa, ierr)
163    call reduce_sum_mpi( mean_total_local, mean_total )
164    ! min_total : min on all procs of sum of all timers on all threads
165    min_total_local = MINVAL( SUM(elapsed(:,1:nb_id), 2) )
166    !call MPI_Reduce(min_total_local, min_total, 1, MPI_DOUBLE, MPI_MIN, 0, comm_icosa, ierr)
167    call reduce_min_mpi(min_total_local, min_total)
168    ! max_total : max on all procs of sum of all timers on all threads
169    max_total_local = MAXVAL( SUM(elapsed(:,1:nb_id), 2) )
170    !call MPI_Reduce(max_total_local, max_total, 1, MPI_DOUBLE, MPI_MAX, 0, comm_icosa, ierr)
171    call reduce_max_mpi(max_total_local, max_total)
172
173    if(is_mpi_root) PRINT *, '---------------------- Profiling (Physic) --------------'
174    if(is_mpi_root) PRINT ('(A15, F12.3, A, F12.3, A, F12.3, A)'), 'Total (Physic) (s) : ', mean_total, "  [", min_total, ",", max_total, "]"
175    if(is_mpi_root) PRINT ('(A15, A7, A47, A47)'), " ","  %  ","  Process mean(s)[ min, max, diff/mean(%) ]  ", "Process 1 threads(s)[ min, max, diff/mean(%) ]"
176    DO i=1,nb_id
177       omp_mean_time_local = SUM(elapsed(:,i))/omp_size ! Mean time of timer on local threads
178       omp_min_time_local = MINVAL(elapsed(:,i)) ! Max time of timer on local threads
179       omp_max_time_local = MAXVAL(elapsed(:,i)) ! Min time of timer on local threads
180       !omp_inbalance_local = 100*(max_time-min_time)/min_time
181       
182       ! mean_time : mean of timer on all threads and all procs
183       !call MPI_Reduce(omp_mean_time_local, mean_time, 1, MPI_DOUBLE, MPI_SUM, 0, comm_icosa, ierr);
184       call reduce_sum_mpi( omp_mean_time_local, mean_time )
185       mean_time = mean_time/mpi_size
186       ! max_process_time : max on procs of mean omp time
187       !call MPI_Reduce(omp_mean_time_local, max_process, 1, MPI_DOUBLE, MPI_MAX, 0, comm_icosa, ierr)
188       call reduce_max_mpi(omp_mean_time_local, max_process)
189       ! min_process_time : min on procs of mean omp time
190       !call MPI_Reduce(omp_mean_time_local, min_process, 1, MPI_DOUBLE, MPI_MIN, 0, comm_icosa, ierr)
191       call reduce_min_mpi(omp_mean_time_local, min_process)
192
193       percent = 100*mean_time/mean_total
194
195       if(is_mpi_root) PRINT ('(A15, F5.1, F12.3, A, F12.3, A, F12.3, A, F6.1, A, F12.3, A, F12.3, A, F12.3, A, F6.1, A)'), name(i), percent, mean_time, " [", min_process, ",", max_process, "," , inbalance(mean_time, min_process, max_process),  "]", omp_mean_time_local, " [", omp_min_time_local, ",", omp_max_time_local, "," , inbalance(omp_mean_time_local, omp_min_time_local, omp_max_time_local) , "]"
196    END DO
197    if(is_mpi_root) PRINT *, '---------------------- Profiling -----------------------'
198    !$OMP END MASTER
199    call exit_profile("print_profile")
200  contains
201    function inbalance( mean, min, max )
202      real :: mean, min, max
203      real :: inbalance
204      inbalance = 100*(max-min)/(mean+1E-10)
205    end function
206  END SUBROUTINE print_profile
207
208  subroutine checked_omp_barrier(message)
209    use omp_lib, only : omp_get_team_size
210    character(*), intent(in) :: message
211    integer, save :: count
212    count = 0
213    !$omp barrier
214    !$omp atomic update
215    count = count + 1
216    !$omp barrier
217    if( count /= omp_get_team_size(1) ) then
218       print *, "Error in barrier in physic timers : ", trim(message)
219       print *, count,"/",omp_get_team_size(1), "threads participated"
220       call abort_physic(modname, "Timer barrier",1)
221     end if
222   end subroutine checked_omp_barrier
223   
224END MODULE profiling_physic_mod
Note: See TracBrowser for help on using the repository browser.