Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fld_file_data.f90
Go to the documentation of this file.
1
9 use num_types, only : rp
10 use math, only : cmult, add2
11 use vector, only : vector_t, vector_ptr_t
13 use field, only : field_t
14 use field_list, only : field_list_t
15 use logger, only : neko_log, log_size
16 use device, only : host_to_device
17 use dofmap, only : dofmap_t
18 use space, only : space_t, gll
21 use mesh, only : mesh_t
22 implicit none
23 private
24
25 type, public :: fld_file_data_t
26 type(vector_t) :: x
27 type(vector_t) :: y
28 type(vector_t) :: z
29 type(vector_t) :: u
30 type(vector_t) :: v
31 type(vector_t) :: w
32 type(vector_t) :: p
33 type(vector_t) :: t
34 integer, allocatable :: idx(:)
35 type(vector_t), allocatable :: s(:)
36 integer :: gdim
37 integer :: n_scalars = 0
38 real(kind=rp) :: time = 0.0
39 integer :: glb_nelv = 0
40 integer :: nelv = 0
41 integer :: offset_el = 0
42 integer :: lx = 0
43 integer :: ly = 0
44 integer :: lz = 0
45 integer :: t_counter = 0
46 ! meta file information (if any)
48 integer :: meta_nsamples = 0
49 integer :: meta_start_counter = 0
51 character(len=1024) :: fld_series_fname
52 contains
53 procedure, pass(this) :: init => fld_file_data_init
54 procedure, pass(this) :: free => fld_file_data_free
55 procedure, pass(this) :: scale => fld_file_data_scale
56 procedure, pass(this) :: add => fld_file_data_add
57 procedure, pass(this) :: size => fld_file_data_size
58 procedure, pass(this) :: get_list => fld_file_data_get_list
59 procedure, pass(this) :: init_same => fld_file_data_init_same
60 procedure, pass(this) :: init_n_fields => fld_file_data_init_n_fields
62 procedure, pass(this) :: generate_interpolator => &
65 procedure, pass(this) :: import_fields => fld_file_data_import_fields
66 end type fld_file_data_t
67
68contains
69
99 subroutine fld_file_data_import_fields(this, u, v, w, p, t, &
100 s_target_list, s_index_list, interpolate, tolerance)
101 class(fld_file_data_t), intent(inout) :: this
102 type(field_t), pointer, intent(inout), optional :: u,v,w,p,t
103 type(field_list_t), intent(inout), optional :: s_target_list
104 integer, intent(in), optional :: s_index_list(:)
105 logical, intent(in), optional :: interpolate
106 real(kind=rp), intent(in), optional :: tolerance
107
108 integer :: i
109
110 ! ---- For the mesh to mesh interpolation
111 logical :: mesh_mismatch
112 logical :: interpolate_
113 type(global_interpolation_t) :: global_interp
114 type(dofmap_t), pointer :: dof
115 type(mesh_t) , pointer :: msh
116 ! -----
117
118 ! ---- For space to space interpolation
119 type(space_t) :: prev_Xh
120 type(space_t) , pointer :: Xh
121 type(interpolator_t) :: space_interp
122 ! ----
123
124 character(len=LOG_SIZE) :: log_buf
125
126 ! ---- Default values
127 interpolate_ = .false.
128 if (present(interpolate)) interpolate_ = interpolate
129 ! ----
130
131 !
132 ! Handle the passing of arguments and pointers
133 !
134 dof => null()
135 msh => null()
136 xh => null()
137
138 if (present(u)) then
139 dof => u%dof; msh => u%msh; xh => u%Xh
140 else if (present(v)) then
141 dof => v%dof; msh => v%msh; xh => v%Xh
142 else if (present(w)) then
143 dof => w%dof; msh => w%msh; xh => w%Xh
144 else if (present(p)) then
145 dof => p%dof; msh => p%msh; xh => p%Xh
146 else if (present(t)) then
147 dof => t%dof; msh => t%msh; xh => t%Xh
148 else if (present(s_target_list)) then
149 if (s_target_list%size() .eq. 0) then
150 call neko_error("Scalar target list is empty")
151 else
152 dof => s_target_list%items(1)%ptr%dof
153 msh => s_target_list%items(1)%ptr%msh
154 xh => s_target_list%items(1)%ptr%Xh
155 end if
156 else
157 call neko_error("At least one field must be passed")
158 end if
159
160 !
161 ! Check that the data in the fld file matches the current case.
162 ! Note that this is a safeguard and there are corner cases where
163 ! two different meshes have the same dimension and same # of elements
164 ! but this should be enough to cover obvious cases.
165 !
166 mesh_mismatch = (this%glb_nelv .ne. msh%glb_nelv .or. &
167 this%gdim .ne. msh%gdim)
168
169 if (mesh_mismatch .and. .not. interpolate_) then
170 call neko_error("The fld file must match the current mesh! &
171 &Use 'interpolate': 'true' to enable interpolation.")
172 else if (.not. mesh_mismatch .and. interpolate_) then
173 call neko_log%warning("You have activated interpolation but you might &
174 &still be using the same mesh.")
175 end if
176
177 !
178 ! Handling of interpolation and I/O
179 !
180 if (interpolate_) then
181
182 ! Throw error if dof or msh are not specified
183 ! This should never be thrown, but just in case.
184 if (.not. associated(dof) .or. .not. associated(msh)) then
185 call neko_error("both dof and msh must be associated")
186 end if
187
188 ! Generates an interpolator object and performs the point search
189 call this%generate_interpolator(global_interp, dof, msh, &
190 tolerance = tolerance)
191
192 ! Evaluate all the fields
193 if (present(u)) call global_interp%evaluate(u%x(:,1,1,1), this%u%x, &
194 on_host = .false.)
195 if (present(v)) call global_interp%evaluate(v%x(:,1,1,1), this%v%x, &
196 on_host = .false.)
197 if (present(w)) call global_interp%evaluate(w%x(:,1,1,1), this%w%x, &
198 on_host = .false.)
199 if (present(p)) call global_interp%evaluate(p%x(:,1,1,1), this%p%x, &
200 on_host = .false.)
201 if (present(t)) call global_interp%evaluate(t%x(:,1,1,1), this%t%x, &
202 on_host = .false.)
203 if (present(s_target_list)) then
204
205 ! If the index list exists, use it as a "mask"
206 if (present(s_index_list)) then
207 do i = 1, size(s_index_list)
208 ! Take care that if we set i=0 we want temperature
209 if (s_index_list(i) .eq. 0) then
210 call global_interp%evaluate(s_target_list%x(i), &
211 this%t%x, on_host = .false.)
212 else
213 ! For scalar fields, require indices in 1:this%n_scalars
214 if (s_index_list(i) < 1 .or. &
215 s_index_list(i) > this%n_scalars) then
216 call neko_error("s_index_list entry out of bounds")
217 end if
218 call global_interp%evaluate(s_target_list%x(i), &
219 this%s(s_index_list(i))%x, on_host = .false.)
220 end if
221 end do
222
223 ! otherwise, just copy element-to-element
224 else
225 do i = 1, s_target_list%size()
226 call global_interp%evaluate(s_target_list%x(i), this%s(i)%x, &
227 on_host = .false.)
228 end do
229 end if ! present s_index_list
230 end if ! present s_tgt
231
232 call global_interp%free()
233
234 else ! No interpolation, but potentially just from different spaces
235
236 ! throw an error is the space is not passed
237 if (.not. associated(xh)) call neko_error("Xh is not associated")
238
239 ! Build a space_t object from the data in the fld file
240 call prev_xh%init(gll, this%lx, this%ly, this%lz)
241 call space_interp%init(xh, prev_xh)
242
243 ! Do the space-to-space interpolation
244 if (present(u)) call space_interp%map(u%x, this%u%x, this%nelv, xh)
245 if (present(v)) call space_interp%map(v%x, this%v%x, this%nelv, xh)
246 if (present(w)) call space_interp%map(w%x, this%w%x, this%nelv, xh)
247 if (present(p)) call space_interp%map(p%x, this%p%x, this%nelv, xh)
248 if (present(t)) call space_interp%map(t%x, this%t%x, this%nelv, xh)
249 if (present(s_target_list)) then
250
251 ! If the index list exists, use it as a "mask"
252 if (present(s_index_list)) then
253 do i = 1, size(s_index_list)
254
255 ! 0 means we want temperature
256 if (s_index_list(i) .eq. 0) then
257 call space_interp%map(s_target_list%x(i), &
258 this%t%x, this%nelv, xh)
259 else
260 call space_interp%map(s_target_list%x(i), &
261 this%s(s_index_list(i))%x, this%nelv, xh)
262 end if
263 end do
264
265 ! otherwise, just copy element-to-element
266 else
267 do i = 1, s_target_list%size()
268 call space_interp%map(s_target_list%x(i), this%s(i)%x, &
269 this%nelv, xh)
270 end do
271 end if ! present s_index_list
272 end if ! present s_tgt
273
274 call space_interp%free
275 call prev_xh%free
276
277 end if ! if interpolate
278
279 nullify(dof)
280 nullify(xh)
281 nullify(msh)
282
283 end subroutine fld_file_data_import_fields
284
288 subroutine fld_file_data_init(this, nelv, offset_el)
289 class(fld_file_data_t), intent(inout) :: this
290 integer, intent(in), optional :: nelv, offset_el
291
292 call this%free()
293 if (present(nelv)) this%nelv = nelv
294 if (present(offset_el)) this%offset_el = offset_el
295
296 end subroutine fld_file_data_init
297
299 function fld_file_data_size(this) result(i)
300 class(fld_file_data_t) :: this
301 integer :: i
302 i = 0
303 if (this%u%size() .gt. 0) i = i + 1
304 if (this%v%size() .gt. 0) i = i + 1
305 if (this%w%size() .gt. 0) i = i + 1
306 if (this%p%size() .gt. 0) i = i + 1
307 if (this%t%size() .gt. 0) i = i + 1
308 i = i + this%n_scalars
309
310 end function fld_file_data_size
311
313 subroutine fld_file_data_init_same(this, fld_file, n)
314 class(fld_file_data_t), target, intent(inout) :: this
315 class(fld_file_data_t), target, intent(in) :: fld_file
316 integer, intent(in) :: n
317 integer :: i, j
318
319 if (fld_file%u%size() .gt. 0) then
320 call this%u%init(n)
321 end if
322 if (fld_file%v%size() .gt. 0) then
323 call this%v%init(n)
324 end if
325 if (fld_file%w%size() .gt. 0) then
326 call this%w%init(n)
327 end if
328 if (fld_file%p%size() .gt. 0) then
329 call this%p%init(n)
330 end if
331 if (fld_file%t%size() .gt. 0) then
332 call this%t%init(n)
333 end if
334 this%n_scalars = fld_file%n_scalars
335 allocate(this%s(fld_file%n_scalars))
336 do j = 1, fld_file%n_scalars
337 call this%s(j)%init(n)
338 end do
339
340 end subroutine fld_file_data_init_same
341
343 subroutine fld_file_data_init_n_fields(this, n_fields, n)
344 class(fld_file_data_t), target, intent(inout) :: this
345 integer, intent(in) :: n, n_fields
346 integer :: i, j
347
348
349 if (n_fields .gt. 0) then
350 call this%u%init(n)
351 end if
352 if (n_fields .gt. 1) then
353 call this%v%init(n)
354 end if
355 if (n_fields .gt. 2) then
356 call this%w%init(n)
357 end if
358 if (n_fields .gt. 3) then
359 call this%p%init(n)
360 end if
361 if (n_fields .gt. 4) then
362 call this%t%init(n)
363 end if
364 if (n_fields .gt. 5) then
365 this%n_scalars = n_fields-5
366 allocate(this%s(this%n_scalars))
367 do j = 1, this%n_scalars
368 call this%s(j)%init(n)
369 end do
370 end if
371
372 end subroutine fld_file_data_init_n_fields
373
375 subroutine fld_file_data_get_list(this, ptr_list, n)
376 class(fld_file_data_t), target, intent(in) :: this
377 integer, intent(in) :: n
378 integer :: i, j
379 type(vector_ptr_t), intent(inout) :: ptr_list(n)
380 i = 1
381 if (this%u%size() .gt. 0) then
382 ptr_list(i)%ptr => this%u
383 i = i + 1
384 end if
385 if (this%v%size() .gt. 0) then
386 ptr_list(i)%ptr => this%v
387 i = i + 1
388 end if
389 if (this%w%size() .gt. 0) then
390 ptr_list(i)%ptr => this%w
391 i = i + 1
392 end if
393 if (this%p%size() .gt. 0) then
394 ptr_list(i)%ptr => this%p
395 i = i + 1
396 end if
397 if (this%t%size() .gt. 0) then
398 ptr_list(i)%ptr => this%t
399 i = i + 1
400 end if
401 do j = 1, this%n_scalars
402 ptr_list(i)%ptr => this%s(j)
403 i = i +1
404 end do
405
406 end subroutine fld_file_data_get_list
407
408
409
411 subroutine fld_file_data_scale(this, c)
412 class(fld_file_data_t), intent(inout) :: this
413 real(kind=rp), intent(in) :: c
414 integer :: i
415
416 if (this%u%size() .gt. 0) call cmult(this%u%x, c, this%u%size())
417 if (this%v%size() .gt. 0) call cmult(this%v%x, c, this%v%size())
418 if (this%w%size() .gt. 0) call cmult(this%w%x, c, this%w%size())
419 if (this%p%size() .gt. 0) call cmult(this%p%x, c, this%p%size())
420 if (this%t%size() .gt. 0) call cmult(this%t%x, c, this%t%size())
421
422 do i = 1, this%n_scalars
423 if (this%s(i)%size() .gt. 0) call cmult(this%s(i)%x, c, this%s(i)%size())
424 end do
425
426 end subroutine fld_file_data_scale
427
429 subroutine fld_file_data_add(this, other)
430 class(fld_file_data_t), intent(inout) :: this
431 class(fld_file_data_t), intent(in) :: other
432 integer :: i, n
433
434 if (this%u%size() .gt. 0) call add2(this%u%x, other%u%x, this%u%size())
435 if (this%v%size() .gt. 0) call add2(this%v%x, other%v%x, this%v%size())
436 if (this%w%size() .gt. 0) call add2(this%w%x, other%w%x, this%w%size())
437 if (this%p%size() .gt. 0) call add2(this%p%x, other%p%x, this%p%size())
438 if (this%t%size() .gt. 0) call add2(this%t%x, other%t%x, this%t%size())
439
440 do i = 1, this%n_scalars
441 if (this%s(i)%size() .gt. 0) call add2(this%s(i)%x, other%s(i)%x, &
442 this%s(i)%size())
443 end do
444 end subroutine fld_file_data_add
445
447 subroutine fld_file_data_free(this)
448 class(fld_file_data_t), intent(inout) :: this
449 integer :: i
450 call this%x%free()
451 call this%y%free()
452 call this%z%free()
453 call this%u%free()
454 call this%v%free()
455 call this%w%free()
456 call this%p%free()
457 call this%t%free()
458 if (allocated(this%s)) then
459 do i = 1, this%n_scalars
460 call this%s(i)%free()
461 end do
462 deallocate(this%s)
463 end if
464 this%n_scalars = 0
465 this%time = 0.0
466 this%glb_nelv = 0
467 this%nelv = 0
468 this%offset_el = 0
469 this%lx = 0
470 this%ly = 0
471 this%lz = 0
472 this%t_counter = 0
473 this%meta_nsamples = 0
474 this%meta_start_counter = 0
475 if (allocated(this%idx)) deallocate(this%idx)
476 end subroutine fld_file_data_free
477
484 subroutine fld_file_data_generate_interpolator(this, global_interp, to_dof, &
485 to_msh, tolerance)
486 class(fld_file_data_t), intent(in) :: this
487 type(global_interpolation_t), intent(inout) :: global_interp
488 type(dofmap_t), intent(in), target :: to_dof
489 type(mesh_t), intent(in), target :: to_msh
490 real(kind=rp), intent(in) :: tolerance
491
492 ! --- variables for interpolation
493 type(space_t) :: fld_xh
494 real(kind=rp), allocatable :: x_coords(:,:,:,:), y_coords(:,:,:,:), &
495 z_coords(:,:,:,:)
496 real(kind=rp) :: center_x, center_y, center_z
497 integer :: e, i
498 ! ---
499
500 type(space_t), pointer :: to_Xh
501 to_xh => to_dof%Xh
502
503 ! Safeguard in case we didn't read mesh information
504 if (.not. allocated(this%x%x) .or. &
505 .not. allocated(this%y%x) .or. &
506 .not. allocated(this%z%x)) call neko_error("Unable to retrieve &
507 &mesh information from fld data.")
508
509 ! Create a space based on the fld data
510 call fld_xh%init(gll, this%lx, this%ly, this%lz)
511
512 ! These are the coordinates of our current dofmap
513 ! that we use for the interpolation
514 allocate(x_coords(to_xh%lx, to_xh%ly, to_xh%lz, to_msh%nelv))
515 allocate(y_coords(to_xh%lx, to_xh%ly, to_xh%lz, to_msh%nelv))
516 allocate(z_coords(to_xh%lx, to_xh%ly, to_xh%lz, to_msh%nelv))
517
522 do e = 1, to_msh%nelv
523 center_x = 0d0
524 center_y = 0d0
525 center_z = 0d0
526 do i = 1, to_xh%lxyz
527 center_x = center_x + to_dof%x(i, 1, 1, e)
528 center_y = center_y + to_dof%y(i, 1, 1, e)
529 center_z = center_z + to_dof%z(i, 1, 1, e)
530 end do
531 center_x = center_x / to_xh%lxyz
532 center_y = center_y / to_xh%lxyz
533 center_z = center_z / to_xh%lxyz
534 do i = 1, to_xh%lxyz
535 x_coords(i, 1, 1, e) = to_dof%x(i, 1, 1, e) - &
536 tolerance * (to_dof%x(i, 1, 1, e) - center_x)
537 y_coords(i, 1, 1, e) = to_dof%y(i, 1, 1, e) - &
538 tolerance * (to_dof%y(i, 1, 1, e) - center_y)
539 z_coords(i, 1, 1, e) = to_dof%z(i, 1, 1, e) - &
540 tolerance * (to_dof%z(i, 1, 1, e) - center_z)
541 end do
542 end do
543
544 ! The initialization is done based on the variables created from
545 ! fld data
546 call global_interp%init(this%x%x, this%y%x, this%z%x, this%gdim, &
547 this%nelv, fld_xh, tol = tolerance)
548 call global_interp%find_points(x_coords, y_coords, z_coords, &
549 to_dof%size())
550
551 deallocate(x_coords)
552 deallocate(y_coords)
553 deallocate(z_coords)
554 call fld_xh%free()
555
557
558end module fld_file_data
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
integer function fld_file_data_size(this)
Get the number of initialized fields in this fld file.
subroutine fld_file_data_init_same(this, fld_file, n)
Genereate same fields as in another fld_file.
subroutine fld_file_data_get_list(this, ptr_list, n)
Get a list with pointers to the fields in the fld file.
subroutine fld_file_data_free(this)
Deallocate fld file data type.
subroutine fld_file_data_init(this, nelv, offset_el)
Initializes a fld_file_data object.
subroutine fld_file_data_init_n_fields(this, n_fields, n)
Generate same fields as in another fld_file.
subroutine fld_file_data_scale(this, c)
Scale the values stored in this fld_file_data.
subroutine fld_file_data_generate_interpolator(this, global_interp, to_dof, to_msh, tolerance)
Generates a global_interpolation object to interpolate the fld data.
subroutine fld_file_data_add(this, other)
Add the values in another fld file to this.
subroutine fld_file_data_import_fields(this, u, v, w, p, t, s_target_list, s_index_list, interpolate, tolerance)
Imports fields from an fld_file_data object, potentially with interpolation.
NEKTON fld file format.
Definition fld_file.f90:35
Implements global_interpolation given a dofmap.
Routines to interpolate between different spaces.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:411
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
Utilities.
Definition utils.f90:35
integer function, public extract_fld_file_index(fld_filename, default_index)
Extracts the index of a field file. For example, "myfield.f00045" will return 45. If the suffix of th...
Definition utils.f90:157
integer, parameter, public neko_fname_len
Definition utils.f90:40
Defines a vector.
Definition vector.f90:34
field_list_t, To be able to group fields together
Implements global interpolation for arbitrary points in the domain.
Interpolation between two space::space_t.
The function space for the SEM solution fields.
Definition space.f90:63