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