SCALE-RM
scale_rm_statistics.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 #include "inc_openmp.h"
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use mpi
20  use scale_precision
21  use scale_stdio
22  use scale_prof
24  use scale_process, only: &
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: stat_setup
34  public :: stat_total
35  public :: stat_detail
36 
37  interface stat_total
38  module procedure stat_total_2d
39  module procedure stat_total_3d
40  end interface stat_total
41 
42  !-----------------------------------------------------------------------------
43  !
44  !++ Public parameters & variables
45  !
46  logical, public :: statistics_checktotal = .false.
47 
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private procedure
51  !
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private parameters & variables
55  !
56  logical, private :: statistics_use_globalcomm = .false.
57 
58  !-----------------------------------------------------------------------------
59 contains
60  !-----------------------------------------------------------------------------
62  subroutine stat_setup
63  use scale_process, only: &
65  implicit none
66 
67  namelist / param_statistics / &
69  statistics_use_globalcomm
70 
71  integer :: ierr
72  !---------------------------------------------------------------------------
73 
74  if( io_l ) write(io_fid_log,*)
75  if( io_l ) write(io_fid_log,*) '++++++ Module[STATISTICS] / Categ[ATMOS-RM COMM] / Origin[SCALElib]'
76 
77  !--- read namelist
78  rewind(io_fid_conf)
79  read(io_fid_conf,nml=param_statistics,iostat=ierr)
80  if( ierr < 0 ) then !--- missing
81  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
82  elseif( ierr > 0 ) then !--- fatal error
83  write(*,*) 'xxx Not appropriate names in namelist PARAM_STATISTICS. Check!'
84  call prc_mpistop
85  endif
86  if( io_lnml ) write(io_fid_log,nml=param_statistics)
87 
88  if( io_l ) write(io_fid_log,*)
89  if( io_l ) write(io_fid_log,*) '*** Caluculate statistics? : ', statistics_checktotal
90  if( io_l ) write(io_fid_log,*) '*** Allow global communication for statistics? : ', statistics_use_globalcomm
91  if ( statistics_use_globalcomm ) then
92  if( io_l ) write(io_fid_log,*) '*** Global total is calculated using MPI_ALLreduce.'
93  else
94  if( io_l ) write(io_fid_log,*) '*** Local total is calculated in each process.'
95  endif
96 
97  return
98  end subroutine stat_setup
99 
100  !-----------------------------------------------------------------------------
102  subroutine stat_total_2d( allstatval, var, varname )
103  use scale_process, only: &
104  prc_myrank, &
106  use scale_comm, only: &
108  use scale_grid_real, only: &
109  area => real_area
110  implicit none
111 
112  real(RP), intent(out) :: allstatval
113  real(RP), intent(in) :: var(ia,ja)
114  character(len=*), intent(in) :: varname
115 
116  character(len=H_SHORT) :: varname_trim
117  real(RP) :: statval
118 
119  integer :: ierr
120  integer :: i, j
121  !---------------------------------------------------------------------------
122 
123  varname_trim = trim(varname)
124 
125  statval = 0.0_rp
126  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2) reduction(+:statval)
127  do j = js, je
128  do i = is, ie
129  statval = statval + var(i,j) * area(i,j)
130  enddo
131  enddo
132 
133  if ( .NOT. ( statval > -1.0_rp .OR. statval < 1.0_rp ) ) then ! must be NaN
134  write(*,*) 'xxx [STAT_total] NaN is detected for ', varname_trim, ' in rank ', prc_myrank
135  call prc_mpistop
136  endif
137 
138  if ( statistics_use_globalcomm ) then
139  call prof_rapstart('COMM_Allreduce', 2)
140  ! All reduce
141  call mpi_allreduce( statval, &
142  allstatval, &
143  1, &
144  comm_datatype, &
145  mpi_sum, &
146  prc_local_comm_world, &
147  ierr )
148 
149  call prof_rapend ('COMM_Allreduce', 2)
150 
151  ! statistics over the all node
152  if ( varname_trim /= "" ) then ! if varname is empty, suppress output
153  if( io_l ) write(io_fid_log,'(1x,A,A,A,ES24.17)') &
154  '[', varname_trim, '] SUM(global) =', allstatval
155  endif
156  else
157  allstatval = statval
158 
159  ! statistics on each node
160  if ( varname_trim /= "" ) then ! if varname is empty, suppress output
161  if( io_l ) write(io_fid_log,'(1x,A,A,A,ES24.17)') &
162  '[', varname_trim, '] SUM(local) =', statval
163  endif
164  endif
165 
166  return
167  end subroutine stat_total_2d
168 
169  !-----------------------------------------------------------------------------
171  subroutine stat_total_3d( allstatval, var, varname )
172  use scale_process, only: &
173  prc_myrank, &
175  use scale_comm, only: &
177  use scale_grid_real, only: &
178  vol => real_vol
179  implicit none
180 
181  real(RP), intent(out) :: allstatval
182  real(RP), intent(in) :: var(ka,ia,ja)
183  character(len=*), intent(in) :: varname
184 
185  character(len=H_SHORT) :: varname_trim
186  real(RP) :: statval
187 
188  integer :: ierr
189  integer :: k, i, j
190  !---------------------------------------------------------------------------
191 
192  varname_trim = trim(varname)
193 
194  statval = 0.0_rp
195  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2) reduction(+:statval)
196  do j = js, je
197  do i = is, ie
198  do k = ks, ke
199  statval = statval + var(k,i,j) * vol(k,i,j)
200  enddo
201  enddo
202  enddo
203 
204  if ( .NOT. ( statval > -1.0_rp .OR. statval < 1.0_rp ) ) then ! must be NaN
205  write(*,*) 'xxx [STAT_total] NaN is detected for ', varname_trim, ' in rank ', prc_myrank
206  call prc_mpistop
207  endif
208 
209  if ( statistics_use_globalcomm ) then
210  call prof_rapstart('COMM_Allreduce', 2)
211  ! All reduce
212  call mpi_allreduce( statval, &
213  allstatval, &
214  1, &
215  comm_datatype, &
216  mpi_sum, &
217  prc_local_comm_world, &
218  ierr )
219 
220  call prof_rapend ('COMM_Allreduce', 2)
221 
222  ! statistics over the all node
223  if ( varname_trim /= "" ) then ! if varname is empty, suppress output
224  if( io_l ) write(io_fid_log,'(1x,A,A,A,ES24.17)') &
225  '[', varname_trim, '] SUM(global) =', allstatval
226  endif
227  else
228  allstatval = statval
229 
230  ! statistics on each node
231  if ( varname_trim /= "" ) then ! if varname is empty, suppress output
232  if( io_l ) write(io_fid_log,'(1x,A,A,A,ES24.17)') &
233  '[', varname_trim, '] SUM(local) =', statval
234  endif
235  endif
236 
237  return
238  end subroutine stat_total_3d
239 
240  !-----------------------------------------------------------------------------
242  subroutine stat_detail(var, varname, supress_globalcomm)
243  use scale_process, only: &
244  prc_nprocs, &
245  prc_myrank
246  use scale_const, only: &
247  const_undef8, &
249  use scale_comm, only: &
251  implicit none
252 
253  real(RP), intent(inout) :: var(:,:,:,:)
254  character(len=*), intent(in) :: varname(:)
255  logical, intent(in), optional :: supress_globalcomm
256 
257  logical , allocatable :: halomask (:,:,:)
258  real(RP), allocatable :: statval (:,:,:)
259  integer, allocatable :: statidx (:,:,:,:)
260  real(RP), allocatable :: allstatval(:,:)
261  integer, allocatable :: allstatidx(:,:,:)
262  integer :: ksize, vsize
263  logical :: do_globalcomm
264 
265  integer :: ierr
266  integer :: v, p
267  !---------------------------------------------------------------------------
268 
269  do_globalcomm = statistics_use_globalcomm
270  if ( present(supress_globalcomm) ) then
271  if ( supress_globalcomm ) then
272  do_globalcomm = .false.
273  endif
274  endif
275 
276  ksize = size(var(:,:,:,:),1)
277  vsize = size(var(:,:,:,:),4)
278 
279  allocate( halomask(ksize,ia,ja) ); halomask(:,:,:) = .false.
280 
281  if ( ksize == ka ) then
282  halomask(ks:ke,is:ie,js:je) = .true.
283  else
284  halomask(:,is:ie,js:je) = .true.
285  endif
286 
287  allocate( statval( vsize,2,0:prc_nprocs-1) ); statval(:,:,:) = const_undef8
288  allocate( statidx(3,vsize,2,0:prc_nprocs-1) ); statidx(:,:,:,:) = const_undef2
289 
290  allocate( allstatval( vsize,2) ); allstatval(:,:) = const_undef8
291  allocate( allstatidx(1,vsize,2) ); allstatidx(:,:,:) = const_undef2
292 
293  if( io_l ) write(io_fid_log,*)
294  if( io_l ) write(io_fid_log,*) '*** Variable Statistics ***'
295  do v = 1, vsize
296  statval( v,1,prc_myrank) = maxval(var(:,:,:,v),mask=halomask)
297  statval( v,2,prc_myrank) = minval(var(:,:,:,v),mask=halomask)
298  statidx(:,v,1,prc_myrank) = maxloc(var(:,:,:,v),mask=halomask)
299  statidx(:,v,2,prc_myrank) = minloc(var(:,:,:,v),mask=halomask)
300  enddo
301 
302  if ( do_globalcomm ) then
303  call prof_rapstart('COMM_Bcast', 2)
304  do p = 0, prc_nprocs-1
305 
306  call mpi_bcast( statval(1,1,p), &
307  vsize*2, &
308  comm_datatype, &
309  p, &
310  prc_local_comm_world, &
311  ierr )
312 
313  call mpi_bcast( statidx(1,1,1,p), &
314  3*vsize*2, &
315  mpi_integer, &
316  p, &
317  prc_local_comm_world, &
318  ierr )
319 
320  enddo
321  call prof_rapend ('COMM_Bcast', 2)
322 
323  do v = 1, vsize
324  allstatval(v,1) = maxval(statval(v,1,:))
325  allstatval(v,2) = minval(statval(v,2,:))
326  allstatidx(:,v,1) = maxloc(statval(v,1,:))-1
327  allstatidx(:,v,2) = minloc(statval(v,2,:))-1
328  if( io_l ) write(io_fid_log,*) '[', trim(varname(v)), ']'
329  if( io_l ) write(io_fid_log,'(1x,A,ES17.10,A,4(I5,A))') ' MAX =', &
330  allstatval( v,1), '(', &
331  allstatidx(1,v,1), ',', &
332  statidx(1,v,1,allstatidx(1,v,1)),',', &
333  statidx(2,v,1,allstatidx(1,v,1)),',', &
334  statidx(3,v,1,allstatidx(1,v,1)),')'
335  if( io_l ) write(io_fid_log,'(1x,A,ES17.10,A,4(I5,A))') ' MIN =', &
336  allstatval( v,2), '(', &
337  allstatidx(1,v,2), ',', &
338  statidx(1,v,2,allstatidx(1,v,2)),',', &
339  statidx(2,v,2,allstatidx(1,v,2)),',', &
340  statidx(3,v,2,allstatidx(1,v,2)),')'
341  enddo
342  else
343  ! statistics on each node
344  do v = 1, vsize
345  if( io_l ) write(io_fid_log,*) '*** [', trim(varname(v)), ']'
346  if( io_l ) write(io_fid_log,'(1x,A,ES17.10,A,3(I5,A))') '*** MAX = ', &
347  statval( v,1,prc_myrank),'(', &
348  statidx(1,v,1,prc_myrank),',', &
349  statidx(2,v,1,prc_myrank),',', &
350  statidx(3,v,1,prc_myrank),')'
351  if( io_l ) write(io_fid_log,'(1x,A,ES17.10,A,3(I5,A))') '*** MIN = ', &
352  statval( v,2,prc_myrank),'(', &
353  statidx(1,v,2,prc_myrank),',', &
354  statidx(2,v,2,prc_myrank),',', &
355  statidx(3,v,2,prc_myrank),')'
356  enddo
357  endif
358 
359  if( io_l ) write(io_fid_log,*)
360 
361  deallocate( halomask )
362 
363  deallocate( statval )
364  deallocate( statidx )
365 
366  deallocate( allstatval )
367  deallocate( allstatidx )
368 
369  return
370  end subroutine stat_detail
371 
372 end module scale_rm_statistics
373 !-------------------------------------------------------------------------------
integer, public is
start point of inner domain: x, local
logical, public statistics_checktotal
calc&report variable totals to logfile?
integer, public comm_datatype
datatype of variable
Definition: scale_comm.F90:117
integer, public je
end point of inner domain: y, local
integer, public prc_local_comm_world
local communicator
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
module Statistics
module grid index
integer, public ia
of x whole cells (local, with HALO)
module GRID (real space)
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:,:), allocatable, public real_area
horizontal area [m2]
module COMMUNICATION
Definition: scale_comm.F90:23
subroutine, public stat_setup
Setup.
integer, public js
start point of inner domain: y, local
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module PROCESS
subroutine stat_total_2d(allstatval, var, varname)
Calc volume/area-weighted global sum.
real(rp), dimension(:,:,:), allocatable, public real_vol
control volume [m3]
subroutine, public stat_detail(var, varname, supress_globalcomm)
Search global maximum & minimum value.
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public prc_myrank
process num in local communicator
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:42
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
module PRECISION
subroutine stat_total_3d(allstatval, var, varname)
Calc volume/area-weighted global sum.
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
integer, public prc_nprocs
myrank in local communicator
integer, public ja
of y whole cells (local, with HALO)