Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
rea_file.f90
Go to the documentation of this file.
1! Copyright (c) 2019-2022, 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!
37 use num_types, only : rp, dp
38 use mesh, only : mesh_t, neko_msh_max_zlbls
39 use point, only : point_t
40 use map, only : map_t
41 use rea, only : rea_t, rea_free
42 use re2_file, only : re2_file_t
43 use map_file, only : map_file_t
44 use comm
45 use datadist, only : linear_dist_t
46 use htable, only : htable_pt_t
49 implicit none
50 private
51
52 ! Defines the conventions for conversion of rea labels to labeled zones.
53 integer, public :: neko_w_bc_label = -1
54 integer, public :: neko_v_bc_label = -1
55 integer, public :: neko_o_bc_label = -1
56 integer, public :: neko_sym_bc_label = -1
57 integer, public :: neko_on_bc_label = -1
58 integer, public :: neko_shl_bc_label = -1
59
61 type, public, extends(generic_file_t) :: rea_file_t
62 contains
63 procedure :: read => rea_file_read
64 procedure :: write => rea_file_write
65 end type rea_file_t
66
67contains
68
70 subroutine rea_file_read(this, data)
71 class(rea_file_t) :: this
72 class(*), target, intent(inout) :: data
73 type(mesh_t), pointer :: msh
74 real(kind=dp), pointer :: params(:)
75 character(len=3), pointer :: cbc(:, :)
76 integer, allocatable :: curve_type(:, :)
77 logical, allocatable :: curve_element(:)
78 character(len=1) :: chtemp
79 integer :: ndim, nparam, nskip, nlogic, ncurve
80 integer :: nelgs, nelgv, i, j, ierr, l
81 integer :: el_idx, pt_idx
82 logical :: read_param, read_bcs, read_map
83 real(kind=dp) :: xc(8), yc(8), zc(8), curve(5)
84 real(kind=dp), allocatable :: bc_data(:, :, :)
85 real(kind=dp), allocatable :: curve_data(:, :, :)
86 type(point_t) :: p(8)
87 type(re2_file_t) :: re2_file
88 type(map_file_t) :: map_file
89 character(len=1024) :: re2_fname, map_fname, fname
90 integer :: start_el, end_el, nel, edge
91 type(linear_dist_t) :: dist
92 type(map_t) :: nm
93 type(htable_pt_t) :: htp
94 integer :: sym_facet, pids(4), p_el_idx, p_facet
95 integer :: off
96 integer, parameter, dimension(6) :: facet_map = [3, 2, 4, 1, 5, 6]
97 logical :: curve_skip = .false.
98 character(len=LOG_SIZE) :: log_buf
99 integer :: file_unit
100 ! ---- Offsets for conversion from internal zones to labeled zones
101 integer :: user_labeled_zones(NEKO_MSH_MAX_ZLBLS)
102 integer :: labeled_zone_offsets(NEKO_MSH_MAX_ZLBLS)
103 integer :: total_labeled_zone_offset, current_internal_zone
104 total_labeled_zone_offset = 0
105 labeled_zone_offsets = 0
106 user_labeled_zones = 0
107 current_internal_zone = 1
108 ! ----
109
110 select type (data)
111 type is (rea_t)
112 call rea_free(data)
113 msh => data%msh
114 params => data%params
115 cbc => data%cbc
116 read_param = .true.
117 read_bcs = .true.
118 type is (mesh_t)
119 msh => data
120 read_param = .false.
121 read_bcs = .false.
122 class default
123 call neko_error('Invalid output data')
124 end select
125
126 if (read_param .and. read_bcs .and. pe_size .gt. 1) then
127 call neko_error('Reading NEKTON session data only implemented in serial')
128 end if
129
130 fname = this%get_fname()
131 open(newunit = file_unit, file = trim(fname), status = 'old', &
132 iostat = ierr)
133 call neko_log%message('Reading NEKTON file ' // fname)
134
135 read(file_unit, *)
136 read(file_unit, *)
137 read(file_unit, *) ndim
138 read(file_unit, *) nparam
139
140 if (.not. read_param) then
141 ! Skip parameters
142 do i = 1, nparam
143 read(file_unit, *)
144 end do
145 else
146 allocate(params(nparam))
147 do i = 1, nparam
148 read(file_unit, *) params(i)
149 end do
150 end if
151
152 ! Skip passive scalars
153 read(file_unit, *) nskip
154 do i = 1, nskip
155 read(file_unit, *)
156 end do
157
158 ! Skip logic switches
159 read(file_unit, *) nlogic
160 do i = 1, nlogic
161 read(file_unit, *)
162 end do
163
164 ! Read mesh info
165 read(file_unit, *)
166 read(file_unit, *)
167 read(file_unit, *) nelgs, ndim, nelgv
168 if (nelgs .lt. 0) then
169 re2_fname = trim(fname(1:scan(trim(fname), '.', back = .true.))) // 're2'
170 call re2_file%init(re2_fname)
171 call re2_file%read(msh)
172 else
173 write(log_buf, 1) ndim, nelgv
1741 format('gdim = ', i1, ', nelements = ', i7)
175 call neko_log%message(log_buf)
176
177 call filename_chsuffix(fname, map_fname, 'map')
178 inquire(file = map_fname, exist = read_map)
179 if (read_map) then
180 call nm%init(nelgv, 2**ndim)
181 call map_file%init(map_fname)
182 call map_file%read(nm)
183 else
184 call neko_log%warning('No NEKTON map file found')
185 end if
186
187 ! Use a load-balanced linear distribution
188 dist = linear_dist_t(nelgv, pe_rank, pe_size, neko_comm)
189 nel = dist%num_local()
190 start_el = dist%start_idx() + 1
191 end_el = dist%end_idx() + 1
192
193 call msh%init(ndim, dist)
194
195 call htp%init(2**ndim * nel, ndim)
196
197 el_idx = 1
198 pt_idx = 0
199 do i = 1, nelgv
200 read(file_unit, *)
201 if (ndim .eq. 2) then
202 read(file_unit, *) (xc(j), j = 1, 4)
203 read(file_unit, *) (yc(j), j = 1, 4)
204 if (i .ge. start_el .and. i .le. end_el) then
205 do j = 1, 4
206 call p(j)%init(real(xc(j), dp), real(yc(j), dp), &
207 real(0d0, dp))
208 call rea_file_add_point(htp, p(j), pt_idx)
209 end do
210 ! swap vertices to keep symmetric vertex numbering in neko
211 call msh%add_element(el_idx, el_idx, p(1), p(2), p(4), p(3))
212 end if
213 else if (ndim .eq. 3) then
214 read(file_unit, *) (xc(j), j = 1, 4)
215 read(file_unit, *) (yc(j), j = 1, 4)
216 read(file_unit, *) (zc(j), j = 1, 4)
217 read(file_unit, *) (xc(j), j = 5, 8)
218 read(file_unit, *) (yc(j), j = 5, 8)
219 read(file_unit, *) (zc(j), j = 5, 8)
220 if (i .ge. start_el .and. i .le. end_el) then
221 do j = 1, 8
222 call p(j)%init(real(xc(j), dp), real(yc(j), dp), &
223 real(zc(j), dp))
224 call rea_file_add_point(htp, p(j), pt_idx)
225 end do
226 ! swap vertices to keep symmetric vertex numbering in neko
227 call msh%add_element(el_idx, el_idx, &
228 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
229 end if
230 end if
231 if (i .ge. start_el .and. i .le. end_el) then
232 el_idx = el_idx + 1
233 end if
234 end do
235
236 call htp%free()
237
238 read(file_unit, *)
239 read(file_unit, *) ncurve
240 allocate(curve_data(5, 8, nelgv))
241 allocate(curve_element(nelgv))
242 allocate(curve_type(8, nelgv))
243 do i = 1, nelgv
244 curve_element(i) = .false.
245 do j = 1, 8
246 curve_type(j, i) = 0
247 do l = 1, 5
248 curve_data(l, j, i) = 0d0
249 end do
250 end do
251 end do
252 do i = 1, ncurve
253 read(file_unit, *) edge, el_idx, (curve(j), j = 1, 5), chtemp
254 do j = 1, 5
255 curve_data(j, edge, el_idx) = curve(j)
256 end do
257 curve_element(el_idx) = .true.
258 select case (trim(chtemp))
259 case ('s')
260 curve_type(edge, el_idx) = 1
261 curve_skip = .true.
262 case ('e')
263 curve_type(edge, el_idx) = 2
264 curve_skip = .true.
265 case ('C')
266 curve_type(edge, el_idx) = 3
267 case ('m')
268 curve_type(edge, el_idx) = 4
269 end select
270 end do
271 if (curve_skip) then
272 call neko_log%warning('Curve type: s, e are not supported, ' // &
273 'treating mesh as non-curved.')
274 else
275 do el_idx = 1, nelgv
276 if (curve_element(el_idx)) then
277 call msh%mark_curve_element(el_idx, &
278 curve_data(1, 1, el_idx), curve_type(1, el_idx))
279 end if
280 end do
281 end if
282 deallocate(curve_data)
283 deallocate(curve_element)
284 deallocate(curve_type)
285
286 ! Read fluid boundary conditions
287 read(file_unit, *)
288 read(file_unit, *)
289 if (.not. read_bcs) then ! Mark zones in the mesh
290 call neko_log%message("Reading boundary conditions", neko_log_debug)
291 allocate(cbc(6, nelgv))
292 allocate(bc_data(6, 2*ndim, nelgv))
293 off = 0
294 !Fix for different horrible .rea periodic bc formats.
295 if (nelgv .lt. 1000) off = 1
296 do i = 1, nelgv
297 if (i .ge. start_el .and. i .le. end_el) then
298 el_idx = i - start_el + 1
299 do j = 1, 2*ndim
300 read(file_unit, *) cbc(j, i), (bc_data(l, j, i), l = 1, 6)
301 sym_facet = facet_map(j)
302
303 select case (trim(cbc(j, i)))
304 case ('W')
305 if (neko_w_bc_label .eq. -1) then
306 neko_w_bc_label = current_internal_zone
307 current_internal_zone = current_internal_zone + 1
308 end if
309
310 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
311 cbc(j, i), &
312 neko_w_bc_label, total_labeled_zone_offset, &
313 labeled_zone_offsets(neko_w_bc_label) .eq. 0)
314
315 labeled_zone_offsets(neko_w_bc_label) = 1
316
317 case ('v', 'V')
318 if (neko_v_bc_label .eq. -1) then
319 neko_v_bc_label = current_internal_zone
320 current_internal_zone = current_internal_zone + 1
321 end if
322
323 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
324 cbc(j, i), &
325 neko_v_bc_label, total_labeled_zone_offset, &
326 labeled_zone_offsets(neko_v_bc_label) .eq. 0)
327
328 labeled_zone_offsets(neko_v_bc_label) = 1
329 case ('O', 'o')
330 if (neko_o_bc_label .eq. -1) then
331 neko_o_bc_label = current_internal_zone
332 current_internal_zone = current_internal_zone + 1
333 end if
334
335 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
336 cbc(j, i), &
337 neko_o_bc_label, total_labeled_zone_offset, &
338 labeled_zone_offsets(neko_o_bc_label) .eq. 0)
339
340 labeled_zone_offsets(neko_o_bc_label) = 1
341 case ('SYM')
342 if (neko_sym_bc_label .eq. -1) then
343 neko_sym_bc_label = current_internal_zone
344 current_internal_zone = current_internal_zone + 1
345 end if
346
347 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
348 cbc(j, i), &
349 neko_sym_bc_label, total_labeled_zone_offset, &
350 labeled_zone_offsets(neko_sym_bc_label) .eq. 0)
351
352 labeled_zone_offsets(neko_sym_bc_label) = 1
353 case ('ON', 'on')
354 if (neko_on_bc_label .eq. -1) then
355 neko_on_bc_label = current_internal_zone
356 current_internal_zone = current_internal_zone + 1
357 end if
358
359 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
360 cbc(j, i), &
361 neko_on_bc_label, total_labeled_zone_offset, &
362 labeled_zone_offsets(neko_on_bc_label) .eq. 0)
363
364 labeled_zone_offsets(neko_on_bc_label) = 1
365 case ('s', 'sl', 'sh', 'shl', 'S', 'SL', 'SH', 'SHL')
366 if (neko_shl_bc_label .eq. -1) then
367 neko_shl_bc_label = current_internal_zone
368 current_internal_zone = current_internal_zone + 1
369 end if
370
371 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
372 cbc(j, i), &
373 neko_shl_bc_label, total_labeled_zone_offset, &
374 labeled_zone_offsets(neko_shl_bc_label) .eq. 0)
375
376 labeled_zone_offsets(neko_shl_bc_label) = 1
377 case ('P')
378 p_el_idx = int(bc_data(2+off, j, i))
379 p_facet = facet_map(int(bc_data(3+off, j, i)))
380 call msh%get_facet_ids(sym_facet, el_idx, pids)
381 call msh%mark_periodic_facet(sym_facet, el_idx, &
382 p_facet, p_el_idx, pids)
383 end select
384 end do
385 end if
386 end do
387
388 do i = 1, nelgv
389 if (i .ge. start_el .and. i .le. end_el) then
390 el_idx = i - start_el + 1
391 do j = 1, 2*ndim
392 sym_facet = facet_map(j)
393 select case (trim(cbc(j, i)))
394 case ('P')
395 p_el_idx = int(bc_data(2+off, j, i))
396 p_facet = facet_map(int(bc_data(3+off, j, i)))
397 call msh%create_periodic_ids(sym_facet, el_idx, &
398 p_facet, p_el_idx)
399 end select
400 end do
401 end if
402 end do
403 do i = 1, nelgv
404 if (i .ge. start_el .and. i .le. end_el) then
405 el_idx = i - start_el + 1
406 do j = 1, 2*ndim
407 sym_facet = facet_map(j)
408 select case (trim(cbc(j, i)))
409 case ('P')
410 p_el_idx = int(bc_data(2+off, j, i))
411 p_facet = facet_map(int(bc_data(3+off, j, i)))
412 call msh%create_periodic_ids(sym_facet, el_idx, &
413 p_facet, p_el_idx)
414 end select
415 end do
416 end if
417 end do
418 do i = 1, nelgv
419 if (i .ge. start_el .and. i .le. end_el) then
420 el_idx = i - start_el + 1
421 do j = 1, 2*ndim
422 sym_facet = facet_map(j)
423 select case (trim(cbc(j, i)))
424 case ('P')
425 p_el_idx = int(bc_data(2+off, j, i))
426 p_facet = facet_map(int(bc_data(3+off, j, i)))
427 call msh%create_periodic_ids(sym_facet, el_idx, &
428 p_facet, p_el_idx)
429 end select
430 end do
431 end if
432 end do
433 deallocate(cbc)
434 deallocate(bc_data)
435 else ! Store bcs in a NEKTON session structure
436 allocate(cbc(6, nelgv))
437 do i = 1, nelgv
438 do j = 1, 2*ndim
439 read(file_unit, '(a1, a3)') chtemp, cbc(j, i)
440 end do
441 end do
442 end if
443
444 call msh%finalize()
445
446 call neko_log%message('Done')
447 close(file_unit)
448 end if
449
450 end subroutine rea_file_read
451
452 subroutine rea_file_write(this, data, t)
453 class(rea_file_t), intent(inout) :: this
454 class(*), target, intent(in) :: data
455 real(kind=rp), intent(in), optional :: t
456 end subroutine rea_file_write
457
458 subroutine rea_file_add_point(htp, p, idx)
459 type(htable_pt_t), intent(inout) :: htp
460 type(point_t), intent(inout) :: p
461 integer, intent(inout) :: idx
462 integer :: tmp
463
464 if (htp%get(p, tmp) .gt. 0) then
465 idx = idx + 1
466 call htp%set(p, idx)
467 call p%set_id(idx)
468 else
469 call p%set_id(tmp)
470 end if
471
472 end subroutine rea_file_add_point
473
484 subroutine rea_file_mark_labeled_bc(msh, el_idx, facet, type, label, &
485 offset, print_info)
486 type(mesh_t), intent(inout) :: msh
487 integer, intent(in) :: el_idx
488 integer, intent(in) :: facet
489 character(len=*), intent(in) :: type
490 integer, intent(in) :: label
491 integer, intent(in) :: offset
492 logical, intent(in) :: print_info
493
494 integer :: mark_label
495 character(len=LOG_SIZE) :: log_buf
496
497 mark_label = offset + label
498
499 if (mark_label .lt. 1 .or. mark_label .gt. neko_msh_max_zlbls) then
500 call neko_error( &
501 "You have reached the maximum amount of allowed labeled " // &
502 "zones (max allowed: 20). This happened when " // &
503 "converting re2 internal labels like e.g. " // &
504 "'w', 'V' or 'o' to labeled zones. Please " // &
505 "reduce the number of labeled zones that you " // &
506 "have defined or make sure that they are " // &
507 "labeled from [1, ..., 20].")
508 end if
509
510 if (print_info) then
511 write(log_buf, '(A3, A, I2)') trim(type), &
512 ' => Labeled index ', mark_label
513 call neko_log%message(log_buf)
514 end if
515 call msh%mark_labeled_facet(facet, el_idx, mark_label)
516
517 end subroutine rea_file_mark_labeled_bc
518
519end module rea_file
double real
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 a domain as a subset of facets in a mesh.
Definition curve.f90:2
Defines practical data distributions.
Definition datadist.f90:34
Module for file I/O operations.
Definition file.f90:34
Implements a hash table ADT.
Definition htable.f90:36
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:87
type(log_t), public neko_log
Global log stream.
Definition log.f90:77
integer, parameter, public log_size
Definition log.f90:46
NEKTON map file.
Definition map_file.f90:35
NEKTON map.
Definition map.f90:3
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:63
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements a point.
Definition point.f90:35
NEKTON mesh data in re2 format.
Definition re2_file.f90:35
integer, public neko_on_bc_label
Definition re2_file.f90:61
integer, public neko_w_bc_label
Definition re2_file.f90:57
integer, public neko_o_bc_label
Definition re2_file.f90:59
integer, public neko_v_bc_label
Definition re2_file.f90:58
integer, public neko_shl_bc_label
Definition re2_file.f90:62
integer, public neko_sym_bc_label
Definition re2_file.f90:60
NEKTON session data reader.
Definition rea_file.f90:35
subroutine rea_file_read(this, data)
Load NEKTON session data from an ascii file.
Definition rea_file.f90:71
subroutine rea_file_add_point(htp, p, idx)
Definition rea_file.f90:459
subroutine rea_file_write(this, data, t)
Definition rea_file.f90:453
subroutine rea_file_mark_labeled_bc(msh, el_idx, facet, type, label, offset, print_info)
Mark a labeled zone based on an offset.
Definition rea_file.f90:486
NEKTON session data.
Definition rea.f90:4
subroutine, public rea_free(r)
Free a NEKTON session data.
Definition rea.f90:26
Utilities.
Definition utils.f90:35
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:140
Load-balanced linear distribution .
Definition datadist.f90:50
A generic file handler.
Point based hash table.
Definition htable.f90:112
NEKTON vertex mapping.
Definition map.f90:8
Interface for NEKTON map files.
Definition map_file.f90:45
A point in with coordinates .
Definition point.f90:43
Interface for NEKTON re2 files.
Definition re2_file.f90:65
NEKTON session data struct.
Definition rea.f90:13
Interface for NEKTON ascii files.
Definition rea_file.f90:61