56 real(kind=
rp) :: start_time = 0.0_rp
58 integer :: n_fields = 0
64 integer :: output_dim = 0
86 coef, avg_dir, name, path)
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
98 this%start_time = start_time
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'&
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'
112 fname =
'user_stats.fld'
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)
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'
131 fname =
'user_stats.csv'
133 call this%map_1d%init_char(coef, avg_dir, 1e-7_rp)
137 call this%init_base(fname)
139 call this%fields%init(n_fields)
140 this%n_fields = n_fields
141 this%mean_fields => mean_fields
143 this%fields%items(i)%ptr => this%mean_fields(i)%mf
152 call this%free_base()
154 nullify(this%mean_fields)
155 call this%map_1d%free()
156 call this%map_2d%free()
157 call this%fields%free()
164 real(kind=
rp),
intent(in) :: t
168 real(kind=
rp) :: u, v, w, p
170 associate(out_fields => this%fields%items)
171 if (t .ge. this%start_time)
then
173 do i = 1,
size(out_fields)
174 call device_memcpy(out_fields(i)%ptr%x, out_fields(i)%ptr%x_d,&
176 sync = (i .eq.
size(out_fields)))
179 if (this%output_dim .eq. 1)
then
180 call this%map_1d%average_planes(avg_output_1d, &
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)
187 call this%file_%write(this%fields, t)
189 do i = 1, this%n_fields
190 call this%mean_fields(i)%reset()
Copy data between host and device (or device and device)
Device abstraction, common interface for various accelerators.
integer, parameter, public device_to_host
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.
Maps a 3D dofmap to a 2D spectral element grid.
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.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
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...
Computes the temporal mean of a field.
Output for a list of mean fields.
Abstract type defining an output type.