Neko 1.99.1
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
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 ! Check if position 5 is a temperature field by name
222 if (trim(data%name(5)) .eq. 'temperature') then
223 ! Position 5 is temperature, remaining fields are scalars
224 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
225 n_scalar_fields = data%size() - 5
226 allocate(scalar_fields(n_scalar_fields))
227 do i = 1, n_scalar_fields
228 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
229 end do
230 write_temperature = .true.
231 else
232 ! All remaining fields are scalars (no temperature field)
233 n_scalar_fields = data%size() - 4
234 allocate(scalar_fields(n_scalar_fields))
235 do i = 1, n_scalar_fields
236 scalar_fields(i)%ptr => data%items(i+4)%ptr%x(:,1,1,1)
237 end do
238 write_temperature = .false.
239 end if
240 write_pressure = .true.
241 write_velocity = .true.
242 case default
243 call neko_error('This many fields not supported yet, fld_file')
244 end select
245 dof => data%dof(1)
246 class default
247 call neko_error('Invalid data')
248 end select
249 ! Fix things for pointers that do not exist in all data types...
250 if (associated(dof)) then
251 x%ptr => dof%x(:,1,1,1)
252 y%ptr => dof%y(:,1,1,1)
253 z%ptr => dof%z(:,1,1,1)
254 msh => dof%msh
255 xh => dof%Xh
256 end if
257
258 if (associated(msh)) then
259 nelv = msh%nelv
260 glb_nelv = msh%glb_nelv
261 offset_el = msh%offset_el
262 gdim = msh%gdim
263 ! Store global idx of each element
264 allocate(idx(msh%nelv))
265 do i = 1, msh%nelv
266 idx(i) = msh%elements(i)%e%id()
267 end do
268 end if
269
270 if (associated(xh)) then
271 lx = xh%lx
272 ly = xh%ly
273 lz = xh%lz
274 end if
275
276 lxyz = lx*ly*lz
277 n = nelv*lxyz
278
279 if (this%dp_precision) then
280 fld_data_size = mpi_double_precision_size
281 else
282 fld_data_size = mpi_real_size
283 end if
284 if (this%dp_precision) then
285 allocate(tmp_dp(gdim*n))
286 else
287 allocate(tmp_sp(gdim*n))
288 end if
289
290
291 !
292 ! Create fld header for NEKTON's multifile output
293 !
294
295 call this%increment_counter()
296 write_mesh = (this%get_counter() .eq. this%get_start_counter())
297 call mpi_allreduce(mpi_in_place, write_mesh, 1, &
298 mpi_logical, mpi_lor, neko_comm)
299 call mpi_allreduce(mpi_in_place, write_velocity, 1, &
300 mpi_logical, mpi_lor, neko_comm)
301 call mpi_allreduce(mpi_in_place, write_pressure, 1, &
302 mpi_logical, mpi_lor, neko_comm)
303 call mpi_allreduce(mpi_in_place, write_temperature, 1, &
304 mpi_logical, mpi_lor, neko_comm)
305 call mpi_allreduce(mpi_in_place, n_scalar_fields, 1, &
306 mpi_integer, mpi_max, neko_comm)
307
308 ! Build rdcode note that for field_t, we only support scalar
309 ! fields at the moment
310 rdcode = ' '
311 i = 1
312 if (write_mesh) then
313 rdcode(i) = 'X'
314 i = i + 1
315 end if
316 if (write_velocity) then
317 rdcode(i) = 'U'
318 i = i + 1
319 end if
320 if (write_pressure) then
321 rdcode(i) = 'P'
322 i = i + 1
323 end if
324 if (write_temperature) then
325 rdcode(i) = 'T'
326 i = i + 1
327 end if
328 if (n_scalar_fields .gt. 0 ) then
329 rdcode(i) = 'S'
330 i = i + 1
331 write(rdcode(i), '(i1)') (n_scalar_fields)/10
332 i = i + 1
333 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
334 i = i + 1
335 end if
336
338 write(hdr, 1) fld_data_size, lx, ly, lz, glb_nelv, glb_nelv,&
339 time, this%get_counter(), 1, 1, (rdcode(i), i = 1, 10)
3401 format('#std', 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
341 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
342
343 ! Change to NEKTON's fld file format
344 fname = this%get_fld_fname()
345
346 call mpi_file_open(neko_comm, trim(fname), &
347 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, &
348 ierr)
349
350 call mpi_file_write_all(fh, hdr, 132, mpi_character, status, ierr)
351 mpi_offset = 132 * mpi_character_size
352
353 call mpi_file_write_all(fh, test_pattern, 1, mpi_real, status, ierr)
354 mpi_offset = mpi_offset + mpi_real_size
355
356 byte_offset = mpi_offset + &
357 int(offset_el, i8) * int(mpi_integer_size, i8)
358 call mpi_file_write_at_all(fh, byte_offset, idx, nelv, &
359 mpi_integer, status, ierr)
360 mpi_offset = mpi_offset + int(glb_nelv, i8) * int(mpi_integer_size, i8)
361 deallocate(idx)
362 if (write_mesh) then
363
364 byte_offset = mpi_offset + int(offset_el, i8) * &
365 (int(gdim*lxyz, i8) * &
366 int(fld_data_size, i8))
367 call fld_file_write_vector_field(this, fh, byte_offset, &
368 x%ptr, y%ptr, z%ptr, &
369 n, gdim, lxyz, nelv)
370 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
371 (int(gdim *lxyz, i8) * &
372 int(fld_data_size, i8))
373 end if
374 if (write_velocity) then
375 byte_offset = mpi_offset + int(offset_el, i8) * &
376 (int(gdim * (lxyz), i8) * int(fld_data_size, i8))
377 call fld_file_write_vector_field(this, fh, byte_offset, &
378 u%ptr, v%ptr, w%ptr, n, gdim, lxyz, nelv)
379
380 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
381 (int(gdim * (lxyz), i8) * &
382 int(fld_data_size, i8))
383
384 end if
385
386 if (write_pressure) then
387 byte_offset = mpi_offset + int(offset_el, i8) * &
388 (int((lxyz), i8) * int(fld_data_size, i8))
389 call fld_file_write_field(this, fh, byte_offset, p%ptr, n)
390 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
391 (int((lxyz), i8) * int(fld_data_size, i8))
392 end if
393
394 if (write_temperature) then
395 byte_offset = mpi_offset + int(offset_el, i8) * &
396 (int((lxyz), i8) * &
397 int(fld_data_size, i8))
398 call fld_file_write_field(this, fh, byte_offset, tem%ptr, n)
399 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
400 (int((lxyz), i8) * &
401 int(fld_data_size, i8))
402 end if
403
404 temp_offset = mpi_offset
405
406 do i = 1, n_scalar_fields
407 ! Without this redundant if statement Cray optimizes this loop to
408 ! Oblivion
409 if (i .eq. 2) then
410 mpi_offset = int(temp_offset, i8) + int(1_i8*glb_nelv, i8) * &
411 (int(lxyz, i8) * int(fld_data_size, i8))
412 end if
413 byte_offset = int(mpi_offset, i8) + int(offset_el, i8) * &
414 (int((lxyz), i8) * &
415 int(fld_data_size, i8))
416 call fld_file_write_field(this, fh, byte_offset, scalar_fields(i)%ptr, n)
417 mpi_offset = int(mpi_offset, i8) + int(glb_nelv, i8) * &
418 (int(lxyz, i8) * &
419 int(fld_data_size, i8))
420 end do
421
422 if (gdim .eq. 3) then
423
425 if (write_mesh) then
426 !The offset is:
427 ! mpioff + element_off * 2(min max value)
428 ! * 4(single precision) * gdim(dimensions)
429 byte_offset = int(mpi_offset, i8) + &
430 int(offset_el, i8) * &
431 int(2, i8) * &
432 int(mpi_real_size, i8) * &
433 int(gdim, i8)
434 call fld_file_write_metadata_vector(this, fh, byte_offset, &
435 x%ptr, y%ptr, z%ptr, gdim, lxyz, nelv)
436 mpi_offset = int(mpi_offset, i8) + &
437 int(glb_nelv, i8) * &
438 int(2, i8) * &
439 int(mpi_real_size, i8) * &
440 int(gdim, i8)
441 end if
442
443 if (write_velocity) then
444 byte_offset = int(mpi_offset, i8) + &
445 int(offset_el, i8) * &
446 int(2, i8) * &
447 int(mpi_real_size, i8) * &
448 int(gdim, i8)
449 call fld_file_write_metadata_vector(this, fh, byte_offset, &
450 u%ptr, v%ptr, w%ptr, gdim, lxyz, nelv)
451 mpi_offset = int(mpi_offset, i8) + &
452 int(glb_nelv, i8) * &
453 int(2, i8) * &
454 int(mpi_real_size, i8) * &
455 int(gdim, i8)
456
457 end if
458
459 if (write_pressure) then
460 byte_offset = int(mpi_offset, i8) + &
461 int(offset_el, i8) * &
462 int(2, i8) * &
463 int(mpi_real_size, i8)
464 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
465 p%ptr, lxyz, nelv)
466 mpi_offset = int(mpi_offset, i8) + &
467 int(glb_nelv, i8) * &
468 int(2, i8) * &
469 int(mpi_real_size, i8)
470
471 end if
472
473 if (write_temperature) then
474 byte_offset = int(mpi_offset, i8) + &
475 int(offset_el, i8) * &
476 int(2, i8) * &
477 int(mpi_real_size, i8)
478 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
479 tem%ptr, lxyz, nelv)
480 mpi_offset = int(mpi_offset, i8) + &
481 int(glb_nelv, i8) * &
482 int(2, i8) * &
483 int(mpi_real_size, i8)
484
485 end if
486
487
488
489 temp_offset = mpi_offset
490
491 do i = 1, n_scalar_fields
492 ! Without this redundant if statement, Cray optimizes this loop to
493 ! Oblivion
494 if (i .eq. 2) then
495 mpi_offset = int(temp_offset, i8) + &
496 int(1_i8*glb_nelv, i8) * &
497 int(2, i8) * &
498 int(mpi_real_size, i8)
499 end if
500
501 byte_offset = int(mpi_offset, i8) + &
502 int(offset_el, i8) * &
503 int(2, i8) * &
504 int(mpi_real_size, i8)
505 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
506 scalar_fields(i)%ptr, lxyz, nelv)
507 mpi_offset = int(mpi_offset, i8) + &
508 int(glb_nelv, i8) * &
509 int(2, i8) * &
510 int(mpi_real_size, i8)
511 end do
512 end if
513
514
515 call mpi_file_sync(fh, ierr)
516 call mpi_file_close(fh, ierr)
517 ! Write metadata file
518 if (pe_rank .eq. 0) then
519 call filename_name(this%get_base_fname(), name)
520
521 open(newunit = file_unit, &
522 file = this%get_meta_fname(), status = 'replace')
523 ! The following string will specify that the files in the file series
524 ! are defined by the filename followed by a 0.
525 ! This 0 is necessary as it specifies the index of number of files
526 ! the output file is split across.
527 ! In the past, many .f files were generated for each write.
528 ! To be consistent with this the trailing 0 is still necessary today.
529 write(file_unit, fmt = '(A,A,A)') 'filetemplate: ', &
530 trim(name), '%01d.f%05d'
531 write(file_unit, fmt = '(A,i5)') 'firsttimestep: ', &
532 this%get_start_counter()
533 write(file_unit, fmt = '(A,i5)') 'numtimesteps: ', &
534 (this%get_counter() + 1) - this%get_start_counter()
535 close(file_unit)
536 end if
537
538 if (allocated(tmp_dp)) deallocate(tmp_dp)
539 if (allocated(tmp_sp)) deallocate(tmp_sp)
540 end subroutine fld_file_write
541
542 subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, &
543 gdim, lxyz, nelv)
544 class(fld_file_t), intent(inout) :: this
545 type(mpi_file), intent(inout) :: fh
546 integer, intent(in) :: gdim, lxyz, nelv
547 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
548 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
549 integer :: el, j, ierr, nout
550 type(mpi_status) :: status
551 real(kind=sp) :: buffer(2*gdim*nelv)
552
553 j = 1
554 do el = 1, nelv
555 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
556 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
557 buffer(j+2) = real(vlmin(y(1, el), lxyz), sp)
558 buffer(j+3) = real(vlmax(y(1, el), lxyz), sp)
559 j = j + 4
560 if (gdim .eq. 3) then
561 buffer(j+0) = real(vlmin(z(1, el), lxyz), sp)
562 buffer(j+1) = real(vlmax(z(1, el), lxyz), sp)
563 j = j + 2
564 end if
565 end do
566
567 ! write out data
568 nout = 2*gdim*nelv
569
570 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
571 mpi_real, status, ierr)
572
573 end subroutine fld_file_write_metadata_vector
574
575 subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, &
576 nelv)
577 class(fld_file_t), intent(inout) :: this
578 type(mpi_file), intent(inout) :: fh
579 integer, intent(in) :: lxyz, nelv
580 real(kind=rp), intent(in) :: x(lxyz, nelv)
581 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
582 integer :: el, j, ierr, nout
583 type(mpi_status) :: status
584 real(kind=sp) :: buffer(2*nelv)
585
586 j = 1
587 do el = 1, nelv
588 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
589 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
590 j = j + 2
591 end do
592
593 ! write out data
594 nout = 2*nelv
595
596 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
597 mpi_real, status, ierr)
598
599 end subroutine fld_file_write_metadata_scalar
600
601 subroutine fld_file_write_field(this, fh, byte_offset, p, n)
602 class(fld_file_t), intent(inout) :: this
603 type(mpi_file), intent(inout) :: fh
604 integer, intent(inout) :: n
605 real(kind=rp), intent(inout) :: p(n)
606 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
607 integer :: i, ierr
608 type(mpi_status) :: status
609
610 if ( this%dp_precision) then
611 do i = 1, n
612 tmp_dp(i) = real(p(i), dp)
613 end do
614
615 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, n, &
616 mpi_double_precision, status, ierr)
617 else
618 do i = 1, n
619 tmp_sp(i) = real(p(i), sp)
620 end do
621 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, n, &
622 mpi_real, status, ierr)
623 end if
624
625 end subroutine fld_file_write_field
626
627 subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, &
628 gdim, lxyz, nelv)
629 class(fld_file_t), intent(inout) :: this
630 type(mpi_file), intent(inout) :: fh
631 integer, intent(in) :: n, gdim, lxyz, nelv
632 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
633 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
634 integer :: i, el, j, ierr
635 type(mpi_status) :: status
636
637 if (this%dp_precision) then
638 i = 1
639 do el = 1, nelv
640 do j = 1, lxyz
641 tmp_dp(i) = real(x(j, el), dp)
642 i = i +1
643 end do
644 do j = 1, lxyz
645 tmp_dp(i) = real(y(j, el), dp)
646 i = i +1
647 end do
648 if (gdim .eq. 3) then
649 do j = 1, lxyz
650 tmp_dp(i) = real(z(j, el), dp)
651 i = i +1
652 end do
653 end if
654 end do
655 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, gdim*n, &
656 mpi_double_precision, status, ierr)
657 else
658 i = 1
659 do el = 1, nelv
660 do j = 1, lxyz
661 tmp_sp(i) = real(x(j, el), sp)
662 i = i +1
663 end do
664 do j = 1, lxyz
665 tmp_sp(i) = real(y(j, el), sp)
666 i = i +1
667 end do
668 if (gdim .eq. 3) then
669 do j = 1, lxyz
670 tmp_sp(i) = real(z(j, el), sp)
671 i = i +1
672 end do
673 end if
674 end do
675 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, gdim*n, &
676 mpi_real, status, ierr)
677 end if
678
679
680 end subroutine fld_file_write_vector_field
681
683 subroutine fld_file_read(this, data)
684 class(fld_file_t) :: this
685 class(*), target, intent(inout) :: data
686 character(len= 132) :: hdr
687 integer :: ierr, suffix_pos, i, j
688 type(mpi_file) :: fh
689 type(mpi_status) :: status
690 character(len= 1024) :: fname, base_fname, meta_fname, string, path
691 logical :: meta_file, read_mesh, read_velocity, read_pressure
692 logical :: read_temp
693 character(len=6) :: suffix
694 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
695 integer :: lx, ly, lz, glb_nelv, counter, lxyz
696 integer :: FLD_DATA_SIZE, n_scalars, n
697 integer :: file_unit
698 real(kind=rp) :: time
699 real(kind=sp) :: temp
700 type(linear_dist_t) :: dist
701 real(kind=sp), parameter :: test_pattern = 6.54321
702 character :: rdcode(10), temp_str(4)
703 character(len=LOG_SIZE) :: log_buf
704
705 select type (data)
706 type is (fld_file_data_t)
707 call filename_chsuffix(this%get_base_fname(), meta_fname, 'nek5000')
708
709 inquire(file = trim(meta_fname), exist = meta_file)
710 if (meta_file .and. data%meta_nsamples .eq. 0) then
711 if (pe_rank .eq. 0) then
712 open(newunit = file_unit, file = trim(meta_fname))
713 read(file_unit, fmt = '(A)') string
714 read(string(14:), fmt = '(A)') string
715 string = trim(string)
716
717 data%fld_series_fname = string(:scan(trim(string), '%')-1)
718 data%fld_series_fname = adjustl(data%fld_series_fname)
719 data%fld_series_fname = trim(data%fld_series_fname)//'0'
720
721 read(file_unit, fmt = '(A)') string
722 read(string(scan(string, ':')+1:), *) data%meta_start_counter
723 read(file_unit, fmt = '(A)') string
724 read(string(scan(string, ':')+1:), *) data%meta_nsamples
725 close(file_unit)
726
727 write(log_buf,*) 'Reading meta file for fld series'
728 call neko_log%message(log_buf)
729 write(log_buf,*) 'Name: ', trim(data%fld_series_fname)
730 call neko_log%message(log_buf)
731 write(log_buf,*) 'Start counter: ', data%meta_start_counter
732 call neko_log%message(log_buf)
733 write(log_buf,*) 'Nsamples: ', data%meta_nsamples
734 call neko_log%message(log_buf)
735
736 end if
737 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
738 neko_comm, ierr)
739 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
740 neko_comm, ierr)
741 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
742 neko_comm, ierr)
743
744 if (this%get_counter() .eq. -1) then
745 call this%set_start_counter(data%meta_start_counter)
746 call this%set_counter(data%meta_start_counter)
747 end if
748 end if
749
750 if (meta_file) then
751 call filename_path(this%get_base_fname(), path)
752 write(suffix, '(a,i5.5)') 'f', this%get_counter()
753 fname = trim(path) // trim(data%fld_series_fname) // '.' // suffix
754 if (this%get_counter() .ge. &
755 data%meta_nsamples+data%meta_start_counter) then
756 call neko_error('Trying to read more fld files than exist')
757 end if
758 else
759 write(suffix, '(a,i5.5)') 'f', this%get_counter()
760 call filename_chsuffix(trim(this%get_base_fname()), fname, suffix)
761 end if
762 call mpi_file_open(neko_comm, trim(fname), &
763 mpi_mode_rdonly, mpi_info_null, fh, ierr)
764
765 if (ierr .ne. 0) call neko_error("Could not read "//trim(fname))
766
767 call neko_log%message('Reading fld file ' // trim(fname))
768
769 call mpi_file_read_all(fh, hdr, 132, mpi_character, status, ierr)
770 ! This read can prorbably be done wihtout the temp variables,
771 ! temp_str, i, j
772
773 read(hdr, 1) temp_str, fld_data_size, lx, ly, lz, glb_nelv, glb_nelv, &
774 time, counter, i, j, (rdcode(i), i = 1, 10)
7751 format(4a, 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
776 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
777 if (data%nelv .eq. 0) then
778 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
779 data%nelv = dist%num_local()
780 data%offset_el = dist%start_idx()
781 end if
782 data%lx = lx
783 data%ly = ly
784 data%lz = lz
785 data%glb_nelv = glb_nelv
786 data%t_counter = counter
787 data%time = time
788 lxyz = lx * ly * lz
789 n = lxyz * data%nelv
790
791 if (lz .eq. 1) then
792 data%gdim = 2
793 else
794 data%gdim = 3
795 end if
796
797
798 if (fld_data_size .eq. mpi_double_precision_size) then
799 this%dp_precision = .true.
800 else
801 this%dp_precision = .false.
802 end if
803 if (this%dp_precision) then
804 allocate(tmp_dp(data%gdim*n))
805 else
806 allocate(tmp_sp(data%gdim*n))
807 end if
808
809
810 i = 1
811 read_mesh = .false.
812 read_velocity = .false.
813 read_pressure = .false.
814 read_temp = .false.
815 if (rdcode(i) .eq. 'X') then
816 read_mesh = .true.
817 call data%x%init(n)
818 call data%y%init(n)
819 call data%z%init(n)
820 i = i + 1
821 end if
822 if (rdcode(i) .eq. 'U') then
823 read_velocity = .true.
824 call data%u%init(n)
825 call data%v%init(n)
826 call data%w%init(n)
827 i = i + 1
828 end if
829 if (rdcode(i) .eq. 'P') then
830 read_pressure = .true.
831 call data%p%init(n)
832 i = i + 1
833 end if
834 if (rdcode(i) .eq. 'T') then
835 read_temp = .true.
836 call data%t%init(n)
837 i = i + 1
838 end if
839 n_scalars = 0
840 if (rdcode(i) .eq. 'S') then
841 i = i + 1
842 read(rdcode(i),*) n_scalars
843 n_scalars = n_scalars*10
844 i = i + 1
845 read(rdcode(i),*) j
846 n_scalars = n_scalars+j
847 i = i + 1
848 if (allocated(data%s)) then
849 if (data%n_scalars .ne. n_scalars) then
850 do j = 1, data%n_scalars
851 call data%s(j)%free()
852 end do
853 deallocate(data%s)
854 data%n_scalars = n_scalars
855 allocate(data%s(n_scalars))
856 do j = 1, data%n_scalars
857 call data%s(j)%init(n)
858 end do
859 end if
860 else
861 data%n_scalars = n_scalars
862 allocate(data%s(data%n_scalars))
863 do j = 1, data%n_scalars
864 call data%s(j)%init(n)
865 end do
866 end if
867 i = i + 1
868 end if
869
870 mpi_offset = 132 * mpi_character_size
871 call mpi_file_read_at_all(fh, mpi_offset, temp, 1, &
872 mpi_real, status, ierr)
873 if (temp .ne. test_pattern) then
874 call neko_error('Incorrect format for fld file, &
875 &test pattern does not match.')
876 end if
877 mpi_offset = mpi_offset + mpi_real_size
878
879
880 if (allocated(data%idx)) then
881 if (size(data%idx) .ne. data%nelv) then
882 deallocate(data%idx)
883 allocate(data%idx(data%nelv))
884 end if
885 else
886 allocate(data%idx(data%nelv))
887 end if
888
889 byte_offset = mpi_offset + &
890 int(data%offset_el, i8) * int(mpi_integer_size, i8)
891
892 call mpi_file_read_at_all(fh, byte_offset, data%idx, data%nelv, &
893 mpi_integer, status, ierr)
894
895 mpi_offset = mpi_offset + &
896 int(data%glb_nelv, i8) * int(mpi_integer_size, i8)
897
898 if (read_mesh) then
899 byte_offset = mpi_offset + int(data%offset_el, i8) * &
900 (int(data%gdim*lxyz, i8) * &
901 int(fld_data_size, i8))
902 call fld_file_read_vector_field(this, fh, byte_offset, &
903 data%x, data%y, data%z, data)
904 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
905 (int(data%gdim *lxyz, i8) * &
906 int(fld_data_size, i8))
907 end if
908
909 if (read_velocity) then
910 byte_offset = mpi_offset + int(data%offset_el, i8) * &
911 (int(data%gdim*lxyz, i8) * &
912 int(fld_data_size, i8))
913 call fld_file_read_vector_field(this, fh, byte_offset, &
914 data%u, data%v, data%w, data)
915 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
916 (int(data%gdim *lxyz, i8) * &
917 int(fld_data_size, i8))
918 end if
919
920 if (read_pressure) then
921 byte_offset = mpi_offset + int(data%offset_el, i8) * &
922 (int(lxyz, i8) * &
923 int(fld_data_size, i8))
924 call fld_file_read_field(this, fh, byte_offset, data%p, data)
925 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
926 (int(lxyz, i8) * &
927 int(fld_data_size, i8))
928 end if
929
930 if (read_temp) then
931 byte_offset = mpi_offset + int(data%offset_el, i8) * &
932 (int(lxyz, i8) * &
933 int(fld_data_size, i8))
934 call fld_file_read_field(this, fh, byte_offset, data%t, data)
935 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
936 (int(lxyz, i8) * &
937 int(fld_data_size, i8))
938 end if
939
940 do i = 1, n_scalars
941 byte_offset = mpi_offset + int(data%offset_el, i8) * &
942 (int(lxyz, i8) * &
943 int(fld_data_size, i8))
944 call fld_file_read_field(this, fh, byte_offset, data%s(i), data)
945 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
946 (int(lxyz, i8) * &
947 int(fld_data_size, i8))
948 end do
949
950 call this%increment_counter()
951
952 if (allocated(tmp_dp)) deallocate(tmp_dp)
953 if (allocated(tmp_sp)) deallocate(tmp_sp)
954 class default
955 call neko_error('Currently we only read into fld_file_data_t, &
956 &please use that data structure instead.')
957 end select
958
959 end subroutine fld_file_read
960
961 subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
962 class(fld_file_t), intent(inout) :: this
963 type(vector_t), intent(inout) :: x
964 type(fld_file_data_t) :: fld_data
965 integer(kind=MPI_OFFSET_KIND) :: byte_offset
966 type(mpi_file) :: fh
967 type(mpi_status) :: status
968 integer :: n, ierr, lxyz, i
969
970 n = x%size()
971 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
972
973 if (this%dp_precision) then
974 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, n, &
975 mpi_double_precision, status, ierr)
976 else
977 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, n, &
978 mpi_real, status, ierr)
979 end if
980
981 if (this%dp_precision) then
982 do i = 1, n
983 x%x(i) = tmp_dp(i)
984 end do
985 else
986 do i = 1, n
987 x%x(i) = tmp_sp(i)
988 end do
989 end if
990
991
992 end subroutine fld_file_read_field
993
994
995 subroutine fld_file_read_vector_field(this, fh, byte_offset, &
996 x, y, z, fld_data)
997 class(fld_file_t), intent(inout) :: this
998 type(vector_t), intent(inout) :: x, y, z
999 type(fld_file_data_t) :: fld_data
1000 integer(kind=MPI_OFFSET_KIND) :: byte_offset
1001 type(mpi_file) :: fh
1002 type(mpi_status) :: status
1003 integer :: n, ierr, lxyz, i, j, e, nd
1004
1005 n = x%size()
1006 nd = n*fld_data%gdim
1007 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
1008
1009 if (this%dp_precision) then
1010 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, nd, &
1011 mpi_double_precision, status, ierr)
1012 else
1013 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, nd, &
1014 mpi_real, status, ierr)
1015 end if
1016
1017
1018 if (this%dp_precision) then
1019 i = 1
1020 do e = 1, fld_data%nelv
1021 do j = 1, lxyz
1022 x%x((e-1)*lxyz+j) = tmp_dp(i)
1023 i = i +1
1024 end do
1025 do j = 1, lxyz
1026 y%x((e-1)*lxyz+j) = tmp_dp(i)
1027 i = i +1
1028 end do
1029 if (fld_data%gdim .eq. 3) then
1030 do j = 1, lxyz
1031 z%x((e-1)*lxyz+j) = tmp_dp(i)
1032 i = i +1
1033 end do
1034 end if
1035 end do
1036 else
1037 i = 1
1038 do e = 1, fld_data%nelv
1039 do j = 1, lxyz
1040 x%x((e-1)*lxyz+j) = tmp_sp(i)
1041 i = i +1
1042 end do
1043 do j = 1, lxyz
1044 y%x((e-1)*lxyz+j) = tmp_sp(i)
1045 i = i +1
1046 end do
1047 if (fld_data%gdim .eq. 3) then
1048 do j = 1, lxyz
1049 z%x((e-1)*lxyz+j) = tmp_sp(i)
1050 i = i +1
1051 end do
1052 end if
1053 end do
1054 end if
1055
1056 end subroutine fld_file_read_vector_field
1057
1058 subroutine fld_file_set_precision(this, precision)
1059 class(fld_file_t) :: this
1060 integer, intent(in) :: precision
1061
1062 if (precision .eq. dp) then
1063 this%dp_precision = .true.
1064 else if (precision .eq. sp) then
1065 this%dp_precision = .false.
1066 else
1067 call neko_error('Invalid precision')
1068 end if
1069
1070 end subroutine fld_file_set_precision
1071
1072 function fld_file_get_fld_fname(this) result(fname)
1073 class(fld_file_t), intent(in) :: this
1074 character(len=1024) :: fname
1075 character(len=1024) :: path, name, id_str
1076 integer :: suffix_pos
1077
1078 call filename_path(this%get_base_fname(), path)
1079 call filename_name(this%get_base_fname(), name)
1080
1081 write(fname, '(a,a,a,i5.5)') trim(path), trim(name), &
1082 '0.f', this%get_counter()
1083
1084 end function fld_file_get_fld_fname
1085
1086 function fld_file_get_meta_fname(this) result(fname)
1087 class(fld_file_t), intent(in) :: this
1088 character(len=1024) :: fname
1089 character(len=1024) :: path, name, id_str
1090
1091 call filename_path(this%get_base_fname(), path)
1092 call filename_name(this%get_base_fname(), name)
1093
1094 write(id_str, '(i5,a)') this%get_start_counter(), '.nek5000'
1095 write(fname, '(a,a,a)') trim(path), trim(name), trim(adjustl(id_str))
1096
1097 end function fld_file_get_meta_fname
1098
1099end 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: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
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:997
subroutine fld_file_read(this, data)
Load a field from a NEKTON fld file.
Definition fld_file.f90:684
subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, gdim, lxyz, nelv)
Definition fld_file.f90:629
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:544
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:962
subroutine fld_file_write_field(this, fh, byte_offset, p, n)
Definition fld_file.f90:602
subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, nelv)
Definition fld_file.f90:577
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:70
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
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:63
integer, public mpi_character_size
Size of MPI type character.
Definition mpi_types.f90:64
integer, public mpi_real_size
Size of MPI type real.
Definition mpi_types.f90:62
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:65
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