Neko 1.99.1
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
37 use num_types, only : dp, sp, i8, rp
39 use field_list, only : field_list_t
40 use dofmap, only : dofmap_t
41 use vector, only : vector_t
42 use space, only : space_t
43 use mesh, only : mesh_t
44 use structs, only : array_ptr_t
46 use utils, only : neko_error
47 use datadist, only : linear_dist_t
48 use comm, only : pe_size, pe_rank, neko_comm
49#ifdef HAVE_ADIOS2_FORTRAN
50 use adios2
51 use mpi_f08, only : mpi_bcast, mpi_character, mpi_integer
52#endif
53 use buffer_1d, only : buffer_1d_t
54 use buffer_4d, only : buffer_4d_t
56 use buffer, only : buffer_t
57 implicit none
58 private
59
60 class(buffer_t), private, allocatable :: outbuf_points
61 class(buffer_t), private, allocatable :: outbuf_npar
62
63#ifdef HAVE_ADIOS2_FORTRAN
64 type(adios2_adios) :: adios
65 type(adios2_io) :: iowriter
66 type(adios2_io) :: ioreader
67#endif
68
70 type, public, extends(generic_file_t) :: bp_file_t
71 logical :: dp_precision = .false.
72 integer :: layout = 1
73 contains
74 procedure :: read => bp_file_read
75 procedure :: write => bp_file_write
76 procedure :: set_precision => bp_file_set_precision
77 procedure :: set_layout => bp_file_set_layout
78 end type bp_file_t
79
80contains
81
82#ifdef HAVE_ADIOS2_FORTRAN
83
84 subroutine bp_file_write(this, data, t)
85 class(bp_file_t), intent(inout) :: this
86 class(*), target, intent(in) :: data
87 real(kind=rp), intent(in), optional :: t
88 type(array_ptr_t) :: x, y, z, u, v, w, p, tem
89 type(mesh_t), pointer :: msh
90 type(dofmap_t), pointer :: dof
91 type(space_t), pointer :: xh
92 real(kind=dp) :: time
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
100 integer :: npar
101 integer, allocatable :: idx(:)
102 type(array_ptr_t), allocatable :: scalar_fields(:)
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
110 integer :: file_unit
111
112 if (present(t)) then
113 time = real(t, dp)
114 else
115 time = 0.0_rp
116 end if
117
118 nullify(msh)
119 nullify(dof)
120 nullify(xh)
121 n_scalar_fields = 0
122 write_velocity = .false.
123 write_pressure = .false.
124 write_temperature = .false.
125
127 select type (data)
128 type is (fld_file_data_t)
129 npar = data%size()
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
134 u%ptr => data%u%x
135 write_velocity = .true.
136 end if
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
140 p%ptr => data%p%x
141 write_pressure = .true.
142 end if
143 if (data%t%n .gt. 0) then
144 write_temperature = .true.
145 tem%ptr => data%t%x
146 end if
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
151 end do
152 nelv = data%nelv
153 lx = data%lx
154 ly = data%ly
155 lz = data%lz
156 gdim = data%gdim
157 glb_nelv = data%glb_nelv
158 offset_el = data%offset_el
159
160 allocate(idx(nelv))
161 do i = 1, nelv
162 idx(i) = data%idx(i)
163 end do
164 type is (field_list_t)
165 npar = data%size()
166 select case (data%size())
167 case (1)
168 p%ptr => data%items(1)%ptr%x(:,1,1,1)
169 write_pressure = .true.
170 case (2)
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.
175 case (3)
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.
180 case (4)
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.
187 case (5:99)
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)
192 ! Check if position 5 is a temperature field by name
193 if (trim(data%name(5)) .eq. 'temperature') then
194 ! Position 5 is temperature, remaining fields are scalars
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)
200 end do
201 write_temperature = .true.
202 else
203 ! All remaining fields are scalars (no temperature field)
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)
208 end do
209 write_temperature = .false.
210 end if
211 write_pressure = .true.
212 write_velocity = .true.
213 case default
214 call neko_error('This many fields not supported yet, bp_file')
215 end select
216 dof => data%dof(1)
217 class default
218 call neko_error('Invalid data')
219 end select
220
221 ! Fix things for pointers that do not exist in all data types...
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)
226 msh => dof%msh
227 xh => dof%Xh
228 end if
229
230 if (associated(msh)) then
231 nelv = msh%nelv
232 glb_nelv = msh%glb_nelv
233 offset_el = msh%offset_el
234 gdim = msh%gdim
235 ! Store global idx of each element
236 allocate(idx(msh%nelv))
237 do i = 1, msh%nelv
238 idx(i) = msh%elements(i)%e%id()
239 end do
240 end if
241
242 if (associated(xh)) then
243 lx = xh%lx
244 ly = xh%ly
245 lz = xh%lz
246 end if
247
248 lxyz = lx*ly*lz
249 n = nelv*lxyz
250
251 if (this%dp_precision) then
252 adios2_type = adios2_type_dp
253 else
254 adios2_type = adios2_type_real
255 end if
256
257 if (.not. allocated(outbuf_points)) allocate(buffer_1d_t::outbuf_points)
258 call outbuf_points%init(this%dp_precision, gdim, glb_nelv, offset_el, &
259 nelv, lx, ly, lz)
260
261 write(*,*) "writing layout ", this%layout
262 if (.not. allocated(outbuf_npar)) then
263 if (this%layout .eq. 1) then
264 allocate(buffer_1d_t::outbuf_npar)
265 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
266 nelv, lx, ly, lz)
267 else if (this%layout .eq. 2) then
268 allocate(buffer_4d_t::outbuf_npar)
269 call outbuf_npar%init(this%dp_precision, gdim, glb_nelv, offset_el, &
270 nelv, lx, ly, lz)
271 else if (this%layout .eq. 4) then
273 call outbuf_npar%init(this%dp_precision, npar, glb_nelv, offset_el, &
274 nelv, lx, ly, lz)
275 else
276 call neko_error('Invalid buffer')
277 end if
278 end if
279
280 !
281 ! Create fld header for NEKTON's multifile output
282 !
283
284 write_mesh = (this%counter .eq. this%start_counter)
285
286 ! Build rdcode note that for field_t, we only support scalar
287 ! fields at the moment
288 rdcode = ' '
289 i = 1
290 if (write_mesh) then
291 rdcode(i) = 'X'
292 i = i + 1
293 end if
294 if (write_velocity) then
295 rdcode(i) = 'U'
296 i = i + 1
297 end if
298 if (write_pressure) then
299 rdcode(i) = 'P'
300 i = i + 1
301 end if
302 if (write_temperature) then
303 rdcode(i) = 'T'
304 i = i + 1
305 end if
306 if (n_scalar_fields .gt. 0 ) then
307 rdcode(i) = 'S'
308 i = i + 1
309 write(rdcode(i), '(i1)') (n_scalar_fields)/10
310 i = i + 1
311 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
312 i = i + 1
313 end if
314
315 ! Create binary header information
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)
320
321 ! Adapt filename with counter
324 suffix_pos = filename_suffix_pos(this%fname)
325 write(id_str, '(i5.5,a)') this%counter, '.bp'
326 fname = trim(this%fname(1:suffix_pos-1)) // "0." // id_str
327
328 if (.not. adios%valid) then
330 call adios2_init(adios, 'adios2.xml', neko_comm%mpi_val, ierr)
331 end if
332
333 if (.not. iowriter%valid) then
334 call adios2_declare_io(iowriter, adios, 'writer', ierr)
335 call adios2_set_engine(iowriter, 'BP5', ierr)
336 end if
337
338 call adios2_open(bpwriter, iowriter, trim(fname), adios2_mode_write, &
339 neko_comm%mpi_val, ierr)
340 call adios2_begin_step(bpwriter, ierr)
341
342 ! Write header
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)
347 end if
348 call adios2_put(bpwriter, variable_hdr, hdr, adios2_mode_sync, ierr)
349
350 ! Write element idxs
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)
359 else
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)
363 end if
364 call adios2_put(bpwriter, variable_idx, idx, adios2_mode_sync, ierr)
365
366 deallocate(idx)
367
368 if (write_mesh) then
369 call outbuf_points%fill(x%ptr, n)
370 call outbuf_points%define(variable, iowriter, 'points-x', ierr)
371 call outbuf_points%write(bpwriter, variable, ierr)
372 call outbuf_points%fill(y%ptr, n)
373 call outbuf_points%define(variable, iowriter, 'points-y', ierr)
374 call outbuf_points%write(bpwriter, variable, ierr)
375 call outbuf_points%fill(z%ptr, n)
376 call outbuf_points%define(variable, iowriter, 'points-z', ierr)
377 call outbuf_points%write(bpwriter, variable, ierr)
378 end if
379
380 if (write_velocity) then
381 call outbuf_npar%fill(u%ptr, n)
382 if (this%layout .le. 3) then
383 call outbuf_npar%define(variable, iowriter, 'velocity-u', ierr)
384 call outbuf_npar%write(bpwriter, variable, ierr)
385 end if
386 call outbuf_npar%fill(v%ptr, n)
387 if (this%layout .le. 3) then
388 call outbuf_npar%define(variable, iowriter, 'velocity-v', ierr)
389 call outbuf_npar%write(bpwriter, variable, ierr)
390 end if
391 call outbuf_npar%fill(w%ptr, n)
392 if (this%layout .le. 3) then
393 call outbuf_npar%define(variable, iowriter, 'velocity-w', ierr)
394 call outbuf_npar%write(bpwriter, variable, ierr)
395 end if
396 end if
397
398 if (write_pressure) then
399 call outbuf_npar%fill(p%ptr, n)
400 if (this%layout .le. 3) then
401 call outbuf_npar%define(variable, iowriter, 'pressure', ierr)
402 call outbuf_npar%write(bpwriter, variable, ierr)
403 end if
404 end if
405
406 if (write_temperature) then
407 call outbuf_npar%fill(tem%ptr, n)
408 if (this%layout .le. 3) then
409 call outbuf_npar%define(variable, iowriter, 'temperature', ierr)
410 call outbuf_npar%write(bpwriter, variable, ierr)
411 end if
412 end if
413
414 do i = 1, n_scalar_fields
415 call outbuf_npar%fill(scalar_fields(i)%ptr, n)
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)
419 call outbuf_npar%write(bpwriter, variable, ierr)
420 end if
421 end do
422
423 if (this%layout .gt. 3) then
424 call outbuf_npar%define(variable, iowriter, 'fields', ierr)
425 call outbuf_npar%write(bpwriter, variable, ierr)
426 end if
427
428 call adios2_end_step(bpwriter, ierr)
429 call adios2_close(bpwriter, ierr)
430
431 ! Write metadata file
432 if (pe_rank .eq. 0) then
433 tslash_pos = filename_tslash_pos(this%fname)
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'
444 close(file_unit)
445 end if
446
447 this%counter = this%counter + 1
448
449 if (allocated(outbuf_points)) deallocate(outbuf_points)
450 if (allocated(outbuf_npar)) deallocate(outbuf_npar)
451 end subroutine bp_file_write
452
454 subroutine bp_file_read(this, data)
455 class(bp_file_t) :: this
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
461 logical :: read_temp
462 character(len=8) :: id_str
463 integer :: lx, ly, lz, glb_nelv, counter, lxyz
464 integer :: npar
465 integer :: adios2_type, n_scalars, n
466 real(kind=rp) :: time
467 type(linear_dist_t) :: dist
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
473 integer :: file_unit
474
475 select type (data)
476 type is (fld_file_data_t)
477 suffix_pos = filename_suffix_pos(this%fname)
478 meta_fname = trim(this%fname(1:suffix_pos-1))
479 call filename_chsuffix(meta_fname, meta_fname,'adios2')
480
482 inquire(file = trim(meta_fname), exist = meta_file)
483 if (meta_file .and. data%meta_nsamples .eq. 0) then
484 if (pe_rank .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
495
496 close(file_unit)
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
501 end if
502 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
503 neko_comm, ierr)
504 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
505 neko_comm, ierr)
506 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
507 neko_comm, ierr)
508 if (this%counter .eq. 0) this%counter = data%meta_start_counter
509 end if
510
511 if (meta_file) then
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')
516 end if
517 else
518 ! @todo write into function because of code duplication
519 !suffix_pos = filename_suffix_pos(this%fname)
520 !write(id_str, '(i5.5,a)') this%counter, '.bp'
521 !fname = trim(this%fname(1:suffix_pos-1))//'.'//id_str
522 fname = this%fname
523 end if
524
525 if (.not.adios%valid) then
526 ! @todo enable parsing XML filename
527 call adios2_init(adios, 'adios2.xml', neko_comm%mpi_val, ierr)
528 end if
529 if (.not.ioreader%valid) then
530 call adios2_declare_io(ioreader, adios, 'reader', ierr)
531 call adios2_set_engine(ioreader, 'BP5', ierr)
532 end if
533
534 ! @todo check if engines and variables should be brought in to local
535 ! subroutine scope
536 call adios2_open(bpreader, ioreader, trim(fname), adios2_mode_read, &
537 neko_comm%mpi_val, ierr)
538 call adios2_begin_step(bpreader, ierr)
539
540 ! Read header and adjust data accordingly
541 if (.not.variable_hdr%valid) then
542 call adios2_inquire_variable(variable_hdr, ioreader, 'header', ierr)
543 end if
544 call adios2_get(bpreader, variable_hdr, hdr, adios2_mode_sync, ierr)
545
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
551 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
552 data%nelv = dist%num_local()
553 data%offset_el = dist%start_idx()
554 end if
555 data%lx = lx
556 data%ly = ly
557 data%lz = lz
558 data%glb_nelv = glb_nelv
559 data%t_counter = counter
560 data%time = time
561 lxyz = lx * ly * lz
562 n = lxyz * data%nelv
563
564 if (lz .eq. 1) then
565 data%gdim = 2
566 else
567 data%gdim = 3
568 end if
569
570 if (adios2_type .eq. adios2_type_dp) then
571 this%dp_precision = .true.
572 else
573 this%dp_precision = .false.
574 end if
575
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)
579
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
586 if (.not. allocated(inpbuf)) allocate(buffer_4d_npar_t::inpbuf)
587 end if
588
589 select type (inpbuf)
590 type is (buffer_1d_t)
591 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
592 data%offset_el, data%nelv, lx, ly, lz)
593 type is (buffer_4d_t)
594 call inpbuf%init(this%dp_precision, data%gdim, data%glb_nelv, &
595 data%offset_el, data%nelv, lx, ly, lz)
596 type is (buffer_4d_npar_t)
597 call inpbuf%init(this%dp_precision, npar, data%glb_nelv, &
598 data%offset_el, data%nelv, lx, ly, lz)
599 class default
600 call neko_error('Invalid buffer')
601 end select
602
603 i = 1
604 read_mesh = .false.
605 read_velocity = .false.
606 read_pressure = .false.
607 read_temp = .false.
608 if (rdcode(i) .eq. 'X') then
609 read_mesh = .true.
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)
613 i = i + 1
614 end if
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)
620 i = i + 1
621 end if
622 if (rdcode(i) .eq. 'P') then
623 read_pressure = .true.
624 if (data%p%n .ne. n) call data%p%init(n)
625 i = i + 1
626 end if
627 if (rdcode(i) .eq. 'T') then
628 read_temp = .true.
629 if (data%t%n .ne. n) call data%t%init(n)
630 i = i + 1
631 end if
632 n_scalars = 0
633 if (rdcode(i) .eq. 'S') then
634 i = i + 1
635 read(rdcode(i),*) n_scalars
636 n_scalars = n_scalars*10
637 i = i + 1
638 read(rdcode(i),*) j
639 n_scalars = n_scalars+j
640 i = i + 1
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()
645 end do
646 deallocate(data%s)
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)
651 end do
652 end if
653 else
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)
658 end do
659 end if
660 i = i + 1
661 end if
662
663 if (allocated(data%idx)) then
664 if (size(data%idx) .ne. data%nelv) then
665 deallocate(data%idx)
666 allocate(data%idx(data%nelv))
667 end if
668 else
669 allocate(data%idx(data%nelv))
670 end if
671
672 ! Read element idxs
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)
679 end if
680 call adios2_get(bpreader, variable_idx, data%idx, adios2_mode_sync, ierr)
681
682 if (read_mesh) then
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)
692 end if
693
694 if (this%layout .eq. 3) then
695 call inpbuf%inquire(variable, ioreader, 'fields', ierr)
696 call inpbuf%read(bpreader, variable, ierr)
697 end if
698
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)
703 end if
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)
708 end if
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)
713 end if
714 call inpbuf%copy(data%w)
715 end if
716
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)
721 end if
722 call inpbuf%copy(data%p)
723 end if
724
725 if (read_temp) then
726 if (this%layout .le. 3) then
727 call inpbuf%inquire(variable, ioreader, 'temperature', ierr)
728 call inpbuf%read(bpreader, variable, ierr)
729 end if
730 call inpbuf%copy(data%t)
731 end if
732
733 do i = 1, n_scalars
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)
738 end if
739 call inpbuf%copy(data%s(i))
740 end do
741
742 call adios2_end_step(bpreader, ierr)
743 call adios2_close(bpreader, ierr)
744
745 this%counter = this%counter + 1
746
747 if (allocated(inpbuf_points)) deallocate(inpbuf_points)
748 if (allocated(inpbuf)) deallocate(inpbuf)
749 class default
750 call neko_error('Currently we only read into fld_file_data_t,&
751 please use that data structure instead.&
752 (output_format.adios2)')
753 end select
754
755 end subroutine bp_file_read
756
757#else
758
759 subroutine bp_file_write(this, data, t)
760 class(bp_file_t), intent(inout) :: this
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')
764 end subroutine bp_file_write
765
766 subroutine bp_file_read(this, data)
767 class(bp_file_t) :: this
768 class(*), target, intent(inout) :: data
769 call neko_error('Neko needs to be built with ADIOS2 Fortran support')
770 end subroutine bp_file_read
771
772#endif
773
774 subroutine bp_file_set_precision(this, precision)
775 class(bp_file_t) :: this
776 integer, intent(in) :: precision
777
778 if (precision .eq. dp) then
779 this%dp_precision = .true.
780 else if (precision .eq. sp) then
781 this%dp_precision = .false.
782 else
783 call neko_error('Invalid precision')
784 end if
785
786 end subroutine bp_file_set_precision
787
788 subroutine bp_file_set_layout(this, layout)
789 class(bp_file_t) :: this
790 integer, intent(in) :: layout
791
793 if (layout .ge. 1 .and. layout .le. 3) then
794 this%layout = layout
795 else
796 call neko_error('Invalid data layout')
797 end if
798
799 end subroutine bp_file_set_layout
800
801end module bp_file
double real
ADIOS2 bp file format.
Definition bp_file.F90:36
subroutine bp_file_read(this, data)
Definition bp_file.F90:767
class(buffer_t), allocatable, private outbuf_npar
Definition bp_file.F90:61
class(buffer_t), allocatable, private outbuf_points
Definition bp_file.F90:60
subroutine bp_file_write(this, data, t)
Definition bp_file.F90:760
subroutine bp_file_set_layout(this, layout)
Definition bp_file.F90:789
subroutine bp_file_set_precision(this, precision)
Definition bp_file.F90:775
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:58
integer, public pe_rank
MPI rank.
Definition comm.F90:55
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
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:78
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Definition utils.f90:64
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:57
Defines a vector.
Definition vector.f90:34
Interface for ADIOS2 bp files.
Definition bp_file.F90:70
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:62
Pointer to array.
Definition structs.f90:14