50#ifdef HAVE_ADIOS2_FORTRAN
52 use mpi_f08,
only : mpi_bcast, mpi_character, mpi_integer
64#ifdef HAVE_ADIOS2_FORTRAN
65 type(adios2_adios) :: adios
66 type(adios2_io) :: iowriter
67 type(adios2_io) :: ioreader
72 logical :: dp_precision = .false.
83#ifdef HAVE_ADIOS2_FORTRAN
87 class(*),
target,
intent(in) :: data
88 real(kind=
rp),
intent(in),
optional :: t
90 type(
mesh_t),
pointer :: msh
94 character(len=132) :: hdr
95 character :: rdcode(10)
96 character(len=8) :: id_str
97 character(len=1024) :: fname, base_fname
98 character(len=1024) :: start_field
99 integer :: i, j, ierr, n, suffix_pos, tslash_pos
100 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
102 integer,
allocatable :: idx(:)
104 integer :: n_scalar_fields
105 logical :: write_mesh, write_velocity, write_pressure, write_temperature
106 integer :: adios2_type
107 type(adios2_engine) :: bpwriter
108 type(adios2_variable) :: variable_idx, variable_hdr, variable, variable_msh
109 type(adios2_variable) :: variable_v, variable_p, variable_temp
110 integer(kind=8),
dimension(1) :: shape_dims, start_dims, count_dims
123 write_velocity = .false.
124 write_pressure = .false.
125 write_temperature = .false.
131 if (data%x%size() .gt. 0) x%ptr => data%x%x
132 if (data%y%size() .gt. 0) y%ptr => data%y%x
133 if (data%z%size() .gt. 0) z%ptr => data%z%x
134 if (data%u%size() .gt. 0)
then
136 write_velocity = .true.
138 if (data%v%size() .gt. 0) v%ptr => data%v%x
139 if (data%w%size() .gt. 0) w%ptr => data%w%x
140 if (data%p%size() .gt. 0)
then
142 write_pressure = .true.
144 if (data%t%size() .gt. 0)
then
145 write_temperature = .true.
148 n_scalar_fields = data%n_scalars
149 allocate(scalar_fields(n_scalar_fields))
150 do i = 1, n_scalar_fields
151 scalar_fields(i)%ptr => data%s(i)%x
158 glb_nelv = data%glb_nelv
159 offset_el = data%offset_el
167 select case (data%size())
169 p%ptr => data%items(1)%ptr%x(:,1,1,1)
170 write_pressure = .true.
172 p%ptr => data%items(1)%ptr%x(:,1,1,1)
173 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
174 write_pressure = .true.
175 write_temperature = .true.
177 u%ptr => data%items(1)%ptr%x(:,1,1,1)
178 v%ptr => data%items(2)%ptr%x(:,1,1,1)
179 w%ptr => data%items(3)%ptr%x(:,1,1,1)
180 write_velocity = .true.
182 p%ptr => data%items(1)%ptr%x(:,1,1,1)
183 u%ptr => data%items(2)%ptr%x(:,1,1,1)
184 v%ptr => data%items(3)%ptr%x(:,1,1,1)
185 w%ptr => data%items(4)%ptr%x(:,1,1,1)
186 write_pressure = .true.
187 write_velocity = .true.
189 p%ptr => data%items(1)%ptr%x(:,1,1,1)
190 u%ptr => data%items(2)%ptr%x(:,1,1,1)
191 v%ptr => data%items(3)%ptr%x(:,1,1,1)
192 w%ptr => data%items(4)%ptr%x(:,1,1,1)
194 if (trim(data%name(5)) .eq.
'temperature')
then
196 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
197 n_scalar_fields = data%size() - 5
198 allocate(scalar_fields(n_scalar_fields))
199 do i = 1, n_scalar_fields
200 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
202 write_temperature = .true.
205 n_scalar_fields = data%size() - 4
206 allocate(scalar_fields(n_scalar_fields))
207 do i = 1, n_scalar_fields
208 scalar_fields(i)%ptr => data%items(i+4)%ptr%x(:,1,1,1)
210 write_temperature = .false.
212 write_pressure = .true.
213 write_velocity = .true.
215 call neko_error(
'This many fields not supported yet, bp_file')
223 if (
associated(dof))
then
224 x%ptr => dof%x(:,1,1,1)
225 y%ptr => dof%y(:,1,1,1)
226 z%ptr => dof%z(:,1,1,1)
231 if (
associated(msh))
then
233 glb_nelv = msh%glb_nelv
234 offset_el = msh%offset_el
237 allocate(idx(msh%nelv))
239 idx(i) = msh%elements(i)%e%id()
243 if (
associated(xh))
then
252 if (this%dp_precision)
then
253 adios2_type = adios2_type_dp
255 adios2_type = adios2_type_real
259 call outbuf_points%init(this%dp_precision, gdim, glb_nelv, offset_el, &
262 write(*,*)
"writing layout ", this%layout
264 if (this%layout .eq. 1)
then
266 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
268 else if (this%layout .eq. 2)
then
270 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
272 else if (this%layout .eq. 4)
then
285 write_mesh = (this%get_counter() .eq. this%get_start_counter())
295 if (write_velocity)
then
299 if (write_pressure)
then
303 if (write_temperature)
then
307 if (n_scalar_fields .gt. 0 )
then
310 write(rdcode(i),
'(i1)') (n_scalar_fields)/10
312 write(rdcode(i),
'(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
317 write(hdr, 1) adios2_type, lx, ly, lz, this%layout, glb_nelv, &
318 time, this%get_counter(),
npar, (rdcode(i), i = 1, 10)
3191
format(
'#std',1x, i1, 1x,
i2, 1x,
i2, 1x,
i2, 1x, i10, 1x, i10, 1x, &
320 e20.13, 1x, i9, 1x, i6, 1x, 10a)
325 base_fname = this%get_base_fname()
327 write(id_str,
'(i5.5,a)') this%get_counter(),
'.bp'
328 fname = trim(base_fname(1:suffix_pos-1)) //
"0." // id_str
330 if (.not. adios%valid)
then
332 call adios2_init(adios,
'adios2.xml',
neko_comm%mpi_val, ierr)
335 if (.not. iowriter%valid)
then
336 call adios2_declare_io(iowriter, adios,
'writer', ierr)
337 call adios2_set_engine(iowriter,
'BP5', ierr)
340 call adios2_open(bpwriter, iowriter, trim(fname), adios2_mode_write, &
342 call adios2_begin_step(bpwriter, ierr)
345 call adios2_inquire_variable(variable_hdr, iowriter,
'header', ierr)
346 if (.not.variable_hdr%valid)
then
347 call adios2_define_variable(variable_hdr, iowriter,
'header', &
348 adios2_type_character, ierr)
350 call adios2_put(bpwriter, variable_hdr, hdr, adios2_mode_sync, ierr)
353 shape_dims = (int(glb_nelv,
i8))
354 start_dims = (int(offset_el,
i8))
355 count_dims = (int(nelv,
i8))
356 call adios2_inquire_variable(variable_idx, iowriter,
'idx', ierr)
357 if (.not.variable_idx%valid)
then
358 call adios2_define_variable(variable_idx, iowriter,
'idx', &
359 adios2_type_integer4,
size(shape_dims), shape_dims, start_dims, &
360 count_dims, .false., ierr)
362 call adios2_set_shape(variable_idx,
size(shape_dims), shape_dims, ierr)
363 call adios2_set_selection(variable_idx,
size(start_dims), &
364 start_dims, count_dims, ierr)
366 call adios2_put(bpwriter, variable_idx, idx, adios2_mode_sync, ierr)
372 call outbuf_points%define(variable, iowriter,
'points-x', ierr)
375 call outbuf_points%define(variable, iowriter,
'points-y', ierr)
378 call outbuf_points%define(variable, iowriter,
'points-z', ierr)
382 if (write_velocity)
then
384 if (this%layout .le. 3)
then
385 call outbuf_npar%define(variable, iowriter,
'velocity-u', ierr)
389 if (this%layout .le. 3)
then
390 call outbuf_npar%define(variable, iowriter,
'velocity-v', ierr)
394 if (this%layout .le. 3)
then
395 call outbuf_npar%define(variable, iowriter,
'velocity-w', ierr)
400 if (write_pressure)
then
402 if (this%layout .le. 3)
then
403 call outbuf_npar%define(variable, iowriter,
'pressure', ierr)
408 if (write_temperature)
then
410 if (this%layout .le. 3)
then
411 call outbuf_npar%define(variable, iowriter,
'temperature', ierr)
416 do i = 1, n_scalar_fields
418 if (this%layout .le. 3)
then
419 write(id_str,
'(a,i1,i1)')
's', i / 10, i - 10*(i / 10)
420 call outbuf_npar%define(variable, iowriter, trim(id_str), ierr)
425 if (this%layout .gt. 3)
then
426 call outbuf_npar%define(variable, iowriter,
'fields', ierr)
430 call adios2_end_step(bpwriter, ierr)
431 call adios2_close(bpwriter, ierr)
436 write(start_field,
"(I5,A7)") this%get_start_counter(),
'.adios2'
437 open(newunit = file_unit, &
438 file = trim(base_fname(1:suffix_pos - 1)) &
439 // trim(adjustl(start_field)), status =
'replace')
440 write(file_unit, fmt =
'(A,A,A)')
'filetemplate: ', &
441 base_fname(tslash_pos+1:suffix_pos-1),
'%01d.%05d.bp'
442 write(file_unit, fmt =
'(A,i5)')
'firsttimestep: ', &
443 this%get_start_counter()
444 write(file_unit, fmt =
'(A,i5)')
'numtimesteps: ', &
445 (this%get_counter() + 1) - this%get_start_counter()
446 write(file_unit, fmt =
'(A)')
'type: adios2-bp'
450 call this%increment_counter()
459 class(*),
target,
intent(inout) :: data
460 character(len=132) :: hdr
461 integer :: ierr, suffix_pos, i, j
462 character(len=1024) :: fname, meta_fname, string, base_fname
463 logical :: meta_file, read_mesh, read_velocity, read_pressure
465 character(len=8) :: id_str
466 integer :: lx, ly, lz, glb_nelv, counter, lxyz
468 integer :: adios2_type, n_scalars, n
469 real(kind=
rp) :: time
471 character :: rdcode(10), temp_str(4)
472 class(
buffer_t),
allocatable :: inpbuf_points, inpbuf
473 type(adios2_engine) :: bpreader
474 type(adios2_variable) :: variable_hdr, variable_idx, variable
475 integer(kind=8),
dimension(1) :: start_dims, count_dims
480 base_fname = this%get_base_fname()
482 meta_fname = trim(base_fname(1:suffix_pos-1))
486 inquire(
file = trim(meta_fname), exist = meta_file)
487 if (meta_file .and. data%meta_nsamples .eq. 0)
then
489 open(newunit = file_unit,
file = trim(meta_fname))
490 read(file_unit, fmt =
'(A)') string
491 read(string(14:), fmt =
'(A)') string
492 string = trim(string)
493 data%fld_series_fname = string(:scan(trim(string),
'%') - 1)
494 data%fld_series_fname = trim(data%fld_series_fname) //
'0'
495 read(file_unit, fmt =
'(A)') string
496 read(string(scan(string,
':')+1:), *) data%meta_start_counter
497 read(file_unit, fmt =
'(A)') string
498 read(string(scan(string,
':')+1:), *) data%meta_nsamples
501 write(*,*)
'Reading meta file for bp series'
502 write(*,*)
'Name: ', trim(data%fld_series_fname)
503 write(*,*)
'Start counter: ', data%meta_start_counter, &
504 'Nsamples: ', data%meta_nsamples
506 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
508 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
510 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
512 if (this%get_counter() .eq. 0) &
513 call this%set_counter(data%meta_start_counter)
517 write(id_str,
'(i5.5,a)') this%get_counter(),
'.bp'
518 fname = trim(data%fld_series_fname) //
'.' // id_str
519 if (this%get_counter() .ge. data%meta_nsamples + &
520 data%meta_start_counter)
then
521 call neko_error(
'Trying to read more bp files than exist')
531 if (.not.adios%valid)
then
533 call adios2_init(adios,
'adios2.xml',
neko_comm%mpi_val, ierr)
535 if (.not.ioreader%valid)
then
536 call adios2_declare_io(ioreader, adios,
'reader', ierr)
537 call adios2_set_engine(ioreader,
'BP5', ierr)
542 call adios2_open(bpreader, ioreader, trim(fname), adios2_mode_read, &
544 call adios2_begin_step(bpreader, ierr)
547 if (.not.variable_hdr%valid)
then
548 call adios2_inquire_variable(variable_hdr, ioreader,
'header', ierr)
550 call adios2_get(bpreader, variable_hdr, hdr, adios2_mode_sync, ierr)
552 read(hdr, 1) temp_str, adios2_type, lx, ly, lz, this%layout, glb_nelv,&
553 time, counter,
npar, (rdcode(i),i = 1,10)
5541
format(4a, 1x, i1, 1x,
i2, 1x,
i2, 1x,
i2, 1x, i10, 1x, i10, 1x, &
555 e20.13, 1x, i9, 1x, i6, 1x, 10a)
556 if (data%nelv .eq. 0)
then
558 data%nelv = dist%num_local()
559 data%offset_el = dist%start_idx()
564 data%glb_nelv = glb_nelv
565 data%t_counter = counter
576 if (adios2_type .eq. adios2_type_dp)
then
577 this%dp_precision = .true.
579 this%dp_precision = .false.
582 if (.not.
allocated(inpbuf_points))
allocate(
buffer_1d_t::inpbuf_points)
583 call inpbuf_points%init(this%dp_precision, data%gdim, data%glb_nelv, &
584 data%offset_el, data%nelv, lx, ly, lz)
586 write(*,*)
"layout ", this%layout
587 if (this%layout .eq. 1)
then
588 if (.not.
allocated(inpbuf))
allocate(
buffer_1d_t::inpbuf)
589 else if (this%layout .eq. 2)
then
590 if (.not.
allocated(inpbuf))
allocate(
buffer_4d_t::inpbuf)
591 else if (this%layout .eq. 3)
then
597 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
598 data%offset_el, data%nelv, lx, ly, lz)
600 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
601 data%offset_el, data%nelv, lx, ly, lz)
603 call inpbuf%init(this%dp_precision,
npar, data%glb_nelv, &
604 data%offset_el, data%nelv, lx, ly, lz)
611 read_velocity = .false.
612 read_pressure = .false.
614 if (rdcode(i) .eq.
'X')
then
616 if (data%x%size() .ne. n)
call data%x%init(n)
617 if (data%y%size() .ne. n)
call data%y%init(n)
618 if (data%z%size() .ne. n)
call data%z%init(n)
621 if (rdcode(i) .eq.
'U')
then
622 read_velocity = .true.
623 if (data%u%size() .ne. n)
call data%u%init(n)
624 if (data%v%size() .ne. n)
call data%v%init(n)
625 if (data%w%size() .ne. n)
call data%w%init(n)
628 if (rdcode(i) .eq.
'P')
then
629 read_pressure = .true.
630 if (data%p%size() .ne. n)
call data%p%init(n)
633 if (rdcode(i) .eq.
'T')
then
635 if (data%t%size() .ne. n)
call data%t%init(n)
639 if (rdcode(i) .eq.
'S')
then
641 read(rdcode(i),*) n_scalars
642 n_scalars = n_scalars*10
645 n_scalars = n_scalars+j
647 if (
allocated(data%s))
then
648 if (data%n_scalars .ne. n_scalars)
then
649 do j = 1, data%n_scalars
650 call data%s(j)%free()
653 data%n_scalars = n_scalars
654 allocate(data%s(n_scalars))
655 do j = 1, data%n_scalars
656 call data%s(j)%init(n)
660 data%n_scalars = n_scalars
661 allocate(data%s(data%n_scalars))
662 do j = 1, data%n_scalars
663 call data%s(j)%init(n)
669 if (
allocated(data%idx))
then
670 if (
size(data%idx) .ne. data%nelv)
then
672 allocate(data%idx(data%nelv))
675 allocate(data%idx(data%nelv))
679 start_dims = (int(data%offset_el,
i8))
680 count_dims = (int(data%nelv,
i8))
681 call adios2_inquire_variable(variable_idx, ioreader,
'idx', ierr)
682 if (variable_idx%valid)
then
683 call adios2_set_selection(variable_idx,
size(start_dims), &
684 start_dims, count_dims, ierr)
686 call adios2_get(bpreader, variable_idx, data%idx, adios2_mode_sync, ierr)
689 call inpbuf_points%inquire(variable, ioreader,
'points-x', ierr)
690 call inpbuf_points%read(bpreader, variable, ierr)
691 call inpbuf_points%copy(data%x)
692 call inpbuf_points%inquire(variable, ioreader,
'points-y', ierr)
693 call inpbuf_points%read(bpreader, variable, ierr)
694 call inpbuf_points%copy(data%y)
695 call inpbuf_points%inquire(variable, ioreader,
'points-z', ierr)
696 call inpbuf_points%read(bpreader, variable, ierr)
697 call inpbuf_points%copy(data%z)
700 if (this%layout .eq. 3)
then
701 call inpbuf%inquire(variable, ioreader,
'fields', ierr)
702 call inpbuf%read(bpreader, variable, ierr)
705 if (read_velocity)
then
706 if (this%layout .le. 3)
then
707 call inpbuf%inquire(variable, ioreader,
'velocity-u', ierr)
708 call inpbuf%read(bpreader, variable, ierr)
710 call inpbuf%copy(data%u)
711 if (this%layout .le. 3)
then
712 call inpbuf%inquire(variable, ioreader,
'velocity-v', ierr)
713 call inpbuf%read(bpreader, variable, ierr)
715 call inpbuf%copy(data%v)
716 if (this%layout .le. 3)
then
717 call inpbuf%inquire(variable, ioreader,
'velocity-w', ierr)
718 call inpbuf%read(bpreader, variable, ierr)
720 call inpbuf%copy(data%w)
723 if (read_pressure)
then
724 if (this%layout .le. 3)
then
725 call inpbuf%inquire(variable, ioreader,
'pressure', ierr)
726 call inpbuf%read(bpreader, variable, ierr)
728 call inpbuf%copy(data%p)
732 if (this%layout .le. 3)
then
733 call inpbuf%inquire(variable, ioreader,
'temperature', ierr)
734 call inpbuf%read(bpreader, variable, ierr)
736 call inpbuf%copy(data%t)
740 if (this%layout .le. 3)
then
741 write(id_str,
'(a,i1,i1)')
's', i/10, i-10*(i/10)
742 call inpbuf%inquire(variable, ioreader, trim(id_str), ierr)
743 call inpbuf%read(bpreader, variable, ierr)
745 call inpbuf%copy(data%s(i))
748 call adios2_end_step(bpreader, ierr)
749 call adios2_close(bpreader, ierr)
751 call this%increment_counter()
753 if (
allocated(inpbuf_points))
deallocate(inpbuf_points)
754 if (
allocated(inpbuf))
deallocate(inpbuf)
756 call neko_error(
'Currently we only read into fld_file_data_t,&
757 please use that data structure instead.&
758 (output_format.adios2)')
767 class(*),
target,
intent(in) :: data
768 real(kind=
rp),
intent(in),
optional :: t
769 call neko_error(
'Neko needs to be built with ADIOS2 Fortran support')
774 class(*),
target,
intent(inout) :: data
775 call neko_error(
'Neko needs to be built with ADIOS2 Fortran support')
782 integer,
intent(in) :: precision
784 if (precision .eq.
dp)
then
785 this%dp_precision = .true.
786 else if (precision .eq.
sp)
then
787 this%dp_precision = .false.
796 integer,
intent(in) :: layout
799 if (layout .ge. 1 .and. layout .le. 3)
then
subroutine bp_file_read(this, data)
class(buffer_t), allocatable, private outbuf_npar
class(buffer_t), allocatable, private outbuf_points
subroutine bp_file_write(this, data, t)
subroutine bp_file_set_layout(this, layout)
subroutine bp_file_set_precision(this, precision)
Generic buffer that is extended with buffers of varying rank.
Generic buffer that is extended with buffers of varying rank.
Generic buffer that is extended with buffers of varying rank.
Generic buffer that is extended with buffers of varying rank.
integer, public pe_size
MPI size of communicator.
integer, public pe_rank
MPI rank.
type(mpi_comm), public neko_comm
MPI communicator.
Defines practical data distributions.
Defines a mapping of the degrees of freedom.
Module for file I/O operations.
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
integer, parameter, public i2
integer, parameter, public i8
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
Defines a function space.
Defines structs that are used... Dont know if we should keep it though.
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Interface for ADIOS2 bp files.
Load-balanced linear distribution .
field_list_t, To be able to group fields together
The function space for the SEM solution fields.