79 class(*),
target,
intent(in) :: data
80 real(kind=
rp),
intent(in),
optional :: t
82 real(kind=
rp),
allocatable,
target :: tempo(:)
83 type(
mesh_t),
pointer :: msh
87 character(len= 132) :: hdr
88 character :: rdcode(10)
89 character(len=6) :: id_str
90 character(len= 1024) :: fname
91 character(len= 1024) :: start_field
92 integer :: i, ierr, n, suffix_pos, tslash_pos
93 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
94 integer,
allocatable :: idx(:)
95 type(mpi_status) :: status
97 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset, temp_offset
98 real(kind=
sp),
parameter :: test_pattern = 6.54321
100 logical :: write_mesh, write_velocity, write_pressure, write_temperature
101 integer :: fld_data_size, n_scalar_fields
113 write_pressure = .false.
114 write_velocity = .false.
115 write_temperature = .false.
124 glb_nelv = data%glb_nelv
125 offset_el = data%offset_el
127 if (data%x%n .gt. 0) x%ptr => data%x%x
128 if (data%y%n .gt. 0) y%ptr => data%y%x
129 if (data%z%n .gt. 0) z%ptr => data%z%x
130 if (gdim .eq. 2) z%ptr => data%y%x
131 if (data%u%n .gt. 0)
then
133 write_velocity = .true.
135 if (data%v%n .gt. 0) v%ptr => data%v%x
136 if (data%w%n .gt. 0) w%ptr => data%w%x
137 if (data%p%n .gt. 0)
then
139 write_pressure = .true.
141 if (data%t%n .gt. 0)
then
142 write_temperature = .true.
147 if (gdim .eq. 2 .and. data%w%n .gt. 0)
then
148 n_scalar_fields = data%n_scalars + 1
149 allocate(scalar_fields(n_scalar_fields))
150 do i = 1, n_scalar_fields -1
151 scalar_fields(i)%ptr => data%s(i)%x
153 scalar_fields(n_scalar_fields)%ptr => data%w%x
155 n_scalar_fields = data%n_scalars
156 allocate(scalar_fields(n_scalar_fields+1))
157 do i = 1, n_scalar_fields
158 scalar_fields(i)%ptr => data%s(i)%x
160 scalar_fields(n_scalar_fields+1)%ptr => data%w%x
165 if (nelv .eq. 0)
then
182 p%ptr => data%x(:,1,1,1)
184 write_pressure = .true.
185 write_velocity = .false.
187 select case (data%size())
189 p%ptr => data%items(1)%ptr%x(:,1,1,1)
190 write_pressure = .true.
191 write_velocity = .false.
193 p%ptr => data%items(1)%ptr%x(:,1,1,1)
194 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
195 write_pressure = .true.
196 write_temperature = .true.
198 u%ptr => data%items(1)%ptr%x(:,1,1,1)
199 v%ptr => data%items(2)%ptr%x(:,1,1,1)
200 w%ptr => data%items(3)%ptr%x(:,1,1,1)
201 write_velocity = .true.
203 p%ptr => data%items(1)%ptr%x(:,1,1,1)
204 u%ptr => data%items(2)%ptr%x(:,1,1,1)
205 v%ptr => data%items(3)%ptr%x(:,1,1,1)
206 w%ptr => data%items(4)%ptr%x(:,1,1,1)
207 write_pressure = .true.
208 write_velocity = .true.
210 p%ptr => data%items(1)%ptr%x(:,1,1,1)
211 u%ptr => data%items(2)%ptr%x(:,1,1,1)
212 v%ptr => data%items(3)%ptr%x(:,1,1,1)
213 w%ptr => data%items(4)%ptr%x(:,1,1,1)
214 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
215 n_scalar_fields = data%size() - 5
216 allocate(scalar_fields(n_scalar_fields))
217 do i = 1, n_scalar_fields
218 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
220 write_pressure = .true.
221 write_velocity = .true.
222 write_temperature = .true.
224 call neko_error(
'This many fields not supported yet, fld_file')
229 u%ptr => data%u%mf%x(:,1,1,1)
230 v%ptr => data%v%mf%x(:,1,1,1)
231 w%ptr => data%w%mf%x(:,1,1,1)
232 p%ptr => data%p%mf%x(:,1,1,1)
234 write_pressure = .true.
235 write_velocity = .true.
237 u%ptr => data%uu%mf%x(:,1,1,1)
238 v%ptr => data%vv%mf%x(:,1,1,1)
239 w%ptr => data%ww%mf%x(:,1,1,1)
240 p%ptr => data%pp%mf%x(:,1,1,1)
241 dof => data%pp%mf%dof
242 write_pressure = .true.
243 write_velocity = .true.
248 if (
associated(dof))
then
249 x%ptr => dof%x(:,1,1,1)
250 y%ptr => dof%y(:,1,1,1)
251 z%ptr => dof%z(:,1,1,1)
256 if (
associated(msh))
then
258 glb_nelv = msh%glb_nelv
259 offset_el = msh%offset_el
262 allocate(idx(msh%nelv))
264 idx(i) = msh%elements(i)%e%id()
268 if (
associated(xh))
then
277 if (this%dp_precision)
then
282 if (this%dp_precision)
then
293 write_mesh = (this%counter .eq. this%start_counter)
294 call mpi_allreduce(mpi_in_place, write_mesh, 1, &
296 call mpi_allreduce(mpi_in_place, write_velocity, 1, &
298 call mpi_allreduce(mpi_in_place, write_pressure, 1, &
300 call mpi_allreduce(mpi_in_place, write_temperature, 1, &
302 call mpi_allreduce(mpi_in_place, n_scalar_fields, 1, &
313 if (write_velocity)
then
317 if (write_pressure)
then
321 if (write_temperature)
then
325 if (n_scalar_fields .gt. 0 )
then
328 write(rdcode(i),
'(i1)') (n_scalar_fields)/10
330 write(rdcode(i),
'(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
335 write(hdr, 1) fld_data_size, lx, ly, lz, glb_nelv, glb_nelv,&
336 time, this%counter, 1, 1, (rdcode(i),i = 1, 10)
3371
format(
'#std', 1x, i1, 1x,
i2, 1x,
i2, 1x,
i2, 1x, i10, 1x, i10, &
338 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
342 write(id_str,
'(a,i5.5)')
'f', this%counter
343 fname = trim(this%fname(1:suffix_pos-1))//
'0.'//id_str
345 call mpi_file_open(
neko_comm, trim(fname), &
346 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, &
349 call mpi_file_write_all(fh, hdr, 132, mpi_character, status, ierr)
352 call mpi_file_write_all(fh, test_pattern, 1, mpi_real, status, ierr)
355 byte_offset = mpi_offset + &
357 call mpi_file_write_at_all(fh, byte_offset, idx, nelv, &
358 mpi_integer, status, ierr)
363 byte_offset = mpi_offset + int(offset_el,
i8) * &
364 (int(gdim*lxyz,
i8) * &
365 int(fld_data_size,
i8))
367 x%ptr, y%ptr, z%ptr, &
369 mpi_offset = mpi_offset + int(glb_nelv,
i8) * &
370 (int(gdim *lxyz,
i8) * &
371 int(fld_data_size,
i8))
373 if (write_velocity)
then
374 byte_offset = mpi_offset + int(offset_el,
i8) * &
375 (int(gdim * (lxyz),
i8) * int(fld_data_size,
i8))
377 u%ptr, v%ptr, w%ptr, n, gdim, lxyz, nelv)
379 mpi_offset = mpi_offset + int(glb_nelv,
i8) * &
380 (int(gdim * (lxyz),
i8) * &
381 int(fld_data_size,
i8))
385 if (write_pressure)
then
386 byte_offset = mpi_offset + int(offset_el,
i8) * &
387 (int((lxyz),
i8) * int(fld_data_size,
i8))
389 mpi_offset = mpi_offset + int(glb_nelv,
i8) * &
390 (int((lxyz),
i8) * int(fld_data_size,
i8))
393 if (write_temperature)
then
394 byte_offset = mpi_offset + int(offset_el,
i8) * &
396 int(fld_data_size,
i8))
398 mpi_offset = mpi_offset + int(glb_nelv,
i8) * &
400 int(fld_data_size,
i8))
403 temp_offset = mpi_offset
405 do i = 1, n_scalar_fields
409 mpi_offset = int(temp_offset,
i8) + int(1_i8*glb_nelv,
i8) * &
410 (int(lxyz,
i8) * int(fld_data_size,
i8))
412 byte_offset = int(mpi_offset,
i8) + int(offset_el,
i8) * &
414 int(fld_data_size,
i8))
416 mpi_offset = int(mpi_offset,
i8) + int(glb_nelv,
i8) * &
418 int(fld_data_size,
i8))
421 if (gdim .eq. 3)
then
428 byte_offset = int(mpi_offset,
i8) + &
429 int(offset_el,
i8) * &
434 x%ptr, y%ptr, z%ptr, gdim, lxyz, nelv)
435 mpi_offset = int(mpi_offset,
i8) + &
436 int(glb_nelv,
i8) * &
442 if (write_velocity)
then
443 byte_offset = int(mpi_offset,
i8) + &
444 int(offset_el,
i8) * &
449 u%ptr, v%ptr, w%ptr, gdim, lxyz, nelv)
450 mpi_offset = int(mpi_offset,
i8) + &
451 int(glb_nelv,
i8) * &
458 if (write_pressure)
then
459 byte_offset = int(mpi_offset,
i8) + &
460 int(offset_el,
i8) * &
465 mpi_offset = int(mpi_offset,
i8) + &
466 int(glb_nelv,
i8) * &
472 if (write_temperature)
then
473 byte_offset = int(mpi_offset,
i8) + &
474 int(offset_el,
i8) * &
479 mpi_offset = int(mpi_offset,
i8) + &
480 int(glb_nelv,
i8) * &
488 temp_offset = mpi_offset
490 do i = 1, n_scalar_fields
494 mpi_offset = int(temp_offset,
i8) + &
495 int(1_i8*glb_nelv,
i8) * &
500 byte_offset = int(mpi_offset,
i8) + &
501 int(offset_el,
i8) * &
505 scalar_fields(i)%ptr, lxyz, nelv)
506 mpi_offset = int(mpi_offset,
i8) + &
507 int(glb_nelv,
i8) * &
514 call mpi_file_sync(fh, ierr)
515 call mpi_file_close(fh, ierr)
519 write(start_field,
"(I5,A8)") this%start_counter,
'.nek5000'
521 file = trim(this%fname(1:suffix_pos-1)) // &
522 trim(adjustl(start_field)), status =
'replace')
523 write(9, fmt =
'(A,A,A)')
'filetemplate: ', &
524 this%fname(tslash_pos+1:suffix_pos-1),
'%01d.f%05d'
525 write(9, fmt =
'(A,i5)')
'firsttimestep: ', this%start_counter
526 write(9, fmt =
'(A,i5)')
'numtimesteps: ', &
527 (this%counter + 1) - this%start_counter
528 write(9, fmt =
'(A)')
'type: binary'
532 this%counter = this%counter + 1
681 class(*),
target,
intent(inout) :: data
682 character(len= 132) :: hdr
683 integer :: ierr, suffix_pos, i, j
685 type(mpi_status) :: status
686 character(len= 1024) :: fname, meta_fname, string, path
687 logical :: meta_file, read_mesh, read_velocity, read_pressure
689 character(len=6) :: id_str
690 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
691 integer :: lx, ly, lz, glb_nelv, counter, lxyz
692 integer :: FLD_DATA_SIZE, n_scalars, n
693 real(kind=
rp) :: time
694 real(kind=
sp) :: temp
696 real(kind=
sp),
parameter :: test_pattern = 6.54321
697 character :: rdcode(10), temp_str(4)
703 inquire(
file = trim(meta_fname), exist = meta_file)
704 if (meta_file .and. data%meta_nsamples .eq. 0)
then
706 open(unit = 9,
file = trim(meta_fname))
707 read(9, fmt =
'(A)') string
708 read(string(14:), fmt =
'(A)') string
709 string = trim(string)
710 data%fld_series_fname = string(:scan(trim(string),
'%')-1)
711 data%fld_series_fname = adjustl(data%fld_series_fname)
712 data%fld_series_fname = trim(data%fld_series_fname)//
'0'
713 read(9, fmt =
'(A)') string
714 read(string(scan(string,
':')+1:), *) data%meta_start_counter
715 read(9, fmt =
'(A)') string
716 read(string(scan(string,
':')+1:), *) data%meta_nsamples
719 write(*,*)
'Reading meta file for fld series'
720 write(*,*)
'Name: ', trim(data%fld_series_fname)
721 write(*,*)
'Start counter: ', data%meta_start_counter, &
722 'Nsamples: ', data%meta_nsamples
724 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
726 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
728 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
730 if (this%counter .eq. 0) this%counter = data%meta_start_counter
734 write(id_str,
'(a,i5.5)')
'f', this%counter
735 path = trim(meta_fname(1:scan(meta_fname,
'/', .true. )))
736 fname = trim(path)//trim(data%fld_series_fname)//
'.'//id_str
737 if (this%counter .ge. data%meta_nsamples+data%meta_start_counter)
then
738 call neko_error(
'Trying to read more fld files than exist')
742 write(id_str,
'(a,i5.5)')
'f', this%counter
743 fname = trim(this%fname(1:suffix_pos-1))//
'.'//id_str
745 call mpi_file_open(
neko_comm, trim(fname), &
746 mpi_mode_rdonly, mpi_info_null, fh, ierr)
748 if (ierr .ne. 0)
call neko_error(
"Could not read "//trim(fname))
750 call mpi_file_read_all(fh, hdr, 132, mpi_character, status, ierr)
754 read(hdr, 1) temp_str, fld_data_size, lx, ly, lz, glb_nelv, glb_nelv, &
755 time, counter, i, j, (rdcode(i), i = 1, 10)
7561
format(4a, 1x, i1, 1x,
i2, 1x,
i2, 1x,
i2, 1x, i10, 1x, i10, &
757 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
758 if (data%nelv .eq. 0)
then
760 data%nelv = dist%num_local()
761 data%offset_el = dist%start_idx()
766 data%glb_nelv = glb_nelv
767 data%t_counter = counter
780 this%dp_precision = .true.
782 this%dp_precision = .false.
784 if (this%dp_precision)
then
785 allocate(
tmp_dp(data%gdim*n))
787 allocate(
tmp_sp(data%gdim*n))
793 read_velocity = .false.
794 read_pressure = .false.
796 if (rdcode(i) .eq.
'X')
then
798 if (data%x%n .ne. n)
call data%x%init(n)
799 if (data%y%n .ne. n)
call data%y%init(n)
800 if (data%z%n .ne. n)
call data%z%init(n)
803 if (rdcode(i) .eq.
'U')
then
804 read_velocity = .true.
805 if (data%u%n .ne. n)
call data%u%init(n)
806 if (data%v%n .ne. n)
call data%v%init(n)
807 if (data%w%n .ne. n)
call data%w%init(n)
810 if (rdcode(i) .eq.
'P')
then
811 read_pressure = .true.
812 if (data%p%n .ne. n)
call data%p%init(n)
815 if (rdcode(i) .eq.
'T')
then
817 if (data%t%n .ne. n)
call data%t%init(n)
821 if (rdcode(i) .eq.
'S')
then
823 read(rdcode(i),*) n_scalars
824 n_scalars = n_scalars*10
827 n_scalars = n_scalars+j
829 if (
allocated(data%s))
then
830 if (data%n_scalars .ne. n_scalars)
then
831 do j = 1, data%n_scalars
832 call data%s(j)%free()
835 data%n_scalars = n_scalars
836 allocate(data%s(n_scalars))
837 do j = 1, data%n_scalars
838 call data%s(j)%init(n)
842 data%n_scalars = n_scalars
843 allocate(data%s(data%n_scalars))
844 do j = 1, data%n_scalars
845 call data%s(j)%init(n)
852 call mpi_file_read_at_all(fh, mpi_offset, temp, 1, &
853 mpi_real, status, ierr)
854 if (temp .ne. test_pattern)
then
855 call neko_error(
'Incorrect format for fld file, &
856 &test pattern does not match.')
861 if (
allocated(data%idx))
then
862 if (
size(data%idx) .ne. data%nelv)
then
864 allocate(data%idx(data%nelv))
867 allocate(data%idx(data%nelv))
870 byte_offset = mpi_offset + &
873 call mpi_file_read_at_all(fh, byte_offset, data%idx, data%nelv, &
874 mpi_integer, status, ierr)
876 mpi_offset = mpi_offset + &
880 byte_offset = mpi_offset + int(data%offset_el,
i8) * &
881 (int(data%gdim*lxyz,
i8) * &
882 int(fld_data_size,
i8))
884 data%x, data%y, data%z, data)
885 mpi_offset = mpi_offset + int(data%glb_nelv,
i8) * &
886 (int(data%gdim *lxyz,
i8) * &
887 int(fld_data_size,
i8))
890 if (read_velocity)
then
891 byte_offset = mpi_offset + int(data%offset_el,
i8) * &
892 (int(data%gdim*lxyz,
i8) * &
893 int(fld_data_size,
i8))
895 data%u, data%v, data%w, data)
896 mpi_offset = mpi_offset + int(data%glb_nelv,
i8) * &
897 (int(data%gdim *lxyz,
i8) * &
898 int(fld_data_size,
i8))
901 if (read_pressure)
then
902 byte_offset = mpi_offset + int(data%offset_el,
i8) * &
904 int(fld_data_size,
i8))
906 mpi_offset = mpi_offset + int(data%glb_nelv,
i8) * &
908 int(fld_data_size,
i8))
912 byte_offset = mpi_offset + int(data%offset_el,
i8) * &
914 int(fld_data_size,
i8))
916 mpi_offset = mpi_offset + int(data%glb_nelv,
i8) * &
918 int(fld_data_size,
i8))
922 byte_offset = mpi_offset + int(data%offset_el,
i8) * &
924 int(fld_data_size,
i8))
926 mpi_offset = mpi_offset + int(data%glb_nelv,
i8) * &
928 int(fld_data_size,
i8))
931 this%counter = this%counter + 1
936 call neko_error(
'Currently we only read into fld_file_data_t, &
937 &please use that data structure instead.')