Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
bp_file.F90
Go to the documentation of this file.
1! Copyright (c) 2024, Gregor Weiss (HLRS)
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!
36module bp_file
38 use num_types, only : dp, sp, i8, rp
40 use field_list, only : field_list_t
41 use dofmap, only : dofmap_t
42 use vector, only : vector_t
43 use space, only : space_t
44 use mesh, only : mesh_t
45 use structs, only : array_ptr_t
47 use utils, only : neko_error
48 use datadist, only : linear_dist_t
49 use comm, only : pe_size, pe_rank, neko_comm
50#ifdef HAVE_ADIOS2_FORTRAN
51 use adios2
52 use mpi_f08, only : mpi_bcast, mpi_character, mpi_integer
53#endif
54 use buffer_1d, only : buffer_1d_t
55 use buffer_4d, only : buffer_4d_t
57 use buffer, only : buffer_t
58 implicit none
59 private
60
61 class(buffer_t), private, allocatable :: outbuf_points
62 class(buffer_t), private, allocatable :: outbuf_npar
63
64#ifdef HAVE_ADIOS2_FORTRAN
65 type(adios2_adios) :: adios
66 type(adios2_io) :: iowriter
67 type(adios2_io) :: ioreader
68#endif
69
71 type, public, extends(generic_file_t) :: bp_file_t
72 logical :: dp_precision = .false.
73 integer :: layout = 1
74 contains
75 procedure :: read => bp_file_read
76 procedure :: write => bp_file_write
77 procedure :: set_precision => bp_file_set_precision
78 procedure :: set_layout => bp_file_set_layout
79 end type bp_file_t
80
81contains
82
83#ifdef HAVE_ADIOS2_FORTRAN
84
85 subroutine bp_file_write(this, data, t)
86 class(bp_file_t), intent(inout) :: this
87 class(*), target, intent(in) :: data
88 real(kind=rp), intent(in), optional :: t
89 type(array_ptr_t) :: x, y, z, u, v, w, p, tem
90 type(mesh_t), pointer :: msh
91 type(dofmap_t), pointer :: dof
92 type(space_t), pointer :: xh
93 real(kind=dp) :: time
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
101 integer :: npar
102 integer, allocatable :: idx(:)
103 type(array_ptr_t), allocatable :: scalar_fields(:)
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
111 integer :: file_unit
112
113 if (present(t)) then
114 time = real(t, dp)
115 else
116 time = 0.0_rp
117 end if
118
119 nullify(msh)
120 nullify(dof)
121 nullify(xh)
122 n_scalar_fields = 0
123 write_velocity = .false.
124 write_pressure = .false.
125 write_temperature = .false.
126
128 select type (data)
129 type is (fld_file_data_t)
130 npar = data%size()
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
135 u%ptr => data%u%x
136 write_velocity = .true.
137 end if
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
141 p%ptr => data%p%x
142 write_pressure = .true.
143 end if
144 if (data%t%size() .gt. 0) then
145 write_temperature = .true.
146 tem%ptr => data%t%x
147 end if
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
152 end do
153 nelv = data%nelv
154 lx = data%lx
155 ly = data%ly
156 lz = data%lz
157 gdim = data%gdim
158 glb_nelv = data%glb_nelv
159 offset_el = data%offset_el
160
161 allocate(idx(nelv))
162 do i = 1, nelv
163 idx(i) = data%idx(i)
164 end do
165 type is (field_list_t)
166 npar = data%size()
167 select case (data%size())
168 case (1)
169 p%ptr => data%items(1)%ptr%x(:,1,1,1)
170 write_pressure = .true.
171 case (2)
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.
176 case (3)
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.
181 case (4)
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.
188 case (5:99)
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)
193 ! Check if position 5 is a temperature field by name
194 if (trim(data%name(5)) .eq. 'temperature') then
195 ! Position 5 is temperature, remaining fields are scalars
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)
201 end do
202 write_temperature = .true.
203 else
204 ! All remaining fields are scalars (no temperature field)
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)
209 end do
210 write_temperature = .false.
211 end if
212 write_pressure = .true.
213 write_velocity = .true.
214 case default
215 call neko_error('This many fields not supported yet, bp_file')
216 end select
217 dof => data%dof(1)
218 class default
219 call neko_error('Invalid data')
220 end select
221
222 ! Fix things for pointers that do not exist in all data types...
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)
227 msh => dof%msh
228 xh => dof%Xh
229 end if
230
231 if (associated(msh)) then
232 nelv = msh%nelv
233 glb_nelv = msh%glb_nelv
234 offset_el = msh%offset_el
235 gdim = msh%gdim
236 ! Store global idx of each element
237 allocate(idx(msh%nelv))
238 do i = 1, msh%nelv
239 idx(i) = msh%elements(i)%e%id()
240 end do
241 end if
242
243 if (associated(xh)) then
244 lx = xh%lx
245 ly = xh%ly
246 lz = xh%lz
247 end if
248
249 lxyz = lx*ly*lz
250 n = nelv*lxyz
251
252 if (this%dp_precision) then
253 adios2_type = adios2_type_dp
254 else
255 adios2_type = adios2_type_real
256 end if
257
258 if (.not. allocated(outbuf_points)) allocate(buffer_1d_t::outbuf_points)
259 call outbuf_points%init(this%dp_precision, gdim, glb_nelv, offset_el, &
260 nelv, lx, ly, lz)
261
262 write(*,*) "writing layout ", this%layout
263 if (.not. allocated(outbuf_npar)) then
264 if (this%layout .eq. 1) then
265 allocate(buffer_1d_t::outbuf_npar)
266 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
267 nelv, lx, ly, lz)
268 else if (this%layout .eq. 2) then
269 allocate(buffer_4d_t::outbuf_npar)
270 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
271 nelv, lx, ly, lz)
272 else if (this%layout .eq. 4) then
274 call outbuf_npar%init(this%dp_precision, npar, glb_nelv, offset_el, &
275 nelv, lx, ly, lz)
276 else
277 call neko_error('Invalid buffer')
278 end if
279 end if
280
281 !
282 ! Create fld header for NEKTON's multifile output
283 !
284
285 write_mesh = (this%get_counter() .eq. this%get_start_counter())
286
287 ! Build rdcode note that for field_t, we only support scalar
288 ! fields at the moment
289 rdcode = ' '
290 i = 1
291 if (write_mesh) then
292 rdcode(i) = 'X'
293 i = i + 1
294 end if
295 if (write_velocity) then
296 rdcode(i) = 'U'
297 i = i + 1
298 end if
299 if (write_pressure) then
300 rdcode(i) = 'P'
301 i = i + 1
302 end if
303 if (write_temperature) then
304 rdcode(i) = 'T'
305 i = i + 1
306 end if
307 if (n_scalar_fields .gt. 0 ) then
308 rdcode(i) = 'S'
309 i = i + 1
310 write(rdcode(i), '(i1)') (n_scalar_fields)/10
311 i = i + 1
312 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
313 i = i + 1
314 end if
315
316 ! Create binary header information
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)
321
322 ! Adapt filename with counter
325 base_fname = this%get_base_fname()
326 suffix_pos = filename_suffix_pos(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
329
330 if (.not. adios%valid) then
332 call adios2_init(adios, 'adios2.xml', neko_comm%mpi_val, ierr)
333 end if
334
335 if (.not. iowriter%valid) then
336 call adios2_declare_io(iowriter, adios, 'writer', ierr)
337 call adios2_set_engine(iowriter, 'BP5', ierr)
338 end if
339
340 call adios2_open(bpwriter, iowriter, trim(fname), adios2_mode_write, &
341 neko_comm%mpi_val, ierr)
342 call adios2_begin_step(bpwriter, ierr)
343
344 ! Write header
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)
349 end if
350 call adios2_put(bpwriter, variable_hdr, hdr, adios2_mode_sync, ierr)
351
352 ! Write element idxs
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)
361 else
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)
365 end if
366 call adios2_put(bpwriter, variable_idx, idx, adios2_mode_sync, ierr)
367
368 deallocate(idx)
369
370 if (write_mesh) then
371 call outbuf_points%fill(x%ptr, n)
372 call outbuf_points%define(variable, iowriter, 'points-x', ierr)
373 call outbuf_points%write(bpwriter, variable, ierr)
374 call outbuf_points%fill(y%ptr, n)
375 call outbuf_points%define(variable, iowriter, 'points-y', ierr)
376 call outbuf_points%write(bpwriter, variable, ierr)
377 call outbuf_points%fill(z%ptr, n)
378 call outbuf_points%define(variable, iowriter, 'points-z', ierr)
379 call outbuf_points%write(bpwriter, variable, ierr)
380 end if
381
382 if (write_velocity) then
383 call outbuf_npar%fill(u%ptr, n)
384 if (this%layout .le. 3) then
385 call outbuf_npar%define(variable, iowriter, 'velocity-u', ierr)
386 call outbuf_npar%write(bpwriter, variable, ierr)
387 end if
388 call outbuf_npar%fill(v%ptr, n)
389 if (this%layout .le. 3) then
390 call outbuf_npar%define(variable, iowriter, 'velocity-v', ierr)
391 call outbuf_npar%write(bpwriter, variable, ierr)
392 end if
393 call outbuf_npar%fill(w%ptr, n)
394 if (this%layout .le. 3) then
395 call outbuf_npar%define(variable, iowriter, 'velocity-w', ierr)
396 call outbuf_npar%write(bpwriter, variable, ierr)
397 end if
398 end if
399
400 if (write_pressure) then
401 call outbuf_npar%fill(p%ptr, n)
402 if (this%layout .le. 3) then
403 call outbuf_npar%define(variable, iowriter, 'pressure', ierr)
404 call outbuf_npar%write(bpwriter, variable, ierr)
405 end if
406 end if
407
408 if (write_temperature) then
409 call outbuf_npar%fill(tem%ptr, n)
410 if (this%layout .le. 3) then
411 call outbuf_npar%define(variable, iowriter, 'temperature', ierr)
412 call outbuf_npar%write(bpwriter, variable, ierr)
413 end if
414 end if
415
416 do i = 1, n_scalar_fields
417 call outbuf_npar%fill(scalar_fields(i)%ptr, n)
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)
421 call outbuf_npar%write(bpwriter, variable, ierr)
422 end if
423 end do
424
425 if (this%layout .gt. 3) then
426 call outbuf_npar%define(variable, iowriter, 'fields', ierr)
427 call outbuf_npar%write(bpwriter, variable, ierr)
428 end if
429
430 call adios2_end_step(bpwriter, ierr)
431 call adios2_close(bpwriter, ierr)
432
433 ! Write metadata file
434 if (pe_rank .eq. 0) then
435 tslash_pos = filename_tslash_pos(base_fname)
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'
447 close(file_unit)
448 end if
449
450 call this%increment_counter()
451
452 if (allocated(outbuf_points)) deallocate(outbuf_points)
453 if (allocated(outbuf_npar)) deallocate(outbuf_npar)
454 end subroutine bp_file_write
455
457 subroutine bp_file_read(this, data)
458 class(bp_file_t) :: this
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
464 logical :: read_temp
465 character(len=8) :: id_str
466 integer :: lx, ly, lz, glb_nelv, counter, lxyz
467 integer :: npar
468 integer :: adios2_type, n_scalars, n
469 real(kind=rp) :: time
470 type(linear_dist_t) :: dist
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
476 integer :: file_unit
477
478 select type (data)
479 type is (fld_file_data_t)
480 base_fname = this%get_base_fname()
481 suffix_pos = filename_suffix_pos(base_fname)
482 meta_fname = trim(base_fname(1:suffix_pos-1))
483 call filename_chsuffix(meta_fname, meta_fname, 'adios2')
484
486 inquire(file = trim(meta_fname), exist = meta_file)
487 if (meta_file .and. data%meta_nsamples .eq. 0) then
488 if (pe_rank .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
499
500 close(file_unit)
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
505 end if
506 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
507 neko_comm, ierr)
508 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
509 neko_comm, ierr)
510 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
511 neko_comm, ierr)
512 if (this%get_counter() .eq. 0) &
513 call this%set_counter(data%meta_start_counter)
514 end if
515
516 if (meta_file) then
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')
522 end if
523 else
524 ! @todo write into function because of code duplication
525 !suffix_pos = filename_suffix_pos(base_fname)
526 !write(id_str, '(i5.5,a)') this%get_counter(), '.bp'
527 !fname = trim(base_fname(1:suffix_pos-1))//'.'//id_str
528 fname = base_fname
529 end if
530
531 if (.not.adios%valid) then
532 ! @todo enable parsing XML filename
533 call adios2_init(adios, 'adios2.xml', neko_comm%mpi_val, ierr)
534 end if
535 if (.not.ioreader%valid) then
536 call adios2_declare_io(ioreader, adios, 'reader', ierr)
537 call adios2_set_engine(ioreader, 'BP5', ierr)
538 end if
539
540 ! @todo check if engines and variables should be brought in to local
541 ! subroutine scope
542 call adios2_open(bpreader, ioreader, trim(fname), adios2_mode_read, &
543 neko_comm%mpi_val, ierr)
544 call adios2_begin_step(bpreader, ierr)
545
546 ! Read header and adjust data accordingly
547 if (.not.variable_hdr%valid) then
548 call adios2_inquire_variable(variable_hdr, ioreader, 'header', ierr)
549 end if
550 call adios2_get(bpreader, variable_hdr, hdr, adios2_mode_sync, ierr)
551
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
557 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
558 data%nelv = dist%num_local()
559 data%offset_el = dist%start_idx()
560 end if
561 data%lx = lx
562 data%ly = ly
563 data%lz = lz
564 data%glb_nelv = glb_nelv
565 data%t_counter = counter
566 data%time = time
567 lxyz = lx * ly * lz
568 n = lxyz * data%nelv
569
570 if (lz .eq. 1) then
571 data%gdim = 2
572 else
573 data%gdim = 3
574 end if
575
576 if (adios2_type .eq. adios2_type_dp) then
577 this%dp_precision = .true.
578 else
579 this%dp_precision = .false.
580 end if
581
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)
585
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
592 if (.not. allocated(inpbuf)) allocate(buffer_4d_npar_t::inpbuf)
593 end if
594
595 select type (inpbuf)
596 type is (buffer_1d_t)
597 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
598 data%offset_el, data%nelv, lx, ly, lz)
599 type is (buffer_4d_t)
600 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
601 data%offset_el, data%nelv, lx, ly, lz)
602 type is (buffer_4d_npar_t)
603 call inpbuf%init(this%dp_precision, npar, data%glb_nelv, &
604 data%offset_el, data%nelv, lx, ly, lz)
605 class default
606 call neko_error('Invalid buffer')
607 end select
608
609 i = 1
610 read_mesh = .false.
611 read_velocity = .false.
612 read_pressure = .false.
613 read_temp = .false.
614 if (rdcode(i) .eq. 'X') then
615 read_mesh = .true.
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)
619 i = i + 1
620 end if
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)
626 i = i + 1
627 end if
628 if (rdcode(i) .eq. 'P') then
629 read_pressure = .true.
630 if (data%p%size() .ne. n) call data%p%init(n)
631 i = i + 1
632 end if
633 if (rdcode(i) .eq. 'T') then
634 read_temp = .true.
635 if (data%t%size() .ne. n) call data%t%init(n)
636 i = i + 1
637 end if
638 n_scalars = 0
639 if (rdcode(i) .eq. 'S') then
640 i = i + 1
641 read(rdcode(i),*) n_scalars
642 n_scalars = n_scalars*10
643 i = i + 1
644 read(rdcode(i),*) j
645 n_scalars = n_scalars+j
646 i = i + 1
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()
651 end do
652 deallocate(data%s)
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)
657 end do
658 end if
659 else
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)
664 end do
665 end if
666 i = i + 1
667 end if
668
669 if (allocated(data%idx)) then
670 if (size(data%idx) .ne. data%nelv) then
671 deallocate(data%idx)
672 allocate(data%idx(data%nelv))
673 end if
674 else
675 allocate(data%idx(data%nelv))
676 end if
677
678 ! Read element idxs
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)
685 end if
686 call adios2_get(bpreader, variable_idx, data%idx, adios2_mode_sync, ierr)
687
688 if (read_mesh) then
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)
698 end if
699
700 if (this%layout .eq. 3) then
701 call inpbuf%inquire(variable, ioreader, 'fields', ierr)
702 call inpbuf%read(bpreader, variable, ierr)
703 end if
704
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)
709 end if
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)
714 end if
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)
719 end if
720 call inpbuf%copy(data%w)
721 end if
722
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)
727 end if
728 call inpbuf%copy(data%p)
729 end if
730
731 if (read_temp) then
732 if (this%layout .le. 3) then
733 call inpbuf%inquire(variable, ioreader, 'temperature', ierr)
734 call inpbuf%read(bpreader, variable, ierr)
735 end if
736 call inpbuf%copy(data%t)
737 end if
738
739 do i = 1, n_scalars
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)
744 end if
745 call inpbuf%copy(data%s(i))
746 end do
747
748 call adios2_end_step(bpreader, ierr)
749 call adios2_close(bpreader, ierr)
750
751 call this%increment_counter()
752
753 if (allocated(inpbuf_points)) deallocate(inpbuf_points)
754 if (allocated(inpbuf)) deallocate(inpbuf)
755 class default
756 call neko_error('Currently we only read into fld_file_data_t,&
757 please use that data structure instead.&
758 (output_format.adios2)')
759 end select
760
761 end subroutine bp_file_read
762
763#else
764
765 subroutine bp_file_write(this, data, t)
766 class(bp_file_t), intent(inout) :: this
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')
770 end subroutine bp_file_write
771
772 subroutine bp_file_read(this, data)
773 class(bp_file_t) :: this
774 class(*), target, intent(inout) :: data
775 call neko_error('Neko needs to be built with ADIOS2 Fortran support')
776 end subroutine bp_file_read
777
778#endif
779
780 subroutine bp_file_set_precision(this, precision)
781 class(bp_file_t) :: this
782 integer, intent(in) :: precision
783
784 if (precision .eq. dp) then
785 this%dp_precision = .true.
786 else if (precision .eq. sp) then
787 this%dp_precision = .false.
788 else
789 call neko_error('Invalid precision')
790 end if
791
792 end subroutine bp_file_set_precision
793
794 subroutine bp_file_set_layout(this, layout)
795 class(bp_file_t) :: this
796 integer, intent(in) :: layout
797
799 if (layout .ge. 1 .and. layout .le. 3) then
800 this%layout = layout
801 else
802 call neko_error('Invalid data layout')
803 end if
804
805 end subroutine bp_file_set_layout
806
807end module bp_file
double real
ADIOS2 bp file format.
Definition bp_file.F90:36
subroutine bp_file_read(this, data)
Definition bp_file.F90:773
class(buffer_t), allocatable, private outbuf_npar
Definition bp_file.F90:62
class(buffer_t), allocatable, private outbuf_points
Definition bp_file.F90:61
subroutine bp_file_write(this, data, t)
Definition bp_file.F90:766
subroutine bp_file_set_layout(this, layout)
Definition bp_file.F90:795
subroutine bp_file_set_precision(this, precision)
Definition bp_file.F90:781
Generic buffer that is extended with buffers of varying rank.
Definition buffer_1d.F90:34
Generic buffer that is extended with buffers of varying rank.
integer, private npar
Generic buffer that is extended with buffers of varying rank.
Definition buffer_4d.F90:34
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:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Defines practical data distributions.
Definition datadist.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
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...
Defines a mesh.
Definition mesh.f90:34
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:140
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Definition utils.f90:65
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:58
Defines a vector.
Definition vector.f90:34
Interface for ADIOS2 bp files.
Definition bp_file.F90:71
Load-balanced linear distribution .
Definition datadist.f90:50
field_list_t, To be able to group fields together
A generic file handler.
The function space for the SEM solution fields.
Definition space.f90:63
Pointer to array.
Definition structs.f90:14