Neko  0.8.99
A portable framework for high-order spectral element flow simulations
statistics.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module stats
35  use num_types, only : rp, dp
36  use stats_quant, only : stats_quant_t
37  use logger, only : log_size, neko_log
38  use comm
39  implicit none
40  private
41 
43  type, private :: quantp_t
44  class(stats_quant_t), pointer :: quantp
45  end type quantp_t
46 
48  type, public :: stats_t
49  type(quantp_t), allocatable :: quant_list(:)
50  integer :: n
51  integer :: size
52  real(kind=rp) :: t_begin
53  real(kind=rp) :: t_diff
54  integer :: samp_interval
55  contains
56  procedure, pass(this) :: init => stats_init
57  procedure, pass(this) :: free => stats_free
58  procedure, pass(this) :: add => stats_add
59  procedure, pass(this) :: eval => stats_eval
60  end type stats_t
61 
62 contains
63 
65  subroutine stats_init(this, T_begin, samp_interval, size)
66  class(stats_t), intent(inout) :: this
67  real(kind=rp), intent(in) :: t_begin
68  integer, intent(in) :: samp_interval
69  integer, intent(inout), optional ::size
70  integer :: n, i
71 
72  call this%free()
73 
74  if (present(size)) then
75  n = size
76  else
77  n = 1
78  end if
79 
80  allocate(this%quant_list(n))
81 
82  do i = 1, n
83  this%quant_list(i)%quantp => null()
84  end do
85 
86  this%n = 0
87  this%size = n
88  this%T_begin = t_begin
89  this%samp_interval = samp_interval
90  this%t_diff = 0.0
91 
92  end subroutine stats_init
93 
95  subroutine stats_free(this)
96  class(stats_t), intent(inout) :: this
97 
98  if (allocated(this%quant_list)) then
99  deallocate(this%quant_list)
100  end if
101 
102  this%n = 0
103  this%size = 0
104  end subroutine stats_free
105 
107  subroutine stats_add(this, quant)
108  class(stats_t), intent(inout) :: this
109  class(stats_quant_t), intent(inout), target :: quant
110  type(quantp_t), allocatable :: tmp(:)
111 
112  if (this%n .ge. this%size) then
113  allocate(tmp(this%size * 2))
114  tmp(1:this%size) = this%quant_list
115  call move_alloc(tmp, this%quant_list)
116  this%size = this%size * 2
117  end if
118 
119  this%n = this%n + 1
120  this%quant_list(this%n)%quantp => quant
121  end subroutine stats_add
122 
124  subroutine stats_eval(this, t, dt, tstep)
125  class(stats_t), intent(inout) :: this
126  real(kind=rp), intent(in) :: t
127  real(kind=rp), intent(in) :: dt
128  integer, intent(in) :: tstep
129  integer :: i, ierr
130  character(len=LOG_SIZE) :: log_buf
131  real(kind=rp) :: sample_start_time, sample_end_time
132  real(kind=dp) :: sample_time
133 
134  if (t .ge. this%T_begin .and. this%n .gt. 0) then
135  this%t_diff = this%t_diff + dt
136  ! There is technically an issue here for the last sample if we
137  ! reset the stats If the reset is not on a multiple of
138  ! samp_interval the weight of the last sample is wrong.
139  if (mod(tstep,this%samp_interval) .eq. 0) then
140  call neko_log%section('Statistics')
141  call mpi_barrier(neko_comm, ierr)
142  sample_start_time = mpi_wtime()
143  do i = 1, this%n
144  call this%quant_list(i)%quantp%update(this%t_diff)
145  end do
146  this%t_diff = 0.0_rp
147  call mpi_barrier(neko_comm, ierr)
148  sample_end_time = mpi_wtime()
149  sample_time = sample_end_time - sample_start_time
150  write(log_buf,'(A17,1x,F10.6,A,F9.6)') 'Sampling at time:', t, &
151  ' Sampling time (s): ', sample_time
152  call neko_log%message(log_buf)
153  call neko_log%end_section()
154  end if
155  end if
156  end subroutine stats_eval
157 
158 end module stats
Definition: comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
Logging routines.
Definition: log.f90:34
type(log_t), public neko_log
Global log stream.
Definition: log.f90:61
integer, parameter, public log_size
Definition: log.f90:40
integer, parameter, public dp
Definition: num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a statistical quantity.
Definition: stats_quant.f90:34
Defines a container for all statistics.
Definition: statistics.f90:34
subroutine stats_init(this, T_begin, samp_interval, size)
Initialize statistics, computed after T_begin.
Definition: statistics.f90:66
subroutine stats_eval(this, t, dt, tstep)
Evaluated all statistical quantities.
Definition: statistics.f90:125
subroutine stats_free(this)
Deallocate.
Definition: statistics.f90:96
subroutine stats_add(this, quant)
Add a statistic quantitiy quant to the backend.
Definition: statistics.f90:108
Pointer to an arbitrary quantitiy.
Definition: statistics.f90:43
Statistics backend.
Definition: statistics.f90:48
Abstract type defining a statistical quantity.
Definition: stats_quant.f90:40