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 if (allocated(tempo)) deallocate(tempo)
529 if (allocated(scalar_fields)) deallocate(scalar_fields)
530
531 end subroutine fld_file_write
532
533 subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, &
534 gdim, lxyz, nelv)
535 class(fld_file_t), intent(inout) :: this
536 type(mpi_file), intent(inout) :: fh
537 integer, intent(in) :: gdim, lxyz, nelv
538 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
539 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
540 integer :: el, j, ierr, nout
541 type(mpi_status) :: status
542 real(kind=sp) :: buffer(2*gdim*nelv)
543
544 j = 1
545 do el = 1, nelv
546 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
547 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
548 buffer(j+2) = real(vlmin(y(1, el), lxyz), sp)
549 buffer(j+3) = real(vlmax(y(1, el), lxyz), sp)
550 j = j + 4
551 if (gdim .eq. 3) then
552 buffer(j+0) = real(vlmin(z(1, el), lxyz), sp)
553 buffer(j+1) = real(vlmax(z(1, el), lxyz), sp)
554 j = j + 2
555 end if
556 end do
557
558 ! write out data
559 nout = 2*gdim*nelv
560
561 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
562 mpi_real, status, ierr)
563
564 end subroutine fld_file_write_metadata_vector
565
566 subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, &
567 nelv)
568 class(fld_file_t), intent(inout) :: this
569 type(mpi_file), intent(inout) :: fh
570 integer, intent(in) :: lxyz, nelv
571 real(kind=rp), intent(in) :: x(lxyz, nelv)
572 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
573 integer :: el, j, ierr, nout
574 type(mpi_status) :: status
575 real(kind=sp) :: buffer(2*nelv)
576
577 j = 1
578 do el = 1, nelv
579 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
580 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
581 j = j + 2
582 end do
583
584 ! write out data
585 nout = 2*nelv
586
587 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
588 mpi_real, status, ierr)
589
590 end subroutine fld_file_write_metadata_scalar
591
592 subroutine fld_file_write_field(this, fh, byte_offset, p, n)
593 class(fld_file_t), intent(inout) :: this
594 type(mpi_file), intent(inout) :: fh
595 integer, intent(inout) :: n
596 real(kind=rp), intent(inout) :: p(n)
597 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
598 integer :: i, ierr
599 type(mpi_status) :: status
600
601 if ( this%dp_precision) then
602 do i = 1, n
603 tmp_dp(i) = real(p(i), dp)
604 end do
605
606 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, n, &
607 mpi_double_precision, status, ierr)
608 else
609 do i = 1, n
610 tmp_sp(i) = real(p(i), sp)
611 end do
612 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, n, &
613 mpi_real, status, ierr)
614 end if
615
616 end subroutine fld_file_write_field
617
618 subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, &
619 gdim, lxyz, nelv)
620 class(fld_file_t), intent(inout) :: this
621 type(mpi_file), intent(inout) :: fh
622 integer, intent(in) :: n, gdim, lxyz, nelv
623 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
624 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
625 integer :: i, el, j, ierr
626 type(mpi_status) :: status
627
628 if (this%dp_precision) then
629 i = 1
630 do el = 1, nelv
631 do j = 1, lxyz
632 tmp_dp(i) = real(x(j, el), dp)
633 i = i +1
634 end do
635 do j = 1, lxyz
636 tmp_dp(i) = real(y(j, el), dp)
637 i = i +1
638 end do
639 if (gdim .eq. 3) then
640 do j = 1, lxyz
641 tmp_dp(i) = real(z(j, el), dp)
642 i = i +1
643 end do
644 end if
645 end do
646 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, gdim*n, &
647 mpi_double_precision, status, ierr)
648 else
649 i = 1
650 do el = 1, nelv
651 do j = 1, lxyz
652 tmp_sp(i) = real(x(j, el), sp)
653 i = i +1
654 end do
655 do j = 1, lxyz
656 tmp_sp(i) = real(y(j, el), sp)
657 i = i +1
658 end do
659 if (gdim .eq. 3) then
660 do j = 1, lxyz
661 tmp_sp(i) = real(z(j, el), sp)
662 i = i +1
663 end do
664 end if
665 end do
666 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, gdim*n, &
667 mpi_real, status, ierr)
668 end if
669
670
671 end subroutine fld_file_write_vector_field
672
674 subroutine fld_file_read(this, data)
675 class(fld_file_t) :: this
676 class(*), target, intent(inout) :: data
677 character(len= 132) :: hdr
678 integer :: ierr, suffix_pos, i, j
679 type(mpi_file) :: fh
680 type(mpi_status) :: status
681 character(len= 1024) :: fname, base_fname, meta_fname, string, path
682 logical :: meta_file, read_mesh, read_velocity, read_pressure
683 logical :: read_temp
684 character(len=6) :: suffix
685 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
686 integer :: lx, ly, lz, glb_nelv, counter, lxyz
687 integer :: FLD_DATA_SIZE, n_scalars, n
688 integer :: file_unit
689 real(kind=rp) :: time
690 real(kind=sp) :: temp
691 type(linear_dist_t) :: dist
692 real(kind=sp), parameter :: test_pattern = 6.54321
693 character :: rdcode(10), temp_str(4)
694 character(len=LOG_SIZE) :: log_buf
695
696 select type (data)
697 type is (fld_file_data_t)
698 call filename_chsuffix(this%get_base_fname(), meta_fname, 'nek5000')
699
700 inquire(file = trim(meta_fname), exist = meta_file)
701 if (meta_file .and. data%meta_nsamples .eq. 0) then
702 if (pe_rank .eq. 0) then
703 open(newunit = file_unit, file = trim(meta_fname))
704 read(file_unit, fmt = '(A)') string
705 read(string(14:), fmt = '(A)') string
706 string = trim(string)
707
708 data%fld_series_fname = string(:scan(trim(string), '%')-1)
709 data%fld_series_fname = adjustl(data%fld_series_fname)
710 data%fld_series_fname = trim(data%fld_series_fname)//'0'
711
712 read(file_unit, fmt = '(A)') string
713 read(string(scan(string, ':')+1:), *) data%meta_start_counter
714 read(file_unit, fmt = '(A)') string
715 read(string(scan(string, ':')+1:), *) data%meta_nsamples
716 close(file_unit)
717
718 write(log_buf,*) 'Reading meta file for fld series'
719 call neko_log%message(log_buf)
720 write(log_buf,*) 'Name: ', trim(data%fld_series_fname)
721 call neko_log%message(log_buf)
722 write(log_buf,*) 'Start counter: ', data%meta_start_counter
723 call neko_log%message(log_buf)
724 write(log_buf,*) 'Nsamples: ', data%meta_nsamples
725 call neko_log%message(log_buf)
726
727 end if
728 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
729 neko_comm, ierr)
730 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
731 neko_comm, ierr)
732 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
733 neko_comm, ierr)
734
735 if (this%get_counter() .eq. -1) then
736 call this%set_start_counter(data%meta_start_counter)
737 call this%set_counter(data%meta_start_counter)
738 end if
739 end if
740
741 if (meta_file) then
742 call filename_path(this%get_base_fname(), path)
743 write(suffix, '(a,i5.5)') 'f', this%get_counter()
744 fname = trim(path) // trim(data%fld_series_fname) // '.' // suffix
745 if (this%get_counter() .ge. &
746 data%meta_nsamples+data%meta_start_counter) then
747 call neko_error('Trying to read more fld files than exist')
748 end if
749 else
750 write(suffix, '(a,i5.5)') 'f', this%get_counter()
751 call filename_chsuffix(trim(this%get_base_fname()), fname, suffix)
752 end if
753 call mpi_file_open(neko_comm, trim(fname), &
754 mpi_mode_rdonly, mpi_info_null, fh, ierr)
755
756 if (ierr .ne. 0) call neko_error("Could not read "//trim(fname))
757
758 call neko_log%message('Reading fld file ' // trim(fname))
759
760 call mpi_file_read_all(fh, hdr, 132, mpi_character, status, ierr)
761 ! This read can prorbably be done wihtout the temp variables,
762 ! temp_str, i, j
763
764 read(hdr, 1) temp_str, fld_data_size, lx, ly, lz, glb_nelv, glb_nelv, &
765 time, counter, i, j, (rdcode(i), i = 1, 10)
7661 format(4a, 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
767 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
768 if (data%nelv .eq. 0) then
769 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
770 data%nelv = dist%num_local()
771 data%offset_el = dist%start_idx()
772 end if
773 data%lx = lx
774 data%ly = ly
775 data%lz = lz
776 data%glb_nelv = glb_nelv
777 data%t_counter = counter
778 data%time = time
779 lxyz = lx * ly * lz
780 n = lxyz * data%nelv
781
782 if (lz .eq. 1) then
783 data%gdim = 2
784 else
785 data%gdim = 3
786 end if
787
788
789 if (fld_data_size .eq. mpi_double_precision_size) then
790 this%dp_precision = .true.
791 else
792 this%dp_precision = .false.
793 end if
794 if (this%dp_precision) then
795 allocate(tmp_dp(data%gdim*n))
796 else
797 allocate(tmp_sp(data%gdim*n))
798 end if
799
800
801 i = 1
802 read_mesh = .false.
803 read_velocity = .false.
804 read_pressure = .false.
805 read_temp = .false.
806 if (rdcode(i) .eq. 'X') then
807 read_mesh = .true.
808 call data%x%init(n)
809 call data%y%init(n)
810 call data%z%init(n)
811 i = i + 1
812 end if
813 if (rdcode(i) .eq. 'U') then
814 read_velocity = .true.
815 call data%u%init(n)
816 call data%v%init(n)
817 call data%w%init(n)
818 i = i + 1
819 end if
820 if (rdcode(i) .eq. 'P') then
821 read_pressure = .true.
822 call data%p%init(n)
823 i = i + 1
824 end if
825 if (rdcode(i) .eq. 'T') then
826 read_temp = .true.
827 call data%t%init(n)
828 i = i + 1
829 end if
830 n_scalars = 0
831 if (rdcode(i) .eq. 'S') then
832 i = i + 1
833 read(rdcode(i),*) n_scalars
834 n_scalars = n_scalars*10
835 i = i + 1
836 read(rdcode(i),*) j
837 n_scalars = n_scalars+j
838 i = i + 1
839 if (allocated(data%s)) then
840 if (data%n_scalars .ne. n_scalars) then
841 do j = 1, data%n_scalars
842 call data%s(j)%free()
843 end do
844 deallocate(data%s)
845 data%n_scalars = n_scalars
846 allocate(data%s(n_scalars))
847 do j = 1, data%n_scalars
848 call data%s(j)%init(n)
849 end do
850 end if
851 else
852 data%n_scalars = n_scalars
853 allocate(data%s(data%n_scalars))
854 do j = 1, data%n_scalars
855 call data%s(j)%init(n)
856 end do
857 end if
858 i = i + 1
859 end if
860
861 mpi_offset = 132 * mpi_character_size
862 call mpi_file_read_at_all(fh, mpi_offset, temp, 1, &
863 mpi_real, status, ierr)
864 if (.not. sabscmp(temp, test_pattern, epsilon(1.0_sp))) then
865 call neko_error('Incorrect format for fld file, &
866 &test pattern does not match.')
867 end if
868 mpi_offset = mpi_offset + mpi_real_size
869
870
871 if (allocated(data%idx)) then
872 if (size(data%idx) .ne. data%nelv) then
873 deallocate(data%idx)
874 allocate(data%idx(data%nelv))
875 end if
876 else
877 allocate(data%idx(data%nelv))
878 end if
879
880 byte_offset = mpi_offset + &
881 int(data%offset_el, i8) * int(mpi_integer_size, i8)
882
883 call mpi_file_read_at_all(fh, byte_offset, data%idx, data%nelv, &
884 mpi_integer, status, ierr)
885
886 mpi_offset = mpi_offset + &
887 int(data%glb_nelv, i8) * int(mpi_integer_size, i8)
888
889 if (read_mesh) then
890 byte_offset = mpi_offset + int(data%offset_el, i8) * &
891 (int(data%gdim*lxyz, i8) * &
892 int(fld_data_size, i8))
893 call fld_file_read_vector_field(this, fh, byte_offset, &
894 data%x, data%y, data%z, data)
895 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
896 (int(data%gdim *lxyz, i8) * &
897 int(fld_data_size, i8))
898 end if
899
900 if (read_velocity) then
901 byte_offset = mpi_offset + int(data%offset_el, i8) * &
902 (int(data%gdim*lxyz, i8) * &
903 int(fld_data_size, i8))
904 call fld_file_read_vector_field(this, fh, byte_offset, &
905 data%u, data%v, data%w, data)
906 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
907 (int(data%gdim *lxyz, i8) * &
908 int(fld_data_size, i8))
909 end if
910
911 if (read_pressure) then
912 byte_offset = mpi_offset + int(data%offset_el, i8) * &
913 (int(lxyz, i8) * &
914 int(fld_data_size, i8))
915 call fld_file_read_field(this, fh, byte_offset, data%p, data)
916 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
917 (int(lxyz, i8) * &
918 int(fld_data_size, i8))
919 end if
920
921 if (read_temp) then
922 byte_offset = mpi_offset + int(data%offset_el, i8) * &
923 (int(lxyz, i8) * &
924 int(fld_data_size, i8))
925 call fld_file_read_field(this, fh, byte_offset, data%t, data)
926 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
927 (int(lxyz, i8) * &
928 int(fld_data_size, i8))
929 end if
930
931 do i = 1, n_scalars
932 byte_offset = mpi_offset + int(data%offset_el, i8) * &
933 (int(lxyz, i8) * &
934 int(fld_data_size, i8))
935 call fld_file_read_field(this, fh, byte_offset, data%s(i), data)
936 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
937 (int(lxyz, i8) * &
938 int(fld_data_size, i8))
939 end do
940
941 call this%increment_counter()
942
943 if (allocated(tmp_dp)) deallocate(tmp_dp)
944 if (allocated(tmp_sp)) deallocate(tmp_sp)
945 class default
946 call neko_error('Currently we only read into fld_file_data_t, &
947 &please use that data structure instead.')
948 end select
949
950 end subroutine fld_file_read
951
952 subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
953 class(fld_file_t), intent(inout) :: this
954 type(vector_t), intent(inout) :: x
955 type(fld_file_data_t) :: fld_data
956 integer(kind=MPI_OFFSET_KIND) :: byte_offset
957 type(mpi_file) :: fh
958 type(mpi_status) :: status
959 integer :: n, ierr, lxyz, i
960
961 n = x%size()
962 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
963
964 if (this%dp_precision) then
965 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, n, &
966 mpi_double_precision, status, ierr)
967 else
968 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, n, &
969 mpi_real, status, ierr)
970 end if
971
972 if (this%dp_precision) then
973 do i = 1, n
974 x%x(i) = tmp_dp(i)
975 end do
976 else
977 do i = 1, n
978 x%x(i) = tmp_sp(i)
979 end do
980 end if
981
982
983 end subroutine fld_file_read_field
984
985
986 subroutine fld_file_read_vector_field(this, fh, byte_offset, &
987 x, y, z, fld_data)
988 class(fld_file_t), intent(inout) :: this
989 type(vector_t), intent(inout) :: x, y, z
990 type(fld_file_data_t) :: fld_data
991 integer(kind=MPI_OFFSET_KIND) :: byte_offset
992 type(mpi_file) :: fh
993 type(mpi_status) :: status
994 integer :: n, ierr, lxyz, i, j, e, nd
995
996 n = x%size()
997 nd = n*fld_data%gdim
998 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
999
1000 if (this%dp_precision) then
1001 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, nd, &
1002 mpi_double_precision, status, ierr)
1003 else
1004 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, nd, &
1005 mpi_real, status, ierr)
1006 end if
1007
1008
1009 if (this%dp_precision) then
1010 i = 1
1011 do e = 1, fld_data%nelv
1012 do j = 1, lxyz
1013 x%x((e-1)*lxyz+j) = tmp_dp(i)
1014 i = i +1
1015 end do
1016 do j = 1, lxyz
1017 y%x((e-1)*lxyz+j) = tmp_dp(i)
1018 i = i +1
1019 end do
1020 if (fld_data%gdim .eq. 3) then
1021 do j = 1, lxyz
1022 z%x((e-1)*lxyz+j) = tmp_dp(i)
1023 i = i +1
1024 end do
1025 end if
1026 end do
1027 else
1028 i = 1
1029 do e = 1, fld_data%nelv
1030 do j = 1, lxyz
1031 x%x((e-1)*lxyz+j) = tmp_sp(i)
1032 i = i +1
1033 end do
1034 do j = 1, lxyz
1035 y%x((e-1)*lxyz+j) = tmp_sp(i)
1036 i = i +1
1037 end do
1038 if (fld_data%gdim .eq. 3) then
1039 do j = 1, lxyz
1040 z%x((e-1)*lxyz+j) = tmp_sp(i)
1041 i = i +1
1042 end do
1043 end if
1044 end do
1045 end if
1046
1047 end subroutine fld_file_read_vector_field
1048
1049 subroutine fld_file_set_precision(this, precision)
1050 class(fld_file_t) :: this
1051 integer, intent(in) :: precision
1052
1053 if (precision .eq. dp) then
1054 this%dp_precision = .true.
1055 else if (precision .eq. sp) then
1056 this%dp_precision = .false.
1057 else
1058 call neko_error('Invalid precision')
1059 end if
1060
1061 end subroutine fld_file_set_precision
1062
1063 function fld_file_get_fld_fname(this) result(fname)
1064 class(fld_file_t), intent(in) :: this
1065 character(len=1024) :: fname
1066 character(len=1024) :: path, name, id_str
1067 integer :: suffix_pos
1068
1069 call filename_path(this%get_base_fname(), path)
1070 call filename_name(this%get_base_fname(), name)
1071
1072 write(fname, '(a,a,a,i5.5)') trim(path), trim(name), &
1073 '0.f', this%get_counter()
1074
1075 end function fld_file_get_fld_fname
1076
1077 function fld_file_get_meta_fname(this) result(fname)
1078 class(fld_file_t), intent(in) :: this
1079 character(len=1024) :: fname
1080 character(len=1024) :: path, name, id_str
1081
1082 call filename_path(this%get_base_fname(), path)
1083 call filename_name(this%get_base_fname(), name)
1084
1085 write(id_str, '(i5,a)') this%get_start_counter(), '.nek5000'
1086 write(fname, '(a,a,a)') trim(path), trim(name), trim(adjustl(id_str))
1087
1088 end function fld_file_get_meta_fname
1089
1090end 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:988
subroutine fld_file_read(this, data)
Load a field from a NEKTON fld file.
Definition fld_file.f90:675
subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, gdim, lxyz, nelv)
Definition fld_file.f90:620
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:535
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:953
subroutine fld_file_write_field(this, fh, byte_offset, p, n)
Definition fld_file.f90:593
subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, nelv)
Definition fld_file.f90:568
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