49#ifdef HAVE_ADIOS2_FORTRAN
51 use mpi_f08,
only : mpi_bcast, mpi_character, mpi_integer
63#ifdef HAVE_ADIOS2_FORTRAN
64 type(adios2_adios) :: adios
65 type(adios2_io) :: iowriter
66 type(adios2_io) :: ioreader
71 logical :: dp_precision = .false.
82#ifdef HAVE_ADIOS2_FORTRAN
86 class(*),
target,
intent(in) :: data
87 real(kind=
rp),
intent(in),
optional :: t
89 type(
mesh_t),
pointer :: msh
93 character(len=132) :: hdr
94 character :: rdcode(10)
95 character(len=8) :: id_str
96 character(len=1024) :: fname
97 character(len=1024) :: start_field
98 integer :: i, j, ierr, n, suffix_pos, tslash_pos
99 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
101 integer,
allocatable :: idx(:)
103 integer :: n_scalar_fields
104 logical :: write_mesh, write_velocity, write_pressure, write_temperature
105 integer :: adios2_type
106 type(adios2_engine) :: bpwriter
107 type(adios2_variable) :: variable_idx, variable_hdr, variable, variable_msh
108 type(adios2_variable) :: variable_v, variable_p, variable_temp
109 integer(kind=8),
dimension(1) :: shape_dims, start_dims, count_dims
122 write_velocity = .false.
123 write_pressure = .false.
124 write_temperature = .false.
130 if (data%x%n .gt. 0) x%ptr => data%x%x
131 if (data%y%n .gt. 0) y%ptr => data%y%x
132 if (data%z%n .gt. 0) z%ptr => data%z%x
133 if (data%u%n .gt. 0)
then
135 write_velocity = .true.
137 if (data%v%n .gt. 0) v%ptr => data%v%x
138 if (data%w%n .gt. 0) w%ptr => data%w%x
139 if (data%p%n .gt. 0)
then
141 write_pressure = .true.
143 if (data%t%n .gt. 0)
then
144 write_temperature = .true.
147 n_scalar_fields = data%n_scalars
148 allocate(scalar_fields(n_scalar_fields))
149 do i = 1, n_scalar_fields
150 scalar_fields(i)%ptr => data%s(i)%x
157 glb_nelv = data%glb_nelv
158 offset_el = data%offset_el
166 select case (data%size())
168 p%ptr => data%items(1)%ptr%x(:,1,1,1)
169 write_pressure = .true.
171 p%ptr => data%items(1)%ptr%x(:,1,1,1)
172 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
173 write_pressure = .true.
174 write_temperature = .true.
176 u%ptr => data%items(1)%ptr%x(:,1,1,1)
177 v%ptr => data%items(2)%ptr%x(:,1,1,1)
178 w%ptr => data%items(3)%ptr%x(:,1,1,1)
179 write_velocity = .true.
181 p%ptr => data%items(1)%ptr%x(:,1,1,1)
182 u%ptr => data%items(2)%ptr%x(:,1,1,1)
183 v%ptr => data%items(3)%ptr%x(:,1,1,1)
184 w%ptr => data%items(4)%ptr%x(:,1,1,1)
185 write_pressure = .true.
186 write_velocity = .true.
188 p%ptr => data%items(1)%ptr%x(:,1,1,1)
189 u%ptr => data%items(2)%ptr%x(:,1,1,1)
190 v%ptr => data%items(3)%ptr%x(:,1,1,1)
191 w%ptr => data%items(4)%ptr%x(:,1,1,1)
193 if (trim(data%name(5)) .eq.
'temperature')
then
195 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
196 n_scalar_fields = data%size() - 5
197 allocate(scalar_fields(n_scalar_fields))
198 do i = 1, n_scalar_fields
199 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
201 write_temperature = .true.
204 n_scalar_fields = data%size() - 4
205 allocate(scalar_fields(n_scalar_fields))
206 do i = 1, n_scalar_fields
207 scalar_fields(i)%ptr => data%items(i+4)%ptr%x(:,1,1,1)
209 write_temperature = .false.
211 write_pressure = .true.
212 write_velocity = .true.
214 call neko_error(
'This many fields not supported yet, bp_file')
222 if (
associated(dof))
then
223 x%ptr => dof%x(:,1,1,1)
224 y%ptr => dof%y(:,1,1,1)
225 z%ptr => dof%z(:,1,1,1)
230 if (
associated(msh))
then
232 glb_nelv = msh%glb_nelv
233 offset_el = msh%offset_el
236 allocate(idx(msh%nelv))
238 idx(i) = msh%elements(i)%e%id()
242 if (
associated(xh))
then
251 if (this%dp_precision)
then
252 adios2_type = adios2_type_dp
254 adios2_type = adios2_type_real
258 call outbuf_points%init(this%dp_precision, gdim, glb_nelv, offset_el, &
261 write(*,*)
"writing layout ", this%layout
263 if (this%layout .eq. 1)
then
265 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
267 else if (this%layout .eq. 2)
then
269 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
271 else if (this%layout .eq. 4)
then
284 write_mesh = (this%counter .eq. this%start_counter)
294 if (write_velocity)
then
298 if (write_pressure)
then
302 if (write_temperature)
then
306 if (n_scalar_fields .gt. 0 )
then
309 write(rdcode(i),
'(i1)') (n_scalar_fields)/10
311 write(rdcode(i),
'(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
316 write(hdr, 1) adios2_type, lx, ly, lz, this%layout, glb_nelv, &
317 time, this%counter,
npar, (rdcode(i), i = 1, 10)
3181
format(
'#std',1x, i1, 1x,
i2, 1x,
i2, 1x,
i2, 1x, i10, 1x, i10, 1x, &
319 e20.13, 1x, i9, 1x, i6, 1x, 10a)
325 write(id_str,
'(i5.5,a)') this%counter,
'.bp'
326 fname = trim(this%fname(1:suffix_pos-1)) //
"0." // id_str
328 if (.not. adios%valid)
then
330 call adios2_init(adios,
'adios2.xml',
neko_comm%mpi_val, ierr)
333 if (.not. iowriter%valid)
then
334 call adios2_declare_io(iowriter, adios,
'writer', ierr)
335 call adios2_set_engine(iowriter,
'BP5', ierr)
338 call adios2_open(bpwriter, iowriter, trim(fname), adios2_mode_write, &
340 call adios2_begin_step(bpwriter, ierr)
343 call adios2_inquire_variable(variable_hdr, iowriter,
'header', ierr)
344 if (.not.variable_hdr%valid)
then
345 call adios2_define_variable(variable_hdr, iowriter,
'header', &
346 adios2_type_character, ierr)
348 call adios2_put(bpwriter, variable_hdr, hdr, adios2_mode_sync, ierr)
351 shape_dims = (int(glb_nelv,
i8))
352 start_dims = (int(offset_el,
i8))
353 count_dims = (int(nelv,
i8))
354 call adios2_inquire_variable(variable_idx, iowriter,
'idx', ierr)
355 if (.not.variable_idx%valid)
then
356 call adios2_define_variable(variable_idx, iowriter,
'idx', &
357 adios2_type_integer4,
size(shape_dims), shape_dims, start_dims, &
358 count_dims, .false., ierr)
360 call adios2_set_shape(variable_idx,
size(shape_dims), shape_dims, ierr)
361 call adios2_set_selection(variable_idx,
size(start_dims), &
362 start_dims, count_dims, ierr)
364 call adios2_put(bpwriter, variable_idx, idx, adios2_mode_sync, ierr)
370 call outbuf_points%define(variable, iowriter,
'points-x', ierr)
373 call outbuf_points%define(variable, iowriter,
'points-y', ierr)
376 call outbuf_points%define(variable, iowriter,
'points-z', ierr)
380 if (write_velocity)
then
382 if (this%layout .le. 3)
then
383 call outbuf_npar%define(variable, iowriter,
'velocity-u', ierr)
387 if (this%layout .le. 3)
then
388 call outbuf_npar%define(variable, iowriter,
'velocity-v', ierr)
392 if (this%layout .le. 3)
then
393 call outbuf_npar%define(variable, iowriter,
'velocity-w', ierr)
398 if (write_pressure)
then
400 if (this%layout .le. 3)
then
401 call outbuf_npar%define(variable, iowriter,
'pressure', ierr)
406 if (write_temperature)
then
408 if (this%layout .le. 3)
then
409 call outbuf_npar%define(variable, iowriter,
'temperature', ierr)
414 do i = 1, n_scalar_fields
416 if (this%layout .le. 3)
then
417 write(id_str,
'(a,i1,i1)')
's', i / 10, i - 10*(i / 10)
418 call outbuf_npar%define(variable, iowriter, trim(id_str), ierr)
423 if (this%layout .gt. 3)
then
424 call outbuf_npar%define(variable, iowriter,
'fields', ierr)
428 call adios2_end_step(bpwriter, ierr)
429 call adios2_close(bpwriter, ierr)
434 write(start_field,
"(I5,A7)") this%start_counter,
'.adios2'
435 open(unit = newunit(file_unit), &
436 file = trim(this%fname(1:suffix_pos - 1)) &
437 // trim(adjustl(start_field)), status=
'replace')
438 write(file_unit, fmt =
'(A,A,A)')
'filetemplate: ', &
439 this%fname(tslash_pos+1:suffix_pos-1),
'%01d.%05d.bp'
440 write(file_unit, fmt =
'(A,i5)')
'firsttimestep: ', this%start_counter
441 write(file_unit, fmt =
'(A,i5)')
'numtimesteps: ', &
442 (this%counter + 1) - this%start_counter
443 write(file_unit, fmt =
'(A)')
'type: adios2-bp'
447 this%counter = this%counter + 1
456 class(*),
target,
intent(inout) :: data
457 character(len=132) :: hdr
458 integer :: ierr, suffix_pos, i, j
459 character(len=1024) :: fname, meta_fname, string
460 logical :: meta_file, read_mesh, read_velocity, read_pressure
462 character(len=8) :: id_str
463 integer :: lx, ly, lz, glb_nelv, counter, lxyz
465 integer :: adios2_type, n_scalars, n
466 real(kind=
rp) :: time
468 character :: rdcode(10), temp_str(4)
469 class(
buffer_t),
allocatable :: inpbuf_points, inpbuf
470 type(adios2_engine) :: bpreader
471 type(adios2_variable) :: variable_hdr, variable_idx, variable
472 integer(kind=8),
dimension(1) :: start_dims, count_dims
478 meta_fname = trim(this%fname(1:suffix_pos-1))
482 inquire(
file = trim(meta_fname), exist = meta_file)
483 if (meta_file .and. data%meta_nsamples .eq. 0)
then
485 open(unit = newunit(file_unit),
file = trim(meta_fname))
486 read(file_unit, fmt =
'(A)') string
487 read(string(14:), fmt =
'(A)') string
488 string = trim(string)
489 data%fld_series_fname = string(:scan(trim(string),
'%') - 1)
490 data%fld_series_fname = trim(data%fld_series_fname) //
'0'
491 read(file_unit, fmt =
'(A)') string
492 read(string(scan(string,
':')+1:), *) data%meta_start_counter
493 read(file_unit, fmt =
'(A)') string
494 read(string(scan(string,
':')+1:), *) data%meta_nsamples
497 write(*,*)
'Reading meta file for bp series'
498 write(*,*)
'Name: ', trim(data%fld_series_fname)
499 write(*,*)
'Start counter: ', data%meta_start_counter, &
500 'Nsamples: ', data%meta_nsamples
502 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
504 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
506 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
508 if (this%counter .eq. 0) this%counter = data%meta_start_counter
512 write(id_str,
'(i5.5,a)') this%counter,
'.bp'
513 fname = trim(data%fld_series_fname) //
'.' // id_str
514 if (this%counter .ge. data%meta_nsamples+data%meta_start_counter)
then
515 call neko_error(
'Trying to read more bp files than exist')
525 if (.not.adios%valid)
then
527 call adios2_init(adios,
'adios2.xml',
neko_comm%mpi_val, ierr)
529 if (.not.ioreader%valid)
then
530 call adios2_declare_io(ioreader, adios,
'reader', ierr)
531 call adios2_set_engine(ioreader,
'BP5', ierr)
536 call adios2_open(bpreader, ioreader, trim(fname), adios2_mode_read, &
538 call adios2_begin_step(bpreader, ierr)
541 if (.not.variable_hdr%valid)
then
542 call adios2_inquire_variable(variable_hdr, ioreader,
'header', ierr)
544 call adios2_get(bpreader, variable_hdr, hdr, adios2_mode_sync, ierr)
546 read(hdr, 1) temp_str, adios2_type, lx, ly, lz, this%layout, glb_nelv,&
547 time, counter,
npar, (rdcode(i),i = 1,10)
5481
format(4a, 1x, i1, 1x,
i2, 1x,
i2, 1x,
i2, 1x, i10, 1x, i10, 1x, &
549 e20.13, 1x, i9, 1x, i6, 1x, 10a)
550 if (data%nelv .eq. 0)
then
552 data%nelv = dist%num_local()
553 data%offset_el = dist%start_idx()
558 data%glb_nelv = glb_nelv
559 data%t_counter = counter
570 if (adios2_type .eq. adios2_type_dp)
then
571 this%dp_precision = .true.
573 this%dp_precision = .false.
576 if (.not.
allocated(inpbuf_points))
allocate(
buffer_1d_t::inpbuf_points)
577 call inpbuf_points%init(this%dp_precision, data%gdim, data%glb_nelv, &
578 data%offset_el, data%nelv, lx, ly, lz)
580 write(*,*)
"layout ", this%layout
581 if (this%layout .eq. 1)
then
582 if (.not.
allocated(inpbuf))
allocate(
buffer_1d_t::inpbuf)
583 else if (this%layout .eq. 2)
then
584 if (.not.
allocated(inpbuf))
allocate(
buffer_4d_t::inpbuf)
585 else if (this%layout .eq. 3)
then
591 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
592 data%offset_el, data%nelv, lx, ly, lz)
594 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
595 data%offset_el, data%nelv, lx, ly, lz)
597 call inpbuf%init(this%dp_precision,
npar, data%glb_nelv, &
598 data%offset_el, data%nelv, lx, ly, lz)
605 read_velocity = .false.
606 read_pressure = .false.
608 if (rdcode(i) .eq.
'X')
then
610 if (data%x%n .ne. n)
call data%x%init(n)
611 if (data%y%n .ne. n)
call data%y%init(n)
612 if (data%z%n .ne. n)
call data%z%init(n)
615 if (rdcode(i) .eq.
'U')
then
616 read_velocity = .true.
617 if (data%u%n .ne. n)
call data%u%init(n)
618 if (data%v%n .ne. n)
call data%v%init(n)
619 if (data%w%n .ne. n)
call data%w%init(n)
622 if (rdcode(i) .eq.
'P')
then
623 read_pressure = .true.
624 if (data%p%n .ne. n)
call data%p%init(n)
627 if (rdcode(i) .eq.
'T')
then
629 if (data%t%n .ne. n)
call data%t%init(n)
633 if (rdcode(i) .eq.
'S')
then
635 read(rdcode(i),*) n_scalars
636 n_scalars = n_scalars*10
639 n_scalars = n_scalars+j
641 if (
allocated(data%s))
then
642 if (data%n_scalars .ne. n_scalars)
then
643 do j = 1, data%n_scalars
644 call data%s(j)%free()
647 data%n_scalars = n_scalars
648 allocate(data%s(n_scalars))
649 do j = 1, data%n_scalars
650 call data%s(j)%init(n)
654 data%n_scalars = n_scalars
655 allocate(data%s(data%n_scalars))
656 do j = 1, data%n_scalars
657 call data%s(j)%init(n)
663 if (
allocated(data%idx))
then
664 if (
size(data%idx) .ne. data%nelv)
then
666 allocate(data%idx(data%nelv))
669 allocate(data%idx(data%nelv))
673 start_dims = (int(data%offset_el,
i8))
674 count_dims = (int(data%nelv,
i8))
675 call adios2_inquire_variable(variable_idx, ioreader,
'idx', ierr)
676 if (variable_idx%valid)
then
677 call adios2_set_selection(variable_idx,
size(start_dims), &
678 start_dims, count_dims, ierr)
680 call adios2_get(bpreader, variable_idx, data%idx, adios2_mode_sync, ierr)
683 call inpbuf_points%inquire(variable, ioreader,
'points-x', ierr)
684 call inpbuf_points%read(bpreader, variable, ierr)
685 call inpbuf_points%copy(data%x)
686 call inpbuf_points%inquire(variable, ioreader,
'points-y', ierr)
687 call inpbuf_points%read(bpreader, variable, ierr)
688 call inpbuf_points%copy(data%y)
689 call inpbuf_points%inquire(variable, ioreader,
'points-z', ierr)
690 call inpbuf_points%read(bpreader, variable, ierr)
691 call inpbuf_points%copy(data%z)
694 if (this%layout .eq. 3)
then
695 call inpbuf%inquire(variable, ioreader,
'fields', ierr)
696 call inpbuf%read(bpreader, variable, ierr)
699 if (read_velocity)
then
700 if (this%layout .le. 3)
then
701 call inpbuf%inquire(variable, ioreader,
'velocity-u', ierr)
702 call inpbuf%read(bpreader, variable, ierr)
704 call inpbuf%copy(data%u)
705 if (this%layout .le. 3)
then
706 call inpbuf%inquire(variable, ioreader,
'velocity-v', ierr)
707 call inpbuf%read(bpreader, variable, ierr)
709 call inpbuf%copy(data%v)
710 if (this%layout .le. 3)
then
711 call inpbuf%inquire(variable, ioreader,
'velocity-w', ierr)
712 call inpbuf%read(bpreader, variable, ierr)
714 call inpbuf%copy(data%w)
717 if (read_pressure)
then
718 if (this%layout .le. 3)
then
719 call inpbuf%inquire(variable, ioreader,
'pressure', ierr)
720 call inpbuf%read(bpreader, variable, ierr)
722 call inpbuf%copy(data%p)
726 if (this%layout .le. 3)
then
727 call inpbuf%inquire(variable, ioreader,
'temperature', ierr)
728 call inpbuf%read(bpreader, variable, ierr)
730 call inpbuf%copy(data%t)
734 if (this%layout .le. 3)
then
735 write(id_str,
'(a,i1,i1)')
's', i/10, i-10*(i/10)
736 call inpbuf%inquire(variable, ioreader, trim(id_str), ierr)
737 call inpbuf%read(bpreader, variable, ierr)
739 call inpbuf%copy(data%s(i))
742 call adios2_end_step(bpreader, ierr)
743 call adios2_close(bpreader, ierr)
745 this%counter = this%counter + 1
747 if (
allocated(inpbuf_points))
deallocate(inpbuf_points)
748 if (
allocated(inpbuf))
deallocate(inpbuf)
750 call neko_error(
'Currently we only read into fld_file_data_t,&
751 please use that data structure instead.&
752 (output_format.adios2)')
761 class(*),
target,
intent(in) :: data
762 real(kind=
rp),
intent(in),
optional :: t
763 call neko_error(
'Neko needs to be built with ADIOS2 Fortran support')
768 class(*),
target,
intent(inout) :: data
769 call neko_error(
'Neko needs to be built with ADIOS2 Fortran support')
776 integer,
intent(in) :: precision
778 if (precision .eq.
dp)
then
779 this%dp_precision = .true.
780 else if (precision .eq.
sp)
then
781 this%dp_precision = .false.
790 integer,
intent(in) :: layout
793 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.