Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_output.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2023, 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
40 use field_list, only : field_list_t
42 use device
43 use output, only : output_t
44 use scalars, only : scalars_t
45 use registry, only : neko_registry
46 use field, only : field_t
47 use fld_file, only : fld_file_t
48 implicit none
49 private
50
52 type, public, extends(output_t) :: fluid_output_t
53 type(field_list_t) :: fluid
54 logical :: always_write_mesh = .false.
55 contains
56 procedure, pass(this) :: init => fluid_output_init
57 procedure, pass(this) :: sample => fluid_output_sample
58 procedure, pass(this) :: free => fluid_output_free
59 end type fluid_output_t
60
61contains
62
63 subroutine fluid_output_init(this, precision, fluid, scalar_fields, name, &
64 path, fmt, layout, always_write_mesh)
65 class(fluid_output_t), intent(inout) :: this
66 integer, intent(inout) :: precision
67 class(fluid_scheme_base_t), intent(in), target :: fluid
68 class(scalars_t), intent(in), optional, target :: scalar_fields
69 character(len=*), intent(in), optional :: name
70 character(len=*), intent(in), optional :: path
71 character(len=*), intent(in), optional :: fmt
72 logical, intent(in), optional :: always_write_mesh
73 integer, intent(in), optional :: layout
74 character(len=1024) :: fname
75 integer :: i, j, n_scalars
76 character(len=10) :: suffix
77 logical :: has_max_wave_speed, has_density
78 type(field_t), pointer :: max_wave_speed_field
79
80 suffix = '.fld'
81 if (present(fmt)) then
82 if (fmt .eq. 'adios2') then
83 suffix = '.bp'
84 end if
85 end if
86
87 if (present(always_write_mesh)) then
88 this%always_write_mesh = always_write_mesh
89 end if
90
91 if (present(name) .and. present(path)) then
92 fname = trim(path) // trim(name) // trim(suffix)
93 else if (present(name)) then
94 fname = trim(name) // trim(suffix)
95 else if (present(path)) then
96 fname = trim(path) // 'field' // trim(suffix)
97 else
98 fname = 'field' // trim(suffix)
99 end if
100
101 if (present(layout)) then
102 call this%init_base(fname, precision, layout)
103 else
104 call this%init_base(fname, precision)
105 end if
106
107 ! Calculate total number of fields
108 n_scalars = 0
109 if (present(scalar_fields)) then
110 n_scalars = size(scalar_fields%scalar_fields)
111 end if
112
113 ! Check if max_wave_speed field exists (for compressible flows)
114 has_max_wave_speed = neko_registry%field_exists("max_wave_speed")
115
116 ! Check if density field exists (for compressible flows)
117 ! We need to check the solver type here since the incompressible
118 ! solver also has a rho field due to the material properties
119 select type (fluid)
121 has_density = associated(fluid%rho)
122 class default
123 has_density = .false.
124 end select
125
126 ! Initialize field list with appropriate size
127 ! Standard fields: p, u, v, w (4)
128 ! Scalar fields: n_scalars
129 ! Compressible fields: density + max_wave_speed (2 additional)
130 i = 4
131
132 if (has_density) then
133 i = i + 1
134 end if
135
136 if (has_max_wave_speed) then
137 i = i + 1
138 end if
139
140 call this%fluid%init(i + n_scalars)
141
142 call this%fluid%assign(1, fluid%p)
143 call this%fluid%assign(2, fluid%u)
144 call this%fluid%assign(3, fluid%v)
145 call this%fluid%assign(4, fluid%w)
146
147 ! Assign all scalar fields first
148 i = 4
149 if (present(scalar_fields)) then
150 do j = 1, n_scalars
151 i = i + 1
152 call this%fluid%assign(i, scalar_fields%scalar_fields(j)%s)
153 end do
154 end if
155
156 ! Add density field if it exists (for compressible flows)
157 if (has_density) then
158 i = i + 1
159 call this%fluid%assign(i, fluid%rho)
160 end if
161
162 ! Add max_wave_speed field if it exists (for compressible flows)
163 if (has_max_wave_speed) then
164 i = i + 1
165 max_wave_speed_field => neko_registry%get_field("max_wave_speed")
166 call this%fluid%assign(i, max_wave_speed_field)
167 end if
168
169 end subroutine fluid_output_init
170
172 subroutine fluid_output_free(this)
173 class(fluid_output_t), intent(inout) :: this
174
175 call this%free_base()
176 call this%fluid%free()
177
178 end subroutine fluid_output_free
179
181 subroutine fluid_output_sample(this, t)
182 class(fluid_output_t), intent(inout) :: this
183 real(kind=rp), intent(in) :: t
184 integer :: i
185
186 if (neko_bcknd_device .eq. 1) then
187
188 associate(fields => this%fluid%items)
189 do i = 1, size(fields)
190 call device_memcpy(fields(i)%ptr%x, fields(i)%ptr%x_d, &
191 fields(i)%ptr%dof%size(), device_to_host, &
192 sync = (i .eq. size(fields))) ! Sync on the last field
193 end do
194 end associate
195
196 end if
197
198 select type (ft => this%file_%file_type)
199 ! Only fld files have the option to write the mesh at command
200 type is (fld_file_t)
201 ft%write_mesh = this%always_write_mesh
202 call ft%write(this%fluid, t)
203 class default
204 call ft%write(this%fluid, t)
205 end select
206
207 end subroutine fluid_output_sample
208
209end module fluid_output
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
Defines a field.
Definition field.f90:34
NEKTON fld file format.
Definition fld_file.f90:35
Defines an output for a fluid.
subroutine fluid_output_sample(this, t)
Sample a fluid solution at time t.
subroutine fluid_output_free(this)
Destroy a fluid output list.
subroutine fluid_output_init(this, precision, fluid, scalar_fields, name, path, fmt, layout, always_write_mesh)
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
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:149
Contains the scalar_scheme_t type.
Contains the scalars_t type that manages multiple scalar fields.
Definition scalars.f90:35
field_list_t, To be able to group fields together
Interface for NEKTON fld files.
Definition fld_file.f90:64
Base type of all fluid formulations.
Abstract type defining an output type.
Definition output.f90:41
Base type for a scalar advection-diffusion solver.
Type to manage multiple scalar transport equations.
Definition scalars.f90:61