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