Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
mean_field_output.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2025, 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 num_types, only : rp
36 use field_list, only : field_list_t
39 use map_2d, only : map_2d_t
40 use map_1d, only : map_1d_t
41 use coefs, only : coef_t
43 use mean_field, only : mean_field_t
44 use output, only : output_t
45 use matrix, only : matrix_t
46 implicit none
47 private
48
50 type, public, extends(output_t) :: mean_field_output_t
52 type(mean_field_t), pointer :: mean_fields(:)
54 type(field_list_t) :: fields
56 real(kind=rp) :: start_time = 0.0_rp
58 integer :: n_fields = 0
64 integer :: output_dim = 0
65 contains
67 procedure, pass(this) :: init => mean_field_output_init
69 procedure, pass(this) :: free => mean_field_output_free
71 procedure, pass(this) :: sample => mean_field_output_sample
72 end type mean_field_output_t
73
74contains
75
85 subroutine mean_field_output_init(this, mean_fields, n_fields, start_time, &
86 coef, avg_dir, name, path)
87 class(mean_field_output_t), intent(inout):: this
88 integer, intent(in) :: n_fields
89 type(mean_field_t), intent(inout), target :: mean_fields(n_fields)
90 type(coef_t), intent(inout) :: coef
91 character(len=*), intent(in) :: avg_dir
92 character(len=*), intent(in), optional :: name
93 character(len=*), intent(in), optional :: path
94 real(kind=rp), intent(in) :: start_time
95 character(len=1024) :: fname
96 integer :: i
97
98 this%start_time = start_time
99
100 if (trim(avg_dir) .eq. 'none' .or. &
101 trim(avg_dir) .eq. 'x' .or.&
102 trim(avg_dir) .eq. 'y' .or.&
103 trim(avg_dir) .eq. 'z'&
104 ) then
105 if (present(name) .and. present(path)) then
106 fname = trim(path) // trim(name) // '.fld'
107 else if (present(name)) then
108 fname = trim(name) // '.fld'
109 else if (present(path)) then
110 fname = trim(path) // 'user_stats.fld'
111 else
112 fname = 'user_stats.fld'
113 end if
114
115 this%output_dim = 3
116
117 if (trim(avg_dir) .eq. 'x' .or.&
118 trim(avg_dir) .eq. 'y' .or.&
119 trim(avg_dir) .eq. 'z' ) then
120 call this%map_2d%init_char(coef, avg_dir, 1e-7_rp)
121 this%output_dim = 2
122 end if
123 else
124 if (present(name) .and. present(path)) then
125 fname = trim(path) // trim(name) // '.csv'
126 else if (present(name)) then
127 fname = trim(name) // '.csv'
128 else if (present(path)) then
129 fname = trim(path) // 'user_stats.csv'
130 else
131 fname = 'user_stats.csv'
132 end if
133 call this%map_1d%init_char(coef, avg_dir, 1e-7_rp)
134 this%output_dim = 1
135 end if
136
137 call this%init_base(fname)
138
139 call this%fields%init(n_fields)
140 this%n_fields = n_fields
141 this%mean_fields => mean_fields
142 do i = 1, n_fields
143 this%fields%items(i)%ptr => this%mean_fields(i)%mf
144 end do
145
146 end subroutine mean_field_output_init
147
149 subroutine mean_field_output_free(this)
150 class(mean_field_output_t), intent(inout) :: this
151
152 call this%free_base()
153
154 nullify(this%mean_fields)
155 call this%map_1d%free()
156 call this%map_2d%free()
157 call this%fields%free()
158
159 end subroutine mean_field_output_free
160
162 subroutine mean_field_output_sample(this, t)
163 class(mean_field_output_t), intent(inout) :: this
164 real(kind=rp), intent(in) :: t
165 integer :: i
166 type(fld_file_data_t) :: output_2d
167 type(matrix_t) :: avg_output_1d
168 real(kind=rp) :: u, v, w, p
169
170 associate(out_fields => this%fields%items)
171 if (t .ge. this%start_time) then
172 if ( neko_bcknd_device .eq. 1) then
173 do i = 1, size(out_fields)
174 call device_memcpy(out_fields(i)%ptr%x, out_fields(i)%ptr%x_d,&
175 out_fields(i)%ptr%dof%size(), device_to_host, &
176 sync = (i .eq. size(out_fields))) ! Sync on last field
177 end do
178 end if
179 if (this%output_dim .eq. 1) then
180 call this%map_1d%average_planes(avg_output_1d, &
181 this%fields)
182 call this%file_%write(avg_output_1d, t)
183 else if (this%output_dim .eq. 2) then
184 call this%map_2d%average(output_2d, this%fields)
185 call this%file_%write(output_2d, t)
186 else
187 call this%file_%write(this%fields, t)
188 end if
189 do i = 1, this%n_fields
190 call this%mean_fields(i)%reset()
191 end do
192 end if
193 end associate
194
195 end subroutine mean_field_output_sample
196
197end module mean_field_output
Copy data between host and device (or device and device)
Definition device.F90:71
Coefficients.
Definition coef.f90:34
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
Defines an output for a list of mean fields.
subroutine mean_field_output_free(this)
Destructor.
subroutine mean_field_output_init(this, mean_fields, n_fields, start_time, coef, avg_dir, name, path)
Constructor.
subroutine mean_field_output_sample(this, t)
Sample the mean solution at time t and reset.
Implements mean_field_t.
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
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
field_list_t, To be able to group fields together
Type that encapsulates a mapping from each gll point in the mesh to its corresponding (global) GLL po...
Definition map_1d.f90:27
Computes the temporal mean of a field.
Output for a list of mean fields.
Abstract type defining an output type.
Definition output.f90:41