Neko 1.99.1
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
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
131 open(newunit=file_unit,file=trim(this%fname), status='old', iostat=ierr)
132 call neko_log%message('Reading NEKTON file ' // this%fname)
133
134 read(file_unit, *)
135 read(file_unit, *)
136 read(file_unit, *) ndim
137 read(file_unit, *) nparam
138
139 if (.not. read_param) then
140 ! Skip parameters
141 do i = 1, nparam
142 read(file_unit, *)
143 end do
144 else
145 allocate(params(nparam))
146 do i = 1, nparam
147 read(file_unit, *) params(i)
148 end do
149 end if
150
151 ! Skip passive scalars
152 read(file_unit, *) nskip
153 do i = 1, nskip
154 read(file_unit, *)
155 end do
156
157 ! Skip logic switches
158 read(file_unit, *) nlogic
159 do i = 1, nlogic
160 read(file_unit, *)
161 end do
162
163 ! Read mesh info
164 read(file_unit, *)
165 read(file_unit, *)
166 read(file_unit, *) nelgs,ndim, nelgv
167 if (nelgs .lt. 0) then
168 re2_fname = trim(this%fname(1:scan(trim(this%fname), &
169 '.', 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(this%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 p(j) = point_t(real(xc(j),dp), real(yc(j),dp),real(0d0,dp))
207 call rea_file_add_point(htp, p(j), pt_idx)
208 end do
209 ! swap vertices to keep symmetric vertex numbering in neko
210 call msh%add_element(el_idx, el_idx, p(1), p(2), p(4), p(3))
211 end if
212 else if (ndim .eq. 3) then
213 read(file_unit, *) (xc(j),j=1,4)
214 read(file_unit, *) (yc(j),j=1,4)
215 read(file_unit, *) (zc(j),j=1,4)
216 read(file_unit, *) (xc(j),j=5,8)
217 read(file_unit, *) (yc(j),j=5,8)
218 read(file_unit, *) (zc(j),j=5,8)
219 if (i .ge. start_el .and. i .le. end_el) then
220 do j = 1, 8
221 p(j) = point_t(real(xc(j),dp), real(yc(j),dp), real(zc(j),dp))
222 call rea_file_add_point(htp, p(j), pt_idx)
223 end do
224 ! swap vertices to keep symmetric vertex numbering in neko
225 call msh%add_element(el_idx, el_idx, &
226 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
227 end if
228 end if
229 if (i .ge. start_el .and. i .le. end_el) then
230 el_idx = el_idx + 1
231 end if
232 end do
233
234 call htp%free()
235
236 read(file_unit, *)
237 read(file_unit, *) ncurve
238 allocate(curve_data(5,8,nelgv))
239 allocate(curve_element(nelgv))
240 allocate(curve_type(8,nelgv))
241 do i = 1, nelgv
242 curve_element(i) = .false.
243 do j = 1, 8
244 curve_type(j,i) = 0
245 do l = 1, 5
246 curve_data(l,j,i) = 0d0
247 end do
248 end do
249 end do
250 do i = 1, ncurve
251 read(file_unit, *) edge, el_idx, (curve(j),j=1,5), chtemp
252 do j = 1, 5
253 curve_data(j,edge,el_idx) = curve(j)
254 end do
255 curve_element(el_idx) = .true.
256 select case(trim(chtemp))
257 case ('s')
258 curve_type(edge,el_idx) = 1
259 curve_skip = .true.
260 case ('e')
261 curve_type(edge,el_idx) = 2
262 curve_skip = .true.
263 case ('C')
264 curve_type(edge,el_idx) = 3
265 case ('m')
266 curve_type(edge,el_idx) = 4
267 end select
268 end do
269 if (curve_skip) then
270 call neko_log%warning('Curve type: s, e are not supported, treating mesh as non-curved.')
271 else
272 do el_idx = 1, nelgv
273 if (curve_element(el_idx)) then
274 call msh%mark_curve_element(el_idx, &
275 curve_data(1,1,el_idx), curve_type(1,el_idx))
276 end if
277 end do
278 end if
279 deallocate(curve_data)
280 deallocate(curve_element)
281 deallocate(curve_type)
282
283 ! Read fluid boundary conditions
284 read(file_unit,*)
285 read(file_unit,*)
286 if (.not. read_bcs) then ! Mark zones in the mesh
287 call neko_log%message("Reading boundary conditions", neko_log_debug)
288 allocate(cbc(6,nelgv))
289 allocate(bc_data(6,2*ndim,nelgv))
290 off = 0
291 !Fix for different horrible .rea periodic bc formats.
292 if (nelgv .lt. 1000) off = 1
293 do i = 1, nelgv
294 if (i .ge. start_el .and. i .le. end_el) then
295 el_idx = i - start_el + 1
296 do j = 1, 2*ndim
297 read(file_unit, *) cbc(j, i), (bc_data(l,j,i),l=1,6)
298 sym_facet = facet_map(j)
299
300 select case(trim(cbc(j,i)))
301 case ('W')
302 if (neko_w_bc_label .eq. -1) then
303 neko_w_bc_label = current_internal_zone
304 current_internal_zone = current_internal_zone + 1
305 end if
306
307 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
308 cbc(j,i), &
309 neko_w_bc_label, total_labeled_zone_offset, &
310 labeled_zone_offsets(neko_w_bc_label) .eq. 0)
311
312 labeled_zone_offsets(neko_w_bc_label) = 1
313
314 case ('v', 'V')
315 if (neko_v_bc_label .eq. -1) then
316 neko_v_bc_label = current_internal_zone
317 current_internal_zone = current_internal_zone + 1
318 end if
319
320 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
321 cbc(j,i), &
322 neko_v_bc_label, total_labeled_zone_offset, &
323 labeled_zone_offsets(neko_v_bc_label) .eq. 0)
324
325 labeled_zone_offsets(neko_v_bc_label) = 1
326 case ('O', 'o')
327 if (neko_o_bc_label .eq. -1) then
328 neko_o_bc_label = current_internal_zone
329 current_internal_zone = current_internal_zone + 1
330 end if
331
332 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
333 cbc(j,i), &
334 neko_o_bc_label, total_labeled_zone_offset, &
335 labeled_zone_offsets(neko_o_bc_label) .eq. 0)
336
337 labeled_zone_offsets(neko_o_bc_label) = 1
338 case ('SYM')
339 if (neko_sym_bc_label .eq. -1) then
340 neko_sym_bc_label = current_internal_zone
341 current_internal_zone = current_internal_zone + 1
342 end if
343
344 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
345 cbc(j,i), &
346 neko_sym_bc_label, total_labeled_zone_offset, &
347 labeled_zone_offsets(neko_sym_bc_label) .eq. 0)
348
349 labeled_zone_offsets(neko_sym_bc_label) = 1
350 case ('ON', 'on')
351 if (neko_on_bc_label .eq. -1) then
352 neko_on_bc_label = current_internal_zone
353 current_internal_zone = current_internal_zone + 1
354 end if
355
356 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
357 cbc(j,i), &
358 neko_on_bc_label, total_labeled_zone_offset, &
359 labeled_zone_offsets(neko_on_bc_label) .eq. 0)
360
361 labeled_zone_offsets(neko_on_bc_label) = 1
362 case ('s', 'sl', 'sh', 'shl', 'S', 'SL', 'SH', 'SHL')
363 if (neko_shl_bc_label .eq. -1) then
364 neko_shl_bc_label = current_internal_zone
365 current_internal_zone = current_internal_zone + 1
366 end if
367
368 call rea_file_mark_labeled_bc(msh, el_idx, sym_facet, &
369 cbc(j,i), &
370 neko_shl_bc_label, total_labeled_zone_offset, &
371 labeled_zone_offsets(neko_shl_bc_label) .eq. 0)
372
373 labeled_zone_offsets(neko_shl_bc_label) = 1
374 case ('P')
375 p_el_idx = int(bc_data(2+off,j,i))
376 p_facet = facet_map(int(bc_data(3+off,j,i)))
377 call msh%get_facet_ids(sym_facet, el_idx, pids)
378 call msh%mark_periodic_facet(sym_facet, el_idx, &
379 p_facet, p_el_idx, pids)
380 end select
381 end do
382 end if
383 end do
384
385 do i = 1, nelgv
386 if (i .ge. start_el .and. i .le. end_el) then
387 el_idx = i - start_el + 1
388 do j = 1, 2*ndim
389 sym_facet = facet_map(j)
390 select case(trim(cbc(j,i)))
391 case ('P')
392 p_el_idx = int(bc_data(2+off,j,i))
393 p_facet = facet_map(int(bc_data(3+off,j,i)))
394 call msh%create_periodic_ids(sym_facet, el_idx, &
395 p_facet, p_el_idx)
396 end select
397 end do
398 end if
399 end do
400 do i = 1, nelgv
401 if (i .ge. start_el .and. i .le. end_el) then
402 el_idx = i - start_el + 1
403 do j = 1, 2*ndim
404 sym_facet = facet_map(j)
405 select case(trim(cbc(j,i)))
406 case ('P')
407 p_el_idx = int(bc_data(2+off,j,i))
408 p_facet = facet_map(int(bc_data(3+off,j,i)))
409 call msh%create_periodic_ids(sym_facet, el_idx, &
410 p_facet, p_el_idx)
411 end select
412 end do
413 end if
414 end do
415 do i = 1, nelgv
416 if (i .ge. start_el .and. i .le. end_el) then
417 el_idx = i - start_el + 1
418 do j = 1, 2*ndim
419 sym_facet = facet_map(j)
420 select case(trim(cbc(j,i)))
421 case ('P')
422 p_el_idx = int(bc_data(2+off,j,i))
423 p_facet = facet_map(int(bc_data(3+off,j,i)))
424 call msh%create_periodic_ids(sym_facet, el_idx, &
425 p_facet, p_el_idx)
426 end select
427 end do
428 end if
429 end do
430 deallocate(cbc)
431 deallocate(bc_data)
432 else ! Store bcs in a NEKTON session structure
433 allocate(cbc(6,nelgv))
434 do i = 1, nelgv
435 do j = 1, 2*ndim
436 read(file_unit,'(a1, a3)') chtemp, cbc(j, i)
437 end do
438 end do
439 end if
440
441 call msh%finalize()
442
443 call neko_log%message('Done')
444 close(file_unit)
445 endif
446
447 end subroutine rea_file_read
448
449 subroutine rea_file_write(this, data, t)
450 class(rea_file_t), intent(inout) :: this
451 class(*), target, intent(in) :: data
452 real(kind=rp), intent(in), optional :: t
453 end subroutine rea_file_write
454
455 subroutine rea_file_add_point(htp, p, idx)
456 type(htable_pt_t), intent(inout) :: htp
457 type(point_t), intent(inout) :: p
458 integer, intent(inout) :: idx
459 integer :: tmp
460
461 if (htp%get(p, tmp) .gt. 0) then
462 idx = idx + 1
463 call htp%set(p, idx)
464 call p%set_id(idx)
465 else
466 call p%set_id(tmp)
467 end if
468
469 end subroutine rea_file_add_point
470
481 subroutine rea_file_mark_labeled_bc(msh, el_idx, facet, type, label, offset, print_info)
482 type(mesh_t), intent(inout) :: msh
483 integer, intent(in) :: el_idx
484 integer, intent(in) :: facet
485 character(len=*), intent(in) :: type
486 integer, intent(in) :: label
487 integer, intent(in) :: offset
488 logical, intent(in) :: print_info
489
490 integer :: mark_label
491 character(len=LOG_SIZE) :: log_buf
492
493 mark_label = offset + label
494
495 if (mark_label .lt. 1 .or. mark_label .gt. neko_msh_max_zlbls) then
496 call neko_error("You have reached the maximum amount of allowed labeled&
497 & zones (max allowed: 20). This happened when converting re2 internal labels&
498 & like e.g. 'w', 'V' or 'o' to labeled zones. Please reduce the number of&
499 & labeled zones that you have defined or make sure that they are labeled&
500 & from [1,...,20].")
501 end if
502
503 if (print_info) then
504 write (log_buf, "(A3,A,I2)") trim(type), " => Labeled index ", mark_label
505 call neko_log%message(log_buf)
506 end if
507 call msh%mark_labeled_facet(facet, el_idx, mark_label)
508
509 end subroutine rea_file_mark_labeled_bc
510
511end module rea_file
double real
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:58
integer, public pe_rank
MPI rank.
Definition comm.F90:55
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
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:78
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
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:62
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:456
subroutine rea_file_write(this, data, t)
Definition rea_file.f90:450
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:482
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:78
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