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