Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fld_file.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2023, The Neko Authors
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!
36 use num_types, only : rp, dp, sp, i8
38 use field, only : field_t
39 use field_list, only : field_list_t
40 use dofmap, only : dofmap_t
41 use space, only : space_t
42 use structs, only : array_ptr_t
43 use vector, only : vector_t
45 use vector, only : vector_t
46 use space, only : space_t
47 use logger, only: neko_log, log_size
48 use mesh, only : mesh_t
51 use comm
52 use datadist, only : linear_dist_t
53 use math, only : vlmin, vlmax, sabscmp
56 use mpi_f08
57 implicit none
58 private
59
60 real(kind=dp), private, allocatable :: tmp_dp(:)
61 real(kind=sp), private, allocatable :: tmp_sp(:)
62
64 type, public, extends(generic_file_t) :: fld_file_t
65 logical :: dp_precision = .false.
66 contains
67 procedure :: read => fld_file_read
68 procedure :: write => fld_file_write
69 procedure :: set_precision => fld_file_set_precision
70 procedure :: get_fld_fname => fld_file_get_fld_fname
71 procedure :: get_meta_fname => fld_file_get_meta_fname
72 end type fld_file_t
73
74
75contains
76
79 subroutine fld_file_write(this, data, t)
80 class(fld_file_t), intent(inout) :: this
81 class(*), target, intent(in) :: data
82 real(kind=rp), intent(in), optional :: t
83 type(array_ptr_t) :: x, y, z, u, v, w, p, tem
84 real(kind=rp), allocatable, target :: tempo(:)
85 type(mesh_t), pointer :: msh
86 type(dofmap_t), pointer :: dof
87 type(space_t), pointer :: Xh
88 real(kind=dp) :: time
89 character(len= 132) :: hdr
90 character :: rdcode(10)
91 character(len=6) :: id_str
92 character(len= 1024) :: fname
93 character(len= 1024) :: name
94 integer :: file_unit
95 integer :: i, ierr, n, suffix_pos, tslash_pos
96 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
97 integer, allocatable :: idx(:)
98 type(mpi_status) :: status
99 type(mpi_file) :: fh
100 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset, temp_offset
101 real(kind=sp), parameter :: test_pattern = 6.54321
102 type(array_ptr_t), allocatable :: scalar_fields(:)
103 logical :: write_mesh, write_velocity, write_pressure, write_temperature
104 integer :: fld_data_size, n_scalar_fields
105
106 if (present(t)) then
107 time = real(t, dp)
108 else
109 time = 0d0
110 end if
111
112 nullify(msh)
113 nullify(dof)
114 nullify(xh)
115 n_scalar_fields = 0
116 write_pressure = .false.
117 write_velocity = .false.
118 write_temperature = .false.
119
120 select type (data)
121 type is (fld_file_data_t)
122 nelv = data%nelv
123 lx = data%lx
124 ly = data%ly
125 lz = data%lz
126 gdim = data%gdim
127 glb_nelv = data%glb_nelv
128 offset_el = data%offset_el
129
130 if (data%x%size() .gt. 0) x%ptr => data%x%x
131 if (data%y%size() .gt. 0) y%ptr => data%y%x
132 if (data%z%size() .gt. 0) z%ptr => data%z%x
133 if (gdim .eq. 2) z%ptr => data%y%x
134 if (data%u%size() .gt. 0) then
135 u%ptr => data%u%x
136 ! In case only u is actually allocated, point the other comps to u
137 ! so that we don't die on trying to write them
138 if (data%v%size() .le. 0) v%ptr => data%u%x
139 if (data%w%size() .le. 0) w%ptr => data%u%x
140 write_velocity = .true.
141 end if
142 if (data%v%size() .gt. 0) v%ptr => data%v%x
143 if (data%w%size() .gt. 0) w%ptr => data%w%x
144 if (data%p%size() .gt. 0) then
145 p%ptr => data%p%x
146 write_pressure = .true.
147 end if
148 if (data%t%size() .gt. 0) then
149 write_temperature = .true.
150 tem%ptr => data%t%x
151 end if
152 ! If gdim = 2 and Z-velocity component exists,
153 ! it is stored in last scalar field
154 if (gdim .eq. 2 .and. data%w%size() .gt. 0) then
155 n_scalar_fields = data%n_scalars + 1
156 allocate(scalar_fields(n_scalar_fields))
157 do i = 1, n_scalar_fields -1
158 scalar_fields(i)%ptr => data%s(i)%x
159 end do
160 scalar_fields(n_scalar_fields)%ptr => data%w%x
161 else
162 n_scalar_fields = data%n_scalars
163 allocate(scalar_fields(n_scalar_fields+1))
164 do i = 1, n_scalar_fields
165 scalar_fields(i)%ptr => data%s(i)%x
166 end do
167 scalar_fields(n_scalar_fields+1)%ptr => data%w%x
168 end if
169 ! This is very stupid...
170 ! Some compilers cannot handle that these pointers dont point to anything
171 ! (although they are not used) this fixes a segfault due to this.
172 if (nelv .eq. 0) then
173 allocate(tempo(1))
174 x%ptr => tempo
175 y%ptr => tempo
176 z%ptr => tempo
177 u%ptr => tempo
178 v%ptr => tempo
179 w%ptr => tempo
180 p%ptr => tempo
181 tem%ptr => tempo
182 end if
183
184 allocate(idx(nelv))
185 do i = 1, nelv
186 idx(i) = data%idx(i)
187 end do
188 type is (field_t)
189 p%ptr => data%x(:,1,1,1)
190 dof => data%dof
191 write_pressure = .true.
192 write_velocity = .false.
193 type is (field_list_t)
194 select case (data%size())
195 case (1)
196 p%ptr => data%items(1)%ptr%x(:,1,1,1)
197 write_pressure = .true.
198 write_velocity = .false.
199 case (2)
200 p%ptr => data%items(1)%ptr%x(:,1,1,1)
201 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
202 write_pressure = .true.
203 write_temperature = .true.
204 case (3)
205 u%ptr => data%items(1)%ptr%x(:,1,1,1)
206 v%ptr => data%items(2)%ptr%x(:,1,1,1)
207 w%ptr => data%items(3)%ptr%x(:,1,1,1)
208 write_velocity = .true.
209 case (4)
210 p%ptr => data%items(1)%ptr%x(:,1,1,1)
211 u%ptr => data%items(2)%ptr%x(:,1,1,1)
212 v%ptr => data%items(3)%ptr%x(:,1,1,1)
213 w%ptr => data%items(4)%ptr%x(:,1,1,1)
214 write_pressure = .true.
215 write_velocity = .true.
216 case (5:99)
217 p%ptr => data%items(1)%ptr%x(:,1,1,1)
218 u%ptr => data%items(2)%ptr%x(:,1,1,1)
219 v%ptr => data%items(3)%ptr%x(:,1,1,1)
220 w%ptr => data%items(4)%ptr%x(:,1,1,1)
221 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
222 n_scalar_fields = data%size() - 5
223 allocate(scalar_fields(n_scalar_fields))
224 do i = 1, n_scalar_fields
225 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
226 end do
227 write_pressure = .true.
228 write_velocity = .true.
229 write_temperature = .true.
230 case default
231 call neko_error('This many fields not supported yet, fld_file')
232 end select
233 dof => data%dof(1)
234 class default
235 call neko_error('Invalid data')
236 end select
237 ! Fix things for pointers that do not exist in all data types...
238 if (associated(dof)) then
239 x%ptr => dof%x(:,1,1,1)
240 y%ptr => dof%y(:,1,1,1)
241 z%ptr => dof%z(:,1,1,1)
242 msh => dof%msh
243 xh => dof%Xh
244 end if
245
246 if (associated(msh)) then
247 nelv = msh%nelv
248 glb_nelv = msh%glb_nelv
249 offset_el = msh%offset_el
250 gdim = msh%gdim
251 ! Store global idx of each element
252 allocate(idx(msh%nelv))
253 do i = 1, msh%nelv
254 idx(i) = msh%elements(i)%e%id()
255 end do
256 end if
257
258 if (associated(xh)) then
259 lx = xh%lx
260 ly = xh%ly
261 lz = xh%lz
262 end if
263
264 lxyz = lx*ly*lz
265 n = nelv*lxyz
266
267 if (this%dp_precision) then
268 fld_data_size = mpi_double_precision_size
269 else
270 fld_data_size = mpi_real_size
271 end if
272 if (this%dp_precision) then
273 allocate(tmp_dp(gdim*n))
274 else
275 allocate(tmp_sp(gdim*n))
276 end if
277
278
279 !
280 ! Create fld header for NEKTON's multifile output
281 !
282
283 call this%increment_counter()
284 write_mesh = (this%get_counter() .eq. this%get_start_counter())
285 call mpi_allreduce(mpi_in_place, write_mesh, 1, &
286 mpi_logical, mpi_lor, neko_comm)
287 call mpi_allreduce(mpi_in_place, write_velocity, 1, &
288 mpi_logical, mpi_lor, neko_comm)
289 call mpi_allreduce(mpi_in_place, write_pressure, 1, &
290 mpi_logical, mpi_lor, neko_comm)
291 call mpi_allreduce(mpi_in_place, write_temperature, 1, &
292 mpi_logical, mpi_lor, neko_comm)
293 call mpi_allreduce(mpi_in_place, n_scalar_fields, 1, &
294 mpi_integer, mpi_max, neko_comm)
295
296 ! Build rdcode note that for field_t, we only support scalar
297 ! fields at the moment
298 rdcode = ' '
299 i = 1
300 if (write_mesh) then
301 rdcode(i) = 'X'
302 i = i + 1
303 end if
304 if (write_velocity) then
305 rdcode(i) = 'U'
306 i = i + 1
307 end if
308 if (write_pressure) then
309 rdcode(i) = 'P'
310 i = i + 1
311 end if
312 if (write_temperature) then
313 rdcode(i) = 'T'
314 i = i + 1
315 end if
316 if (n_scalar_fields .gt. 0 ) then
317 rdcode(i) = 'S'
318 i = i + 1
319 write(rdcode(i), '(i1)') (n_scalar_fields)/10
320 i = i + 1
321 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
322 i = i + 1
323 end if
324
326 write(hdr, 1) fld_data_size, lx, ly, lz, glb_nelv, glb_nelv,&
327 time, this%get_counter(), 1, 1, (rdcode(i), i = 1, 10)
3281 format('#std', 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
329 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
330
331 ! Change to NEKTON's fld file format
332 fname = this%get_fld_fname()
333
334 call mpi_file_open(neko_comm, trim(fname), &
335 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, &
336 ierr)
337
338 call mpi_file_write_all(fh, hdr, 132, mpi_character, status, ierr)
339 mpi_offset = 132 * mpi_character_size
340
341 call mpi_file_write_all(fh, test_pattern, 1, mpi_real, status, ierr)
342 mpi_offset = mpi_offset + mpi_real_size
343
344 byte_offset = mpi_offset + &
345 int(offset_el, i8) * int(mpi_integer_size, i8)
346 call mpi_file_write_at_all(fh, byte_offset, idx, nelv, &
347 mpi_integer, status, ierr)
348 mpi_offset = mpi_offset + int(glb_nelv, i8) * int(mpi_integer_size, i8)
349 deallocate(idx)
350 if (write_mesh) then
351
352 byte_offset = mpi_offset + int(offset_el, i8) * &
353 (int(gdim*lxyz, i8) * &
354 int(fld_data_size, i8))
355 call fld_file_write_vector_field(this, fh, byte_offset, &
356 x%ptr, y%ptr, z%ptr, &
357 n, gdim, lxyz, nelv)
358 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
359 (int(gdim *lxyz, i8) * &
360 int(fld_data_size, i8))
361 end if
362 if (write_velocity) then
363 byte_offset = mpi_offset + int(offset_el, i8) * &
364 (int(gdim * (lxyz), i8) * int(fld_data_size, i8))
365 call fld_file_write_vector_field(this, fh, byte_offset, &
366 u%ptr, v%ptr, w%ptr, n, gdim, lxyz, nelv)
367
368 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
369 (int(gdim * (lxyz), i8) * &
370 int(fld_data_size, i8))
371
372 end if
373
374 if (write_pressure) then
375 byte_offset = mpi_offset + int(offset_el, i8) * &
376 (int((lxyz), i8) * int(fld_data_size, i8))
377 call fld_file_write_field(this, fh, byte_offset, p%ptr, n)
378 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
379 (int((lxyz), i8) * int(fld_data_size, i8))
380 end if
381
382 if (write_temperature) then
383 byte_offset = mpi_offset + int(offset_el, i8) * &
384 (int((lxyz), i8) * &
385 int(fld_data_size, i8))
386 call fld_file_write_field(this, fh, byte_offset, tem%ptr, n)
387 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
388 (int((lxyz), i8) * &
389 int(fld_data_size, i8))
390 end if
391
392 temp_offset = mpi_offset
393
394 do i = 1, n_scalar_fields
395 ! Without this redundant if statement Cray optimizes this loop to
396 ! Oblivion
397 if (i .eq. 2) then
398 mpi_offset = int(temp_offset, i8) + int(1_i8*glb_nelv, i8) * &
399 (int(lxyz, i8) * int(fld_data_size, i8))
400 end if
401 byte_offset = int(mpi_offset, i8) + int(offset_el, i8) * &
402 (int((lxyz), i8) * &
403 int(fld_data_size, i8))
404 call fld_file_write_field(this, fh, byte_offset, scalar_fields(i)%ptr, n)
405 mpi_offset = int(mpi_offset, i8) + int(glb_nelv, i8) * &
406 (int(lxyz, i8) * &
407 int(fld_data_size, i8))
408 end do
409
410 if (gdim .eq. 3) then
411
413 if (write_mesh) then
414 !The offset is:
415 ! mpioff + element_off * 2(min max value)
416 ! * 4(single precision) * gdim(dimensions)
417 byte_offset = int(mpi_offset, i8) + &
418 int(offset_el, i8) * &
419 int(2, i8) * &
420 int(mpi_real_size, i8) * &
421 int(gdim, i8)
422 call fld_file_write_metadata_vector(this, fh, byte_offset, &
423 x%ptr, y%ptr, z%ptr, gdim, lxyz, nelv)
424 mpi_offset = int(mpi_offset, i8) + &
425 int(glb_nelv, i8) * &
426 int(2, i8) * &
427 int(mpi_real_size, i8) * &
428 int(gdim, i8)
429 end if
430
431 if (write_velocity) then
432 byte_offset = int(mpi_offset, i8) + &
433 int(offset_el, i8) * &
434 int(2, i8) * &
435 int(mpi_real_size, i8) * &
436 int(gdim, i8)
437 call fld_file_write_metadata_vector(this, fh, byte_offset, &
438 u%ptr, v%ptr, w%ptr, gdim, lxyz, nelv)
439 mpi_offset = int(mpi_offset, i8) + &
440 int(glb_nelv, i8) * &
441 int(2, i8) * &
442 int(mpi_real_size, i8) * &
443 int(gdim, i8)
444
445 end if
446
447 if (write_pressure) then
448 byte_offset = int(mpi_offset, i8) + &
449 int(offset_el, i8) * &
450 int(2, i8) * &
451 int(mpi_real_size, i8)
452 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
453 p%ptr, lxyz, nelv)
454 mpi_offset = int(mpi_offset, i8) + &
455 int(glb_nelv, i8) * &
456 int(2, i8) * &
457 int(mpi_real_size, i8)
458
459 end if
460
461 if (write_temperature) then
462 byte_offset = int(mpi_offset, i8) + &
463 int(offset_el, i8) * &
464 int(2, i8) * &
465 int(mpi_real_size, i8)
466 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
467 tem%ptr, lxyz, nelv)
468 mpi_offset = int(mpi_offset, i8) + &
469 int(glb_nelv, i8) * &
470 int(2, i8) * &
471 int(mpi_real_size, i8)
472
473 end if
474
475
476
477 temp_offset = mpi_offset
478
479 do i = 1, n_scalar_fields
480 ! Without this redundant if statement, Cray optimizes this loop to
481 ! Oblivion
482 if (i .eq. 2) then
483 mpi_offset = int(temp_offset, i8) + &
484 int(1_i8*glb_nelv, i8) * &
485 int(2, i8) * &
486 int(mpi_real_size, i8)
487 end if
488
489 byte_offset = int(mpi_offset, i8) + &
490 int(offset_el, i8) * &
491 int(2, i8) * &
492 int(mpi_real_size, i8)
493 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
494 scalar_fields(i)%ptr, lxyz, nelv)
495 mpi_offset = int(mpi_offset, i8) + &
496 int(glb_nelv, i8) * &
497 int(2, i8) * &
498 int(mpi_real_size, i8)
499 end do
500 end if
501
502
503 call mpi_file_sync(fh, ierr)
504 call mpi_file_close(fh, ierr)
505 ! Write metadata file
506 if (pe_rank .eq. 0) then
507 call filename_name(this%get_base_fname(), name)
508
509 open(newunit = file_unit, &
510 file = this%get_meta_fname(), status = 'replace')
511 ! The following string will specify that the files in the file series
512 ! are defined by the filename followed by a 0.
513 ! This 0 is necessary as it specifies the index of number of files
514 ! the output file is split across.
515 ! In the past, many .f files were generated for each write.
516 ! To be consistent with this the trailing 0 is still necessary today.
517 write(file_unit, fmt = '(A,A,A)') 'filetemplate: ', &
518 trim(name), '%01d.f%05d'
519 write(file_unit, fmt = '(A,i5)') 'firsttimestep: ', &
520 this%get_start_counter()
521 write(file_unit, fmt = '(A,i5)') 'numtimesteps: ', &
522 (this%get_counter() + 1) - this%get_start_counter()
523 close(file_unit)
524 end if
525
526 if (allocated(tmp_dp)) deallocate(tmp_dp)
527 if (allocated(tmp_sp)) deallocate(tmp_sp)
528 end subroutine fld_file_write
529
530 subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, &
531 gdim, lxyz, nelv)
532 class(fld_file_t), intent(inout) :: this
533 type(mpi_file), intent(inout) :: fh
534 integer, intent(in) :: gdim, lxyz, nelv
535 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
536 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
537 integer :: el, j, ierr, nout
538 type(mpi_status) :: status
539 real(kind=sp) :: buffer(2*gdim*nelv)
540
541 j = 1
542 do el = 1, nelv
543 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
544 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
545 buffer(j+2) = real(vlmin(y(1, el), lxyz), sp)
546 buffer(j+3) = real(vlmax(y(1, el), lxyz), sp)
547 j = j + 4
548 if (gdim .eq. 3) then
549 buffer(j+0) = real(vlmin(z(1, el), lxyz), sp)
550 buffer(j+1) = real(vlmax(z(1, el), lxyz), sp)
551 j = j + 2
552 end if
553 end do
554
555 ! write out data
556 nout = 2*gdim*nelv
557
558 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
559 mpi_real, status, ierr)
560
561 end subroutine fld_file_write_metadata_vector
562
563 subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, &
564 nelv)
565 class(fld_file_t), intent(inout) :: this
566 type(mpi_file), intent(inout) :: fh
567 integer, intent(in) :: lxyz, nelv
568 real(kind=rp), intent(in) :: x(lxyz, nelv)
569 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
570 integer :: el, j, ierr, nout
571 type(mpi_status) :: status
572 real(kind=sp) :: buffer(2*nelv)
573
574 j = 1
575 do el = 1, nelv
576 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
577 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
578 j = j + 2
579 end do
580
581 ! write out data
582 nout = 2*nelv
583
584 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
585 mpi_real, status, ierr)
586
587 end subroutine fld_file_write_metadata_scalar
588
589 subroutine fld_file_write_field(this, fh, byte_offset, p, n)
590 class(fld_file_t), intent(inout) :: this
591 type(mpi_file), intent(inout) :: fh
592 integer, intent(inout) :: n
593 real(kind=rp), intent(inout) :: p(n)
594 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
595 integer :: i, ierr
596 type(mpi_status) :: status
597
598 if ( this%dp_precision) then
599 do i = 1, n
600 tmp_dp(i) = real(p(i), dp)
601 end do
602
603 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, n, &
604 mpi_double_precision, status, ierr)
605 else
606 do i = 1, n
607 tmp_sp(i) = real(p(i), sp)
608 end do
609 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, n, &
610 mpi_real, status, ierr)
611 end if
612
613 end subroutine fld_file_write_field
614
615 subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, &
616 gdim, lxyz, nelv)
617 class(fld_file_t), intent(inout) :: this
618 type(mpi_file), intent(inout) :: fh
619 integer, intent(in) :: n, gdim, lxyz, nelv
620 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
621 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
622 integer :: i, el, j, ierr
623 type(mpi_status) :: status
624
625 if (this%dp_precision) then
626 i = 1
627 do el = 1, nelv
628 do j = 1, lxyz
629 tmp_dp(i) = real(x(j, el), dp)
630 i = i +1
631 end do
632 do j = 1, lxyz
633 tmp_dp(i) = real(y(j, el), dp)
634 i = i +1
635 end do
636 if (gdim .eq. 3) then
637 do j = 1, lxyz
638 tmp_dp(i) = real(z(j, el), dp)
639 i = i +1
640 end do
641 end if
642 end do
643 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, gdim*n, &
644 mpi_double_precision, status, ierr)
645 else
646 i = 1
647 do el = 1, nelv
648 do j = 1, lxyz
649 tmp_sp(i) = real(x(j, el), sp)
650 i = i +1
651 end do
652 do j = 1, lxyz
653 tmp_sp(i) = real(y(j, el), sp)
654 i = i +1
655 end do
656 if (gdim .eq. 3) then
657 do j = 1, lxyz
658 tmp_sp(i) = real(z(j, el), sp)
659 i = i +1
660 end do
661 end if
662 end do
663 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, gdim*n, &
664 mpi_real, status, ierr)
665 end if
666
667
668 end subroutine fld_file_write_vector_field
669
671 subroutine fld_file_read(this, data)
672 class(fld_file_t) :: this
673 class(*), target, intent(inout) :: data
674 character(len= 132) :: hdr
675 integer :: ierr, suffix_pos, i, j
676 type(mpi_file) :: fh
677 type(mpi_status) :: status
678 character(len= 1024) :: fname, base_fname, meta_fname, string, path
679 logical :: meta_file, read_mesh, read_velocity, read_pressure
680 logical :: read_temp
681 character(len=6) :: suffix
682 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
683 integer :: lx, ly, lz, glb_nelv, counter, lxyz
684 integer :: FLD_DATA_SIZE, n_scalars, n
685 integer :: file_unit
686 real(kind=rp) :: time
687 real(kind=sp) :: temp
688 type(linear_dist_t) :: dist
689 real(kind=sp), parameter :: test_pattern = 6.54321
690 character :: rdcode(10), temp_str(4)
691 character(len=LOG_SIZE) :: log_buf
692
693 select type (data)
694 type is (fld_file_data_t)
695 call filename_chsuffix(this%get_base_fname(), meta_fname, 'nek5000')
696
697 inquire(file = trim(meta_fname), exist = meta_file)
698 if (meta_file .and. data%meta_nsamples .eq. 0) then
699 if (pe_rank .eq. 0) then
700 open(newunit = file_unit, file = trim(meta_fname))
701 read(file_unit, fmt = '(A)') string
702 read(string(14:), fmt = '(A)') string
703 string = trim(string)
704
705 data%fld_series_fname = string(:scan(trim(string), '%')-1)
706 data%fld_series_fname = adjustl(data%fld_series_fname)
707 data%fld_series_fname = trim(data%fld_series_fname)//'0'
708
709 read(file_unit, fmt = '(A)') string
710 read(string(scan(string, ':')+1:), *) data%meta_start_counter
711 read(file_unit, fmt = '(A)') string
712 read(string(scan(string, ':')+1:), *) data%meta_nsamples
713 close(file_unit)
714
715 write(log_buf,*) 'Reading meta file for fld series'
716 call neko_log%message(log_buf)
717 write(log_buf,*) 'Name: ', trim(data%fld_series_fname)
718 call neko_log%message(log_buf)
719 write(log_buf,*) 'Start counter: ', data%meta_start_counter
720 call neko_log%message(log_buf)
721 write(log_buf,*) 'Nsamples: ', data%meta_nsamples
722 call neko_log%message(log_buf)
723
724 end if
725 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
726 neko_comm, ierr)
727 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
728 neko_comm, ierr)
729 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
730 neko_comm, ierr)
731
732 if (this%get_counter() .eq. -1) then
733 call this%set_start_counter(data%meta_start_counter)
734 call this%set_counter(data%meta_start_counter)
735 end if
736 end if
737
738 if (meta_file) then
739 call filename_path(this%get_base_fname(), path)
740 write(suffix, '(a,i5.5)') 'f', this%get_counter()
741 fname = trim(path) // trim(data%fld_series_fname) // '.' // suffix
742 if (this%get_counter() .ge. &
743 data%meta_nsamples+data%meta_start_counter) then
744 call neko_error('Trying to read more fld files than exist')
745 end if
746 else
747 write(suffix, '(a,i5.5)') 'f', this%get_counter()
748 call filename_chsuffix(trim(this%get_base_fname()), fname, suffix)
749 end if
750 call mpi_file_open(neko_comm, trim(fname), &
751 mpi_mode_rdonly, mpi_info_null, fh, ierr)
752
753 if (ierr .ne. 0) call neko_error("Could not read "//trim(fname))
754
755 call neko_log%message('Reading fld file ' // trim(fname))
756
757 call mpi_file_read_all(fh, hdr, 132, mpi_character, status, ierr)
758 ! This read can prorbably be done wihtout the temp variables,
759 ! temp_str, i, j
760
761 read(hdr, 1) temp_str, fld_data_size, lx, ly, lz, glb_nelv, glb_nelv, &
762 time, counter, i, j, (rdcode(i), i = 1, 10)
7631 format(4a, 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
764 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
765 if (data%nelv .eq. 0) then
766 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
767 data%nelv = dist%num_local()
768 data%offset_el = dist%start_idx()
769 end if
770 data%lx = lx
771 data%ly = ly
772 data%lz = lz
773 data%glb_nelv = glb_nelv
774 data%t_counter = counter
775 data%time = time
776 lxyz = lx * ly * lz
777 n = lxyz * data%nelv
778
779 if (lz .eq. 1) then
780 data%gdim = 2
781 else
782 data%gdim = 3
783 end if
784
785
786 if (fld_data_size .eq. mpi_double_precision_size) then
787 this%dp_precision = .true.
788 else
789 this%dp_precision = .false.
790 end if
791 if (this%dp_precision) then
792 allocate(tmp_dp(data%gdim*n))
793 else
794 allocate(tmp_sp(data%gdim*n))
795 end if
796
797
798 i = 1
799 read_mesh = .false.
800 read_velocity = .false.
801 read_pressure = .false.
802 read_temp = .false.
803 if (rdcode(i) .eq. 'X') then
804 read_mesh = .true.
805 call data%x%init(n)
806 call data%y%init(n)
807 call data%z%init(n)
808 i = i + 1
809 end if
810 if (rdcode(i) .eq. 'U') then
811 read_velocity = .true.
812 call data%u%init(n)
813 call data%v%init(n)
814 call data%w%init(n)
815 i = i + 1
816 end if
817 if (rdcode(i) .eq. 'P') then
818 read_pressure = .true.
819 call data%p%init(n)
820 i = i + 1
821 end if
822 if (rdcode(i) .eq. 'T') then
823 read_temp = .true.
824 call data%t%init(n)
825 i = i + 1
826 end if
827 n_scalars = 0
828 if (rdcode(i) .eq. 'S') then
829 i = i + 1
830 read(rdcode(i),*) n_scalars
831 n_scalars = n_scalars*10
832 i = i + 1
833 read(rdcode(i),*) j
834 n_scalars = n_scalars+j
835 i = i + 1
836 if (allocated(data%s)) then
837 if (data%n_scalars .ne. n_scalars) then
838 do j = 1, data%n_scalars
839 call data%s(j)%free()
840 end do
841 deallocate(data%s)
842 data%n_scalars = n_scalars
843 allocate(data%s(n_scalars))
844 do j = 1, data%n_scalars
845 call data%s(j)%init(n)
846 end do
847 end if
848 else
849 data%n_scalars = n_scalars
850 allocate(data%s(data%n_scalars))
851 do j = 1, data%n_scalars
852 call data%s(j)%init(n)
853 end do
854 end if
855 i = i + 1
856 end if
857
858 mpi_offset = 132 * mpi_character_size
859 call mpi_file_read_at_all(fh, mpi_offset, temp, 1, &
860 mpi_real, status, ierr)
861 if (.not. sabscmp(temp, test_pattern, epsilon(1.0_sp))) then
862 call neko_error('Incorrect format for fld file, &
863 &test pattern does not match.')
864 end if
865 mpi_offset = mpi_offset + mpi_real_size
866
867
868 if (allocated(data%idx)) then
869 if (size(data%idx) .ne. data%nelv) then
870 deallocate(data%idx)
871 allocate(data%idx(data%nelv))
872 end if
873 else
874 allocate(data%idx(data%nelv))
875 end if
876
877 byte_offset = mpi_offset + &
878 int(data%offset_el, i8) * int(mpi_integer_size, i8)
879
880 call mpi_file_read_at_all(fh, byte_offset, data%idx, data%nelv, &
881 mpi_integer, status, ierr)
882
883 mpi_offset = mpi_offset + &
884 int(data%glb_nelv, i8) * int(mpi_integer_size, i8)
885
886 if (read_mesh) then
887 byte_offset = mpi_offset + int(data%offset_el, i8) * &
888 (int(data%gdim*lxyz, i8) * &
889 int(fld_data_size, i8))
890 call fld_file_read_vector_field(this, fh, byte_offset, &
891 data%x, data%y, data%z, data)
892 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
893 (int(data%gdim *lxyz, i8) * &
894 int(fld_data_size, i8))
895 end if
896
897 if (read_velocity) then
898 byte_offset = mpi_offset + int(data%offset_el, i8) * &
899 (int(data%gdim*lxyz, i8) * &
900 int(fld_data_size, i8))
901 call fld_file_read_vector_field(this, fh, byte_offset, &
902 data%u, data%v, data%w, data)
903 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
904 (int(data%gdim *lxyz, i8) * &
905 int(fld_data_size, i8))
906 end if
907
908 if (read_pressure) then
909 byte_offset = mpi_offset + int(data%offset_el, i8) * &
910 (int(lxyz, i8) * &
911 int(fld_data_size, i8))
912 call fld_file_read_field(this, fh, byte_offset, data%p, data)
913 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
914 (int(lxyz, i8) * &
915 int(fld_data_size, i8))
916 end if
917
918 if (read_temp) then
919 byte_offset = mpi_offset + int(data%offset_el, i8) * &
920 (int(lxyz, i8) * &
921 int(fld_data_size, i8))
922 call fld_file_read_field(this, fh, byte_offset, data%t, data)
923 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
924 (int(lxyz, i8) * &
925 int(fld_data_size, i8))
926 end if
927
928 do i = 1, n_scalars
929 byte_offset = mpi_offset + int(data%offset_el, i8) * &
930 (int(lxyz, i8) * &
931 int(fld_data_size, i8))
932 call fld_file_read_field(this, fh, byte_offset, data%s(i), data)
933 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
934 (int(lxyz, i8) * &
935 int(fld_data_size, i8))
936 end do
937
938 call this%increment_counter()
939
940 if (allocated(tmp_dp)) deallocate(tmp_dp)
941 if (allocated(tmp_sp)) deallocate(tmp_sp)
942 class default
943 call neko_error('Currently we only read into fld_file_data_t, &
944 &please use that data structure instead.')
945 end select
946
947 end subroutine fld_file_read
948
949 subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
950 class(fld_file_t), intent(inout) :: this
951 type(vector_t), intent(inout) :: x
952 type(fld_file_data_t) :: fld_data
953 integer(kind=MPI_OFFSET_KIND) :: byte_offset
954 type(mpi_file) :: fh
955 type(mpi_status) :: status
956 integer :: n, ierr, lxyz, i
957
958 n = x%size()
959 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
960
961 if (this%dp_precision) then
962 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, n, &
963 mpi_double_precision, status, ierr)
964 else
965 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, n, &
966 mpi_real, status, ierr)
967 end if
968
969 if (this%dp_precision) then
970 do i = 1, n
971 x%x(i) = tmp_dp(i)
972 end do
973 else
974 do i = 1, n
975 x%x(i) = tmp_sp(i)
976 end do
977 end if
978
979
980 end subroutine fld_file_read_field
981
982
983 subroutine fld_file_read_vector_field(this, fh, byte_offset, &
984 x, y, z, fld_data)
985 class(fld_file_t), intent(inout) :: this
986 type(vector_t), intent(inout) :: x, y, z
987 type(fld_file_data_t) :: fld_data
988 integer(kind=MPI_OFFSET_KIND) :: byte_offset
989 type(mpi_file) :: fh
990 type(mpi_status) :: status
991 integer :: n, ierr, lxyz, i, j, e, nd
992
993 n = x%size()
994 nd = n*fld_data%gdim
995 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
996
997 if (this%dp_precision) then
998 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, nd, &
999 mpi_double_precision, status, ierr)
1000 else
1001 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, nd, &
1002 mpi_real, status, ierr)
1003 end if
1004
1005
1006 if (this%dp_precision) then
1007 i = 1
1008 do e = 1, fld_data%nelv
1009 do j = 1, lxyz
1010 x%x((e-1)*lxyz+j) = tmp_dp(i)
1011 i = i +1
1012 end do
1013 do j = 1, lxyz
1014 y%x((e-1)*lxyz+j) = tmp_dp(i)
1015 i = i +1
1016 end do
1017 if (fld_data%gdim .eq. 3) then
1018 do j = 1, lxyz
1019 z%x((e-1)*lxyz+j) = tmp_dp(i)
1020 i = i +1
1021 end do
1022 end if
1023 end do
1024 else
1025 i = 1
1026 do e = 1, fld_data%nelv
1027 do j = 1, lxyz
1028 x%x((e-1)*lxyz+j) = tmp_sp(i)
1029 i = i +1
1030 end do
1031 do j = 1, lxyz
1032 y%x((e-1)*lxyz+j) = tmp_sp(i)
1033 i = i +1
1034 end do
1035 if (fld_data%gdim .eq. 3) then
1036 do j = 1, lxyz
1037 z%x((e-1)*lxyz+j) = tmp_sp(i)
1038 i = i +1
1039 end do
1040 end if
1041 end do
1042 end if
1043
1044 end subroutine fld_file_read_vector_field
1045
1046 subroutine fld_file_set_precision(this, precision)
1047 class(fld_file_t) :: this
1048 integer, intent(in) :: precision
1049
1050 if (precision .eq. dp) then
1051 this%dp_precision = .true.
1052 else if (precision .eq. sp) then
1053 this%dp_precision = .false.
1054 else
1055 call neko_error('Invalid precision')
1056 end if
1057
1058 end subroutine fld_file_set_precision
1059
1060 function fld_file_get_fld_fname(this) result(fname)
1061 class(fld_file_t), intent(in) :: this
1062 character(len=1024) :: fname
1063 character(len=1024) :: path, name, id_str
1064 integer :: suffix_pos
1065
1066 call filename_path(this%get_base_fname(), path)
1067 call filename_name(this%get_base_fname(), name)
1068
1069 write(fname, '(a,a,a,i5.5)') trim(path), trim(name), &
1070 '0.f', this%get_counter()
1071
1072 end function fld_file_get_fld_fname
1073
1074 function fld_file_get_meta_fname(this) result(fname)
1075 class(fld_file_t), intent(in) :: this
1076 character(len=1024) :: fname
1077 character(len=1024) :: path, name, id_str
1078
1079 call filename_path(this%get_base_fname(), path)
1080 call filename_name(this%get_base_fname(), name)
1081
1082 write(id_str, '(i5,a)') this%get_start_counter(), '.nek5000'
1083 write(fname, '(a,a,a)') trim(path), trim(name), trim(adjustl(id_str))
1084
1085 end function fld_file_get_meta_fname
1086
1087end module fld_file
double real
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
Defines a field.
Definition field.f90:34
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...
NEKTON fld file format.
Definition fld_file.f90:35
character(len=1024) function fld_file_get_fld_fname(this)
real(kind=dp), dimension(:), allocatable, private tmp_dp
Definition fld_file.f90:60
subroutine fld_file_set_precision(this, precision)
subroutine fld_file_read_vector_field(this, fh, byte_offset, x, y, z, fld_data)
Definition fld_file.f90:985
subroutine fld_file_read(this, data)
Load a field from a NEKTON fld file.
Definition fld_file.f90:672
subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, gdim, lxyz, nelv)
Definition fld_file.f90:617
character(len=1024) function fld_file_get_meta_fname(this)
subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, gdim, lxyz, nelv)
Definition fld_file.f90:532
real(kind=sp), dimension(:), allocatable, private tmp_sp
Definition fld_file.f90:61
subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
Definition fld_file.f90:950
subroutine fld_file_write_field(this, fh, byte_offset, p, n)
Definition fld_file.f90:590
subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, nelv)
Definition fld_file.f90:565
subroutine fld_file_write(this, data, t)
Write fields to a NEKTON fld file.
Definition fld_file.f90:80
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
pure logical function, public sabscmp(x, y, tol)
Return single precision absolute comparison .
Definition math.f90:117
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:602
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:591
Defines a mesh.
Definition mesh.f90:34
MPI derived types.
Definition mpi_types.f90:34
integer, public mpi_double_precision_size
Size of MPI type double precision.
Definition mpi_types.f90:66
integer, public mpi_character_size
Size of MPI type character.
Definition mpi_types.f90:67
integer, public mpi_real_size
Size of MPI type real.
Definition mpi_types.f90:65
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:68
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_name(fname, name)
Extract the base name of a file (without path and suffix)
Definition utils.f90:87
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:140
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:58
subroutine, public filename_path(fname, path)
Extract the path to a file.
Definition utils.f90:72
Defines a vector.
Definition vector.f90:34
Load-balanced linear distribution .
Definition datadist.f90:50
field_list_t, To be able to group fields together
Interface for NEKTON fld files.
Definition fld_file.f90:64
A generic file handler.
The function space for the SEM solution fields.
Definition space.f90:63
Pointer to array.
Definition structs.f90:14