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