Neko  0.8.1
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
36  use stats_quant
37  use logger
38  use comm
39  implicit none
40 
42  type, private :: quantp_t
43  class(stats_quant_t), pointer :: quantp
44  end type quantp_t
45 
47  type :: stats_t
48  type(quantp_t), allocatable :: quant_list(:)
49  integer :: n
50  integer :: size
51  real(kind=rp) :: t_begin
52  real(kind=rp) :: t_diff
53  integer :: samp_interval
54  contains
55  procedure, pass(this) :: init => stats_init
56  procedure, pass(this) :: free => stats_free
57  procedure, pass(this) :: add => stats_add
58  procedure, pass(this) :: eval => stats_eval
59  end type stats_t
60 
61 contains
62 
64  subroutine stats_init(this, T_begin, samp_interval, size)
65  class(stats_t), intent(inout) :: this
66  real(kind=rp), intent(in) :: t_begin
67  integer, intent(in) :: samp_interval
68  integer, intent(inout), optional ::size
69  integer :: n, i
70 
71  call this%free()
72 
73  if (present(size)) then
74  n = size
75  else
76  n = 1
77  end if
78 
79  allocate(this%quant_list(n))
80 
81  do i = 1, n
82  this%quant_list(i)%quantp => null()
83  end do
84 
85  this%n = 0
86  this%size = n
87  this%T_begin = t_begin
88  this%samp_interval = samp_interval
89  this%t_diff = 0.0
90 
91  end subroutine stats_init
92 
94  subroutine stats_free(this)
95  class(stats_t), intent(inout) :: this
96 
97  if (allocated(this%quant_list)) then
98  deallocate(this%quant_list)
99  end if
100 
101  this%n = 0
102  this%size = 0
103  end subroutine stats_free
104 
106  subroutine stats_add(this, quant)
107  class(stats_t), intent(inout) :: this
108  class(stats_quant_t), intent(inout), target :: quant
109  type(quantp_t), allocatable :: tmp(:)
110 
111  if (this%n .ge. this%size) then
112  allocate(tmp(this%size * 2))
113  tmp(1:this%size) = this%quant_list
114  call move_alloc(tmp, this%quant_list)
115  this%size = this%size * 2
116  end if
117 
118  this%n = this%n + 1
119  this%quant_list(this%n)%quantp => quant
120  end subroutine stats_add
121 
123  subroutine stats_eval(this, t, dt, tstep)
124  class(stats_t), intent(inout) :: this
125  real(kind=rp), intent(in) :: t
126  real(kind=rp), intent(in) :: dt
127  integer, intent(in) :: tstep
128  integer :: i, ierr
129  character(len=LOG_SIZE) :: log_buf
130  real(kind=rp) :: sample_start_time, sample_end_time
131  real(kind=dp) :: sample_time
132 
133  if (t .ge. this%T_begin .and. this%n .gt. 0) then
134  this%t_diff = this%t_diff + dt
135  ! There is technically an issue here for the last sample if we
136  ! reset the stats If the reset is not on a multiple of
137  ! samp_interval the weight of the last sample is wrong.
138  if (mod(tstep,this%samp_interval) .eq. 0) then
139  call neko_log%section('Statistics')
140  call mpi_barrier(neko_comm, ierr)
141  sample_start_time = mpi_wtime()
142  do i = 1, this%n
143  call this%quant_list(i)%quantp%update(this%t_diff)
144  end do
145  this%t_diff = 0.0_rp
146  call mpi_barrier(neko_comm, ierr)
147  sample_end_time = mpi_wtime()
148  sample_time = sample_end_time - sample_start_time
149  write(log_buf,'(A17,1x,F10.6,A,F9.6)') 'Sampling at time:', t, &
150  ' Sampling time (s): ', sample_time
151  call neko_log%message(log_buf)
152  call neko_log%end_section()
153  end if
154  end if
155  end subroutine stats_eval
156 
157 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 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:65
subroutine stats_eval(this, t, dt, tstep)
Evaluated all statistical quantities.
Definition: statistics.f90:124
subroutine stats_free(this)
Deallocate.
Definition: statistics.f90:95
subroutine stats_add(this, quant)
Add a statistic quantitiy quant to the backend.
Definition: statistics.f90:107
Pointer to an arbitrary quantitiy.
Definition: statistics.f90:42
Statistics backend.
Definition: statistics.f90:47
Abstract type defining a statistical quantity.
Definition: stats_quant.f90:40