Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
34module 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 mpi_f08, only : mpi_barrier, mpi_wtime
39 use comm, only : neko_comm
40 implicit none
41 private
42
44 type, private :: quantp_t
45 class(stats_quant_t), pointer :: quantp
46 contains
48 procedure, pass(this) :: free => quantp_free
49 end type quantp_t
50
52 type, public :: stats_t
53 type(quantp_t), allocatable :: quant_list(:)
54 integer :: n
55 integer :: size
56 real(kind=rp) :: t_begin
57 real(kind=rp) :: t_diff
58 integer :: samp_interval
59 contains
60 procedure, pass(this) :: init => stats_init
61 procedure, pass(this) :: free => stats_free
62 procedure, pass(this) :: add => stats_add
63 procedure, pass(this) :: eval => stats_eval
64 end type stats_t
65
66contains
68 subroutine quantp_free(this)
69 class(quantp_t), intent(inout) :: this
70
71 nullify(this%quantp)
72 end subroutine quantp_free
73
75 subroutine stats_init(this, T_begin, samp_interval, size)
76 class(stats_t), intent(inout) :: this
77 real(kind=rp), intent(in) :: t_begin
78 integer, intent(in) :: samp_interval
79 integer, intent(inout), optional ::size
80 integer :: n, i
81
82 call this%free()
83
84 if (present(size)) then
85 n = size
86 else
87 n = 1
88 end if
89
90 allocate(this%quant_list(n))
91
92 do i = 1, n
93 this%quant_list(i)%quantp => null()
94 end do
95
96 this%n = 0
97 this%size = n
98 this%T_begin = t_begin
99 this%samp_interval = samp_interval
100 this%t_diff = 0.0
101
102 end subroutine stats_init
103
105 subroutine stats_free(this)
106 class(stats_t), intent(inout) :: this
107 integer :: i
108
109 if (allocated(this%quant_list)) then
110 do i = 1, size(this%quant_list)
111 call this%quant_list(i)%free()
112 end do
113
114 deallocate(this%quant_list)
115 end if
116
117 this%n = 0
118 this%size = 0
119 end subroutine stats_free
120
122 subroutine stats_add(this, quant)
123 class(stats_t), intent(inout) :: this
124 class(stats_quant_t), intent(inout), target :: quant
125 type(quantp_t), allocatable :: tmp(:)
126
127 if (this%n .ge. this%size) then
128 allocate(tmp(this%size * 2))
129 tmp(1:this%size) = this%quant_list
130 call move_alloc(tmp, this%quant_list)
131 this%size = this%size * 2
132 end if
133
134 this%n = this%n + 1
135 this%quant_list(this%n)%quantp => quant
136 end subroutine stats_add
137
139 subroutine stats_eval(this, t, dt, tstep)
140 class(stats_t), intent(inout) :: this
141 real(kind=rp), intent(in) :: t
142 real(kind=rp), intent(in) :: dt
143 integer, intent(in) :: tstep
144 integer :: i, ierr
145 character(len=LOG_SIZE) :: log_buf
146 real(kind=rp) :: sample_start_time, sample_end_time
147 real(kind=dp) :: sample_time
148
149 if (t .ge. this%T_begin .and. this%n .gt. 0) then
150 this%t_diff = this%t_diff + dt
151 ! There is technically an issue here for the last sample if we
152 ! reset the stats If the reset is not on a multiple of
153 ! samp_interval the weight of the last sample is wrong.
154 if (mod(tstep,this%samp_interval) .eq. 0) then
155 call neko_log%section('Statistics')
156 call mpi_barrier(neko_comm, ierr)
157 sample_start_time = mpi_wtime()
158 do i = 1, this%n
159 call this%quant_list(i)%quantp%update(this%t_diff)
160 end do
161 this%t_diff = 0.0_rp
162 call mpi_barrier(neko_comm, ierr)
163 sample_end_time = mpi_wtime()
164 sample_time = sample_end_time - sample_start_time
165 write(log_buf,'(A17,1x,F10.6,A,F9.6)') 'Sampling at time:', t, &
166 ' Sampling time (s): ', sample_time
167 call neko_log%message(log_buf)
168 call neko_log%end_section()
169 end if
170 end if
171 end subroutine stats_eval
172
173end module stats
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
integer, parameter, public log_size
Definition log.f90:46
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.
Defines a container for all statistics.
subroutine stats_eval(this, t, dt, tstep)
Evaluated all statistical quantities.
subroutine stats_free(this)
Deallocate.
subroutine stats_init(this, t_begin, samp_interval, size)
Initialize statistics, computed after T_begin.
subroutine stats_add(this, quant)
Add a statistic quantitiy quant to the backend.
subroutine quantp_free(this)
Destructor for quantp_t.
Pointer to an arbitrary quantitiy.
Statistics backend.
Abstract type defining a statistical quantity.