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  global, &
129  mean, sum )
130  use scale_prc, only: &
131  prc_myrank, &
132  prc_abort
133  use scale_const, only: &
134  eps => const_eps, &
135  undef => const_undef
136  implicit none
137 
138  integer, intent(in) :: IA, IS, IE
139  integer, intent(in) :: JA, JS, JE
140 
141  real(RP), intent(in) :: var(IA,JA)
142  character(len=*), intent(in) :: varname
143  real(RP), intent(in) :: area(IA,JA)
144  real(RP), intent(in) :: total
145 
146  logical, intent(in), optional :: log_suppress
147  logical, intent(in), optional :: global
148  real(RP), intent(out), optional :: mean
149  real(DP), intent(out), optional :: sum
150 
151  real(DP) :: statval
152  real(DP) :: sendbuf(2), recvbuf(2)
153  real(DP) :: sum_, mean_
154 
155  logical :: suppress_, global_
156  integer :: ierr
157  integer :: i, j
158  !---------------------------------------------------------------------------
159 
160  statval = 0.0_dp
161  if ( var(is,js) /= undef ) then
162  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2) reduction(+:statval)
163  do j = js, je
164  do i = is, ie
165  statval = statval + var(i,j) * area(i,j)
166  end do
167  end do
168  end if
169 
170  if ( .NOT. ( statval > -1.0_dp .OR. statval < 1.0_dp ) ) then ! must be NaN
171  log_error("STATISTICS_total_2D",*) 'NaN is detected for ', trim(varname), ' in rank ', prc_myrank
172  call prc_abort
173  endif
174 
175  if ( present(log_suppress) ) then
176  suppress_ = log_suppress
177  else
178  suppress_ = .false.
179  end if
180 
181  if ( present(global) ) then
182  global_ = global
183  else
184  global_ = statistics_use_globalcomm
185  end if
186 
187  if ( global_ ) then
188  call prof_rapstart('COMM_Allreduce', 2)
189  sendbuf(1) = statval
190  sendbuf(2) = total
191  ! All reduce
192  call mpi_allreduce( sendbuf(:), recvbuf(:), &
193  2, &
194  mpi_double_precision, &
195  mpi_sum, &
196  prc_local_comm_world, &
197  ierr )
198  call prof_rapend ('COMM_Allreduce', 2)
199 
200  if ( recvbuf(2) < eps ) then
201  sum_ = undef
202  mean_ = undef
203  else
204  sum_ = recvbuf(1)
205  mean_ = recvbuf(1) / recvbuf(2)
206  end if
207  ! statistics over the all node
208  if ( .not. suppress_ ) then ! if varname is empty, suppress output
209  log_info("STATISTICS_total_2D",'(1x,A,A24,A,ES24.17)') &
210  '[', trim(varname), '] MEAN(global) = ', mean_
211  endif
212  else
213  if ( total < eps ) then
214  sum_ = undef
215  mean_ = undef
216  else
217  sum_ = statval
218  mean_ = statval / total
219  end if
220 
221  ! statistics on each node
222  if ( .not. suppress_ ) then ! if varname is empty, suppress output
223  log_info("STATISTICS_total_2D",'(1x,A,A24,A,ES24.17)') &
224  '[', trim(varname), '] MEAN(local) = ', mean_
225  endif
226  endif
227 
228  if ( present(mean) ) mean = mean_
229  if ( present(sum ) ) sum = sum_
230 
231  return
232  end subroutine statistics_total_2d
233 
234  !-----------------------------------------------------------------------------
236  subroutine statistics_total_3d( &
237  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
238  var, varname, &
239  vol, total, &
240  log_suppress, &
241  global, &
242  mean, sum )
243  use scale_prc, only: &
244  prc_myrank, &
245  prc_abort
246  use scale_const, only: &
247  eps => const_eps , &
248  undef => const_undef
249  implicit none
250 
251  integer, intent(in) :: KA, KS, KE
252  integer, intent(in) :: IA, IS, IE
253  integer, intent(in) :: JA, JS, JE
254 
255  real(RP), intent(in) :: var(KA,IA,JA)
256  character(len=*), intent(in) :: varname
257  real(RP), intent(in) :: vol(KA,IA,JA)
258  real(RP), intent(in) :: total
259 
260  logical, intent(in), optional :: log_suppress
261  logical, intent(in), optional :: global
262  real(RP), intent(out), optional :: mean
263  real(DP), intent(out), optional :: sum
264 
265  real(DP) :: statval
266  real(DP) :: sendbuf(2), recvbuf(2)
267  real(DP) :: mean_, sum_
268 
269  logical :: suppress_, global_
270 
271  real(DP) :: work
272  integer :: ierr
273  integer :: k, i, j
274  !---------------------------------------------------------------------------
275 
276  statval = 0.0_dp
277  if ( var(ks,is,js) /= undef ) then
278  !$omp parallel do OMP_SCHEDULE_ reduction(+:statval) &
279  !$omp private(work)
280  do j = je, js, -1
281  do i = ie, is, -1
282  work = 0.0_rp
283  do k = ke, ks, -1
284  work = work + var(k,i,j) * vol(k,i,j)
285  enddo
286  statval = statval + work
287  enddo
288  enddo
289  end if
290 
291  if ( .NOT. ( statval > -1.0_dp .OR. statval < 1.0_dp ) ) then ! must be NaN
292  log_error("STATISTICS_total_3D",*) 'NaN is detected for ', trim(varname), ' in rank ', prc_myrank
293  call prc_abort
294  endif
295 
296  if ( present(log_suppress) ) then
297  suppress_ = log_suppress
298  else
299  suppress_ = .false.
300  end if
301 
302  if ( present(global) ) then
303  global_ = global
304  else
305  global_ = statistics_use_globalcomm
306  end if
307 
308  if ( global_ ) then
309  call prof_rapstart('COMM_Allreduce', 2)
310  sendbuf(1) = statval
311  sendbuf(2) = total
312  ! All reduce
313  call mpi_allreduce( sendbuf(:), recvbuf(:), &
314  2, &
315  mpi_double_precision, &
316  mpi_sum, &
317  prc_local_comm_world, &
318  ierr )
319  call prof_rapend ('COMM_Allreduce', 2)
320 
321  if ( recvbuf(2) < eps ) then
322  sum_ = undef
323  mean_ = undef
324  else
325  sum_ = recvbuf(1)
326  mean_ = recvbuf(1) / recvbuf(2)
327  end if
328  ! statistics over the all node
329  if ( .not. suppress_ ) then ! if varname is empty, suppress output
330  log_info("STATISTICS_total_3D",'(1x,A,A24,A,ES24.17)') &
331  '[', trim(varname), '] MEAN(global) = ', mean_
332  endif
333  else
334  if ( total < eps ) then
335  sum_ = undef
336  mean_ = undef
337  else
338  sum_ = statval
339  mean_ = statval / total
340  end if
341 
342  ! statistics on each node
343  if ( .not. suppress_ ) then ! if varname is empty, suppress output
344  log_info("STATISTICS_total_3D",'(1x,A,A24,A,ES24.17)') &
345  '[', trim(varname), '] MEAN(local) = ', mean_
346  endif
347  endif
348 
349  if ( present(mean) ) mean = mean_
350  if ( present(sum ) ) sum = sum_
351 
352  return
353  end subroutine statistics_total_3d
354 
355  !-----------------------------------------------------------------------------
357  subroutine statistics_horizontal_mean_2d( &
358  IA, IS, IE, JA, JS, JE, &
359  var, area, &
360  varmean )
361  use scale_const, only: &
362  undef => const_undef
363  integer, intent(in) :: IA, IS, IE
364  integer, intent(in) :: JA, JS, JE
365  real(RP), intent(in) :: var (IA,JA)
366  real(RP), intent(in) :: area(IA,JA)
367  real(RP), intent(out) :: varmean
368 
369  real(DP) :: statval (2)
370  real(DP) :: allstatval(2)
371 
372  integer :: ierr
373  integer :: i, j
374  !---------------------------------------------------------------------------
375 
376  statval(:) = 0.0_dp
377 ! !$omp parallel do reduction(+:statval)
378  do j = js, je
379  do i = is, ie
380  if ( var(i,j) /= undef ) then
381  statval(1) = statval(1) + area(i,j) * var(i,j)
382  statval(2) = statval(2) + area(i,j)
383  endif
384  enddo
385  enddo
386 
387  call prof_rapstart('COMM_Allreduce', 2)
388  ! All reduce
389  call mpi_allreduce( statval(:), &
390  allstatval(:), &
391  2, &
392  mpi_double_precision, &
393  mpi_sum, &
395  ierr )
396  call prof_rapend ('COMM_Allreduce', 2)
397 
398  if ( allstatval(2) > 0.0_dp ) then
399  varmean = allstatval(1) / allstatval(2)
400  else
401  varmean = undef
402  end if
403 
404  return
405  end subroutine statistics_horizontal_mean_2d
406 
407  subroutine statistics_horizontal_mean_3d( &
408  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
409  var, area, &
410  varmean )
411  use scale_const, only: &
412  undef => const_undef
413  integer, intent(in) :: KA, KS, KE
414  integer, intent(in) :: IA, IS, IE
415  integer, intent(in) :: JA, JS, JE
416 
417  real(RP), intent(in) :: var (KA,IA,JA)
418  real(RP), intent(in) :: area( IA,JA)
419  real(RP), intent(out) :: varmean(KA)
420 
421  real(DP) :: statval (2,KA)
422  real(DP) :: allstatval(2,KA)
423 
424  integer :: ierr
425  integer :: k, i, j
426  !---------------------------------------------------------------------------
427 
428  statval(:,:) = 0.0_dp
429 ! !$omp parallel do reduction(+:statval)
430  do j = js, je
431  do i = is, ie
432  do k = ks, ke
433  if ( var(k,i,j) /= undef ) then
434  statval(1,k) = statval(1,k) + area(i,j) * var(k,i,j)
435  statval(2,k) = statval(2,k) + area(i,j)
436  endif
437  enddo
438  enddo
439  enddo
440 
441  call prof_rapstart('COMM_Allreduce', 2)
442  ! All reduce
443  call mpi_allreduce( statval(:,ks:ke), &
444  allstatval(:,ks:ke), &
445  (ke-ks+1)*2, &
446  mpi_double_precision, &
447  mpi_sum, &
449  ierr )
450  call prof_rapend ('COMM_Allreduce', 2)
451 
452  do k = ks, ke
453  if ( allstatval(2,k) > 0.0_dp ) then
454  varmean(k) = allstatval(1,k) / allstatval(2,k)
455  else
456  varmean(k) = undef
457  end if
458  enddo
459  do k = 1, ks-1
460  varmean(k) = undef
461  end do
462  do k = ke+1, ka
463  varmean(k) = undef
464  end do
465 
466  return
467  end subroutine statistics_horizontal_mean_3d
468 
469  !-----------------------------------------------------------------------------
471  subroutine statistics_horizontal_min_2d( &
472  IA, IS, IE, JA, JS, JE, &
473  var, &
474  varmin )
475  use scale_const, only: &
476  undef => const_undef, &
477  huge => const_huge
478  use scale_comm_cartesc, only: &
480  integer, intent(in) :: IA, IS, IE
481  integer, intent(in) :: JA, JS, JE
482 
483  real(RP), intent(in) :: var(IA,JA)
484  real(RP), intent(out) :: varmin
485 
486  real(RP) :: statval
487  real(RP) :: allstatval
488 
489  integer :: ierr
490  integer :: i, j
491  !---------------------------------------------------------------------------
492 
493  statval = huge
494  !$omp parallel do reduction(min:statval)
495  do j = js, je
496  do i = is, ie
497  if ( var(i,j) /= undef .and. var(i,j) < statval ) then
498  statval = var(i,j)
499  endif
500  enddo
501  enddo
502 
503  call prof_rapstart('COMM_Allreduce', 2)
504  ! All reduce
505  call mpi_allreduce( statval, &
506  allstatval, &
507  1, &
508  comm_datatype, &
509  mpi_min, &
511  ierr )
512  call prof_rapend ('COMM_Allreduce', 2)
513 
514  if ( allstatval < huge ) then
515  varmin = allstatval
516  else
517  varmin = undef
518  end if
519 
520  return
521  end subroutine statistics_horizontal_min_2d
522 
523  subroutine statistics_horizontal_min_3d( &
524  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
525  var, &
526  varmin )
527  use scale_const, only: &
528  undef => const_undef, &
529  huge => const_huge
530  use scale_comm_cartesc, only: &
532  integer, intent(in) :: KA, KS, KE
533  integer, intent(in) :: IA, IS, IE
534  integer, intent(in) :: JA, JS, JE
535 
536  real(RP), intent(in) :: var(KA,IA,JA)
537  real(RP), intent(out) :: varmin(KA)
538 
539  real(RP) :: statval (KA)
540  real(RP) :: allstatval(KA)
541 
542  integer :: ierr
543  integer :: k, i, j
544  !---------------------------------------------------------------------------
545 
546  statval(:) = huge
547 ! !$omp parallel do reduction(min:statval)
548  do j = js, je
549  do i = is, ie
550  do k = ks, ke
551  if ( var(k,i,j) /= undef .and. var(k,i,j) < statval(k) ) then
552  statval(k) = var(k,i,j)
553  endif
554  enddo
555  enddo
556  enddo
557 
558  call prof_rapstart('COMM_Allreduce', 2)
559  ! All reduce
560  call mpi_allreduce( statval(ks:ke), &
561  allstatval(ks:ke), &
562  ke-ks+1, &
563  comm_datatype, &
564  mpi_min, &
566  ierr )
567  call prof_rapend ('COMM_Allreduce', 2)
568 
569  do k = ks, ke
570  if ( allstatval(k) < huge ) then
571  varmin(k) = allstatval(k)
572  else
573  varmin(k) = undef
574  end if
575  enddo
576  do k = 1, ks-1
577  varmin(k) = undef
578  end do
579  do k = ke+1, ka
580  varmin(k) = undef
581  end do
582 
583  return
584  end subroutine statistics_horizontal_min_3d
585 
586  !-----------------------------------------------------------------------------
588  subroutine statistics_horizontal_max_2d( &
589  IA, IS, IE, JA, JS, JE, &
590  var, &
591  varmax )
592  use scale_const, only: &
593  undef => const_undef, &
594  huge => const_huge
595  use scale_comm_cartesc, only: &
597  integer, intent(in) :: IA, IS, IE
598  integer, intent(in) :: JA, JS, JE
599 
600  real(RP), intent(in) :: var(IA,JA)
601  real(RP), intent(out) :: varmax
602 
603  real(RP) :: statval
604  real(RP) :: allstatval
605 
606  integer :: ierr
607  integer :: i, j
608  !---------------------------------------------------------------------------
609 
610  statval = - huge
611  !$omp parallel do reduction(max:statval)
612  do j = js, je
613  do i = is, ie
614  if ( var(i,j) /= undef .and. var(i,j) > statval ) then
615  statval = var(i,j)
616  endif
617  enddo
618  enddo
619 
620  call prof_rapstart('COMM_Allreduce', 2)
621  ! All reduce
622  call mpi_allreduce( statval, &
623  allstatval, &
624  1, &
625  comm_datatype, &
626  mpi_max, &
628  ierr )
629  call prof_rapend ('COMM_Allreduce', 2)
630 
631  if ( allstatval > - huge ) then
632  varmax = allstatval
633  else
634  varmax = undef
635  end if
636 
637  return
638  end subroutine statistics_horizontal_max_2d
639 
640  subroutine statistics_horizontal_max_3d( &
641  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
642  var, &
643  varmax )
644  use scale_const, only: &
645  undef => const_undef, &
646  huge => const_huge
647  use scale_comm_cartesc, only: &
649  integer, intent(in) :: KA, KS, KE
650  integer, intent(in) :: IA, IS, IE
651  integer, intent(in) :: JA, JS, JE
652 
653  real(RP), intent(in) :: var(KA,IA,JA)
654  real(RP), intent(out) :: varmax(KA)
655 
656  real(RP) :: statval (KA)
657  real(RP) :: allstatval(KA)
658 
659  integer :: ierr
660  integer :: k, i, j
661  !---------------------------------------------------------------------------
662 
663  statval(:) = - huge
664 ! !$omp parallel do reduction(max:statval)
665  do j = js, je
666  do i = is, ie
667  do k = ks, ke
668  if ( var(k,i,j) /= undef .and. var(k,i,j) > statval(k) ) then
669  statval(k) = var(k,i,j)
670  endif
671  enddo
672  enddo
673  enddo
674 
675  call prof_rapstart('COMM_Allreduce', 2)
676  ! All reduce
677  call mpi_allreduce( statval(ks:ke), &
678  allstatval(ks:ke), &
679  ke-ks+1, &
680  comm_datatype, &
681  mpi_max, &
683  ierr )
684  call prof_rapend ('COMM_Allreduce', 2)
685 
686  do k = ks, ke
687  if ( allstatval(k) > - huge ) then
688  varmax(k) = allstatval(k)
689  else
690  varmax(k) = undef
691  end if
692  enddo
693  do k = 1, ks-1
694  varmax(k) = undef
695  end do
696  do k = ke+1, ka
697  varmax(k) = undef
698  end do
699 
700  return
701  end subroutine statistics_horizontal_max_3d
702 
703  !-----------------------------------------------------------------------------
705  subroutine statistics_detail_3d( &
706  KA, KS, KE, IA, IS, IE, JA, JS, JE, VA, &
707  varname, var, &
708  local )
709  use scale_prc, only: &
710  prc_nprocs
711  use scale_comm_cartesc, only: &
713  implicit none
714 
715  integer, intent(in) :: KA, KS, KE
716  integer, intent(in) :: IA, IS, IE
717  integer, intent(in) :: JA, JS, JE
718  integer, intent(in) :: VA
719 
720  character(len=*), intent(in) :: varname(VA)
721  real(RP), intent(in) :: var(KA,IA,JA,VA)
722 
723  logical, intent(in), optional :: local
724 
725  real(RP) :: statval_l ( VA,2)
726  integer :: statidx_l (3,VA,2)
727  real(RP) :: statval ( VA,2,0:PRC_nprocs-1)
728  integer :: statidx (3,VA,2,0:PRC_nprocs-1)
729  real(RP) :: allstatval(VA,2)
730  integer :: allstatidx(VA,2)
731  logical :: do_globalcomm
732 
733  integer :: k, i, j
734  integer :: ierr
735  integer :: v, p
736  !---------------------------------------------------------------------------
737 
738  do_globalcomm = statistics_use_globalcomm
739  if ( present(local) ) do_globalcomm = ( .not. local )
740 
741  log_newline
742  log_info("STATISTICS_detail_3D",*) 'Variable Statistics '
743  do v = 1, va
744  statval_l( v,:) = var(ks,is,js,v)
745  statidx_l(1,v,:) = ks
746  statidx_l(2,v,:) = is
747  statidx_l(3,v,:) = js
748  do j = js, je
749  do i = is, ie
750  do k = ks, ke
751  if ( var(k,i,j,v) > statval_l(v,1) ) then
752  statval_l( v,1) = var(k,i,j,v)
753  statidx_l(1,v,1) = k
754  statidx_l(2,v,1) = i
755  statidx_l(3,v,1) = j
756  end if
757  if ( var(k,i,j,v) < statval_l(v,2) ) then
758  statval_l( v,2) = var(k,i,j,v)
759  statidx_l(1,v,2) = k
760  statidx_l(2,v,2) = i
761  statidx_l(3,v,2) = j
762  end if
763  end do
764  end do
765  end do
766  enddo
767 
768  if ( do_globalcomm ) then
769  call prof_rapstart('COMM_Bcast', 2)
770 
771  call mpi_allgather( statval_l(:,:), &
772  va*2, &
773  comm_datatype, &
774  statval(:,:,:), &
775  va*2, &
776  comm_datatype, &
777  prc_local_comm_world, &
778  ierr )
779 
780  call mpi_allgather( statidx_l(:,:,:), &
781  3*va*2, &
782  mpi_integer, &
783  statidx(:,:,:,:), &
784  3*va*2, &
785  mpi_integer, &
786  prc_local_comm_world, &
787  ierr )
788 
789  call prof_rapend ('COMM_Bcast', 2)
790 
791  do v = 1, va
792  allstatval(v,1) = statval(v,1,0)
793  allstatval(v,2) = statval(v,2,0)
794  allstatidx(v,:) = 0
795  do p = 1, prc_nprocs-1
796  if ( statval(v,1,p) > allstatval(v,1) ) then
797  allstatval(v,1) = statval(v,1,p)
798  allstatidx(v,1) = p
799  end if
800  if ( statval(v,2,p) < allstatval(v,2) ) then
801  allstatval(v,2) = statval(v,2,p)
802  allstatidx(v,2) = p
803  end if
804  end do
805  log_info_cont(*) '[', trim(varname(v)), ']'
806  log_info_cont('(1x,A,ES17.10,A,4(I5,A))') ' MAX =', &
807  allstatval(v,1), ' (rank=', &
808  allstatidx(v,1), '; ', &
809  statidx(1,v,1,allstatidx(v,1)),',', &
810  statidx(2,v,1,allstatidx(v,1)),',', &
811  statidx(3,v,1,allstatidx(v,1)),')'
812  log_info_cont('(1x,A,ES17.10,A,4(I5,A))') ' MIN =', &
813  allstatval(v,2), ' (rank=', &
814  allstatidx(v,2), '; ', &
815  statidx(1,v,2,allstatidx(v,2)),',', &
816  statidx(2,v,2,allstatidx(v,2)),',', &
817  statidx(3,v,2,allstatidx(v,2)),')'
818  enddo
819  else
820  ! statistics on each node
821  do v = 1, va
822  log_info_cont(*) '[', trim(varname(v)), ']'
823  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') 'MAX = ', &
824  statval_l( v,1),' (', &
825  statidx_l(1,v,1),',', &
826  statidx_l(2,v,1),',', &
827  statidx_l(3,v,1),')'
828  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') 'MIN = ', &
829  statval_l( v,2),' (', &
830  statidx_l(1,v,2),',', &
831  statidx_l(2,v,2),',', &
832  statidx_l(3,v,2),')'
833  enddo
834  endif
835 
836  log_newline
837 
838  return
839  end subroutine statistics_detail_3d
840 
841  subroutine statistics_detail_2d( &
842  IA, IS, IE, JA, JS, JE, VA, &
843  varname, var, &
844  local )
845  use scale_prc, only: &
846  prc_nprocs, &
847  prc_myrank
848  use scale_comm_cartesc, only: &
850  implicit none
851 
852  integer, intent(in) :: IA, IS, IE
853  integer, intent(in) :: JA, JS, JE
854  integer, intent(in) :: VA
855 
856  character(len=*), intent(in) :: varname(VA)
857  real(RP), intent(in) :: var(IA,JA,VA)
858 
859  logical, intent(in), optional :: local ! calc in local node
860 
861  real(RP) :: statval_l ( VA,2)
862  integer :: statidx_l (2,VA,2)
863  real(RP) :: statval ( VA,2,0:PRC_nprocs-1)
864  integer :: statidx (2,VA,2,0:PRC_nprocs-1)
865  real(RP) :: allstatval(VA,2)
866  integer :: allstatidx(VA,2)
867  logical :: do_globalcomm
868 
869  integer :: i, j
870  integer :: ierr
871  integer :: v, p
872  !---------------------------------------------------------------------------
873 
874  do_globalcomm = statistics_use_globalcomm
875  if ( present(local) ) do_globalcomm = ( .not. local )
876 
877  log_newline
878  log_info("STATISTICS_detail_2D",*) 'Variable Statistics '
879  do v = 1, va
880  statval_l( v,:) = var(is,js,v)
881  statidx_l(1,v,:) = is
882  statidx_l(2,v,:) = js
883  do j = js, je
884  do i = is, ie
885  if ( var(i,j,v) > statval_l(v,1) ) then
886  statval_l( v,1) = var(i,j,v)
887  statidx_l(1,v,1) = i
888  statidx_l(2,v,1) = j
889  end if
890  if ( var(i,j,v) < statval_l(v,2) ) then
891  statval_l( v,2) = var(i,j,v)
892  statidx_l(1,v,2) = i
893  statidx_l(2,v,2) = j
894  end if
895  end do
896  end do
897  enddo
898 
899  if ( do_globalcomm ) then
900  call prof_rapstart('COMM_Bcast', 2)
901 
902  call mpi_allgather( statval_l(:,:), &
903  va*2, &
904  comm_datatype, &
905  statval(:,:,:), &
906  va*2, &
907  comm_datatype, &
908  prc_local_comm_world, &
909  ierr )
910 
911  call mpi_allgather( statidx_l(:,:,:), &
912  2*va*2, &
913  mpi_integer, &
914  statidx(:,:,:,:), &
915  2*va*2, &
916  mpi_integer, &
917  prc_local_comm_world, &
918  ierr )
919 
920  call prof_rapend ('COMM_Bcast', 2)
921 
922  do v = 1, va
923  allstatval(v,1) = statval(v,1,0)
924  allstatval(v,2) = statval(v,2,0)
925  allstatidx(v,:) = 0
926  do p = 1, prc_nprocs-1
927  if ( statval(v,1,p) > allstatval(v,1) ) then
928  allstatval(v,1) = statval(v,1,p)
929  allstatidx(v,1) = p
930  end if
931  if ( statval(v,2,p) < allstatval(v,2) ) then
932  allstatval(v,2) = statval(v,2,p)
933  allstatidx(v,2) = p
934  end if
935  end do
936  log_info_cont(*) '[', trim(varname(v)), ']'
937  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') ' MAX =', &
938  allstatval(v,1), ' (rank=', &
939  allstatidx(v,1), '; ', &
940  statidx(1,v,1,allstatidx(v,1)),',', &
941  statidx(2,v,1,allstatidx(v,1)),')'
942  log_info_cont('(1x,A,ES17.10,A,3(I5,A))') ' MIN =', &
943  allstatval(v,2), ' (rank=', &
944  allstatidx(v,2), '; ', &
945  statidx(1,v,2,allstatidx(v,2)),',', &
946  statidx(2,v,2,allstatidx(v,2)),')'
947  enddo
948  else
949  ! statistics on each node
950  do v = 1, va
951  log_info_cont(*) '[', trim(varname(v)), ']'
952  log_info_cont('(1x,A,ES17.10,A,2(I5,A))') 'MAX = ', &
953  statval_l( v,1),' (', &
954  statidx_l(1,v,1),',', &
955  statidx_l(2,v,1),')'
956  log_info_cont('(1x,A,ES17.10,A,2(I5,A))') 'MIN = ', &
957  statval_l( v,2),' (', &
958  statidx_l(1,v,2),',', &
959  statidx_l(2,v,2),')'
960  enddo
961  endif
962 
963  log_newline
964 
965  return
966  end subroutine statistics_detail_2d
967 
968 end module scale_statistics
scale_comm_cartesc::comm_setup
subroutine, public comm_setup
Setup.
Definition: scale_comm_cartesC.F90:145
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
scale_comm_cartesc::comm_datatype
integer, public comm_datatype
datatype of variable
Definition: scale_comm_cartesC.F90:98
scale_statistics::statistics_horizontal_min_2d
subroutine statistics_horizontal_min_2d(IA, IS, IE, JA, JS, JE, var, varmin)
Calc horizontal minimum value.
Definition: scale_statistics.F90:475
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_statistics::statistics_horizontal_max_2d
subroutine statistics_horizontal_max_2d(IA, IS, IE, JA, JS, JE, var, varmax)
Calc horizontal maximum value.
Definition: scale_statistics.F90:592
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_prc::prc_local_comm_world
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:88
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:90
scale_statistics::statistics_total_3d
subroutine statistics_total_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, varname, vol, total, log_suppress, global, mean, sum)
Calc domain sum and volume-weighted mean.
Definition: scale_statistics.F90:243
scale_statistics::statistics_horizontal_mean_2d
subroutine statistics_horizontal_mean_2d(IA, IS, IE, JA, JS, JE, var, area, varmean)
Calc horizontal mean value.
Definition: scale_statistics.F90:361
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_statistics::statistics_horizontal_min_3d
subroutine statistics_horizontal_min_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, varmin)
Definition: scale_statistics.F90:527
scale_io
module STDIO
Definition: scale_io.F90:10
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_huge
real(rp), public const_huge
huge number
Definition: scale_const.F90:35
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_statistics::statistics_horizontal_mean_3d
subroutine statistics_horizontal_mean_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, area, varmean)
Definition: scale_statistics.F90:411
scale_prc::prc_nprocs
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:89
scale_statistics::statistics_detail_2d
subroutine statistics_detail_2d(IA, IS, IE, JA, JS, JE, VA, varname, var, local)
Definition: scale_statistics.F90:845
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:64
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_statistics::statistics_total_2d
subroutine statistics_total_2d(IA, IS, IE, JA, JS, JE, var, varname, area, total, log_suppress, global, mean, sum)
Calc domain sum and area-weighted mean.
Definition: scale_statistics.F90:130
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_statistics::statistics_setup
subroutine, public statistics_setup
Setup.
Definition: scale_statistics.F90:81
scale_statistics::statistics_horizontal_max_3d
subroutine statistics_horizontal_max_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, var, varmax)
Definition: scale_statistics.F90:644
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_statistics::statistics_detail_3d
subroutine statistics_detail_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, VA, varname, var, local)
Search global maximum & minimum value.
Definition: scale_statistics.F90:709
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:56