Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_sgs_stats_output.f90
Go to the documentation of this file.
1! Copyright (c) 2026, 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!
37 use num_types, only : rp
38 use map_1d, only : map_1d_t
39 use map_2d, only : map_2d_t
42 use output, only : output_t
43 use matrix, only : matrix_t
44 implicit none
45 private
46
49 type, public, extends(output_t) :: scalar_sgs_stats_output_t
51 type(scalar_sgs_stats_t), pointer :: stats => null()
56 real(kind=rp) :: t_begin
58 integer :: output_dim
59 contains
61 procedure, pass(this) :: init => scalar_sgs_stats_output_init
63 procedure, pass(this) :: free => scalar_sgs_stats_output_free
65 procedure, pass(this) :: sample => scalar_sgs_stats_output_sample
67
68
69contains
70
72 subroutine scalar_sgs_stats_output_init(this, stats, T_begin, &
73 hom_dir, name, path)
74 class(scalar_sgs_stats_output_t), intent(inout) :: this
75 type(scalar_sgs_stats_t), intent(inout), target :: stats
76 real(kind=rp), intent(in) :: t_begin
77 character(len=*), intent(in) :: hom_dir
78 character(len=*), intent(in), optional :: name
79 character(len=*), intent(in), optional :: path
80 character(len=1024) :: fname
81
82 if (trim(hom_dir) .eq. 'none' .or. &
83 trim(hom_dir) .eq. 'x' .or.&
84 trim(hom_dir) .eq. 'y' .or.&
85 trim(hom_dir) .eq. 'z'&
86 ) then
87 if (present(name) .and. present(path)) then
88 fname = trim(path) // trim(name) // '.fld'
89 else if (present(name)) then
90 fname = trim(name) // '.fld'
91 else if (present(path)) then
92 fname = trim(path) // 'scalar_sgs_stats.fld'
93 else
94 fname = 'scalar_sgs_stats.fld'
95 end if
96
97 this%output_dim = 3
98
99 if (trim(hom_dir) .eq. 'x' .or.&
100 trim(hom_dir) .eq. 'y' .or.&
101 trim(hom_dir) .eq. 'z' ) then
102 call this%map_2d%init_char(stats%coef, hom_dir, 1e-7_rp)
103 this%output_dim = 2
104 end if
105 else
106 if (present(name) .and. present(path)) then
107 fname = trim(path) // trim(name) // '.csv'
108 else if (present(name)) then
109 fname = trim(name) // '.csv'
110 else if (present(path)) then
111 fname = trim(path) // 'scalar_sgs_stats.csv'
112 else
113 fname = 'scalar_sgs_stats.csv'
114 end if
115 call this%map_1d%init_char(stats%coef, hom_dir, 1e-7_rp)
116 this%output_dim = 1
117 end if
118
119 call this%init_base(fname)
120 this%stats => stats
121 this%T_begin = t_begin
122 end subroutine scalar_sgs_stats_output_init
123
126 class(scalar_sgs_stats_output_t), intent(inout) :: this
127
128 call this%free_base()
129
130 nullify(this%stats)
131 call this%map_1d%free()
132 call this%map_2d%free()
133
134 end subroutine scalar_sgs_stats_output_free
135
138 class(scalar_sgs_stats_output_t), intent(inout) :: this
139 real(kind=rp), intent(in) :: t
140 integer :: i
141 type(matrix_t) :: avg_output_1d
142 type(fld_file_data_t) :: output_2d
143 real(kind=rp) :: u, v, w, p
144 associate(out_fields => this%stats%stat_fields%items)
145 if (t .ge. this%T_begin) then
146 if ( neko_bcknd_device .eq. 1) then
147 do i = 1, size(out_fields)
148 call device_memcpy(out_fields(i)%ptr%x, out_fields(i)%ptr%x_d,&
149 out_fields(i)%ptr%dof%size(), device_to_host, &
150 sync = (i .eq. size(out_fields))) ! Sync on last field
151 end do
152 end if
153 if (this%output_dim .eq. 1) then
154 call this%map_1d%average_planes(avg_output_1d, &
155 this%stats%stat_fields)
156 call this%file_%write(avg_output_1d, t)
157 else if (this%output_dim .eq. 2) then
158 call this%map_2d%average(output_2d, this%stats%stat_fields)
159 !Switch around fields to get correct orders
160 do i = 1, this%map_2d%n_2d
161 u = output_2d%v%x(i)
162 v = output_2d%w%x(i)
163 w = output_2d%p%x(i)
164 p = output_2d%u%x(i)
165 output_2d%p%x(i) = p
166 output_2d%u%x(i) = u
167 output_2d%v%x(i) = v
168 output_2d%w%x(i) = w
169 end do
170
171 call this%file_%write(output_2d, t)
172 else
173 call this%file_%write(this%stats%stat_fields, t)
174 end if
175 call this%stats%reset()
176 end if
177 end associate
178 end subroutine scalar_sgs_stats_output_sample
179
Copy data between host and device (or device and device)
Definition device.F90:71
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...
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.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines an output.
Definition output.f90:34
Implements scalar_sgs_stats_ouput_t.
subroutine scalar_sgs_stats_output_free(this)
Destructor.
subroutine scalar_sgs_stats_output_sample(this, t)
Sample scalar_sgs_stats at time t.
subroutine scalar_sgs_stats_output_init(this, stats, t_begin, hom_dir, name, path)
Constructor.
Computes the subgrid-scale contributions for the scalar fluxes. We use the Reynolds decomposition for...
Defines a container for all statistics.
Type that encapsulates a mapping from each gll point in the mesh to its corresponding (global) GLL po...
Definition map_1d.f90:27
Abstract type defining an output type.
Definition output.f90:41
Defines an output for the sgs statistics for scalar computed using the scalar_sgs_stats_t object.