SCALE-RM
scale_statistics.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use mpi ! TODO: to replace functions in scale_comm module
17  use scale_precision
18  use scale_io
19  use scale_prof
20  use scale_prc, only: &
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
29  public :: statistics_setup
30  public :: statistics_total
31  public :: statistics_detail
32  public :: statistics_horizontal_mean
33  public :: statistics_horizontal_min
34  public :: statistics_horizontal_max
35 
36  interface statistics_total
37  module procedure statistics_total_2d
38  module procedure statistics_total_3d
39  end interface statistics_total
40 
41  interface statistics_detail
42  module procedure statistics_detail_2d
43  module procedure statistics_detail_3d
44  end interface statistics_detail
45 
46  interface statistics_horizontal_mean
47  module procedure statistics_horizontal_mean_2d
48  module procedure statistics_horizontal_mean_3d
49  end interface statistics_horizontal_mean
50 
51  interface statistics_horizontal_max
52  module procedure statistics_horizontal_max_2d
53  module procedure statistics_horizontal_max_3d
54  end interface statistics_horizontal_max
55 
56  interface statistics_horizontal_min
57  module procedure statistics_horizontal_min_2d
58  module procedure statistics_horizontal_min_3d
59  end interface statistics_horizontal_min
60  !-----------------------------------------------------------------------------
61  !
62  !++ Public parameters & variables
63  !
64  logical, public :: statistics_checktotal = .false.
65 
66  !-----------------------------------------------------------------------------
67  !
68  !++ Private procedure
69  !
70  !-----------------------------------------------------------------------------
71  !
72  !++ Private parameters & variables
73  !
74  logical, private :: statistics_use_globalcomm = .false.
75 
76  !-----------------------------------------------------------------------------
77 contains
78  !-----------------------------------------------------------------------------
80  subroutine statistics_setup
81  use scale_prc, only: &
82  prc_abort
83  use scale_comm_cartesc, only: &
85  implicit none
86 
87  namelist / param_statistics / &
89  statistics_use_globalcomm
90 
91  integer :: ierr
92  !---------------------------------------------------------------------------
93 
94  call comm_setup
95 
96  log_newline
97  log_info("STATISTICS_setup",*) 'Setup'
98 
99  !--- read namelist
100  rewind(io_fid_conf)
101  read(io_fid_conf,nml=param_statistics,iostat=ierr)
102  if( ierr < 0 ) then !--- missing
103  log_info("STATISTICS_setup",*) 'Not found namelist. Default used.'
104  elseif( ierr > 0 ) then !--- fatal error
105  log_error("STATISTICS_setup",*) 'Not appropriate names in namelist PARAM_STATISTICS. Check!'
106  call prc_abort
107  endif
108  log_nml(param_statistics)
109 
110  log_newline
111  log_info("STATISTICS_setup",*) 'Caluculate total statistics for monitoring? : ', statistics_checktotal
112  if ( statistics_use_globalcomm ) then
113  log_info_cont(*) '=> The total is calculated for the global domain.'
114  else
115  log_info_cont(*) '=> The total is calculated for the local domain.'
116  endif
117 
118  return
119  end subroutine statistics_setup
120 
121  !-----------------------------------------------------------------------------
123  subroutine statistics_total_2d( &
124  IA, IS, IE, JA, JS, JE, &
125  var, varname, &
126  area, total, &
127  log_suppress, &
128  mean, sum )
129  use scale_prc, only: &
130  prc_myrank, &
131  prc_abort
132  use scale_const, only: &
133  undef => const_undef
134  implicit none
135 
136  integer, intent(in) :: IA, IS, IE
137  integer, intent(in) :: JA, JS, JE
138 
139  real(RP), intent(in) :: var(ia,ja)
140  character(len=*), intent(in) :: varname
141  real(RP), intent(in) :: area(ia,ja)
142  real(RP), intent(in) :: total
143 
144  logical, intent(in), optional :: log_suppress
145  real(RP), intent(out), optional :: mean
146  real(DP), intent(out), optional :: sum
147 
148  real(DP) :: statval
149  real(DP) :: sendbuf(2), recvbuf(2)
150  real(DP) :: sum_, mean_
151 
152  logical :: suppress_
153  integer :: ierr
154  integer :: i, j
155  !---------------------------------------------------------------------------
156 
157  statval = 0.0_dp
158  if ( var(is,js) /= undef ) then
159  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2) reduction(+:statval)
160  do j = js, je
161  do i = is, ie
162  statval = statval + var(i,j) * area(i,j)
163  end do
164  end do
165  end if
166 
167  if ( .NOT. ( statval > -1.0_dp .OR. statval < 1.0_dp ) ) then ! must be NaN
168  log_error("STATISTICS_total_2D",*) 'NaN is detected for ', trim(varname), ' in rank ', prc_myrank
169  call prc_abort
170  endif
171 
172  if ( present(log_suppress) ) then
173  suppress_ = log_suppress
174  else
175  suppress_ = .false.
176  end if
177 
178  if ( statistics_use_globalcomm ) then
179  call prof_rapstart('COMM_Allreduce', 2)
180  sendbuf(1) = statval
181  sendbuf(2) = total
182  ! All reduce
183  call mpi_allreduce( sendbuf(:), recvbuf(:), &
184  2, &
185  mpi_double_precision, &
186  mpi_sum, &
187  prc_local_comm_world, &
188  ierr )
189  call prof_rapend ('COMM_Allreduce', 2)
190 
191  sum_ = recvbuf(1)
192  mean_ = recvbuf(1) / recvbuf(2)
193  ! statistics over the all node
194  if ( .not. suppress_ ) then ! if varname is empty, suppress output
195  log_info("STATISTICS_total_2D",'(1x,A,A24,A,ES24.17)') &
196  '[', trim(varname), '] MEAN(global) = ', mean_
197  endif
198  else
199  sum_ = statval
200  mean_ = statval / total
201 
202  ! statistics on each node
203  if ( .not. suppress_ ) then ! if varname is empty, suppress output
204  log_info("STATISTICS_total_2D",'(1x,A,A24,A,ES24.17)') &
205  '[', trim(varname), '] MEAN(local) = ', mean_
206  endif
207  endif
208 
209  if ( present(mean) ) mean = mean_
210  if ( present(sum ) ) sum = sum_
211 
212  return
213  end subroutine statistics_total_2d
214 
215  !-----------------------------------------------------------------------------
217  subroutine statistics_total_3d( &
218  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
219  var, varname, &
220  vol, total, &
221  log_suppress, &
222  mean, sum )
223  use scale_prc, only: &
224  prc_myrank, &
225  prc_abort
226  use scale_const, only: &
227  undef => const_undef
228  implicit none
229 
230  integer, intent(in) :: KA, KS, KE
231  integer, intent(in) :: IA, IS, IE
232  integer, intent(in) :: JA, JS, JE
233 
234  real(RP), intent(in) :: var(ka,ia,ja)
235  character(len=*), intent(in) :: varname
236  real(RP), intent(in) :: vol(ka,ia,ja)
237  real(RP), intent(in) :: total
238 
239  logical, intent(in), optional :: log_suppress
240  real(RP), intent(out), optional :: mean
241  real(DP), intent(out), optional :: sum
242 
243  real(DP) :: statval
244  real(DP) :: sendbuf(2), recvbuf(2)
245  real(DP) :: mean_, sum_
246 
247  logical :: suppress_
248  integer :: ierr
249  integer :: k, i, j
250  !---------------------------------------------------------------------------
251 
252  statval = 0.0_dp
253  if ( var(ks,is,js) /= undef ) then
254  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2) reduction(+:statval)
255  do j = js, je
256  do i = is, ie
257  do k = ks, ke
258  statval = statval + var(k,i,j) * vol(k,i,j)
259  enddo
260  enddo
261  enddo
262  end if
263 
264  if ( .NOT. ( statval > -1.0_dp .OR. statval < 1.0_dp ) ) then ! must be NaN
265  log_error("STATISTICS_total_3D",*) 'NaN is detected for ', trim(varname), ' in rank ', prc_myrank
266  call prc_abort
267  endif
268 
269  if ( present(log_suppress) ) then
270  suppress_ = log_suppress
271  else
272  suppress_ = .false.
273  end if
274 
275  if ( statistics_use_globalcomm ) then
276  call prof_rapstart('COMM_Allreduce', 2)
277  sendbuf(1) = statval
278  sendbuf(2) = total
279  ! All reduce
280  call mpi_allreduce( sendbuf(:), recvbuf(:), &
281  2, &
282  mpi_double_precision, &
283  mpi_sum, &
284  prc_local_comm_world, &
285  ierr )
286  call prof_rapend ('COMM_Allreduce', 2)
287 
288  sum_ = recvbuf(1)
289  mean_ = recvbuf(1) / recvbuf(2)
290  ! statistics over the all node
291  if ( .not. suppress_ ) then ! if varname is empty, suppress output
292  log_info("STATISTICS_total_3D",'(1x,A,A24,A,ES24.17)') &
293  '[', trim(varname), '] MEAN(global) = ', mean_
294  endif
295  else
296  sum_ = statval
297  mean_ = statval / total
298 
299  ! statistics on each node
300  if ( .not. suppress_ ) then ! if varname is empty, suppress output
301  log_info("STATISTICS_total_3D",'(1x,A,A24,A,ES24.17)') &
302  '[', trim(varname), '] MEAN(local) = ', mean_
303  endif
304  endif
305 
306  if ( present(mean) ) mean = mean_
307  if ( present(sum ) ) sum = sum_
308 
309  return
310  end subroutine statistics_total_3d
311 
312  !-----------------------------------------------------------------------------
314  subroutine statistics_horizontal_mean_2d( &
315  IA, IS, IE, JA, JS, JE, &
316  var, area, &
317  varmean )
318  use scale_const, only: &
319  undef => const_undef
320  integer, intent(in) :: IA, IS, IE
321  integer, intent(in) :: JA, JS, JE
322  real(RP), intent(in) :: var (ia,ja)
323  real(RP), intent(in) :: area(ia,ja)
324  real(RP), intent(out) :: varmean
325 
326  real(DP) :: statval (2)
327  real(DP) :: allstatval(2)
328 
329  integer :: ierr
330  integer :: i, j
331  !---------------------------------------------------------------------------
332 
333  statval(:) = 0.0_dp
334  do j = js, je
335  do i = is, ie
336  if ( var(i,j) /= undef ) then
337  statval(1) = statval(1) + area(i,j) * var(i,j)
338  statval(2) = statval(2) + area(i,j)
339  endif
340  enddo
341  enddo
342 
343  call prof_rapstart('COMM_Allreduce', 2)
344  ! All reduce
345  call mpi_allreduce( statval(:), &
346  allstatval(:), &
347  2, &
348  mpi_double_precision, &
349  mpi_sum, &
351  ierr )
352  call prof_rapend ('COMM_Allreduce', 2)
353 
354  if ( allstatval(2) > 0.0_dp ) then
355  varmean = allstatval(1) / allstatval(2)
356  else
357  varmean = undef
358  end if
359 
360  return
361  end subroutine statistics_horizontal_mean_2d
362 
363  subroutine statistics_horizontal_mean_3d( &
364  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
365  var, area, &
366  varmean )
367  use scale_const, only: &
368  undef => const_undef
369  integer, intent(in) :: KA, KS, KE
370  integer, intent(in) :: IA, IS, IE
371  integer, intent(in) :: JA, JS, JE
372 
373  real(RP), intent(in) :: var (ka,ia,ja)
374  real(RP), intent(in) :: area( ia,ja)
375  real(RP), intent(out) :: varmean(ka)
376 
377  real(DP) :: statval (2,ka)
378  real(DP) :: allstatval(2,ka)
379 
380  integer :: ierr
381  integer :: k, i, j
382  !---------------------------------------------------------------------------
383 
384  statval(:,:) = 0.0_dp
385  do j = js, je
386  do i = is, ie
387  do k = ks, ke
388  if ( var(k,i,j) /= undef ) then
389  statval(1,k) = statval(1,k) + area(i,j) * var(k,i,j)
390  statval(2,k) = statval(2,k) + area(i,j)
391  endif
392  enddo
393  enddo
394  enddo
395 
396  call prof_rapstart('COMM_Allreduce', 2)
397  ! All reduce
398  call mpi_allreduce( statval(:,ks:ke), &
399  allstatval(:,ks:ke), &
400  (ke-ks+1)*2, &
401  mpi_double_precision, &
402  mpi_sum, &
404  ierr )
405  call prof_rapend ('COMM_Allreduce', 2)
406 
407  do k = ks, ke
408  if ( allstatval(2,k) > 0.0_dp ) then
409  varmean(k) = allstatval(1,k) / allstatval(2,k)
410  else
411  varmean(k) = undef
412  end if
413  enddo
414  do k = 1, ks-1
415  varmean(k) = undef
416  end do
417  do k = ke+1, ka
418  varmean(k) = undef
419  end do
420 
421  return
422  end subroutine statistics_horizontal_mean_3d
423 
424  !-----------------------------------------------------------------------------
426  subroutine statistics_horizontal_min_2d( &
427  IA, IS, IE, JA, JS, JE, &
428  var, &
429  varmin )
430  use scale_const, only: &
431  undef => const_undef, &
432  huge => const_huge
433  use scale_comm_cartesc, only: &
435  integer, intent(in) :: IA, IS, IE
436  integer, intent(in) :: JA, JS, JE
437 
438  real(RP), intent(in) :: var(ia,ja)
439  real(RP), intent(out) :: varmin
440 
441  real(RP) :: statval
442  real(RP) :: allstatval
443 
444  integer :: ierr
445  integer :: i, j
446  !---------------------------------------------------------------------------
447 
448  statval = huge
449  do j = js, je
450  do i = is, ie
451  if ( var(i,j) /= undef .and. var(i,j) < statval ) then
452  statval = var(i,j)
453  endif
454  enddo
455  enddo
456 
457  call prof_rapstart('COMM_Allreduce', 2)
458  ! All reduce
459  call mpi_allreduce( statval, &
460  allstatval, &
461  1, &
462  comm_datatype, &
463  mpi_min, &
465  ierr )
466  call prof_rapend ('COMM_Allreduce', 2)
467 
468  if ( allstatval < huge ) then
469  varmin = allstatval
470  else
471  varmin = undef
472  end if
473 
474  return
475  end subroutine statistics_horizontal_min_2d
476 
477  subroutine statistics_horizontal_min_3d( &
478  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
479  var, &
480  varmin )
481  use scale_const, only: &
482  undef => const_undef, &
483  huge => const_huge
484  use scale_comm_cartesc, only: &
486  integer, intent(in) :: KA, KS, KE
487  integer, intent(in) :: IA, IS, IE
488  integer, intent(in) :: JA, JS, JE
489 
490  real(RP), intent(in) :: var(ka,ia,ja)
491  real(RP), intent(out) :: varmin(ka)
492 
493  real(RP) :: statval (ka)
494  real(RP) :: allstatval(ka)
495 
496  integer :: ierr
497  integer :: k, i, j
498  !---------------------------------------------------------------------------
499 
500  statval(:) = huge
501  do j = js, je
502  do i = is, ie
503  do k = ks, ke
504  if ( var(k,i,j) /= undef .and. var(k,i,j) < statval(k) ) then
505  statval(k) = var(k,i,j)
506  endif
507  enddo
508  enddo
509  enddo
510 
511  call prof_rapstart('COMM_Allreduce', 2)
512  ! All reduce
513  call mpi_allreduce( statval(ks:ke), &
514  allstatval(ks:ke), &
515  ke-ks+1, &
516  comm_datatype, &
517  mpi_min, &
519  ierr )
520  call prof_rapend ('COMM_Allreduce', 2)
521 
522  do k = ks, ke
523  if ( allstatval(k) < huge ) then
524  varmin(k) = allstatval(k)
525  else
526  varmin(k) = undef
527  end if
528  enddo
529  do k = 1, ks-1
530  varmin(k) = undef
531  end do
532  do k = ke+1, ka
533  varmin(k) = undef
534  end do
535 
536  return
537  end subroutine statistics_horizontal_min_3d
538 
539  !-----------------------------------------------------------------------------
541  subroutine statistics_horizontal_max_2d( &
542  IA, IS, IE, JA, JS, JE, &
543  var, &
544  varmax )
545  use scale_const, only: &
546  undef => const_undef, &
547  huge => const_huge
548  use scale_comm_cartesc, only: &
550  integer, intent(in) :: IA, IS, IE
551  integer, intent(in) :: JA, JS, JE
552 
553  real(RP), intent(in) :: var(ia,ja)
554  real(RP), intent(out) :: varmax
555 
556  real(RP) :: statval
557  real(RP) :: allstatval
558 
559  integer :: ierr
560  integer :: i, j
561  !---------------------------------------------------------------------------
562 
563  statval = - huge
564  do j = js, je
565  do i = is, ie
566  if ( var(i,j) /= undef .and. var(i,j) > statval ) then
567  statval = var(i,j)
568  endif
569  enddo
570  enddo
571 
572  call prof_rapstart('COMM_Allreduce', 2)
573  ! All reduce
574  call mpi_allreduce( statval, &
575  allstatval, &
576  1, &
577  comm_datatype, &
578  mpi_max, &
580  ierr )
581  call prof_rapend ('COMM_Allreduce', 2)
582 
583  if ( allstatval > - huge ) then
584  varmax = allstatval
585  else
586  varmax = undef
587  end if
588 
589  return
590  end subroutine statistics_horizontal_max_2d
591 
592  subroutine statistics_horizontal_max_3d( &
593  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
594  var, &
595  varmax )
596  use scale_const, only: &
597  undef => const_undef, &
598  huge => const_huge
599  use scale_comm_cartesc, only: &
601  integer, intent(in) :: KA, KS, KE
602  integer, intent(in) :: IA, IS, IE
603  integer, intent(in) :: JA, JS, JE
604 
605  real(RP), intent(in) :: var(ka,ia,ja)
606  real(RP), intent(out) :: varmax(ka)
607 
608  real(RP) :: statval (ka)
609  real(RP) :: allstatval(ka)
610 
611  integer :: ierr
612  integer :: k, i, j
613  !---------------------------------------------------------------------------
614 
615  statval(:) = - huge
616  do j = js, je
617  do i = is, ie
618  do k = ks, ke
619  if ( var(k,i,j) /= undef .and. var(k,i,j) > statval(k) ) then
620  statval(k) = var(k,i,j)
621  endif
622  enddo
623  enddo
624  enddo
625 
626  call prof_rapstart('COMM_Allreduce', 2)
627  ! All reduce
628  call mpi_allreduce( statval(ks:ke), &
629  allstatval(ks:ke), &
630  ke-ks+1, &
631  comm_datatype, &
632  mpi_max, &
634  ierr )
635  call prof_rapend ('COMM_Allreduce', 2)
636 
637  do k = ks, ke
638  if ( allstatval(k) > - huge ) then
639  varmax(k) = allstatval(k)
640  else
641  varmax(k) = undef
642  end if
643  enddo
644  do k = 1, ks-1
645  varmax(k) = undef
646  end do
647  do k = ke+1, ka
648  varmax(k) = undef
649  end do
650 
651  return
652  end subroutine statistics_horizontal_max_3d
653 
654  !-----------------------------------------------------------------------------
656  subroutine statistics_detail_3d( &
657  KA, KS, KE, IA, IS, IE, JA, JS, JE, VA, &
658  varname, var, &
659  local )
660  use scale_prc, only: &
661  prc_nprocs
662  use scale_comm_cartesc, only: &
664  implicit none
665 
666  integer, intent(in) :: KA, KS, KE
667  integer, intent(in) :: IA, IS, IE
668  integer, intent(in) :: JA, JS, JE
669  integer, intent(in) :: VA
670 
671  character(len=*), intent(in) :: varname(va)
672  real(RP), intent(in) :: var(ka,ia,ja,va)
673 
674  logical, intent(in), optional :: local
675 
676  real(RP) :: statval_l ( va,2)
677  integer :: statidx_l (3,va,2)
678  real(RP) :: statval ( va,2,0:prc_nprocs-1)
679  integer :: statidx (3,va,2,0:prc_nprocs-1)
680  real(RP) :: allstatval(va,2)
681  integer :: allstatidx(va,2)
682  logical :: do_globalcomm
683 
684  integer :: k, i, j
685  integer :: ierr
686  integer :: v, p
687  !---------------------------------------------------------------------------
688 
689  do_globalcomm = statistics_use_globalcomm
690  if ( present(local) ) do_globalcomm = ( .not. local )
691 
692  log_newline
693  log_info("STATISTICS_detail_3D",*) 'Variable Statistics '
694  do v = 1, va
695  statval_l( v,:) = var(ks,is,js,v)
696  statidx_l(1,v,:) = ks
697  statidx_l(2,v,:) = is
698  statidx_l(3,v,:) = js
699  do j = js, je
700  do i = is, ie
701  do k = ks, ke
702  if ( var(k,i,j,v) > statval_l(v,1) ) then
703  statval_l( v,1) = var(k,i,j,v)
704  statidx_l(1,v,1) = k
705  statidx_l(2,v,1) = i
706  statidx_l(3,v,1) = j
707  end if
708  if ( var(k,i,j,v) < statval_l(v,2) ) then
709  statval_l( v,2) = var(k,i,j,v)
710  statidx_l(1,v,2) = k
711  statidx_l(2,v,2) = i
712  statidx_l(3,v,2) = j
713  end if
714  end do
715  end do
716  end do
717  enddo
718 
719  if ( do_globalcomm ) then
720  call prof_rapstart('COMM_Bcast', 2)
721 
722  call mpi_allgather( statval_l(:,:), &
723  va*2, &
724  comm_datatype, &
725  statval(:,:,:), &
726  va*2, &
727  comm_datatype, &
728  prc_local_comm_world, &
729  ierr )
730 
731  call mpi_allgather( statidx_l(:,:,:), &
732  3*va*2, &
733  mpi_integer, &
734  statidx(:,:,:,:), &
735  3*va*2, &
736  mpi_integer, &
737  prc_local_comm_world, &
738  ierr )
739 
740  call prof_rapend ('COMM_Bcast', 2)
741 
742  do v = 1, va
743  allstatval(v,1) = statval(v,1,0)
744  allstatval(v,2) = statval(v,2,0)
745  allstatidx(v,:) = 0
746  do p = 1, prc_nprocs-1
747  if ( statval(v,1,p) > allstatval(v,1) ) then
748  allstatval(v,1) = statval(v,1,p)
749  allstatidx(v,1) = p
750  end if
751  if ( statval(v,2,p) < allstatval(v,2) ) then
752  allstatval(v,2) = statval(v,2,p)
753  allstatidx(v,2) = p
754  end if
755  end do
756  log_info_cont(*) '[', trim(varname(v)), ']'
757  log_info_cont('(1x,A,ES17.10,A,4(I5,A))') ' MAX =', &
758  allstatval(v,1), ' (rank=', &
759  allstatidx(v,1), '; ', &
760  statidx(1,v,1,allstatidx(v,1)),',', &
761  statidx(2,v,1,allstatidx(v,1)),',', &
762  statidx(3,v,1,allstatidx(v,1)),')'
763  log_info_cont('(1x,A,ES17.10,A,4(I5,A))') ' MIN =', &
764  allstatval(v,2), ' (rank=', &
765  allstatidx(v,2), '; ', &
766  statidx(1,v,2,allstatidx(v,2)),',', &
767  statidx(2,v,2,allstatidx(v,2)),',', &
768  statidx(3,v,2,allstatidx(v,2)),')'
769  enddo
770  else
771  ! statistics on each node
772  do v = 1, va
773  log_info_cont(*) '[', trim(varname(v)), ']'
774  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') 'MAX = ', &
775  statval_l( v,1),' (', &
776  statidx_l(1,v,1),',', &
777  statidx_l(2,v,1),',', &
778  statidx_l(3,v,1),')'
779  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') 'MIN = ', &
780  statval_l( v,2),' (', &
781  statidx_l(1,v,2),',', &
782  statidx_l(2,v,2),',', &
783  statidx_l(3,v,2),')'
784  enddo
785  endif
786 
787  log_newline
788 
789  return
790  end subroutine statistics_detail_3d
791 
792  subroutine statistics_detail_2d( &
793  IA, IS, IE, JA, JS, JE, VA, &
794  varname, var, &
795  local )
796  use scale_prc, only: &
797  prc_nprocs, &
798  prc_myrank
799  use scale_comm_cartesc, only: &
801  implicit none
802 
803  integer, intent(in) :: IA, IS, IE
804  integer, intent(in) :: JA, JS, JE
805  integer, intent(in) :: VA
806 
807  character(len=*), intent(in) :: varname(va)
808  real(RP), intent(in) :: var(ia,ja,va)
809 
810  logical, intent(in), optional :: local ! calc in local node
811 
812  real(RP) :: statval_l ( va,2)
813  integer :: statidx_l (2,va,2)
814  real(RP) :: statval ( va,2,0:prc_nprocs-1)
815  integer :: statidx (2,va,2,0:prc_nprocs-1)
816  real(RP) :: allstatval(va,2)
817  integer :: allstatidx(va,2)
818  logical :: do_globalcomm
819 
820  integer :: i, j
821  integer :: ierr
822  integer :: v, p
823  !---------------------------------------------------------------------------
824 
825  do_globalcomm = statistics_use_globalcomm
826  if ( present(local) ) do_globalcomm = ( .not. local )
827 
828  log_newline
829  log_info("STATISTICS_detail_2D",*) 'Variable Statistics '
830  do v = 1, va
831  statval_l( v,:) = var(is,js,v)
832  statidx_l(1,v,:) = is
833  statidx_l(2,v,:) = js
834  do j = js, je
835  do i = is, ie
836  if ( var(i,j,v) > statval_l(v,1) ) then
837  statval_l( v,1) = var(i,j,v)
838  statidx_l(1,v,1) = i
839  statidx_l(2,v,1) = j
840  end if
841  if ( var(i,j,v) < statval_l(v,2) ) then
842  statval_l( v,2) = var(i,j,v)
843  statidx_l(1,v,2) = i
844  statidx_l(2,v,2) = j
845  end if
846  end do
847  end do
848  enddo
849 
850  if ( do_globalcomm ) then
851  call prof_rapstart('COMM_Bcast', 2)
852 
853  call mpi_allgather( statval_l(:,:), &
854  va*2, &
855  comm_datatype, &
856  statval(:,:,:), &
857  va*2, &
858  comm_datatype, &
859  prc_local_comm_world, &
860  ierr )
861 
862  call mpi_allgather( statidx_l(:,:,:), &
863  2*va*2, &
864  mpi_integer, &
865  statidx(:,:,:,:), &
866  2*va*2, &
867  mpi_integer, &
868  prc_local_comm_world, &
869  ierr )
870 
871  call prof_rapend ('COMM_Bcast', 2)
872 
873  do v = 1, va
874  allstatval(v,1) = statval(v,1,0)
875  allstatval(v,2) = statval(v,2,0)
876  allstatidx(v,:) = 0
877  do p = 1, prc_nprocs-1
878  if ( statval(v,1,p) > allstatval(v,1) ) then
879  allstatval(v,1) = statval(v,1,p)
880  allstatidx(v,1) = p
881  end if
882  if ( statval(v,2,p) < allstatval(v,2) ) then
883  allstatval(v,2) = statval(v,2,p)
884  allstatidx(v,2) = p
885  end if
886  end do
887  log_info_cont(*) '[', trim(varname(v)), ']'
888  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') ' MAX =', &
889  allstatval(v,1), ' (rank=', &
890  allstatidx(v,1), '; ', &
891  statidx(1,v,1,allstatidx(v,1)),',', &
892  statidx(2,v,1,allstatidx(v,1)),')'
893  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') ' MIN =', &
894  allstatval(v,2), ' (rank=', &
895  allstatidx(v,2), '; ', &
896  statidx(1,v,2,allstatidx(v,2)),',', &
897  statidx(2,v,2,allstatidx(v,2)),')'
898  enddo
899  else
900  ! statistics on each node
901  do v = 1, va
902  log_info_cont(*) '[', trim(varname(v)), ']'
903  log_info_cont('(1x,A,ES17.10,A,2(I5,A))') 'MAX = ', &
904  statval_l( v,1),' (', &
905  statidx_l(1,v,1),',', &
906  statidx_l(2,v,1),')'
907  log_info_cont('(1x,A,ES17.10,A,2(I5,A))') 'MIN = ', &
908  statval_l( v,2),' (', &
909  statidx_l(1,v,2),',', &
910  statidx_l(2,v,2),')'
911  enddo
912  endif
913 
914  log_newline
915 
916  return
917  end subroutine statistics_detail_2d
918 
919 end module scale_statistics
real(rp), public const_huge
huge number
Definition: scale_const.F90:35
subroutine statistics_horizontal_mean_2d(IA, IS, IE, JA, JS, JE, var, area, varmean)
Calc horizontal mean value.
integer, public comm_datatype
datatype of variable
subroutine statistics_detail_2d(IA, IS, IE, JA, JS, JE, VA, varname, var, local)
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
logical, public statistics_checktotal
calc&report variable totals to logfile?
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
subroutine statistics_horizontal_mean_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, area, varmean)
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
subroutine statistics_total_2d(IA, IS, IE, JA, JS, JE, var, varname, area, total, log_suppress, mean, sum)
Calc domain sum and area-weighted mean.
module PROCESS
Definition: scale_prc.F90:11
subroutine, public statistics_setup
Setup.
subroutine statistics_horizontal_max_2d(IA, IS, IE, JA, JS, JE, var, varmax)
Calc horizontal maximum value.
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
subroutine statistics_horizontal_min_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, varmin)
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
subroutine statistics_horizontal_min_2d(IA, IS, IE, JA, JS, JE, var, varmin)
Calc horizontal minimum value.
module profiler
Definition: scale_prof.F90:11
subroutine statistics_total_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, varname, vol, total, log_suppress, mean, sum)
Calc domain sum and volume-weighted mean.
module PRECISION
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
module Statistics
subroutine, public comm_setup
Setup.
module STDIO
Definition: scale_io.F90:10
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
subroutine statistics_horizontal_max_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, varmax)
subroutine statistics_detail_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, VA, varname, var, local)
Search global maximum & minimum value.