Neko  0.9.99
A portable framework for high-order spectral element flow simulations
fluid_stats_output.f90
Go to the documentation of this file.
1 ! Copyright (c) 2024, 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 !
35  use fluid_stats, only : fluid_stats_t
36  use neko_config, only : neko_bcknd_device
37  use num_types, only : rp
38  use map_1d, only : map_1d_t
39  use map_2d, only : map_2d_t
40  use fld_file_data, only : fld_file_data_t
41  use device
42  use output, only : output_t
43  use matrix, only : matrix_t
44  implicit none
45  private
46 
49  type, public, extends(output_t) :: fluid_stats_output_t
51  type(fluid_stats_t), pointer :: stats
53  type(map_1d_t) :: map_1d
55  type(map_2d_t) :: map_2d
56  real(kind=rp) :: t_begin
58  integer :: output_dim
59  contains
61  procedure, pass(this) :: init => fluid_stats_output_init
63  procedure, pass(this) :: sample => fluid_stats_output_sample
64  end type fluid_stats_output_t
65 
66 
67 contains
68 
70  subroutine fluid_stats_output_init(this, stats, T_begin, hom_dir, name, path)
71  type(fluid_stats_t), intent(inout), target :: stats
72  real(kind=rp), intent(in) :: t_begin
73  character(len=*), intent(in) :: hom_dir
74  character(len=*), intent(in), optional :: name
75  character(len=*), intent(in), optional :: path
76  class(fluid_stats_output_t), intent(inout) :: this
77  character(len=1024) :: fname
78 
79  if (trim(hom_dir) .eq. 'none' .or. &
80  trim(hom_dir) .eq. 'x' .or.&
81  trim(hom_dir) .eq. 'y' .or.&
82  trim(hom_dir) .eq. 'z'&
83  ) then
84  if (present(name) .and. present(path)) then
85  fname = trim(path) // trim(name) // '.fld'
86  else if (present(name)) then
87  fname = trim(name) // '.fld'
88  else if (present(path)) then
89  fname = trim(path) // 'fluid_stats.fld'
90  else
91  fname = 'fluid_stats.fld'
92  end if
93 
94  this%output_dim = 3
95 
96  if (trim(hom_dir) .eq. 'x' .or.&
97  trim(hom_dir) .eq. 'y' .or.&
98  trim(hom_dir) .eq. 'z' ) then
99  call this%map_2d%init_char(stats%coef, hom_dir, 1e-7_rp)
100  this%output_dim = 2
101  end if
102  else
103  if (present(name) .and. present(path)) then
104  fname = trim(path) // trim(name) // '.csv'
105  else if (present(name)) then
106  fname = trim(name) // '.csv'
107  else if (present(path)) then
108  fname = trim(path) // 'fluid_stats.csv'
109  else
110  fname = 'fluid_stats.csv'
111  end if
112  call this%map_1d%init_char(stats%coef, hom_dir, 1e-7_rp)
113  this%output_dim = 1
114  end if
115 
116  call this%init_base(fname)
117  this%stats => stats
118  this%T_begin = t_begin
119  end subroutine fluid_stats_output_init
120 
122  subroutine fluid_stats_output_sample(this, t)
123  class(fluid_stats_output_t), intent(inout) :: this
124  real(kind=rp), intent(in) :: t
125  integer :: i
126  type(matrix_t) :: avg_output_1d
127  type(fld_file_data_t) :: output_2d
128  real(kind=rp) :: u, v, w, p
129  associate(out_fields => this%stats%stat_fields%items)
130  if (t .ge. this%T_begin) then
131  call this%stats%make_strong_grad()
132  if ( neko_bcknd_device .eq. 1) then
133  do i = 1, size(out_fields)
134  call device_memcpy(out_fields(i)%ptr%x, out_fields(i)%ptr%x_d,&
135  out_fields(i)%ptr%dof%size(), device_to_host, &
136  sync = (i .eq. size(out_fields))) ! Sync on last field
137  end do
138  end if
139  if (this%output_dim .eq. 1) then
140  call this%map_1d%average_planes(avg_output_1d, &
141  this%stats%stat_fields)
142  call this%file_%write(avg_output_1d, t)
143  else if (this%output_dim .eq. 2) then
144  call this%map_2d%average(output_2d, this%stats%stat_fields)
145  !Switch around fields to get correct orders
146  !Put average direction mean_vel in scalar45
147  do i = 1, this%map_2d%n_2d
148  u = output_2d%v%x(i)
149  v = output_2d%w%x(i)
150  w = output_2d%p%x(i)
151  p = output_2d%u%x(i)
152  output_2d%p%x(i) = p
153  output_2d%u%x(i) = u
154  output_2d%v%x(i) = v
155  output_2d%w%x(i) = w
156  end do
157 
158  call this%file_%write(output_2d, t)
159  else
160  call this%file_%write(this%stats%stat_fields, t)
161  end if
162  call this%stats%reset()
163  end if
164  end associate
165  end subroutine fluid_stats_output_sample
166 
167 end module fluid_stats_output
168 
169 
Copy data between host and device (or device and device)
Definition: device.F90:51
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public device_to_host
Definition: device.F90:47
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
Implements fluid_stats_ouput_t.
subroutine fluid_stats_output_sample(this, t)
Sample fluid_stats at time t.
subroutine fluid_stats_output_init(this, stats, T_begin, hom_dir, name, path)
Constructor.
Computes various statistics for the fluid fields. We use the Reynolds decomposition for a field u = ...
Definition: fluid_stats.f90:36
Creates a 1d GLL point map along a specified direction based on the connectivity in the mesh.
Definition: map_1d.f90:3
Maps a 3D dofmap to a 2D spectral element grid.
Definition: map_2d.f90:3
Defines a matrix.
Definition: matrix.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines an output.
Definition: output.f90:34
Defines a container for all statistics.
Definition: statistics.f90:34
Defines an output for the fluid statistics computed using the fluid_stats_t object.
Type that encapsulates a mapping from each gll point in the mesh to its corresponding (global) GLL po...
Definition: map_1d.f90:26
Abstract type defining an output type.
Definition: output.f90:41