Neko 1.99.3
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
49 use mask, only : mask_t
52 use comm
53 use datadist, only : linear_dist_t
54 use math, only : vlmin, vlmax, sabscmp
58 use mpi_f08
59 implicit none
60 private
61
62 real(kind=dp), private, allocatable :: tmp_dp(:)
63 real(kind=sp), private, allocatable :: tmp_sp(:)
64
66 type, public, extends(generic_file_t) :: fld_file_t
67 logical :: dp_precision = .false.
68 logical :: write_mesh = .false.
69 type(mask_t) :: mask
70 contains
71 procedure :: read => fld_file_read
72 procedure :: write => fld_file_write_manager
73 procedure, private :: write_all => fld_file_write
74 procedure, private :: write_masked => fld_file_write_masked
75 procedure :: set_precision => fld_file_set_precision
76 procedure :: set_mask => fld_file_set_mask
77 procedure :: get_fld_fname => fld_file_get_fld_fname
78 procedure :: get_meta_fname => fld_file_get_meta_fname
79 end type fld_file_t
80
81
82contains
83
85 subroutine fld_file_write_manager(this, data, t)
86 class(fld_file_t), intent(inout) :: this
87 class(*), target, intent(in) :: data
88 real(kind=rp), intent(in), optional :: t
89
90 if (this%mask%is_set()) then
91 call this%write_masked(data, this%mask, t)
92 else
93 call this%write_all(data, t)
94 end if
95
96 end subroutine fld_file_write_manager
97
100 subroutine fld_file_write(this, data, t)
101 class(fld_file_t), intent(inout) :: this
102 class(*), target, intent(in) :: data
103 real(kind=rp), intent(in), optional :: t
104 type(array_ptr_t) :: x, y, z, u, v, w, p, tem
105 real(kind=rp), allocatable, target :: tempo(:)
106 type(mesh_t), pointer :: msh
107 type(dofmap_t), pointer :: dof
108 type(space_t), pointer :: Xh
109 real(kind=dp) :: time
110 character(len= 132) :: hdr
111 character :: rdcode(10)
112 character(len=6) :: id_str
113 character(len= 1024) :: fname
114 character(len= 1024) :: name
115 integer :: file_unit
116 integer :: i, ierr, n, suffix_pos, tslash_pos
117 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
118 integer, allocatable :: idx(:)
119 type(mpi_status) :: status
120 type(mpi_file) :: fh
121 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset, temp_offset
122 real(kind=sp), parameter :: test_pattern = 6.54321
123 type(array_ptr_t), allocatable :: scalar_fields(:)
124 logical :: write_mesh, write_velocity, write_pressure, write_temperature
125 integer :: fld_data_size, n_scalar_fields
126
127 if (present(t)) then
128 time = real(t, dp)
129 else
130 time = 0d0
131 end if
132
133 nullify(msh)
134 nullify(dof)
135 nullify(xh)
136 n_scalar_fields = 0
137 write_pressure = .false.
138 write_velocity = .false.
139 write_temperature = .false.
140
141 select type (data)
142 type is (fld_file_data_t)
143 nelv = data%nelv
144 lx = data%lx
145 ly = data%ly
146 lz = data%lz
147 gdim = data%gdim
148 glb_nelv = data%glb_nelv
149 offset_el = data%offset_el
150
151 if (data%x%size() .gt. 0) x%ptr => data%x%x
152 if (data%y%size() .gt. 0) y%ptr => data%y%x
153 if (data%z%size() .gt. 0) z%ptr => data%z%x
154 if (gdim .eq. 2) z%ptr => data%y%x
155 if (data%u%size() .gt. 0) then
156 u%ptr => data%u%x
157 ! In case only u is actually allocated, point the other comps to u
158 ! so that we don't die on trying to write them
159 if (data%v%size() .le. 0) v%ptr => data%u%x
160 if (data%w%size() .le. 0) w%ptr => data%u%x
161 write_velocity = .true.
162 end if
163 if (data%v%size() .gt. 0) v%ptr => data%v%x
164 if (data%w%size() .gt. 0) w%ptr => data%w%x
165 if (data%p%size() .gt. 0) then
166 p%ptr => data%p%x
167 write_pressure = .true.
168 end if
169 if (data%t%size() .gt. 0) then
170 write_temperature = .true.
171 tem%ptr => data%t%x
172 end if
173 ! If gdim = 2 and Z-velocity component exists,
174 ! it is stored in last scalar field
175 if (gdim .eq. 2 .and. data%w%size() .gt. 0) then
176 n_scalar_fields = data%n_scalars + 1
177 allocate(scalar_fields(n_scalar_fields))
178 do i = 1, n_scalar_fields -1
179 scalar_fields(i)%ptr => data%s(i)%x
180 end do
181 scalar_fields(n_scalar_fields)%ptr => data%w%x
182 else
183 n_scalar_fields = data%n_scalars
184 allocate(scalar_fields(n_scalar_fields+1))
185 do i = 1, n_scalar_fields
186 scalar_fields(i)%ptr => data%s(i)%x
187 end do
188 scalar_fields(n_scalar_fields+1)%ptr => data%w%x
189 end if
190 ! This is very stupid...
191 ! Some compilers cannot handle that these pointers dont point to anything
192 ! (although they are not used) this fixes a segfault due to this.
193 if (nelv .eq. 0) then
194 allocate(tempo(1))
195 x%ptr => tempo
196 y%ptr => tempo
197 z%ptr => tempo
198 u%ptr => tempo
199 v%ptr => tempo
200 w%ptr => tempo
201 p%ptr => tempo
202 tem%ptr => tempo
203 end if
204
205 allocate(idx(nelv))
206 do i = 1, nelv
207 idx(i) = data%idx(i)
208 end do
209 type is (field_t)
210 p%ptr => data%x(:,1,1,1)
211 dof => data%dof
212 write_pressure = .true.
213 write_velocity = .false.
214 type is (field_list_t)
215 select case (data%size())
216 case (1)
217 p%ptr => data%items(1)%ptr%x(:,1,1,1)
218 write_pressure = .true.
219 write_velocity = .false.
220 case (2)
221 p%ptr => data%items(1)%ptr%x(:,1,1,1)
222 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
223 write_pressure = .true.
224 write_temperature = .true.
225 case (3)
226 u%ptr => data%items(1)%ptr%x(:,1,1,1)
227 v%ptr => data%items(2)%ptr%x(:,1,1,1)
228 w%ptr => data%items(3)%ptr%x(:,1,1,1)
229 write_velocity = .true.
230 case (4)
231 p%ptr => data%items(1)%ptr%x(:,1,1,1)
232 u%ptr => data%items(2)%ptr%x(:,1,1,1)
233 v%ptr => data%items(3)%ptr%x(:,1,1,1)
234 w%ptr => data%items(4)%ptr%x(:,1,1,1)
235 write_pressure = .true.
236 write_velocity = .true.
237 case (5:99)
238 p%ptr => data%items(1)%ptr%x(:,1,1,1)
239 u%ptr => data%items(2)%ptr%x(:,1,1,1)
240 v%ptr => data%items(3)%ptr%x(:,1,1,1)
241 w%ptr => data%items(4)%ptr%x(:,1,1,1)
242 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
243 n_scalar_fields = data%size() - 5
244 allocate(scalar_fields(n_scalar_fields))
245 do i = 1, n_scalar_fields
246 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
247 end do
248 write_pressure = .true.
249 write_velocity = .true.
250 write_temperature = .true.
251 case default
252 call neko_error('This many fields not supported yet, fld_file')
253 end select
254 dof => data%dof(1)
255 class default
256 call neko_error('Invalid data')
257 end select
258 ! Fix things for pointers that do not exist in all data types...
259 if (associated(dof)) then
260 x%ptr => dof%x(:,1,1,1)
261 y%ptr => dof%y(:,1,1,1)
262 z%ptr => dof%z(:,1,1,1)
263 msh => dof%msh
264 xh => dof%Xh
265 end if
266
267 if (associated(msh)) then
268 nelv = msh%nelv
269 glb_nelv = msh%glb_nelv
270 offset_el = msh%offset_el
271 gdim = msh%gdim
272 ! Store global idx of each element
273 allocate(idx(msh%nelv))
274 do i = 1, msh%nelv
275 idx(i) = msh%elements(i)%e%id()
276 end do
277 end if
278
279 if (associated(xh)) then
280 lx = xh%lx
281 ly = xh%ly
282 lz = xh%lz
283 end if
284
285 lxyz = lx*ly*lz
286 n = nelv*lxyz
287
288 if (this%dp_precision) then
289 fld_data_size = mpi_double_precision_size
290 else
291 fld_data_size = mpi_real_size
292 end if
293 if (this%dp_precision) then
294 allocate(tmp_dp(gdim*n))
295 else
296 allocate(tmp_sp(gdim*n))
297 end if
298
299
300 !
301 ! Create fld header for NEKTON's multifile output
302 !
303
304 call this%increment_counter()
305 ! Check if I should write the mesh. Always override at the start counters
306 if (.not. this%write_mesh) then
307 write_mesh = (this%get_counter() .eq. this%get_start_counter())
308 else
309 write_mesh = this%write_mesh
310 end if
311 call mpi_allreduce(mpi_in_place, write_mesh, 1, &
312 mpi_logical, mpi_lor, neko_comm)
313 call mpi_allreduce(mpi_in_place, write_velocity, 1, &
314 mpi_logical, mpi_lor, neko_comm)
315 call mpi_allreduce(mpi_in_place, write_pressure, 1, &
316 mpi_logical, mpi_lor, neko_comm)
317 call mpi_allreduce(mpi_in_place, write_temperature, 1, &
318 mpi_logical, mpi_lor, neko_comm)
319 call mpi_allreduce(mpi_in_place, n_scalar_fields, 1, &
320 mpi_integer, mpi_max, neko_comm)
321
322 ! Build rdcode note that for field_t, we only support scalar
323 ! fields at the moment
324 rdcode = ' '
325 i = 1
326 if (write_mesh) then
327 rdcode(i) = 'X'
328 i = i + 1
329 end if
330 if (write_velocity) then
331 rdcode(i) = 'U'
332 i = i + 1
333 end if
334 if (write_pressure) then
335 rdcode(i) = 'P'
336 i = i + 1
337 end if
338 if (write_temperature) then
339 rdcode(i) = 'T'
340 i = i + 1
341 end if
342 if (n_scalar_fields .gt. 0 ) then
343 rdcode(i) = 'S'
344 i = i + 1
345 write(rdcode(i), '(i1)') (n_scalar_fields)/10
346 i = i + 1
347 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
348 i = i + 1
349 end if
350
352 write(hdr, 1) fld_data_size, lx, ly, lz, glb_nelv, glb_nelv,&
353 time, this%get_counter(), 1, 1, (rdcode(i), i = 1, 10)
3541 format('#std', 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
355 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
356
357 ! Change to NEKTON's fld file format
358 fname = this%get_fld_fname()
359
360 call mpi_file_open(neko_comm, trim(fname), &
361 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, &
362 ierr)
363
364 call mpi_file_write_all(fh, hdr, 132, mpi_character, status, ierr)
365 mpi_offset = 132 * mpi_character_size
366
367 call mpi_file_write_all(fh, test_pattern, 1, mpi_real, status, ierr)
368 mpi_offset = mpi_offset + mpi_real_size
369
370 byte_offset = mpi_offset + &
371 int(offset_el, i8) * int(mpi_integer_size, i8)
372 call mpi_file_write_at_all(fh, byte_offset, idx, nelv, &
373 mpi_integer, status, ierr)
374 mpi_offset = mpi_offset + int(glb_nelv, i8) * int(mpi_integer_size, i8)
375 deallocate(idx)
376 if (write_mesh) then
377
378 byte_offset = mpi_offset + int(offset_el, i8) * &
379 (int(gdim*lxyz, i8) * &
380 int(fld_data_size, i8))
381 call fld_file_write_vector_field(this, fh, byte_offset, &
382 x%ptr, y%ptr, z%ptr, &
383 n, gdim, lxyz, nelv)
384 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
385 (int(gdim *lxyz, i8) * &
386 int(fld_data_size, i8))
387 end if
388 if (write_velocity) then
389 byte_offset = mpi_offset + int(offset_el, i8) * &
390 (int(gdim * (lxyz), i8) * int(fld_data_size, i8))
391 call fld_file_write_vector_field(this, fh, byte_offset, &
392 u%ptr, v%ptr, w%ptr, n, gdim, lxyz, nelv)
393
394 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
395 (int(gdim * (lxyz), i8) * &
396 int(fld_data_size, i8))
397
398 end if
399
400 if (write_pressure) then
401 byte_offset = mpi_offset + int(offset_el, i8) * &
402 (int((lxyz), i8) * int(fld_data_size, i8))
403 call fld_file_write_field(this, fh, byte_offset, p%ptr, n)
404 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
405 (int((lxyz), i8) * int(fld_data_size, i8))
406 end if
407
408 if (write_temperature) then
409 byte_offset = mpi_offset + int(offset_el, i8) * &
410 (int((lxyz), i8) * &
411 int(fld_data_size, i8))
412 call fld_file_write_field(this, fh, byte_offset, tem%ptr, n)
413 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
414 (int((lxyz), i8) * &
415 int(fld_data_size, i8))
416 end if
417
418 temp_offset = mpi_offset
419
420 do i = 1, n_scalar_fields
421 ! Without this redundant if statement Cray optimizes this loop to
422 ! Oblivion
423 if (i .eq. 2) then
424 mpi_offset = int(temp_offset, i8) + int(1_i8*glb_nelv, i8) * &
425 (int(lxyz, i8) * int(fld_data_size, i8))
426 end if
427 byte_offset = int(mpi_offset, i8) + int(offset_el, i8) * &
428 (int((lxyz), i8) * &
429 int(fld_data_size, i8))
430 call fld_file_write_field(this, fh, byte_offset, scalar_fields(i)%ptr, n)
431 mpi_offset = int(mpi_offset, i8) + int(glb_nelv, i8) * &
432 (int(lxyz, i8) * &
433 int(fld_data_size, i8))
434 end do
435
436 if (gdim .eq. 3) then
437
439 if (write_mesh) then
440 !The offset is:
441 ! mpioff + element_off * 2(min max value)
442 ! * 4(single precision) * gdim(dimensions)
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 x%ptr, y%ptr, z%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 end if
456
457 if (write_velocity) then
458 byte_offset = int(mpi_offset, i8) + &
459 int(offset_el, i8) * &
460 int(2, i8) * &
461 int(mpi_real_size, i8) * &
462 int(gdim, i8)
463 call fld_file_write_metadata_vector(this, fh, byte_offset, &
464 u%ptr, v%ptr, w%ptr, gdim, 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 int(gdim, i8)
470
471 end if
472
473 if (write_pressure) 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 p%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 if (write_temperature) then
488 byte_offset = int(mpi_offset, i8) + &
489 int(offset_el, i8) * &
490 int(2, i8) * &
491 int(mpi_real_size, i8)
492 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
493 tem%ptr, lxyz, nelv)
494 mpi_offset = int(mpi_offset, i8) + &
495 int(glb_nelv, i8) * &
496 int(2, i8) * &
497 int(mpi_real_size, i8)
498
499 end if
500
501
502
503 temp_offset = mpi_offset
504
505 do i = 1, n_scalar_fields
506 ! Without this redundant if statement, Cray optimizes this loop to
507 ! Oblivion
508 if (i .eq. 2) then
509 mpi_offset = int(temp_offset, i8) + &
510 int(1_i8*glb_nelv, i8) * &
511 int(2, i8) * &
512 int(mpi_real_size, i8)
513 end if
514
515 byte_offset = int(mpi_offset, i8) + &
516 int(offset_el, i8) * &
517 int(2, i8) * &
518 int(mpi_real_size, i8)
519 call fld_file_write_metadata_scalar(this, fh, byte_offset, &
520 scalar_fields(i)%ptr, lxyz, nelv)
521 mpi_offset = int(mpi_offset, i8) + &
522 int(glb_nelv, i8) * &
523 int(2, i8) * &
524 int(mpi_real_size, i8)
525 end do
526 end if
527
528
529 call mpi_file_sync(fh, ierr)
530 call mpi_file_close(fh, ierr)
531 ! Write metadata file
532 if (pe_rank .eq. 0) then
533 call filename_name(this%get_base_fname(), name)
534
535 open(newunit = file_unit, &
536 file = this%get_meta_fname(), status = 'replace')
537 ! The following string will specify that the files in the file series
538 ! are defined by the filename followed by a 0.
539 ! This 0 is necessary as it specifies the index of number of files
540 ! the output file is split across.
541 ! In the past, many .f files were generated for each write.
542 ! To be consistent with this the trailing 0 is still necessary today.
543 write(file_unit, fmt = '(A,A,A)') 'filetemplate: ', &
544 trim(name), '%01d.f%05d'
545 write(file_unit, fmt = '(A,i5)') 'firsttimestep: ', &
546 this%get_start_counter()
547 write(file_unit, fmt = '(A,i5)') 'numtimesteps: ', &
548 (this%get_counter() + 1) - this%get_start_counter()
549 close(file_unit)
550 end if
551
552 if (allocated(tmp_dp)) deallocate(tmp_dp)
553 if (allocated(tmp_sp)) deallocate(tmp_sp)
554 if (allocated(tempo)) deallocate(tempo)
555 if (allocated(scalar_fields)) deallocate(scalar_fields)
556
557 end subroutine fld_file_write
558
561 subroutine fld_file_write_masked(this, data, mask, t)
562 class(fld_file_t), intent(inout) :: this
563 class(*), target, intent(in) :: data
564 type(mask_t), intent(in) :: mask
565 real(kind=rp), intent(in), optional :: t
566 type(array_ptr_t) :: x, y, z, u, v, w, p, tem
567 real(kind=rp), allocatable, target :: tempo(:)
568 type(mesh_t), pointer :: msh
569 type(dofmap_t), pointer :: dof
570 type(space_t), pointer :: Xh
571 real(kind=dp) :: time
572 character(len= 132) :: hdr
573 character :: rdcode(10)
574 character(len=6) :: id_str
575 character(len= 1024) :: fname
576 character(len= 1024) :: name
577 integer :: file_unit
578 integer :: i, ierr, n, suffix_pos, tslash_pos
579 integer :: lx, ly, lz, lxyz, gdim, glb_nelv, nelv, offset_el
580 integer, allocatable :: idx(:)
581 type(mpi_status) :: status
582 type(mpi_file) :: fh
583 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset, temp_offset
584 real(kind=sp), parameter :: test_pattern = 6.54321
585 type(array_ptr_t), allocatable :: scalar_fields(:)
586 logical :: write_mesh, write_velocity, write_pressure, write_temperature
587 integer :: fld_data_size, n_scalar_fields
588
589 if (present(t)) then
590 time = real(t, dp)
591 else
592 time = 0d0
593 end if
594
595 nullify(msh)
596 nullify(dof)
597 nullify(xh)
598 n_scalar_fields = 0
599 write_pressure = .false.
600 write_velocity = .false.
601 write_temperature = .false.
602
603 select type (data)
604 type is (fld_file_data_t)
605 nelv = data%nelv
606 lx = data%lx
607 ly = data%ly
608 lz = data%lz
609 gdim = data%gdim
610 glb_nelv = data%glb_nelv
611 offset_el = data%offset_el
612
613 if (data%x%size() .gt. 0) x%ptr => data%x%x
614 if (data%y%size() .gt. 0) y%ptr => data%y%x
615 if (data%z%size() .gt. 0) z%ptr => data%z%x
616 if (gdim .eq. 2) z%ptr => data%y%x
617 if (data%u%size() .gt. 0) then
618 u%ptr => data%u%x
619 ! In case only u is actually allocated, point the other comps to u
620 ! so that we don't die on trying to write them
621 if (data%v%size() .le. 0) v%ptr => data%u%x
622 if (data%w%size() .le. 0) w%ptr => data%u%x
623 write_velocity = .true.
624 end if
625 if (data%v%size() .gt. 0) v%ptr => data%v%x
626 if (data%w%size() .gt. 0) w%ptr => data%w%x
627 if (data%p%size() .gt. 0) then
628 p%ptr => data%p%x
629 write_pressure = .true.
630 end if
631 if (data%t%size() .gt. 0) then
632 write_temperature = .true.
633 tem%ptr => data%t%x
634 end if
635 ! If gdim = 2 and Z-velocity component exists,
636 ! it is stored in last scalar field
637 if (gdim .eq. 2 .and. data%w%size() .gt. 0) then
638 n_scalar_fields = data%n_scalars + 1
639 allocate(scalar_fields(n_scalar_fields))
640 do i = 1, n_scalar_fields -1
641 scalar_fields(i)%ptr => data%s(i)%x
642 end do
643 scalar_fields(n_scalar_fields)%ptr => data%w%x
644 else
645 n_scalar_fields = data%n_scalars
646 allocate(scalar_fields(n_scalar_fields+1))
647 do i = 1, n_scalar_fields
648 scalar_fields(i)%ptr => data%s(i)%x
649 end do
650 scalar_fields(n_scalar_fields+1)%ptr => data%w%x
651 end if
652 ! This is very stupid...
653 ! Some compilers cannot handle that these pointers dont point to anything
654 ! (although they are not used) this fixes a segfault due to this.
655 if (nelv .eq. 0) then
656 allocate(tempo(1))
657 x%ptr => tempo
658 y%ptr => tempo
659 z%ptr => tempo
660 u%ptr => tempo
661 v%ptr => tempo
662 w%ptr => tempo
663 p%ptr => tempo
664 tem%ptr => tempo
665 end if
666
667 allocate(idx(nelv))
668 do i = 1, nelv
669 idx(i) = data%idx(i)
670 end do
671 type is (field_t)
672 p%ptr => data%x(:,1,1,1)
673 dof => data%dof
674 write_pressure = .true.
675 write_velocity = .false.
676 type is (field_list_t)
677 select case (data%size())
678 case (1)
679 p%ptr => data%items(1)%ptr%x(:,1,1,1)
680 write_pressure = .true.
681 write_velocity = .false.
682 case (2)
683 p%ptr => data%items(1)%ptr%x(:,1,1,1)
684 tem%ptr => data%items(2)%ptr%x(:,1,1,1)
685 write_pressure = .true.
686 write_temperature = .true.
687 case (3)
688 u%ptr => data%items(1)%ptr%x(:,1,1,1)
689 v%ptr => data%items(2)%ptr%x(:,1,1,1)
690 w%ptr => data%items(3)%ptr%x(:,1,1,1)
691 write_velocity = .true.
692 case (4)
693 p%ptr => data%items(1)%ptr%x(:,1,1,1)
694 u%ptr => data%items(2)%ptr%x(:,1,1,1)
695 v%ptr => data%items(3)%ptr%x(:,1,1,1)
696 w%ptr => data%items(4)%ptr%x(:,1,1,1)
697 write_pressure = .true.
698 write_velocity = .true.
699 case (5:99)
700 p%ptr => data%items(1)%ptr%x(:,1,1,1)
701 u%ptr => data%items(2)%ptr%x(:,1,1,1)
702 v%ptr => data%items(3)%ptr%x(:,1,1,1)
703 w%ptr => data%items(4)%ptr%x(:,1,1,1)
704 tem%ptr => data%items(5)%ptr%x(:,1,1,1)
705 n_scalar_fields = data%size() - 5
706 allocate(scalar_fields(n_scalar_fields))
707 do i = 1, n_scalar_fields
708 scalar_fields(i)%ptr => data%items(i+5)%ptr%x(:,1,1,1)
709 end do
710 write_pressure = .true.
711 write_velocity = .true.
712 write_temperature = .true.
713 case default
714 call neko_error('This many fields not supported yet, fld_file')
715 end select
716 dof => data%dof(1)
717 class default
718 call neko_error('Invalid data')
719 end select
720 ! Fix things for pointers that do not exist in all data types...
721 if (associated(dof)) then
722 x%ptr => dof%x(:,1,1,1)
723 y%ptr => dof%y(:,1,1,1)
724 z%ptr => dof%z(:,1,1,1)
725 msh => dof%msh
726 xh => dof%Xh
727 end if
728
729 if (associated(msh)) then
730 nelv = msh%nelv
731 glb_nelv = msh%glb_nelv
732 offset_el = msh%offset_el
733 gdim = msh%gdim
734 ! Store global idx of each element
735 allocate(idx(msh%nelv))
736 do i = 1, msh%nelv
737 idx(i) = msh%elements(i)%e%id()
738 end do
739 end if
740
741 if (associated(xh)) then
742 lx = xh%lx
743 ly = xh%ly
744 lz = xh%lz
745 end if
746
747 ! Up to now, all data types have dealt with their stuff.
748 ! Now overwrite with the masked info
749 lxyz = lx*ly*lz
750 nelv = mask%size() / lxyz
751 if (mod(mask%size(), lxyz) /= 0) then
752 call neko_error("Mask size must be a multiple of the number of elements in the mesh.")
753 end if
754 call mpi_allreduce(nelv, glb_nelv, 1, &
755 mpi_integer, mpi_sum, neko_comm)
756 call mpi_scan(nelv, offset_el, 1, &
757 mpi_integer, mpi_sum, neko_comm, ierr)
758 offset_el = offset_el - nelv
759
760 if (allocated(idx)) then
761 deallocate(idx)
762 end if
763
764 allocate(idx(nelv))
765 do i = 1, nelv
766 idx(i) = offset_el + i
767 end do
768 n = nelv*lxyz
769
770 if (this%dp_precision) then
771 fld_data_size = mpi_double_precision_size
772 else
773 fld_data_size = mpi_real_size
774 end if
775 if (this%dp_precision) then
776 allocate(tmp_dp(gdim*n))
777 else
778 allocate(tmp_sp(gdim*n))
779 end if
780
781
782 !
783 ! Create fld header for NEKTON's multifile output
784 !
785
786 call this%increment_counter()
787 ! Check if I should write the mesh. Always override at the start counters
788 if (.not. this%write_mesh) then
789 write_mesh = (this%get_counter() .eq. this%get_start_counter())
790 else
791 write_mesh = this%write_mesh
792 end if
793 call mpi_allreduce(mpi_in_place, write_mesh, 1, &
794 mpi_logical, mpi_lor, neko_comm)
795 call mpi_allreduce(mpi_in_place, write_velocity, 1, &
796 mpi_logical, mpi_lor, neko_comm)
797 call mpi_allreduce(mpi_in_place, write_pressure, 1, &
798 mpi_logical, mpi_lor, neko_comm)
799 call mpi_allreduce(mpi_in_place, write_temperature, 1, &
800 mpi_logical, mpi_lor, neko_comm)
801 call mpi_allreduce(mpi_in_place, n_scalar_fields, 1, &
802 mpi_integer, mpi_max, neko_comm)
803
804 ! Build rdcode note that for field_t, we only support scalar
805 ! fields at the moment
806 rdcode = ' '
807 i = 1
808 if (write_mesh) then
809 rdcode(i) = 'X'
810 i = i + 1
811 end if
812 if (write_velocity) then
813 rdcode(i) = 'U'
814 i = i + 1
815 end if
816 if (write_pressure) then
817 rdcode(i) = 'P'
818 i = i + 1
819 end if
820 if (write_temperature) then
821 rdcode(i) = 'T'
822 i = i + 1
823 end if
824 if (n_scalar_fields .gt. 0 ) then
825 rdcode(i) = 'S'
826 i = i + 1
827 write(rdcode(i), '(i1)') (n_scalar_fields)/10
828 i = i + 1
829 write(rdcode(i), '(i1)') (n_scalar_fields) - 10*((n_scalar_fields)/10)
830 i = i + 1
831 end if
832
834 write(hdr, 1) fld_data_size, lx, ly, lz, glb_nelv, glb_nelv,&
835 time, this%get_counter(), 1, 1, (rdcode(i), i = 1, 10)
8361 format('#std', 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
837 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
838
839 ! Change to NEKTON's fld file format
840 fname = this%get_fld_fname()
841
842 call mpi_file_open(neko_comm, trim(fname), &
843 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, &
844 ierr)
845
846 call mpi_file_write_all(fh, hdr, 132, mpi_character, status, ierr)
847 mpi_offset = 132 * mpi_character_size
848
849 call mpi_file_write_all(fh, test_pattern, 1, mpi_real, status, ierr)
850 mpi_offset = mpi_offset + mpi_real_size
851
852 byte_offset = mpi_offset + &
853 int(offset_el, i8) * int(mpi_integer_size, i8)
854 call mpi_file_write_at_all(fh, byte_offset, idx, nelv, &
855 mpi_integer, status, ierr)
856 mpi_offset = mpi_offset + int(glb_nelv, i8) * int(mpi_integer_size, i8)
857 deallocate(idx)
858 if (write_mesh) then
859
860 byte_offset = mpi_offset + int(offset_el, i8) * &
861 (int(gdim*lxyz, i8) * &
862 int(fld_data_size, i8))
863 call fld_file_write_vector_field_masked(this, fh, byte_offset, &
864 x%ptr, y%ptr, z%ptr, &
865 n, gdim, lxyz, nelv, lx, ly, lz, mask%get())
866 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
867 (int(gdim *lxyz, i8) * &
868 int(fld_data_size, i8))
869 end if
870 if (write_velocity) then
871 byte_offset = mpi_offset + int(offset_el, i8) * &
872 (int(gdim * (lxyz), i8) * int(fld_data_size, i8))
873 call fld_file_write_vector_field_masked(this, fh, byte_offset, &
874 u%ptr, v%ptr, w%ptr, n, gdim, lxyz, nelv, lx, ly, lz, mask%get())
875
876 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
877 (int(gdim * (lxyz), i8) * &
878 int(fld_data_size, i8))
879
880 end if
881
882 if (write_pressure) then
883 byte_offset = mpi_offset + int(offset_el, i8) * &
884 (int((lxyz), i8) * int(fld_data_size, i8))
885 call fld_file_write_field_masked(this, fh, byte_offset, p%ptr, n, mask%get())
886 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
887 (int((lxyz), i8) * int(fld_data_size, i8))
888 end if
889
890 if (write_temperature) then
891 byte_offset = mpi_offset + int(offset_el, i8) * &
892 (int((lxyz), i8) * &
893 int(fld_data_size, i8))
894 call fld_file_write_field_masked(this, fh, byte_offset, tem%ptr, n, mask%get())
895 mpi_offset = mpi_offset + int(glb_nelv, i8) * &
896 (int((lxyz), i8) * &
897 int(fld_data_size, i8))
898 end if
899
900 temp_offset = mpi_offset
901
902 do i = 1, n_scalar_fields
903 ! Without this redundant if statement Cray optimizes this loop to
904 ! Oblivion
905 if (i .eq. 2) then
906 mpi_offset = int(temp_offset, i8) + int(1_i8*glb_nelv, i8) * &
907 (int(lxyz, i8) * int(fld_data_size, i8))
908 end if
909 byte_offset = int(mpi_offset, i8) + int(offset_el, i8) * &
910 (int((lxyz), i8) * &
911 int(fld_data_size, i8))
912 call fld_file_write_field_masked(this, fh, byte_offset, scalar_fields(i)%ptr, n, mask%get())
913 mpi_offset = int(mpi_offset, i8) + int(glb_nelv, i8) * &
914 (int(lxyz, i8) * &
915 int(fld_data_size, i8))
916 end do
917
918 if (gdim .eq. 3) then
919
921 if (write_mesh) then
922 !The offset is:
923 ! mpioff + element_off * 2(min max value)
924 ! * 4(single precision) * gdim(dimensions)
925 byte_offset = int(mpi_offset, i8) + &
926 int(offset_el, i8) * &
927 int(2, i8) * &
928 int(mpi_real_size, i8) * &
929 int(gdim, i8)
930 call fld_file_write_metadata_vector_masked(this, fh, byte_offset, &
931 x%ptr, y%ptr, z%ptr, gdim, lxyz, nelv, lx, ly, lz, n, mask%get())
932 mpi_offset = int(mpi_offset, i8) + &
933 int(glb_nelv, i8) * &
934 int(2, i8) * &
935 int(mpi_real_size, i8) * &
936 int(gdim, i8)
937 end if
938
939 if (write_velocity) then
940 byte_offset = int(mpi_offset, i8) + &
941 int(offset_el, i8) * &
942 int(2, i8) * &
943 int(mpi_real_size, i8) * &
944 int(gdim, i8)
945 call fld_file_write_metadata_vector_masked(this, fh, byte_offset, &
946 u%ptr, v%ptr, w%ptr, gdim, lxyz, nelv, lx, ly, lz, n, mask%get())
947 mpi_offset = int(mpi_offset, i8) + &
948 int(glb_nelv, i8) * &
949 int(2, i8) * &
950 int(mpi_real_size, i8) * &
951 int(gdim, i8)
952
953 end if
954
955 if (write_pressure) then
956 byte_offset = int(mpi_offset, i8) + &
957 int(offset_el, i8) * &
958 int(2, i8) * &
959 int(mpi_real_size, i8)
960 call fld_file_write_metadata_scalar_masked(this, fh, byte_offset, &
961 p%ptr, lxyz, nelv, lx, ly, lz, n, mask%get())
962 mpi_offset = int(mpi_offset, i8) + &
963 int(glb_nelv, i8) * &
964 int(2, i8) * &
965 int(mpi_real_size, i8)
966
967 end if
968
969 if (write_temperature) then
970 byte_offset = int(mpi_offset, i8) + &
971 int(offset_el, i8) * &
972 int(2, i8) * &
973 int(mpi_real_size, i8)
974 call fld_file_write_metadata_scalar_masked(this, fh, byte_offset, &
975 tem%ptr, lxyz, nelv, lx, ly, lz, n, mask%get())
976 mpi_offset = int(mpi_offset, i8) + &
977 int(glb_nelv, i8) * &
978 int(2, i8) * &
979 int(mpi_real_size, i8)
980
981 end if
982
983
984
985 temp_offset = mpi_offset
986
987 do i = 1, n_scalar_fields
988 ! Without this redundant if statement, Cray optimizes this loop to
989 ! Oblivion
990 if (i .eq. 2) then
991 mpi_offset = int(temp_offset, i8) + &
992 int(1_i8*glb_nelv, i8) * &
993 int(2, i8) * &
994 int(mpi_real_size, i8)
995 end if
996
997 byte_offset = int(mpi_offset, i8) + &
998 int(offset_el, i8) * &
999 int(2, i8) * &
1000 int(mpi_real_size, i8)
1001 call fld_file_write_metadata_scalar_masked(this, fh, byte_offset, &
1002 scalar_fields(i)%ptr, lxyz, nelv, lx, ly, lz, n, mask%get())
1003 mpi_offset = int(mpi_offset, i8) + &
1004 int(glb_nelv, i8) * &
1005 int(2, i8) * &
1006 int(mpi_real_size, i8)
1007 end do
1008 end if
1009
1010
1011 call mpi_file_sync(fh, ierr)
1012 call mpi_file_close(fh, ierr)
1013 ! Write metadata file
1014 if (pe_rank .eq. 0) then
1015 call filename_name(this%get_base_fname(), name)
1016
1017 open(newunit = file_unit, &
1018 file = this%get_meta_fname(), status = 'replace')
1019 ! The following string will specify that the files in the file series
1020 ! are defined by the filename followed by a 0.
1021 ! This 0 is necessary as it specifies the index of number of files
1022 ! the output file is split across.
1023 ! In the past, many .f files were generated for each write.
1024 ! To be consistent with this the trailing 0 is still necessary today.
1025 write(file_unit, fmt = '(A,A,A)') 'filetemplate: ', &
1026 trim(name), '%01d.f%05d'
1027 write(file_unit, fmt = '(A,i5)') 'firsttimestep: ', &
1028 this%get_start_counter()
1029 write(file_unit, fmt = '(A,i5)') 'numtimesteps: ', &
1030 (this%get_counter() + 1) - this%get_start_counter()
1031 close(file_unit)
1032 end if
1033
1034 if (allocated(tmp_dp)) deallocate(tmp_dp)
1035 if (allocated(tmp_sp)) deallocate(tmp_sp)
1036 if (allocated(tempo)) deallocate(tempo)
1037 if (allocated(scalar_fields)) deallocate(scalar_fields)
1038
1039 end subroutine fld_file_write_masked
1040
1041 subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, &
1042 gdim, lxyz, nelv)
1043 class(fld_file_t), intent(inout) :: this
1044 type(mpi_file), intent(inout) :: fh
1045 integer, intent(in) :: gdim, lxyz, nelv
1046 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
1047 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1048 integer :: el, j, ierr, nout
1049 type(mpi_status) :: status
1050 real(kind=sp) :: buffer(2*gdim*nelv)
1051
1052 j = 1
1053 do el = 1, nelv
1054 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
1055 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
1056 buffer(j+2) = real(vlmin(y(1, el), lxyz), sp)
1057 buffer(j+3) = real(vlmax(y(1, el), lxyz), sp)
1058 j = j + 4
1059 if (gdim .eq. 3) then
1060 buffer(j+0) = real(vlmin(z(1, el), lxyz), sp)
1061 buffer(j+1) = real(vlmax(z(1, el), lxyz), sp)
1062 j = j + 2
1063 end if
1064 end do
1065
1066 ! write out data
1067 nout = 2*gdim*nelv
1068
1069 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
1070 mpi_real, status, ierr)
1071
1072 end subroutine fld_file_write_metadata_vector
1073
1074 subroutine fld_file_write_metadata_vector_masked(this, fh, byte_offset, x, y, z, &
1075 gdim, lxyz, nelv, lx, ly, lz, n, mask)
1076 class(fld_file_t), intent(inout) :: this
1077 type(mpi_file), intent(inout) :: fh
1078 integer, intent(in) :: gdim, lxyz, nelv, lx, ly, lz, n
1079 real(kind=rp), intent(in) :: x(lxyz, *), y(lxyz, *), z(lxyz, *)
1080 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1081 integer, intent(in) :: mask(n)
1082 integer :: el, j, ierr, nout, i_m, nidx(4), e_m
1083 type(mpi_status) :: status
1084 real(kind=sp) :: buffer(2*gdim*nelv)
1085
1086 j = 1
1087 do el = 1, nelv
1088 i_m = 1 + lxyz * (el - 1) ! Offset in the mask array
1089 nidx = nonlinear_index(mask(i_m), lx, ly, lz)
1090 e_m = nidx(4) ! Actual element in the field array
1091 buffer(j+0) = real(vlmin(x(1, e_m), lxyz), sp)
1092 buffer(j+1) = real(vlmax(x(1, e_m), lxyz), sp)
1093 buffer(j+2) = real(vlmin(y(1, e_m), lxyz), sp)
1094 buffer(j+3) = real(vlmax(y(1, e_m), lxyz), sp)
1095 j = j + 4
1096 if (gdim .eq. 3) then
1097 buffer(j+0) = real(vlmin(z(1, e_m), lxyz), sp)
1098 buffer(j+1) = real(vlmax(z(1, e_m), lxyz), sp)
1099 j = j + 2
1100 end if
1101 end do
1102
1103 ! write out data
1104 nout = 2*gdim*nelv
1105
1106 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
1107 mpi_real, status, ierr)
1108
1110
1111 subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, &
1112 nelv)
1113 class(fld_file_t), intent(inout) :: this
1114 type(mpi_file), intent(inout) :: fh
1115 integer, intent(in) :: lxyz, nelv
1116 real(kind=rp), intent(in) :: x(lxyz, nelv)
1117 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1118 integer :: el, j, ierr, nout
1119 type(mpi_status) :: status
1120 real(kind=sp) :: buffer(2*nelv)
1121
1122 j = 1
1123 do el = 1, nelv
1124 buffer(j+0) = real(vlmin(x(1, el), lxyz), sp)
1125 buffer(j+1) = real(vlmax(x(1, el), lxyz), sp)
1126 j = j + 2
1127 end do
1128
1129 ! write out data
1130 nout = 2*nelv
1131
1132 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
1133 mpi_real, status, ierr)
1134
1135 end subroutine fld_file_write_metadata_scalar
1136
1137 subroutine fld_file_write_metadata_scalar_masked(this, fh, byte_offset, x, lxyz, &
1138 nelv, lx, ly, lz, n, mask)
1139 class(fld_file_t), intent(inout) :: this
1140 type(mpi_file), intent(inout) :: fh
1141 integer, intent(in) :: lxyz, nelv, lx, ly, lz, n
1142 real(kind=rp), intent(in) :: x(lxyz, *)
1143 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1144 integer, intent(in) :: mask(n)
1145 integer :: el, j, ierr, nout, i_m, nidx(4), e_m
1146 type(mpi_status) :: status
1147 real(kind=sp) :: buffer(2*nelv)
1148
1149 j = 1
1150 do el = 1, nelv
1151 i_m = 1 + lxyz * (el - 1) ! Offset in the mask array
1152 nidx = nonlinear_index(mask(i_m), lx, ly, lz)
1153 e_m = nidx(4) ! Actual element in the field array
1154 buffer(j+0) = real(vlmin(x(1, e_m), lxyz), sp)
1155 buffer(j+1) = real(vlmax(x(1, e_m), lxyz), sp)
1156 j = j + 2
1157 end do
1158
1159 ! write out data
1160 nout = 2*nelv
1161
1162 call mpi_file_write_at_all(fh, byte_offset, buffer, nout, &
1163 mpi_real, status, ierr)
1164
1166
1167 subroutine fld_file_write_field(this, fh, byte_offset, p, n)
1168 class(fld_file_t), intent(inout) :: this
1169 type(mpi_file), intent(inout) :: fh
1170 integer, intent(inout) :: n
1171 real(kind=rp), intent(inout) :: p(n)
1172 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1173 integer :: i, ierr
1174 type(mpi_status) :: status
1175
1176 if ( this%dp_precision) then
1177 do i = 1, n
1178 tmp_dp(i) = real(p(i), dp)
1179 end do
1180
1181 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, n, &
1182 mpi_double_precision, status, ierr)
1183 else
1184 do i = 1, n
1185 tmp_sp(i) = real(p(i), sp)
1186 end do
1187 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, n, &
1188 mpi_real, status, ierr)
1189 end if
1190
1191 end subroutine fld_file_write_field
1192
1193 subroutine fld_file_write_field_masked(this, fh, byte_offset, p, n, mask)
1194 class(fld_file_t), intent(inout) :: this
1195 type(mpi_file), intent(inout) :: fh
1196 integer, intent(inout) :: n
1197 real(kind=rp), intent(inout) :: p(:)
1198 integer, intent(in) :: mask(n)
1199 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1200 integer :: i, ierr
1201 type(mpi_status) :: status
1202
1203 if ( this%dp_precision) then
1204 do i = 1, n
1205 tmp_dp(i) = real(p(mask(i)), dp)
1206 end do
1207
1208 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, n, &
1209 mpi_double_precision, status, ierr)
1210 else
1211 do i = 1, n
1212 tmp_sp(i) = real(p(mask(i)), sp)
1213 end do
1214 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, n, &
1215 mpi_real, status, ierr)
1216 end if
1217
1218 end subroutine fld_file_write_field_masked
1219
1220 subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, &
1221 gdim, lxyz, nelv)
1222 class(fld_file_t), intent(inout) :: this
1223 type(mpi_file), intent(inout) :: fh
1224 integer, intent(in) :: n, gdim, lxyz, nelv
1225 real(kind=rp), intent(in) :: x(lxyz, nelv), y(lxyz, nelv), z(lxyz, nelv)
1226 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1227 integer :: i, el, j, ierr
1228 type(mpi_status) :: status
1229
1230 if (this%dp_precision) then
1231 i = 1
1232 do el = 1, nelv
1233 do j = 1, lxyz
1234 tmp_dp(i) = real(x(j, el), dp)
1235 i = i +1
1236 end do
1237 do j = 1, lxyz
1238 tmp_dp(i) = real(y(j, el), dp)
1239 i = i +1
1240 end do
1241 if (gdim .eq. 3) then
1242 do j = 1, lxyz
1243 tmp_dp(i) = real(z(j, el), dp)
1244 i = i +1
1245 end do
1246 end if
1247 end do
1248 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, gdim*n, &
1249 mpi_double_precision, status, ierr)
1250 else
1251 i = 1
1252 do el = 1, nelv
1253 do j = 1, lxyz
1254 tmp_sp(i) = real(x(j, el), sp)
1255 i = i +1
1256 end do
1257 do j = 1, lxyz
1258 tmp_sp(i) = real(y(j, el), sp)
1259 i = i +1
1260 end do
1261 if (gdim .eq. 3) then
1262 do j = 1, lxyz
1263 tmp_sp(i) = real(z(j, el), sp)
1264 i = i +1
1265 end do
1266 end if
1267 end do
1268 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, gdim*n, &
1269 mpi_real, status, ierr)
1270 end if
1271
1272
1273 end subroutine fld_file_write_vector_field
1274
1275 subroutine fld_file_write_vector_field_masked(this, fh, byte_offset, x, y, z, n, &
1276 gdim, lxyz, nelv, lx, ly, lz, mask)
1277 class(fld_file_t), intent(inout) :: this
1278 type(mpi_file), intent(inout) :: fh
1279 integer, intent(in) :: n, gdim, lxyz, lx, ly, lz, nelv
1280 real(kind=rp), intent(in) :: x(lxyz, *), y(lxyz, *), z(lxyz, *)
1281 integer (kind=MPI_OFFSET_KIND), intent(in) :: byte_offset
1282 integer, intent(in) :: mask(n)
1283 integer :: i, el, j, ierr, i_m, e_m, nidx(4)
1284 type(mpi_status) :: status
1285
1286 if (this%dp_precision) then
1287 i = 1
1288 do el = 1, nelv
1289 i_m = 1 + lxyz * (el - 1) ! Offset in the mask array
1290 nidx = nonlinear_index(mask(i_m), lx, ly, lz)
1291 e_m = nidx(4) ! Actual element in the field array
1292 do j = 1, lxyz
1293 tmp_dp(i) = real(x(j, e_m), dp)
1294 i = i +1
1295 end do
1296 do j = 1, lxyz
1297 tmp_dp(i) = real(y(j, e_m), dp)
1298 i = i +1
1299 end do
1300 if (gdim .eq. 3) then
1301 do j = 1, lxyz
1302 tmp_dp(i) = real(z(j, e_m), dp)
1303 i = i +1
1304 end do
1305 end if
1306 end do
1307 call mpi_file_write_at_all(fh, byte_offset, tmp_dp, gdim*n, &
1308 mpi_double_precision, status, ierr)
1309 else
1310 i = 1
1311 do el = 1, nelv
1312 i_m = 1 + lxyz * (el - 1) ! Offset in the mask array
1313 nidx = nonlinear_index(mask(i_m), lx, ly, lz)
1314 e_m = nidx(4) ! Actual element in the field array
1315 do j = 1, lxyz
1316 tmp_sp(i) = real(x(j, e_m), sp)
1317 i = i +1
1318 end do
1319 do j = 1, lxyz
1320 tmp_sp(i) = real(y(j, e_m), sp)
1321 i = i +1
1322 end do
1323 if (gdim .eq. 3) then
1324 do j = 1, lxyz
1325 tmp_sp(i) = real(z(j, e_m), sp)
1326 i = i +1
1327 end do
1328 end if
1329 end do
1330 call mpi_file_write_at_all(fh, byte_offset, tmp_sp, gdim*n, &
1331 mpi_real, status, ierr)
1332 end if
1333
1334
1336
1338 subroutine fld_file_read(this, data)
1339 class(fld_file_t) :: this
1340 class(*), target, intent(inout) :: data
1341 character(len= 132) :: hdr
1342 integer :: ierr, suffix_pos, i, j
1343 type(mpi_file) :: fh
1344 type(mpi_status) :: status
1345 character(len= 1024) :: fname, base_fname, meta_fname, string, path
1346 logical :: meta_file, read_mesh, read_velocity, read_pressure
1347 logical :: read_temp
1348 character(len=6) :: suffix
1349 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
1350 integer :: lx, ly, lz, glb_nelv, counter, lxyz
1351 integer :: FLD_DATA_SIZE, n_scalars, n
1352 integer :: file_unit
1353 real(kind=rp) :: time
1354 real(kind=sp) :: temp
1355 type(linear_dist_t) :: dist
1356 real(kind=sp), parameter :: test_pattern = 6.54321
1357 character :: rdcode(10), temp_str(4)
1358 character(len=LOG_SIZE) :: log_buf
1359
1360 select type (data)
1361 type is (fld_file_data_t)
1362 call filename_chsuffix(this%get_base_fname(), meta_fname, 'nek5000')
1363
1364 inquire(file = trim(meta_fname), exist = meta_file)
1365 if (meta_file .and. data%meta_nsamples .eq. 0) then
1366 if (pe_rank .eq. 0) then
1367 open(newunit = file_unit, file = trim(meta_fname))
1368 read(file_unit, fmt = '(A)') string
1369 read(string(14:), fmt = '(A)') string
1370 string = trim(string)
1371
1372 data%fld_series_fname = string(:scan(trim(string), '%')-1)
1373 data%fld_series_fname = adjustl(data%fld_series_fname)
1374 data%fld_series_fname = trim(data%fld_series_fname)//'0'
1375
1376 read(file_unit, fmt = '(A)') string
1377 read(string(scan(string, ':')+1:), *) data%meta_start_counter
1378 read(file_unit, fmt = '(A)') string
1379 read(string(scan(string, ':')+1:), *) data%meta_nsamples
1380 close(file_unit)
1381
1382 write(log_buf,*) 'Reading meta file for fld series'
1383 call neko_log%message(log_buf)
1384 write(log_buf,*) 'Name: ', trim(data%fld_series_fname)
1385 call neko_log%message(log_buf)
1386 write(log_buf,*) 'Start counter: ', data%meta_start_counter
1387 call neko_log%message(log_buf)
1388 write(log_buf,*) 'Nsamples: ', data%meta_nsamples
1389 call neko_log%message(log_buf)
1390
1391 end if
1392 call mpi_bcast(data%fld_series_fname, 1024, mpi_character, 0, &
1393 neko_comm, ierr)
1394 call mpi_bcast(data%meta_start_counter, 1, mpi_integer, 0, &
1395 neko_comm, ierr)
1396 call mpi_bcast(data%meta_nsamples, 1, mpi_integer, 0, &
1397 neko_comm, ierr)
1398
1399 if (this%get_counter() .eq. -1) then
1400 call this%set_start_counter(data%meta_start_counter)
1401 call this%set_counter(data%meta_start_counter)
1402 end if
1403 end if
1404
1405 if (meta_file) then
1406 call filename_path(this%get_base_fname(), path)
1407 write(suffix, '(a,i5.5)') 'f', this%get_counter()
1408 fname = trim(path) // trim(data%fld_series_fname) // '.' // suffix
1409 if (this%get_counter() .ge. &
1410 data%meta_nsamples+data%meta_start_counter) then
1411 call neko_error('Trying to read more fld files than exist')
1412 end if
1413 else
1414 write(suffix, '(a,i5.5)') 'f', this%get_counter()
1415 call filename_chsuffix(trim(this%get_base_fname()), fname, suffix)
1416 end if
1417 call mpi_file_open(neko_comm, trim(fname), &
1418 mpi_mode_rdonly, mpi_info_null, fh, ierr)
1419
1420 if (ierr .ne. 0) call neko_error("Could not read "//trim(fname))
1421
1422 call neko_log%message('Reading fld file ' // trim(fname))
1423
1424 call mpi_file_read_all(fh, hdr, 132, mpi_character, status, ierr)
1425 ! This read can prorbably be done wihtout the temp variables,
1426 ! temp_str, i, j
1427
1428 read(hdr, 1) temp_str, fld_data_size, lx, ly, lz, glb_nelv, glb_nelv, &
1429 time, counter, i, j, (rdcode(i), i = 1, 10)
14301 format(4a, 1x, i1, 1x, i2, 1x, i2, 1x, i2, 1x, i10, 1x, i10, &
1431 1x, e20.13, 1x, i9, 1x, i6, 1x, i6, 1x, 10a)
1432 if (data%nelv .eq. 0) then
1433 dist = linear_dist_t(glb_nelv, pe_rank, pe_size, neko_comm)
1434 data%nelv = dist%num_local()
1435 data%offset_el = dist%start_idx()
1436 end if
1437 data%lx = lx
1438 data%ly = ly
1439 data%lz = lz
1440 data%glb_nelv = glb_nelv
1441 data%t_counter = counter
1442 data%time = time
1443 lxyz = lx * ly * lz
1444 n = lxyz * data%nelv
1445
1446 if (lz .eq. 1) then
1447 data%gdim = 2
1448 else
1449 data%gdim = 3
1450 end if
1451
1452
1453 if (fld_data_size .eq. mpi_double_precision_size) then
1454 this%dp_precision = .true.
1455 else
1456 this%dp_precision = .false.
1457 end if
1458 if (this%dp_precision) then
1459 allocate(tmp_dp(data%gdim*n))
1460 else
1461 allocate(tmp_sp(data%gdim*n))
1462 end if
1463
1464
1465 i = 1
1466 read_mesh = .false.
1467 read_velocity = .false.
1468 read_pressure = .false.
1469 read_temp = .false.
1470 if (rdcode(i) .eq. 'X') then
1471 read_mesh = .true.
1472 call data%x%init(n)
1473 call data%y%init(n)
1474 call data%z%init(n)
1475 i = i + 1
1476 end if
1477 if (rdcode(i) .eq. 'U') then
1478 read_velocity = .true.
1479 call data%u%init(n)
1480 call data%v%init(n)
1481 call data%w%init(n)
1482 i = i + 1
1483 end if
1484 if (rdcode(i) .eq. 'P') then
1485 read_pressure = .true.
1486 call data%p%init(n)
1487 i = i + 1
1488 end if
1489 if (rdcode(i) .eq. 'T') then
1490 read_temp = .true.
1491 call data%t%init(n)
1492 i = i + 1
1493 end if
1494 n_scalars = 0
1495 if (rdcode(i) .eq. 'S') then
1496 i = i + 1
1497 read(rdcode(i),*) n_scalars
1498 n_scalars = n_scalars*10
1499 i = i + 1
1500 read(rdcode(i),*) j
1501 n_scalars = n_scalars+j
1502 i = i + 1
1503 if (allocated(data%s)) then
1504 if (data%n_scalars .ne. n_scalars) then
1505 do j = 1, data%n_scalars
1506 call data%s(j)%free()
1507 end do
1508 deallocate(data%s)
1509 data%n_scalars = n_scalars
1510 allocate(data%s(n_scalars))
1511 do j = 1, data%n_scalars
1512 call data%s(j)%init(n)
1513 end do
1514 end if
1515 else
1516 data%n_scalars = n_scalars
1517 allocate(data%s(data%n_scalars))
1518 do j = 1, data%n_scalars
1519 call data%s(j)%init(n)
1520 end do
1521 end if
1522 i = i + 1
1523 end if
1524
1525 mpi_offset = 132 * mpi_character_size
1526 call mpi_file_read_at_all(fh, mpi_offset, temp, 1, &
1527 mpi_real, status, ierr)
1528 if (.not. sabscmp(temp, test_pattern, epsilon(1.0_sp))) then
1529 call neko_error('Incorrect format for fld file, &
1530 &test pattern does not match.')
1531 end if
1532 mpi_offset = mpi_offset + mpi_real_size
1533
1534
1535 if (allocated(data%idx)) then
1536 if (size(data%idx) .ne. data%nelv) then
1537 deallocate(data%idx)
1538 allocate(data%idx(data%nelv))
1539 end if
1540 else
1541 allocate(data%idx(data%nelv))
1542 end if
1543
1544 byte_offset = mpi_offset + &
1545 int(data%offset_el, i8) * int(mpi_integer_size, i8)
1546
1547 call mpi_file_read_at_all(fh, byte_offset, data%idx, data%nelv, &
1548 mpi_integer, status, ierr)
1549
1550 mpi_offset = mpi_offset + &
1551 int(data%glb_nelv, i8) * int(mpi_integer_size, i8)
1552
1553 if (read_mesh) then
1554 byte_offset = mpi_offset + int(data%offset_el, i8) * &
1555 (int(data%gdim*lxyz, i8) * &
1556 int(fld_data_size, i8))
1557 call fld_file_read_vector_field(this, fh, byte_offset, &
1558 data%x, data%y, data%z, data)
1559 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
1560 (int(data%gdim *lxyz, i8) * &
1561 int(fld_data_size, i8))
1562 end if
1563
1564 if (read_velocity) then
1565 byte_offset = mpi_offset + int(data%offset_el, i8) * &
1566 (int(data%gdim*lxyz, i8) * &
1567 int(fld_data_size, i8))
1568 call fld_file_read_vector_field(this, fh, byte_offset, &
1569 data%u, data%v, data%w, data)
1570 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
1571 (int(data%gdim *lxyz, i8) * &
1572 int(fld_data_size, i8))
1573 end if
1574
1575 if (read_pressure) then
1576 byte_offset = mpi_offset + int(data%offset_el, i8) * &
1577 (int(lxyz, i8) * &
1578 int(fld_data_size, i8))
1579 call fld_file_read_field(this, fh, byte_offset, data%p, data)
1580 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
1581 (int(lxyz, i8) * &
1582 int(fld_data_size, i8))
1583 end if
1584
1585 if (read_temp) then
1586 byte_offset = mpi_offset + int(data%offset_el, i8) * &
1587 (int(lxyz, i8) * &
1588 int(fld_data_size, i8))
1589 call fld_file_read_field(this, fh, byte_offset, data%t, data)
1590 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
1591 (int(lxyz, i8) * &
1592 int(fld_data_size, i8))
1593 end if
1594
1595 do i = 1, n_scalars
1596 byte_offset = mpi_offset + int(data%offset_el, i8) * &
1597 (int(lxyz, i8) * &
1598 int(fld_data_size, i8))
1599 call fld_file_read_field(this, fh, byte_offset, data%s(i), data)
1600 mpi_offset = mpi_offset + int(data%glb_nelv, i8) * &
1601 (int(lxyz, i8) * &
1602 int(fld_data_size, i8))
1603 end do
1604
1605 call device_sync()
1606 call this%increment_counter()
1607
1608 if (allocated(tmp_dp)) deallocate(tmp_dp)
1609 if (allocated(tmp_sp)) deallocate(tmp_sp)
1610 class default
1611 call neko_error('Currently we only read into fld_file_data_t, &
1612 &please use that data structure instead.')
1613 end select
1614
1615 end subroutine fld_file_read
1616
1617 subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
1618 class(fld_file_t), intent(inout) :: this
1619 type(vector_t), intent(inout) :: x
1620 type(fld_file_data_t) :: fld_data
1621 integer(kind=MPI_OFFSET_KIND) :: byte_offset
1622 type(mpi_file) :: fh
1623 type(mpi_status) :: status
1624 integer :: n, ierr, lxyz, i
1625
1626 n = x%size()
1627 lxyz = fld_data%lx * fld_data%ly * fld_data%lz
1628
1629 if (this%dp_precision) then
1630 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, n, &
1631 mpi_double_precision, status, ierr)
1632 else
1633 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, n, &
1634 mpi_real, status, ierr)
1635 end if
1636
1637 if (this%dp_precision) then
1638 do i = 1, n
1639 x%x(i) = tmp_dp(i)
1640 end do
1641 else
1642 do i = 1, n
1643 x%x(i) = tmp_sp(i)
1644 end do
1645 end if
1646
1647 call x%copy_from(host_to_device, sync = .false.)
1648
1649 end subroutine fld_file_read_field
1650
1651
1652 subroutine fld_file_read_vector_field(this, fh, byte_offset, &
1653 x, y, z, fld_data)
1654 class(fld_file_t), intent(inout) :: this
1655 type(vector_t), intent(inout) :: x, y, z
1656 type(fld_file_data_t) :: fld_data
1657 integer(kind=MPI_OFFSET_KIND) :: byte_offset
1658 type(mpi_file) :: fh
1659 type(mpi_status) :: status
1660 integer :: n, ierr, lxyz, i, j, e, nd
1661
1662 n = x%size()
1663 nd = n*fld_data%gdim
1664 lxyz = fld_data%lx*fld_data%ly*fld_data%lz
1665
1666 if (this%dp_precision) then
1667 call mpi_file_read_at_all(fh, byte_offset, tmp_dp, nd, &
1668 mpi_double_precision, status, ierr)
1669 else
1670 call mpi_file_read_at_all(fh, byte_offset, tmp_sp, nd, &
1671 mpi_real, status, ierr)
1672 end if
1673
1674
1675 if (this%dp_precision) then
1676 i = 1
1677 do e = 1, fld_data%nelv
1678 do j = 1, lxyz
1679 x%x((e-1)*lxyz+j) = tmp_dp(i)
1680 i = i + 1
1681 end do
1682 do j = 1, lxyz
1683 y%x((e-1)*lxyz+j) = tmp_dp(i)
1684 i = i + 1
1685 end do
1686 if (fld_data%gdim .eq. 3) then
1687 do j = 1, lxyz
1688 z%x((e-1)*lxyz+j) = tmp_dp(i)
1689 i = i + 1
1690 end do
1691 end if
1692 end do
1693 else
1694 i = 1
1695 do e = 1, fld_data%nelv
1696 do j = 1, lxyz
1697 x%x((e-1)*lxyz+j) = tmp_sp(i)
1698 i = i +1
1699 end do
1700 do j = 1, lxyz
1701 y%x((e-1)*lxyz+j) = tmp_sp(i)
1702 i = i +1
1703 end do
1704 if (fld_data%gdim .eq. 3) then
1705 do j = 1, lxyz
1706 z%x((e-1)*lxyz+j) = tmp_sp(i)
1707 i = i +1
1708 end do
1709 end if
1710 end do
1711 end if
1712
1713 call x%copy_from(host_to_device, sync = .false.)
1714 call y%copy_from(host_to_device, sync = .false.)
1715 if (fld_data%gdim .eq. 3) then
1716 call z%copy_from(host_to_device, sync = .false.)
1717 end if
1718
1719 end subroutine fld_file_read_vector_field
1720
1721 subroutine fld_file_set_precision(this, precision)
1722 class(fld_file_t) :: this
1723 integer, intent(in) :: precision
1724
1725 if (precision .eq. dp) then
1726 this%dp_precision = .true.
1727 else if (precision .eq. sp) then
1728 this%dp_precision = .false.
1729 else
1730 call neko_error('Invalid precision')
1731 end if
1732
1733 end subroutine fld_file_set_precision
1734
1735 subroutine fld_file_set_mask(this, mask)
1736 class(fld_file_t) :: this
1737 type(mask_t), intent(inout), optional :: mask
1738
1739 if (present(mask)) then
1740 call this%mask%init_from_mask(mask)
1741 else
1742 call this%mask%free()
1743 end if
1744
1745 end subroutine fld_file_set_mask
1746
1747 function fld_file_get_fld_fname(this) result(fname)
1748 class(fld_file_t), intent(in) :: this
1749 character(len=1024) :: fname
1750 character(len=1024) :: path, name, id_str
1751 integer :: suffix_pos
1752
1753 call filename_path(this%get_base_fname(), path)
1754 call filename_name(this%get_base_fname(), name)
1755
1756 write(fname, '(a,a,a,i5.5)') trim(path), trim(name), &
1757 '0.f', this%get_counter()
1758
1759 end function fld_file_get_fld_fname
1760
1761 function fld_file_get_meta_fname(this) result(fname)
1762 class(fld_file_t), intent(in) :: this
1763 character(len=1024) :: fname
1764 character(len=1024) :: path, name, id_str
1765
1766 call filename_path(this%get_base_fname(), path)
1767 call filename_name(this%get_base_fname(), name)
1768
1769 write(id_str, '(i5,a)') this%get_start_counter(), '.nek5000'
1770 write(fname, '(a,a,a)') trim(path), trim(name), trim(adjustl(id_str))
1771
1772 end function fld_file_get_meta_fname
1773
1774end module fld_file
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
double real
Synchronize a device or stream.
Definition device.F90:113
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
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
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
subroutine fld_file_write_metadata_vector_masked(this, fh, byte_offset, x, y, z, gdim, lxyz, nelv, lx, ly, lz, n, mask)
subroutine fld_file_write_field_masked(this, fh, byte_offset, p, n, mask)
character(len=1024) function fld_file_get_fld_fname(this)
real(kind=dp), dimension(:), allocatable, private tmp_dp
Definition fld_file.f90:62
subroutine fld_file_write_vector_field_masked(this, fh, byte_offset, x, y, z, n, gdim, lxyz, nelv, lx, ly, lz, mask)
subroutine fld_file_set_precision(this, precision)
subroutine fld_file_read_vector_field(this, fh, byte_offset, x, y, z, fld_data)
subroutine fld_file_write_metadata_scalar_masked(this, fh, byte_offset, x, lxyz, nelv, lx, ly, lz, n, mask)
subroutine fld_file_write_manager(this, data, t)
Manage writer to use.
Definition fld_file.f90:86
subroutine fld_file_read(this, data)
Load a field from a NEKTON fld file.
subroutine fld_file_write_vector_field(this, fh, byte_offset, x, y, z, n, gdim, lxyz, nelv)
subroutine fld_file_set_mask(this, mask)
character(len=1024) function fld_file_get_meta_fname(this)
subroutine fld_file_write_masked(this, data, mask, t)
Write fields to a NEKTON fld file from a masked array.
Definition fld_file.f90:562
subroutine fld_file_write_metadata_vector(this, fh, byte_offset, x, y, z, gdim, lxyz, nelv)
real(kind=sp), dimension(:), allocatable, private tmp_sp
Definition fld_file.f90:63
subroutine fld_file_read_field(this, fh, byte_offset, x, fld_data)
subroutine fld_file_write_field(this, fh, byte_offset, p, n)
subroutine fld_file_write_metadata_scalar(this, fh, byte_offset, x, lxyz, nelv)
subroutine fld_file_write(this, data, t)
Write fields to a NEKTON fld file.
Definition fld_file.f90:101
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
pure logical function, public sabscmp(x, y, tol)
Return single precision absolute comparison .
Definition math.f90:120
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:652
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:641
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:97
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:150
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition utils.f90:68
subroutine, public filename_path(fname, path)
Extract the path to a file.
Definition utils.f90:82
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:66
A generic file handler.
Type for consistently handling masks in Neko. This type encapsulates the mask array and its associate...
Definition mask.f90:51
The function space for the SEM solution fields.
Definition space.f90:63
Pointer to array.
Definition structs.f90:14