Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fld_file.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!
36 use num_types, only: rp, dp, sp, i8
38 use field, only: field_t
39 use field_list, only: field_list_t
40 use dofmap, only: dofmap_t
41 use space, only: space_t
42 use structs, only: array_ptr_t
43 use vector, only: vector_t
45 use mean_flow, only: mean_flow_t
47 use vector, only : vector_t
48 use space, only : space_t
49 use mesh, only : mesh_t
51 use utils, only: neko_error
52 use comm
53 use datadist, only: linear_dist_t
54 use math, only: vlmin, vlmax
57 use mpi_f08
58 implicit none
59 private
60
61 real(kind=dp), private, allocatable :: tmp_dp(:)
62 real(kind=sp), private, allocatable :: tmp_sp(:)
63
65 type, public, extends(generic_file_t) :: fld_file_t
66 logical :: dp_precision = .false.
67 contains
68 procedure :: read => fld_file_read
69 procedure :: write => fld_file_write
70 procedure :: set_precision => fld_file_set_precision
71 end type fld_file_t
72
73
74contains
75
78 subroutine fld_file_write(this, data, t)
79 class(fld_file_t), intent(inout) :: this
80 class(*), target, intent(in) :: data
81 real(kind=rp), intent(in), optional :: t
82 type(array_ptr_t) :: x, y, z, u, v, w, p, tem
83 real(kind=rp), allocatable, target :: tempo(:)
84 type(mesh_t), pointer :: msh
85 type(dofmap_t), pointer :: dof
86 type(space_t), pointer :: Xh
87 real(kind=dp) :: time
88 character(len= 132) :: hdr
89 character :: rdcode(10)
90 character(len=6) :: id_str
91 character(len= 1024) :: fname
92 character(len= 1024) :: start_field
93 integer :: file_unit
94 integer :: i, ierr, n, suffix_pos, tslash_pos
95 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
96 integer, allocatable :: idx(:)
97 type(mpi_status) :: status
98 type(mpi_file) :: fh
99 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset, temp_offset
100 real(kind=sp), parameter :: test_pattern = 6.54321
101 type(array_ptr_t), allocatable :: scalar_fields(:)
102 logical :: write_mesh, write_velocity, write_pressure, write_temperature
103 integer :: fld_data_size, n_scalar_fields
104
105 if (present(t)) then
106 time = real(t, dp)
107 else
108 time = 0d0
109 end if
110
111 nullify(msh)
112 nullify(dof)
113 nullify(xh)
114 n_scalar_fields = 0
115 write_pressure = .false.
116 write_velocity = .false.
117 write_temperature = .false.
118
119 select type (data)
120 type is (fld_file_data_t)
121 nelv = data%nelv
122 lx = data%lx
123 ly = data%ly
124 lz = data%lz
125 gdim = data%gdim
126 glb_nelv = data%glb_nelv
127 offset_el = data%offset_el
128
129 if (data%x%size() .gt. 0) x%ptr => data%x%x
130 if (data%y%size() .gt. 0) y%ptr => data%y%x
131 if (data%z%size() .gt. 0) z%ptr => data%z%x
132 if (gdim .eq. 2) z%ptr => data%y%x
133 if (data%u%size() .gt. 0) then
134 u%ptr => data%u%x
135 write_velocity = .true.
136 end if
137 if (data%v%size() .gt. 0) v%ptr => data%v%x
138 if (data%w%size() .gt. 0) w%ptr => data%w%x
139 if (data%p%size() .gt. 0) then
140 p%ptr => data%p%x
141 write_pressure = .true.
142 end if
143 if (data%t%size() .gt. 0) then
144 write_temperature = .true.
145 tem%ptr => data%t%x
146 end if
147 ! If gdim = 2 and Z-velocity component exists,
148 ! it is stored in last scalar field
149 if (gdim .eq. 2 .and. data%w%size() .gt. 0) then
150 n_scalar_fields = data%n_scalars + 1
151 allocate(scalar_fields(n_scalar_fields))
152 do i = 1, n_scalar_fields -1
153 scalar_fields(i)%ptr => data%s(i)%x
154 end do
155 scalar_fields(n_scalar_fields)%ptr => data%w%x
156 else
157 n_scalar_fields = data%n_scalars
158 allocate(scalar_fields(n_scalar_fields+1))
159 do i = 1, n_scalar_fields
160 scalar_fields(i)%ptr => data%s(i)%x
161 end do
162 scalar_fields(n_scalar_fields+1)%ptr => data%w%x
163 end if
164 ! This is very stupid...
165 ! Some compilers cannot handle that these pointers dont point to anything
166 ! (although they are not used) this fixes a segfault due to this.
167 if (nelv .eq. 0) then
168 allocate(tempo(1))
169 x%ptr => tempo
170 y%ptr => tempo
171 z%ptr => tempo
172 u%ptr => tempo
173 v%ptr => tempo
174 w%ptr => tempo
175 p%ptr => tempo
176 tem%ptr => tempo
177 end if
178
179 allocate(idx(nelv))
180 do i = 1, nelv
181 idx(i) = data%idx(i)
182 end do
183 type is (field_t)
184 p%ptr => data%x(:,1,1,1)
185 dof => data%dof
186 write_pressure = .true.
187 write_velocity = .false.
188 type is (field_list_t)
189 select case (data%size())
190 case (1)
191 p%ptr => data%items(1)%ptr%x(:,1,1,1)
192 write_pressure = .true.
193 write_velocity = .false.
194 case (2)
195 p%ptr => data%items(1)%ptr%x(:,1,1,1)
196 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
197 write_pressure = .true.
198 write_temperature = .true.
199 case (3)
200 u%ptr => data%items(1)%ptr%x(:,1,1,1)
201 v%ptr => data%items(2)%ptr%x(:,1,1,1)
202 w%ptr => data%items(3)%ptr%x(:,1,1,1)
203 write_velocity = .true.
204 case (4)
205 p%ptr => data%items(1)%ptr%x(:,1,1,1)
206 u%ptr => data%items(2)%ptr%x(:,1,1,1)
207 v%ptr => data%items(3)%ptr%x(:,1,1,1)
208 w%ptr => data%items(4)%ptr%x(:,1,1,1)
209 write_pressure = .true.
210 write_velocity = .true.
211 case (5:99)
212 p%ptr => data%items(1)%ptr%x(:,1,1,1)
213 u%ptr => data%items(2)%ptr%x(:,1,1,1)
214 v%ptr => data%items(3)%ptr%x(:,1,1,1)
215 w%ptr => data%items(4)%ptr%x(:,1,1,1)
216 ! Check if position 5 is a temperature field by name
217 if (trim(data%name(5)) .eq. 'temperature') then
218 ! Position 5 is temperature, remaining fields are scalars
219 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
220 n_scalar_fields = data%size() - 5
221 allocate(scalar_fields(n_scalar_fields))
222 do i = 1, n_scalar_fields
223 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
224 end do
225 write_temperature = .true.
226 else
227 ! All remaining fields are scalars (no temperature field)
228 n_scalar_fields = data%size() - 4
229 allocate(scalar_fields(n_scalar_fields))
230 do i = 1, n_scalar_fields
231 scalar_fields(i)%ptr => data%items(i+4)%ptr%x(:,1,1,1)
232 end do
233 write_temperature = .false.
234 end if
235 write_pressure = .true.
236 write_velocity = .true.
237 case default
238 call neko_error('This many fields not supported yet, fld_file')
239 end select
240 dof => data%dof(1)
241
242 type is (mean_flow_t)
243 u%ptr => data%u%mf%x(:,1,1,1)
244 v%ptr => data%v%mf%x(:,1,1,1)
245 w%ptr => data%w%mf%x(:,1,1,1)
246 p%ptr => data%p%mf%x(:,1,1,1)
247 dof => data%u%mf%dof
248 write_pressure = .true.
249 write_velocity = .true.
250 type is (mean_sqr_flow_t)
251 u%ptr => data%uu%mf%x(:,1,1,1)
252 v%ptr => data%vv%mf%x(:,1,1,1)
253 w%ptr => data%ww%mf%x(:,1,1,1)
254 p%ptr => data%pp%mf%x(:,1,1,1)
255 dof => data%pp%mf%dof
256 write_pressure = .true.
257 write_velocity = .true.
258 class default
259 call neko_error('Invalid data')
260 end select
261 ! Fix things for pointers that do not exist in all data types...
262 if (associated(dof)) then
263 x%ptr => dof%x(:,1,1,1)
264 y%ptr => dof%y(:,1,1,1)
265 z%ptr => dof%z(:,1,1,1)
266 msh => dof%msh
267 xh => dof%Xh
268 end if
269
270 if (associated(msh)) then
271 nelv = msh%nelv
272 glb_nelv = msh%glb_nelv
273 offset_el = msh%offset_el
274 gdim = msh%gdim
275 ! Store global idx of each element
276 allocate(idx(msh%nelv))
277 do i = 1, msh%nelv
278 idx(i) = msh%elements(i)%e%id()
279 end do
280 end if
281
282 if (associated(xh)) then
283 lx = xh%lx
284 ly = xh%ly
285 lz = xh%lz
286 end if
287
288 lxyz = lx*ly*lz
289 n = nelv*lxyz
290
291 if (this%dp_precision) then
292 fld_data_size = mpi_double_precision_size
293 else
294 fld_data_size = mpi_real_size
295 end if
296 if (this%dp_precision) then
297 allocate(tmp_dp(gdim*n))
298 else
299 allocate(tmp_sp(gdim*n))
300 end if
301
302
303 !
304 ! Create fld header for NEKTON's multifile output
305 !
306
307 write_mesh = (this%counter .eq. this%start_counter)
308 call mpi_allreduce(mpi_in_place, write_mesh, 1, &
309 mpi_logical, mpi_lor, neko_comm)
310 call mpi_allreduce(mpi_in_place, write_velocity, 1, &
311 mpi_logical, mpi_lor, neko_comm)
312 call mpi_allreduce(mpi_in_place, write_pressure, 1, &
313 mpi_logical, mpi_lor, neko_comm)
314 call mpi_allreduce(mpi_in_place, write_temperature, 1, &
315 mpi_logical, mpi_lor, neko_comm)
316 call mpi_allreduce(mpi_in_place, n_scalar_fields, 1, &
317 mpi_integer, mpi_max, neko_comm)
318
319 ! Build rdcode note that for field_t, we only support scalar
320 ! fields at the moment
321 rdcode = ' '
322 i = 1
323 if (write_mesh) then
324 rdcode(i) = 'X'
325 i = i + 1
326 end if
327 if (write_velocity) then
328 rdcode(i) = 'U'
329 i = i + 1
330 end if
331 if (write_pressure) then
332 rdcode(i) = 'P'
333 i = i + 1
334 end if
335 if (write_temperature) then
336 rdcode(i) = 'T'
337 i = i + 1
338 end if
339 if (n_scalar_fields .gt. 0 ) then
340 rdcode(i) = 'S'
341 i = i + 1
342 write(rdcode(i), '(i1)') (n_scalar_fields)/10
343 i = i + 1
344 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
345 i = i + 1
346 end if
347
349 write(hdr, 1) fld_data_size, lx, ly, lz, glb_nelv, glb_nelv,&
350 time, this%counter, 1, 1, (rdcode(i),i = 1, 10)
3511 format('#std', 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
352 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
353
354 ! Change to NEKTON's fld file format
355 suffix_pos = filename_suffix_pos(this%fname)
356 write(id_str, '(a,i5.5)') 'f', this%counter
357 fname = trim(this%fname(1:suffix_pos-1))//'0.'//id_str
358
359 call mpi_file_open(neko_comm, trim(fname), &
360 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, &
361 ierr)
362
363 call mpi_file_write_all(fh, hdr, 132, mpi_character, status, ierr)
364 mpi_offset = 132 * mpi_character_size
365
366 call mpi_file_write_all(fh, test_pattern, 1, mpi_real, status, ierr)
367 mpi_offset = mpi_offset + mpi_real_size
368
369 byte_offset = mpi_offset + &
370 int(offset_el, i8) * int(mpi_integer_size, i8)
371 call mpi_file_write_at_all(fh, byte_offset, idx, nelv, &
372 mpi_integer, status, ierr)
373 mpi_offset = mpi_offset + int(glb_nelv, i8) * int(mpi_integer_size, i8)
374 deallocate(idx)
375 if (write_mesh) then
376
377 byte_offset = mpi_offset + int(offset_el, i8) * &
378 (int(gdim*lxyz, i8) * &
379 int(fld_data_size, i8))
380 call fld_file_write_vector_field(this, fh, byte_offset, &
381 x%ptr, y%ptr, z%ptr, &
382 n, gdim, lxyz, nelv)
383 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
384 (int(gdim *lxyz, i8) * &
385 int(fld_data_size, i8))
386 end if
387 if (write_velocity) then
388 byte_offset = mpi_offset + int(offset_el, i8) * &
389 (int(gdim * (lxyz), i8) * int(fld_data_size, i8))
390 call fld_file_write_vector_field(this, fh, byte_offset, &
391 u%ptr, v%ptr, w%ptr, n, gdim, lxyz, nelv)
392
393 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
394 (int(gdim * (lxyz), i8) * &
395 int(fld_data_size, i8))
396
397 end if
398
399 if (write_pressure) then
400 byte_offset = mpi_offset + int(offset_el, i8) * &
401 (int((lxyz), i8) * int(fld_data_size, i8))
402 call fld_file_write_field(this, fh, byte_offset, p%ptr, n)
403 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
404 (int((lxyz), i8) * int(fld_data_size, i8))
405 end if
406
407 if (write_temperature) then
408 byte_offset = mpi_offset + int(offset_el, i8) * &
409 (int((lxyz), i8) * &
410 int(fld_data_size, i8))
411 call fld_file_write_field(this, fh, byte_offset, tem%ptr, n)
412 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
413 (int((lxyz), i8) * &
414 int(fld_data_size, i8))
415 end if
416
417 temp_offset = mpi_offset
418
419 do i = 1, n_scalar_fields
420 ! Without this redundant if statement Cray optimizes this loop to
421 ! Oblivion
422 if (i .eq. 2) then
423 mpi_offset = int(temp_offset, i8) + int(1_i8*glb_nelv, i8) * &
424 (int(lxyz, i8) * int(fld_data_size, i8))
425 end if
426 byte_offset = int(mpi_offset, i8) + int(offset_el, i8) * &
427 (int((lxyz), i8) * &
428 int(fld_data_size, i8))
429 call fld_file_write_field(this, fh, byte_offset, scalar_fields(i)%ptr, n)
430 mpi_offset = int(mpi_offset, i8) + int(glb_nelv, i8) * &
431 (int(lxyz, i8) * &
432 int(fld_data_size, i8))
433 end do
434
435 if (gdim .eq. 3) then
436
438 if (write_mesh) then
439 !The offset is:
440 ! mpioff + element_off * 2(min max value)
441 ! * 4(single precision) * gdim(dimensions)
442 byte_offset = int(mpi_offset, i8) + &
443 int(offset_el, i8) * &
444 int(2, i8) * &
445 int(mpi_real_size, i8) * &
446 int(gdim, i8)
447 call fld_file_write_metadata_vector(this, fh, byte_offset, &
448 x%ptr, y%ptr, z%ptr, gdim, lxyz, nelv)
449 mpi_offset = int(mpi_offset, i8) + &
450 int(glb_nelv, i8) * &
451 int(2, i8) * &
452 int(mpi_real_size, i8) * &
453 int(gdim, i8)
454 end if
455
456 if (write_velocity) then
457 byte_offset = int(mpi_offset, i8) + &
458 int(offset_el, i8) * &
459 int(2, i8) * &
460 int(mpi_real_size, i8) * &
461 int(gdim, i8)
462 call fld_file_write_metadata_vector(this, fh, byte_offset, &
463 u%ptr, v%ptr, w%ptr, gdim, lxyz, nelv)
464 mpi_offset = int(mpi_offset, i8) + &
465 int(glb_nelv, i8) * &
466 int(2, i8) * &
467 int(mpi_real_size, i8) * &
468 int(gdim, i8)
469
470 end if
471
472 if (write_pressure) then
473 byte_offset = int(mpi_offset, i8) + &
474 int(offset_el, i8) * &
475 int(2, i8) * &
476 int(mpi_real_size, i8)
477 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
478 p%ptr, lxyz, nelv)
479 mpi_offset = int(mpi_offset, i8) + &
480 int(glb_nelv, i8) * &
481 int(2, i8) * &
482 int(mpi_real_size, i8)
483
484 end if
485
486 if (write_temperature) then
487 byte_offset = int(mpi_offset, i8) + &
488 int(offset_el, i8) * &
489 int(2, i8) * &
490 int(mpi_real_size, i8)
491 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
492 tem%ptr, lxyz, nelv)
493 mpi_offset = int(mpi_offset, i8) + &
494 int(glb_nelv, i8) * &
495 int(2, i8) * &
496 int(mpi_real_size, i8)
497
498 end if
499
500
501
502 temp_offset = mpi_offset
503
504 do i = 1, n_scalar_fields
505 ! Without this redundant if statement, Cray optimizes this loop to
506 ! Oblivion
507 if (i .eq. 2) then
508 mpi_offset = int(temp_offset, i8) + &
509 int(1_i8*glb_nelv, i8) * &
510 int(2, i8) * &
511 int(mpi_real_size, i8)
512 end if
513
514 byte_offset = int(mpi_offset, i8) + &
515 int(offset_el, i8) * &
516 int(2, i8) * &
517 int(mpi_real_size, i8)
518 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
519 scalar_fields(i)%ptr, lxyz, nelv)
520 mpi_offset = int(mpi_offset, i8) + &
521 int(glb_nelv, i8) * &
522 int(2, i8) * &
523 int(mpi_real_size, i8)
524 end do
525 end if
526
527
528 call mpi_file_sync(fh, ierr)
529 call mpi_file_close(fh, ierr)
530 ! Write metadata file
531 if (pe_rank .eq. 0) then
532 tslash_pos = filename_tslash_pos(this%fname)
533 write(start_field, "(I5,A8)") this%start_counter, '.nek5000'
534 open(newunit = file_unit, &
535 file = trim(this%fname(1:suffix_pos-1)) // &
536 trim(adjustl(start_field)), status = 'replace')
537 ! The following string will specify that the files in the file series
538 ! are defined by the filename followed by a 0.
539 ! This 0 is necessary as it specifies the index of number of files
540 ! the output file is split across.
541 ! In the past, many .f files were generated for each write.
542 ! To be consistent with this the trailing 0 is still necessary today.
543 write(file_unit, fmt = '(A,A,A)') 'filetemplate: ', &
544 this%fname(tslash_pos+1:suffix_pos-1), '%01d.f%05d'
545 write(file_unit, fmt = '(A,i5)') 'firsttimestep: ', this%start_counter
546 write(file_unit, fmt = '(A,i5)') 'numtimesteps: ', &
547 (this%counter + 1) - this%start_counter
548 close(file_unit)
549 end if
550
551 this%counter = this%counter + 1
552
553 if (allocated(tmp_dp)) deallocate(tmp_dp)
554 if (allocated(tmp_sp)) deallocate(tmp_sp)
555 end subroutine fld_file_write
556
557 subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, &
558 gdim, lxyz, nelv)
559 class(fld_file_t), intent(inout) :: this
560 type(mpi_file), intent(inout) :: fh
561 integer, intent(in) :: gdim, lxyz, nelv
562 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
563 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
564 integer :: el, j, ierr, nout
565 type(mpi_status) :: status
566 real(kind=sp) :: buffer(2*gdim*nelv)
567
568 j = 1
569 do el = 1, nelv
570 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
571 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
572 buffer(j+2) = real(vlmin(y(1, el), lxyz), sp)
573 buffer(j+3) = real(vlmax(y(1, el), lxyz), sp)
574 j = j + 4
575 if (gdim .eq. 3) then
576 buffer(j+0) = real(vlmin(z(1, el), lxyz), sp)
577 buffer(j+1) = real(vlmax(z(1, el), lxyz), sp)
578 j = j + 2
579 end if
580 end do
581
582 ! write out data
583 nout = 2*gdim*nelv
584
585 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
586 mpi_real, status, ierr)
587
588 end subroutine fld_file_write_metadata_vector
589
590 subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, &
591 nelv)
592 class(fld_file_t), intent(inout) :: this
593 type(mpi_file), intent(inout) :: fh
594 integer, intent(in) :: lxyz, nelv
595 real(kind=rp), intent(in) :: x(lxyz, nelv)
596 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
597 integer :: el, j, ierr, nout
598 type(mpi_status) :: status
599 real(kind=sp) :: buffer(2*nelv)
600
601 j = 1
602 do el = 1, nelv
603 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
604 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
605 j = j + 2
606 end do
607
608 ! write out data
609 nout = 2*nelv
610
611 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
612 mpi_real, status, ierr)
613
614 end subroutine fld_file_write_metadata_scalar
615
616 subroutine fld_file_write_field(this, fh, byte_offset, p, n)
617 class(fld_file_t), intent(inout) :: this
618 type(mpi_file), intent(inout) :: fh
619 integer, intent(inout) :: n
620 real(kind=rp), intent(inout) :: p(n)
621 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
622 integer :: i, ierr
623 type(mpi_status) :: status
624
625 if ( this%dp_precision) then
626 do i = 1, n
627 tmp_dp(i) = real(p(i), dp)
628 end do
629
630 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, n, &
631 mpi_double_precision, status, ierr)
632 else
633 do i = 1, n
634 tmp_sp(i) = real(p(i), sp)
635 end do
636 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, n, &
637 mpi_real, status, ierr)
638 end if
639
640 end subroutine fld_file_write_field
641
642 subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, &
643 gdim, lxyz, nelv)
644 class(fld_file_t), intent(inout) :: this
645 type(mpi_file), intent(inout) :: fh
646 integer, intent(in) :: n, gdim, lxyz, nelv
647 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
648 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
649 integer :: i, el, j, ierr
650 type(mpi_status) :: status
651
652 if (this%dp_precision) then
653 i = 1
654 do el = 1, nelv
655 do j = 1, lxyz
656 tmp_dp(i) = real(x(j, el), dp)
657 i = i +1
658 end do
659 do j = 1, lxyz
660 tmp_dp(i) = real(y(j, el), dp)
661 i = i +1
662 end do
663 if (gdim .eq. 3) then
664 do j = 1, lxyz
665 tmp_dp(i) = real(z(j, el), dp)
666 i = i +1
667 end do
668 end if
669 end do
670 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, gdim*n, &
671 mpi_double_precision, status, ierr)
672 else
673 i = 1
674 do el = 1, nelv
675 do j = 1, lxyz
676 tmp_sp(i) = real(x(j, el), sp)
677 i = i +1
678 end do
679 do j = 1, lxyz
680 tmp_sp(i) = real(y(j, el), sp)
681 i = i +1
682 end do
683 if (gdim .eq. 3) then
684 do j = 1, lxyz
685 tmp_sp(i) = real(z(j, el), sp)
686 i = i +1
687 end do
688 end if
689 end do
690 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, gdim*n, &
691 mpi_real, status, ierr)
692 end if
693
694
695 end subroutine fld_file_write_vector_field
696
698 subroutine fld_file_read(this, data)
699 class(fld_file_t) :: this
700 class(*), target, intent(inout) :: data
701 character(len= 132) :: hdr
702 integer :: ierr, suffix_pos, i, j
703 type(mpi_file) :: fh
704 type(mpi_status) :: status
705 character(len= 1024) :: fname, meta_fname, string, path
706 logical :: meta_file, read_mesh, read_velocity, read_pressure
707 logical :: read_temp
708 character(len=6) :: id_str
709 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
710 integer :: lx, ly, lz, glb_nelv, counter, lxyz
711 integer :: FLD_DATA_SIZE, n_scalars, n
712 integer :: file_unit
713 real(kind=rp) :: time
714 real(kind=sp) :: temp
715 type(linear_dist_t) :: dist
716 real(kind=sp), parameter :: test_pattern = 6.54321
717 character :: rdcode(10), temp_str(4)
718
719 select type (data)
720 type is (fld_file_data_t)
721 call filename_chsuffix(this%fname, meta_fname, 'nek5000')
722
723 inquire(file = trim(meta_fname), exist = meta_file)
724 if (meta_file .and. data%meta_nsamples .eq. 0) then
725 if (pe_rank .eq. 0) then
726 open(newunit = file_unit, file = trim(meta_fname))
727 read(file_unit, fmt = '(A)') string
728 read(string(14:), fmt = '(A)') string
729 string = trim(string)
730 data%fld_series_fname = string(:scan(trim(string), '%')-1)
731 data%fld_series_fname = adjustl(data%fld_series_fname)
732 data%fld_series_fname = trim(data%fld_series_fname)//'0'
733 read(file_unit, fmt = '(A)') string
734 read(string(scan(string, ':')+1:), *) data%meta_start_counter
735 read(file_unit, fmt = '(A)') string
736 read(string(scan(string, ':')+1:), *) data%meta_nsamples
737
738 close(file_unit)
739 write(*,*) 'Reading meta file for fld series'
740 write(*,*) 'Name: ', trim(data%fld_series_fname)
741 write(*,*) 'Start counter: ', data%meta_start_counter, &
742 'Nsamples: ', data%meta_nsamples
743 end if
744 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
745 neko_comm, ierr)
746 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
747 neko_comm, ierr)
748 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
749 neko_comm, ierr)
750 if (this%counter .eq. 0) this%counter = data%meta_start_counter
751 end if
752
753 if (meta_file) then
754 write(id_str, '(a,i5.5)') 'f', this%counter
755 path = trim(meta_fname(1:scan(meta_fname, '/', .true. )))
756 fname = trim(path)//trim(data%fld_series_fname)//'.'//id_str
757 if (this%counter .ge. data%meta_nsamples+data%meta_start_counter) then
758 call neko_error('Trying to read more fld files than exist')
759 end if
760 else
761 suffix_pos = filename_suffix_pos(this%fname)
762 write(id_str, '(a,i5.5)') 'f', this%counter
763 fname = trim(this%fname(1:suffix_pos-1))//'.'//id_str
764 end if
765 call mpi_file_open(neko_comm, trim(fname), &
766 mpi_mode_rdonly, mpi_info_null, fh, ierr)
767
768 if (ierr .ne. 0) call neko_error("Could not read "//trim(fname))
769
770 call mpi_file_read_all(fh, hdr, 132, mpi_character, status, ierr)
771 ! This read can prorbably be done wihtout the temp variables,
772 ! temp_str, i, j
773
774 read(hdr, 1) temp_str, fld_data_size, lx, ly, lz, glb_nelv, glb_nelv, &
775 time, counter, i, j, (rdcode(i), i = 1, 10)
7761 format(4a, 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
777 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
778 if (data%nelv .eq. 0) then
779 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
780 data%nelv = dist%num_local()
781 data%offset_el = dist%start_idx()
782 end if
783 data%lx = lx
784 data%ly = ly
785 data%lz = lz
786 data%glb_nelv = glb_nelv
787 data%t_counter = counter
788 data%time = time
789 lxyz = lx * ly * lz
790 n = lxyz * data%nelv
791
792 if (lz .eq. 1) then
793 data%gdim = 2
794 else
795 data%gdim = 3
796 end if
797
798
799 if (fld_data_size .eq. mpi_double_precision_size) then
800 this%dp_precision = .true.
801 else
802 this%dp_precision = .false.
803 end if
804 if (this%dp_precision) then
805 allocate(tmp_dp(data%gdim*n))
806 else
807 allocate(tmp_sp(data%gdim*n))
808 end if
809
810
811 i = 1
812 read_mesh = .false.
813 read_velocity = .false.
814 read_pressure = .false.
815 read_temp = .false.
816 if (rdcode(i) .eq. 'X') then
817 read_mesh = .true.
818 call data%x%init(n)
819 call data%y%init(n)
820 call data%z%init(n)
821 i = i + 1
822 end if
823 if (rdcode(i) .eq. 'U') then
824 read_velocity = .true.
825 call data%u%init(n)
826 call data%v%init(n)
827 call data%w%init(n)
828 i = i + 1
829 end if
830 if (rdcode(i) .eq. 'P') then
831 read_pressure = .true.
832 call data%p%init(n)
833 i = i + 1
834 end if
835 if (rdcode(i) .eq. 'T') then
836 read_temp = .true.
837 call data%t%init(n)
838 i = i + 1
839 end if
840 n_scalars = 0
841 if (rdcode(i) .eq. 'S') then
842 i = i + 1
843 read(rdcode(i),*) n_scalars
844 n_scalars = n_scalars*10
845 i = i + 1
846 read(rdcode(i),*) j
847 n_scalars = n_scalars+j
848 i = i + 1
849 if (allocated(data%s)) then
850 if (data%n_scalars .ne. n_scalars) then
851 do j = 1, data%n_scalars
852 call data%s(j)%free()
853 end do
854 deallocate(data%s)
855 data%n_scalars = n_scalars
856 allocate(data%s(n_scalars))
857 do j = 1, data%n_scalars
858 call data%s(j)%init(n)
859 end do
860 end if
861 else
862 data%n_scalars = n_scalars
863 allocate(data%s(data%n_scalars))
864 do j = 1, data%n_scalars
865 call data%s(j)%init(n)
866 end do
867 end if
868 i = i + 1
869 end if
870
871 mpi_offset = 132 * mpi_character_size
872 call mpi_file_read_at_all(fh, mpi_offset, temp, 1, &
873 mpi_real, status, ierr)
874 if (temp .ne. test_pattern) then
875 call neko_error('Incorrect format for fld file, &
876 &test pattern does not match.')
877 end if
878 mpi_offset = mpi_offset + mpi_real_size
879
880
881 if (allocated(data%idx)) then
882 if (size(data%idx) .ne. data%nelv) then
883 deallocate(data%idx)
884 allocate(data%idx(data%nelv))
885 end if
886 else
887 allocate(data%idx(data%nelv))
888 end if
889
890 byte_offset = mpi_offset + &
891 int(data%offset_el, i8) * int(mpi_integer_size, i8)
892
893 call mpi_file_read_at_all(fh, byte_offset, data%idx, data%nelv, &
894 mpi_integer, status, ierr)
895
896 mpi_offset = mpi_offset + &
897 int(data%glb_nelv, i8) * int(mpi_integer_size, i8)
898
899 if (read_mesh) then
900 byte_offset = mpi_offset + int(data%offset_el, i8) * &
901 (int(data%gdim*lxyz, i8) * &
902 int(fld_data_size, i8))
903 call fld_file_read_vector_field(this, fh, byte_offset, &
904 data%x, data%y, data%z, data)
905 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
906 (int(data%gdim *lxyz, i8) * &
907 int(fld_data_size, i8))
908 end if
909
910 if (read_velocity) then
911 byte_offset = mpi_offset + int(data%offset_el, i8) * &
912 (int(data%gdim*lxyz, i8) * &
913 int(fld_data_size, i8))
914 call fld_file_read_vector_field(this, fh, byte_offset, &
915 data%u, data%v, data%w, data)
916 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
917 (int(data%gdim *lxyz, i8) * &
918 int(fld_data_size, i8))
919 end if
920
921 if (read_pressure) then
922 byte_offset = mpi_offset + int(data%offset_el, i8) * &
923 (int(lxyz, i8) * &
924 int(fld_data_size, i8))
925 call fld_file_read_field(this, fh, byte_offset, data%p, data)
926 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
927 (int(lxyz, i8) * &
928 int(fld_data_size, i8))
929 end if
930
931 if (read_temp) then
932 byte_offset = mpi_offset + int(data%offset_el, i8) * &
933 (int(lxyz, i8) * &
934 int(fld_data_size, i8))
935 call fld_file_read_field(this, fh, byte_offset, data%t, data)
936 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
937 (int(lxyz, i8) * &
938 int(fld_data_size, i8))
939 end if
940
941 do i = 1, n_scalars
942 byte_offset = mpi_offset + int(data%offset_el, i8) * &
943 (int(lxyz, i8) * &
944 int(fld_data_size, i8))
945 call fld_file_read_field(this, fh, byte_offset, data%s(i), data)
946 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
947 (int(lxyz, i8) * &
948 int(fld_data_size, i8))
949 end do
950
951 this%counter = this%counter + 1
952
953 if (allocated(tmp_dp)) deallocate(tmp_dp)
954 if (allocated(tmp_sp)) deallocate(tmp_sp)
955 class default
956 call neko_error('Currently we only read into fld_file_data_t, &
957 &please use that data structure instead.')
958 end select
959
960 end subroutine fld_file_read
961
962 subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
963 class(fld_file_t), intent(inout) :: this
964 type(vector_t), intent(inout) :: x
965 type(fld_file_data_t) :: fld_data
966 integer(kind=MPI_OFFSET_KIND) :: byte_offset
967 type(mpi_file) :: fh
968 type(mpi_status) :: status
969 integer :: n, ierr, lxyz, i
970
971 n = x%size()
972 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
973
974 if (this%dp_precision) then
975 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, n, &
976 mpi_double_precision, status, ierr)
977 else
978 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, n, &
979 mpi_real, status, ierr)
980 end if
981
982 if (this%dp_precision) then
983 do i = 1, n
984 x%x(i) = tmp_dp(i)
985 end do
986 else
987 do i = 1, n
988 x%x(i) = tmp_sp(i)
989 end do
990 end if
991
992
993 end subroutine fld_file_read_field
994
995
996 subroutine fld_file_read_vector_field(this, fh, byte_offset, &
997 x, y, z, fld_data)
998 class(fld_file_t), intent(inout) :: this
999 type(vector_t), intent(inout) :: x, y, z
1000 type(fld_file_data_t) :: fld_data
1001 integer(kind=MPI_OFFSET_KIND) :: byte_offset
1002 type(mpi_file) :: fh
1003 type(mpi_status) :: status
1004 integer :: n, ierr, lxyz, i, j, e, nd
1005
1006 n = x%size()
1007 nd = n*fld_data%gdim
1008 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
1009
1010 if (this%dp_precision) then
1011 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, nd, &
1012 mpi_double_precision, status, ierr)
1013 else
1014 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, nd, &
1015 mpi_real, status, ierr)
1016 end if
1017
1018
1019 if (this%dp_precision) then
1020 i = 1
1021 do e = 1, fld_data%nelv
1022 do j = 1, lxyz
1023 x%x((e-1)*lxyz+j) = tmp_dp(i)
1024 i = i +1
1025 end do
1026 do j = 1, lxyz
1027 y%x((e-1)*lxyz+j) = tmp_dp(i)
1028 i = i +1
1029 end do
1030 if (fld_data%gdim .eq. 3) then
1031 do j = 1, lxyz
1032 z%x((e-1)*lxyz+j) = tmp_dp(i)
1033 i = i +1
1034 end do
1035 end if
1036 end do
1037 else
1038 i = 1
1039 do e = 1, fld_data%nelv
1040 do j = 1, lxyz
1041 x%x((e-1)*lxyz+j) = tmp_sp(i)
1042 i = i +1
1043 end do
1044 do j = 1, lxyz
1045 y%x((e-1)*lxyz+j) = tmp_sp(i)
1046 i = i +1
1047 end do
1048 if (fld_data%gdim .eq. 3) then
1049 do j = 1, lxyz
1050 z%x((e-1)*lxyz+j) = tmp_sp(i)
1051 i = i +1
1052 end do
1053 end if
1054 end do
1055 end if
1056
1057 end subroutine fld_file_read_vector_field
1058
1059 subroutine fld_file_set_precision(this, precision)
1060 class(fld_file_t) :: this
1061 integer, intent(in) :: precision
1062
1063 if (precision .eq. dp) then
1064 this%dp_precision = .true.
1065 else if (precision .eq. sp) then
1066 this%dp_precision = .false.
1067 else
1068 call neko_error('Invalid precision')
1069 end if
1070
1071 end subroutine fld_file_set_precision
1072
1073
1074end module fld_file
double real
Generic buffer that is extended with buffers of varying rank.
Definition buffer.F90:34
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:58
integer, public pe_rank
MPI rank.
Definition comm.F90:55
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
Defines practical data distributions.
Definition datadist.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
NEKTON fld file format.
Definition fld_file.f90:35
real(kind=dp), dimension(:), allocatable, private tmp_dp
Definition fld_file.f90:61
subroutine fld_file_set_precision(this, precision)
subroutine fld_file_read_vector_field(this, fh, byte_offset, x, y, z, fld_data)
Definition fld_file.f90:998
subroutine fld_file_read(this, data)
Load a field from a NEKTON fld file.
Definition fld_file.f90:699
subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, gdim, lxyz, nelv)
Definition fld_file.f90:644
subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, gdim, lxyz, nelv)
Definition fld_file.f90:559
real(kind=sp), dimension(:), allocatable, private tmp_sp
Definition fld_file.f90:62
subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
Definition fld_file.f90:963
subroutine fld_file_write_field(this, fh, byte_offset, p, n)
Definition fld_file.f90:617
subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, nelv)
Definition fld_file.f90:592
subroutine fld_file_write(this, data, t)
Write fields to a NEKTON fld file.
Definition fld_file.f90:79
Definition math.f90:60
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:608
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:597
Defines a mean flow field.
Definition mean_flow.f90:34
Defines a mean squared flow field.
Defines a mesh.
Definition mesh.f90:34
MPI derived types.
Definition mpi_types.f90:34
integer, public mpi_double_precision_size
Size of MPI type double precision.
Definition mpi_types.f90:63
integer, public mpi_character_size
Size of MPI type character.
Definition mpi_types.f90:64
integer, public mpi_real_size
Size of MPI type real.
Definition mpi_types.f90:62
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:65
integer, parameter, public i2
Definition num_types.f90:5
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
Defines structs that are used... Dont know if we should keep it though.
Definition structs.f90:2
Utilities.
Definition utils.f90:35
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:78
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Definition utils.f90:64
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:57
Defines a vector.
Definition vector.f90:34
Load-balanced linear distribution .
Definition datadist.f90:50
field_list_t, To be able to group fields together
Interface for NEKTON fld files.
Definition fld_file.f90:65
A generic file handler.
The function space for the SEM solution fields.
Definition space.f90:62
Pointer to array.
Definition structs.f90:14