Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
re2_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, i8, sp, dp
38 use map_file, only : map_file_t
39 use mesh, only: mesh_t
40 use point, only: point_t
42 use datadist, only : linear_dist_t
43 use map, only : map_t
44 use mesh, only : neko_msh_max_zlbls
45 use logger, only: neko_log, log_size
50 use htable, only : htable_pt_t
51 use comm
52 use mpi_f08
53 implicit none
54 private
55
56 ! Defines the conventions for conversion of re2 labels to labeled zones.
57 integer, public :: neko_w_bc_label = -1
58 integer, public :: neko_v_bc_label = -1
59 integer, public :: neko_o_bc_label = -1
60 integer, public :: neko_sym_bc_label = -1
61 integer, public :: neko_on_bc_label = -1
62 integer, public :: neko_shl_bc_label = -1
63
65 type, public, extends(generic_file_t) :: re2_file_t
66 contains
67 procedure :: read => re2_file_read
68 procedure :: write => re2_file_write
69 end type re2_file_t
70
71contains
72
74 subroutine re2_file_read(this, data)
75 class(re2_file_t) :: this
76 class(*), target, intent(inout) :: data
77 type(mesh_t), pointer :: msh
78 integer :: file_unit
79 character(len=5) :: hdr_ver
80 character(len=54) :: hdr_str
81 character(len=80) :: hdr_full
82 integer :: nel, ndim, nelv, ierr, nBCre2
83 type(mpi_status) :: status
84 type(mpi_file) :: fh
85 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
86 real(kind=sp) :: test
87 real(kind=dp) :: t2
88 integer :: ncurv, nbcs
89 type(linear_dist_t) :: dist
90 type(map_t) :: nm
91 type(map_file_t) :: map_file
92 character(len=1024) :: map_fname
93 logical :: read_map
94 integer :: re2_data_xy_size
95 integer :: re2_data_xyz_size
96 integer :: re2_data_cv_size
97 integer :: re2_data_bc_size
98 logical :: v2_format
99 character(len=LOG_SIZE) :: log_buf
100
101 call this%check_exists()
102
103 select type(data)
104 type is (mesh_t)
105 msh => data
106 class default
107 call neko_error('Invalid output data')
108 end select
109
110 v2_format = .false.
111 open(newunit=file_unit,file=trim(this%fname), status='old', iostat=ierr)
112 call neko_log%message('Reading binary NEKTON file ' // this%fname)
113
114 read(file_unit,'(a80)') hdr_full
115 read(hdr_full, '(a5)') hdr_ver
116
117 if (hdr_ver .eq. '#v004') then
118 read(hdr_full, '(a5,i16,i3,i16,i4,a36)') hdr_ver, nel, ndim, nelv, nbcre2, hdr_str
119 v2_format = .true.
120 else if (hdr_ver .eq. '#v002' .or. hdr_ver .eq. '#v003') then
121 read(hdr_full, '(a5,i9,i3,i9,a54)') hdr_ver, nel, ndim, nelv, hdr_str
122 v2_format = .true.
123 else if (hdr_ver .eq. '#v001') then
124 read(hdr_full, '(a5,i9,i3,i9,a54)') hdr_ver, nel, ndim, nelv, hdr_str
125 end if
126
127 if (v2_format) then
128 call mpi_type_size(mpi_re2v2_data_xy, re2_data_xy_size, ierr)
129 call mpi_type_size(mpi_re2v2_data_xyz, re2_data_xyz_size, ierr)
130 call mpi_type_size(mpi_re2v2_data_cv, re2_data_cv_size, ierr)
131 call mpi_type_size(mpi_re2v2_data_bc, re2_data_bc_size, ierr)
132 else
133 call mpi_type_size(mpi_re2v1_data_xy, re2_data_xy_size, ierr)
134 call mpi_type_size(mpi_re2v1_data_xyz, re2_data_xyz_size, ierr)
135 call mpi_type_size(mpi_re2v1_data_cv, re2_data_cv_size, ierr)
136 call mpi_type_size(mpi_re2v1_data_bc, re2_data_bc_size, ierr)
137 end if
138
139 write(log_buf,1) ndim, nelv
1401 format('gdim = ', i1, ', nelements =', i7)
141 call neko_log%message(log_buf)
142 close(file_unit)
143
144 call filename_chsuffix(this%fname, map_fname,'map')
145
146 inquire(file=map_fname, exist=read_map)
147 if (read_map) then
148 call nm%init(nelv, 2**ndim)
149 call map_file%init(map_fname)
150 call map_file%read(nm)
151 else
152 call neko_log%warning('No NEKTON map file found')
153 end if
154
155 call mpi_file_open(neko_comm, trim(this%fname), &
156 mpi_mode_rdonly, mpi_info_null, fh, ierr)
157
158 if (ierr .ne. 0) then
159 call neko_log%error("Can't open binary NEKTON file ")
160 end if
161 dist = linear_dist_t(nelv, pe_rank, pe_size, neko_comm)
162
163
164 call msh%init(ndim, dist)
165
166 ! Set offset (header)
167 mpi_offset = re2_hdr_size * mpi_character_size
168
169 call mpi_file_read_at_all(fh, mpi_offset, test, 1, mpi_real, status, ierr)
170 mpi_offset = mpi_offset + mpi_real_size
171
172 if (abs(re2_endian_test - test) .gt. 1e-4) then
173 call neko_error('Invalid endian of re2 file, byte swap not implemented yet')
174 end if
175
176 call re2_file_read_points(msh, ndim, nel, dist, fh, &
177 mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
178
179
180 ! Set offset to start of curved side data
182 if (ndim .eq. 2) then
183 mpi_offset = mpi_offset + int(dist%num_global(),i8) * int(re2_data_xy_size,i8)
184 else
185 mpi_offset = mpi_offset + int(dist%num_global(),i8) * int(re2_data_xyz_size,i8)
186 end if
187
190 if (v2_format) then
191 call mpi_file_read_at_all(fh, mpi_offset, t2, 1, mpi_double_precision, status, ierr)
192 ncurv = int(t2)
193 mpi_offset = mpi_offset + mpi_double_precision_size
194 call re2_file_read_curve(msh, ncurv, dist, fh, mpi_offset, v2_format)
195 mpi_offset = mpi_offset + int(ncurv,i8) * int(re2_data_cv_size,i8)
196 call mpi_file_read_at_all(fh, mpi_offset, t2, 1, mpi_double_precision, status, ierr)
197 nbcs = int(t2)
198 mpi_offset = mpi_offset + mpi_double_precision_size
199
200 call re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
201 else
202 call mpi_file_read_at_all(fh, mpi_offset, ncurv, 1, mpi_integer, status, ierr)
203 mpi_offset = mpi_offset + mpi_integer_size
204 call re2_file_read_curve(msh, ncurv, dist, fh, mpi_offset, v2_format)
205 mpi_offset = mpi_offset + int(ncurv,i8) * int(re2_data_cv_size,i8)
206 call mpi_file_read_at_all(fh, mpi_offset, nbcs, 1, mpi_integer, status, ierr)
207 mpi_offset = mpi_offset + mpi_integer_size
208
209 call re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
210
211 end if
212
213 call mpi_file_close(fh, ierr)
214 call msh%finalize()
215
216
217 call neko_log%message('Done')
218
219
220 end subroutine re2_file_read
221
222 subroutine re2_file_write(this, data, t)
223 class(re2_file_t), intent(inout) :: this
224 class(*), target, intent(in) :: data
225 real(kind=rp), intent(in), optional :: t
226 type(re2v1_xy_t), allocatable :: re2_data_xy(:)
227 type(re2v1_xyz_t), allocatable :: re2_data_xyz(:)
228 type(mesh_t), pointer :: msh
229 character(len=5), parameter :: RE2_HDR_VER = '#v001'
230 character(len=54), parameter :: RE2_HDR_STR = 'RE2 exported by NEKO'
231 integer :: i, j, ierr, nelgv
232 type(mpi_status) :: status
233 type(mpi_file) :: fh
234 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
235 integer :: element_offset
236 integer :: re2_data_xy_size
237 integer :: re2_data_xyz_size
238
239 select type(data)
240 type is (mesh_t)
241 msh => data
242 class default
243 call neko_error('Invalid output data')
244 end select
245
246 call mpi_type_size(mpi_re2v1_data_xy, re2_data_xy_size, ierr)
247 call mpi_type_size(mpi_re2v1_data_xyz, re2_data_xyz_size, ierr)
248 call mpi_reduce(msh%nelv, nelgv, 1, &
249 mpi_integer, mpi_sum, 0, neko_comm, ierr)
250 element_offset = 0
251 call mpi_exscan(msh%nelv, element_offset, 1, &
252 mpi_integer, mpi_sum, neko_comm, ierr)
253
254 call neko_log%message('Writing data as a binary NEKTON file ' // this%fname)
255
256 if (pe_rank .eq. 0) then
257 open(unit=9,file=trim(this%fname), status='new', iostat=ierr)
258 write(9, '(a5,i9,i3,i9,a54)') re2_hdr_ver, nelgv, msh%gdim,&
259 nelgv, re2_hdr_str
260 close(9)
261 end if
262
263 call mpi_barrier(neko_comm, ierr)
264 call mpi_file_open(neko_comm, trim(this%fname), &
265 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
266 mpi_offset = re2_hdr_size * mpi_character_size
267
268 call mpi_file_write_at(fh, mpi_offset, re2_endian_test, 1, &
269 mpi_real, status, ierr)
270 mpi_offset = mpi_offset + mpi_real_size
271
272 if (msh%gdim .eq. 2) then
273 allocate(re2_data_xy(msh%nelv))
274 do i = 1, msh%nelv
275 re2_data_xy(i)%rgroup = 1.0 ! Not used
276 do j = 1, 4
277 re2_data_xy(i)%x(j) = real(msh%elements(i)%e%pts(j)%p%x(1))
278 re2_data_xy(i)%y(j) = real(msh%elements(i)%e%pts(j)%p%x(2))
279 end do
280 end do
281 mpi_offset = mpi_offset + element_offset * re2_data_xy_size
282 call mpi_file_write_at(fh, mpi_offset, &
283 re2_data_xy, msh%nelv, mpi_re2v1_data_xy, status, ierr)
284
285 deallocate(re2_data_xy)
286
287 else if (msh%gdim .eq. 3) then
288 allocate(re2_data_xyz(msh%nelv))
289 do i = 1, msh%nelv
290 re2_data_xyz(i)%rgroup = 1.0 ! Not used
291 do j = 1, 8
292 re2_data_xyz(i)%x(j) = real(msh%elements(i)%e%pts(j)%p%x(1))
293 re2_data_xyz(i)%y(j) = real(msh%elements(i)%e%pts(j)%p%x(2))
294 re2_data_xyz(i)%z(j) = real(msh%elements(i)%e%pts(j)%p%x(3))
295 end do
296 end do
297 mpi_offset = mpi_offset + element_offset * re2_data_xyz_size
298 call mpi_file_write_at(fh, mpi_offset, &
299 re2_data_xyz, msh%nelv, mpi_re2v1_data_xyz, status, ierr)
300
301 deallocate(re2_data_xyz)
302 else
303 call neko_error("Invalid dimension of mesh")
304 end if
305
306 call mpi_file_close(fh, ierr)
307 call neko_log%message('Done')
308
310
311 end subroutine re2_file_write
312
313 subroutine re2_file_read_points(msh, ndim, nel, dist, fh, &
314 mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
315 type(mesh_t), intent(inout) :: msh
316 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
317 integer, intent(inout) :: ndim
318 integer, intent(inout) :: nel
319 type(mpi_file), intent(inout) :: fh
320 integer, intent(in) :: re2_data_xy_size
321 integer, intent(in) :: re2_data_xyz_size
322 logical, intent(in) :: v2_format
323 type(linear_dist_t) :: dist
324 integer :: element_offset
325 type(re2v1_xy_t), allocatable :: re2v1_data_xy(:)
326 type(re2v1_xyz_t), allocatable :: re2v1_data_xyz(:)
327 type(re2v2_xy_t), allocatable :: re2v2_data_xy(:)
328 type(re2v2_xyz_t), allocatable :: re2v2_data_xyz(:)
329 type(mpi_status) :: status
330 type(htable_pt_t) :: htp
331 type(point_t) :: p(8)
332 integer :: pt_idx, nelv
333 integer :: i, j, ierr
334
335 nelv = dist%num_local()
336 element_offset = dist%start_idx()
337
338 call htp%init(2*nel, ndim)
339 pt_idx = 0
340 if (ndim .eq. 2) then
341 mpi_offset = mpi_offset + element_offset * re2_data_xy_size
342 if (.not. v2_format) then
343 allocate(re2v1_data_xy(nelv))
344 call mpi_file_read_at_all(fh, mpi_offset, &
345 re2v1_data_xy, nelv, mpi_re2v1_data_xy, status, ierr)
346 do i = 1, nelv
347 do j = 1, 4
348 p(j) = point_t(real(re2v1_data_xy(i)%x(j),dp), &
349 real(re2v1_data_xy(i)%y(j),dp), 0.0d0)
350 call re2_file_add_point(htp, p(j), pt_idx)
351 end do
352 if (nelv > 10) then
353 if(mod(i,nelv/10) .eq. 0) write(*,*) i, 'elements read'
354 end if
355 ! swap vertices to keep symmetric vertex numbering in neko
356 call msh%add_element(i, i, p(1), p(2), p(4), p(3))
357 end do
358 deallocate(re2v1_data_xy)
359 else
360 allocate(re2v2_data_xy(nelv))
361 call mpi_file_read_at_all(fh, mpi_offset, &
362 re2v2_data_xy, nelv, mpi_re2v2_data_xy, status, ierr)
363 do i = 1, nelv
364 do j = 1, 4
365 p(j) = point_t(re2v2_data_xy(i)%x(j), &
366 re2v2_data_xy(i)%y(j), 0.0d0)
367 call re2_file_add_point(htp, p(j), pt_idx)
368 end do
369 if (nelv > 10) then
370 if(mod(i,nelv/10) .eq. 0) write(*,*) i, 'elements read'
371 end if
372 ! swap vertices to keep symmetric vertex numbering in neko
373 call msh%add_element(i, i, p(1), p(2), p(4), p(3))
374 end do
375 deallocate(re2v2_data_xy)
376 end if
377 else if (ndim .eq. 3) then
378 mpi_offset = mpi_offset + element_offset * re2_data_xyz_size
379 if (.not. v2_format) then
380 allocate(re2v1_data_xyz(nelv))
381 call mpi_file_read_at_all(fh, mpi_offset, &
382 re2v1_data_xyz, nelv, mpi_re2v1_data_xyz, status, ierr)
383 do i = 1, nelv
384 do j = 1, 8
385 p(j) = point_t(real(re2v1_data_xyz(i)%x(j),dp), &
386 real(re2v1_data_xyz(i)%y(j),dp),&
387 real(re2v1_data_xyz(i)%z(j),dp))
388 call re2_file_add_point(htp, p(j), pt_idx)
389 end do
390 if (nelv > 100) then
391 if(mod(i,nelv/100) .eq. 0) write(*,*) i, 'elements read'
392 end if
393 ! swap vertices to keep symmetric vertex numbering in neko
394 call msh%add_element(i, i, &
395 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
396 end do
397 deallocate(re2v1_data_xyz)
398 else
399 allocate(re2v2_data_xyz(nelv))
400 call mpi_file_read_at_all(fh, mpi_offset, &
401 re2v2_data_xyz, nelv, mpi_re2v2_data_xyz, status, ierr)
402 do i = 1, nelv
403 do j = 1, 8
404 p(j) = point_t(re2v2_data_xyz(i)%x(j), &
405 re2v2_data_xyz(i)%y(j),&
406 re2v2_data_xyz(i)%z(j))
407 call re2_file_add_point(htp, p(j), pt_idx)
408 end do
409 if (nelv > 100) then
410 if(mod(i,nelv/100) .eq. 0) write(*,*) i, 'elements read'
411 end if
412 ! swap vertices to keep symmetric vertex numbering in neko
413 call msh%add_element(i, i, &
414 p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
415 end do
416 deallocate(re2v2_data_xyz)
417 end if
418 end if
419
420 call htp%free()
421 end subroutine re2_file_read_points
422
423 subroutine re2_file_read_curve(msh, ncurve, dist, fh, mpi_offset, v2_format)
424 type(mesh_t), intent(inout) :: msh
425 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
426 integer, intent(inout) :: ncurve
427 type(linear_dist_t) :: dist
428 type(mpi_file), intent(inout) :: fh
429 logical, intent(in) :: v2_format
430 type(mpi_status) :: status
431 integer :: i, j, l, ierr, el_idx, id
432 type(re2v1_curve_t), allocatable :: re2v1_data_curve(:)
433 type(re2v2_curve_t), allocatable :: re2v2_data_curve(:)
434 real(kind=dp), allocatable :: curve_data(:,:,:)
435 integer, allocatable :: curve_type(:,:)
436 logical, allocatable :: curve_element(:)
437 character(len=1) :: chtemp
438 logical :: curve_skip = .false.
439
440 allocate(curve_data(5,12,msh%nelv))
441 allocate(curve_element(msh%nelv))
442 allocate(curve_type(12,msh%nelv))
443 do i = 1, msh%nelv
444 curve_element(i) = .false.
445 do j = 1, 12
446 curve_type(j,i) = 0
447 do l = 1, 5
448 curve_data(l,j,i) = 0d0
449 end do
450 end do
451 end do
452
453 call neko_log%message('Reading curved data')
454 if (.not. v2_format) then
455 allocate(re2v1_data_curve(ncurve))
456 call mpi_file_read_at_all(fh, mpi_offset, re2v1_data_curve, ncurve, &
457 mpi_re2v1_data_cv, status, ierr)
458 else
459 allocate(re2v2_data_curve(ncurve))
460 call mpi_file_read_at_all(fh, mpi_offset, re2v2_data_curve, ncurve, &
461 mpi_re2v2_data_cv, status, ierr)
462 end if
463 !This can probably be made nicer...
464 do i = 1, ncurve
465 if(v2_format) then
466 el_idx = re2v2_data_curve(i)%elem - dist%start_idx()
467 id = re2v2_data_curve(i)%zone
468 chtemp = re2v2_data_curve(i)%type
469 do j = 1, 5
470 curve_data(j,id, el_idx) = re2v2_data_curve(i)%point(j)
471 enddo
472 else
473 el_idx = re2v1_data_curve(i)%elem - dist%start_idx()
474 id = re2v1_data_curve(i)%zone
475 chtemp = re2v1_data_curve(i)%type
476 do j = 1, 5
477 curve_data(j,id, el_idx) = real(re2v1_data_curve(i)%point(j),dp)
478 enddo
479 end if
480
481 curve_element(el_idx) = .true.
482 !This might need to be extended
483 select case(trim(chtemp))
484 case ('s')
485 curve_type(id,el_idx) = 1
486 call neko_log%warning('curve type s not supported, treating mesh as non-curved')
487 curve_skip = .true.
488 exit
489 case ('e')
490 curve_type(id,el_idx) = 2
491 call neko_log%warning('curve type e not supported, treating mesh as non-curved')
492 curve_skip = .true.
493 exit
494 case ('C')
495 curve_type(id,el_idx) = 3
496 case ('m')
497 curve_type(id,el_idx) = 4
498 case default
499 write(*,*) chtemp, 'curve type not supported yet, treating mesh as non-curved',id, el_idx
500
501 curve_skip = .true.
502 end select
503 end do
504
505 if( v2_format) then
506 deallocate(re2v2_data_curve)
507 else
508 deallocate(re2v1_data_curve)
509 end if
510 if (.not. curve_skip) then
511 do el_idx = 1, msh%nelv
512 if (curve_element(el_idx)) then
513 call msh%mark_curve_element(el_idx, &
514 curve_data(1,1,el_idx), curve_type(1,el_idx))
515 end if
516 end do
517 end if
518
519 deallocate(curve_data)
520 deallocate(curve_element)
521 deallocate(curve_type)
522 end subroutine re2_file_read_curve
523
524
525 subroutine re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
526 type(mesh_t), intent(inout) :: msh
527 integer (kind=MPI_OFFSET_KIND) :: mpi_offset
528 integer, intent(inout) :: nbcs
529 type(linear_dist_t) :: dist
530 type(mpi_file), intent(inout) :: fh
531 logical, intent(in) :: v2_format
532 type(mpi_status) :: status
533 integer :: pids(4)
534 integer :: sym_facet, label
535 integer :: p_el_idx, p_facet
536 integer :: i, j, ierr, el_idx
537 integer, parameter, dimension(6) :: facet_map = (/3, 2, 4, 1, 5, 6/)
538 logical :: periodic
539 type(re2v1_bc_t), allocatable :: re2v1_data_bc(:)
540 type(re2v2_bc_t), allocatable :: re2v2_data_bc(:)
541 character(len=LOG_SIZE) :: log_buf
542
543 ! ---- Offsets for conversion from internal zones to labeled zones
544 integer :: user_labeled_zones(NEKO_MSH_MAX_ZLBLS)
545 integer :: labeled_zone_offsets(NEKO_MSH_MAX_ZLBLS)
546 integer :: total_labeled_zone_offset, current_internal_zone
547 total_labeled_zone_offset = 0
548 labeled_zone_offsets = 0
549 user_labeled_zones = 0
550 current_internal_zone = 1
551 ! ----
552
553 call neko_log%message("Reading boundary conditions")
554 if (.not. v2_format) then
555 allocate(re2v1_data_bc(nbcs))
556 call mpi_file_read_at_all(fh, mpi_offset, re2v1_data_bc, nbcs, &
557 mpi_re2v1_data_bc, status, ierr)
558 else
559 allocate(re2v2_data_bc(nbcs))
560 call mpi_file_read_at_all(fh, mpi_offset, re2v2_data_bc, nbcs, &
561 mpi_re2v2_data_bc, status, ierr)
562 end if
563
564 periodic = .false.
565
567 if (v2_format) then ! V2 format
568
569 !
570 ! Search for user labeled zones
571 !
572 do i = 1, nbcs
573
574 el_idx = int(re2v2_data_bc(i)%elem) - dist%start_idx()
575 sym_facet = facet_map(int(re2v2_data_bc(i)%face))
576
577 select case(trim(re2v2_data_bc(i)%type))
578 case ('MSH', 'msh', 'EXO','exo')
579
580 label = int(re2v2_data_bc(i)%bc_data(5))
581
582 if (label .lt. 1 .or. label .gt. neko_msh_max_zlbls) then
583 call neko_error('Invalid label id (valid range [1,...,20])')
584 end if
585
586 if (user_labeled_zones(label) .eq. 0) then
587 write (log_buf, "(A,I2,A)") "Labeled zone ", label, &
588 " found."
589 call neko_log%message(log_buf)
590 user_labeled_zones(label) = 1
591 end if
592
593 call msh%mark_labeled_facet(sym_facet, el_idx, label)
594
595 ! Get the largest label possible in case we need to convert
596 ! re2 labels to labeled zones (see below).
597 total_labeled_zone_offset = max(total_labeled_zone_offset, label)
598
599
600 end select
601 end do
602
603 do i = 1, nbcs
604 el_idx = int(re2v2_data_bc(i)%elem) - dist%start_idx()
605 sym_facet = facet_map(int(re2v2_data_bc(i)%face))
606 select case(trim(re2v2_data_bc(i)%type))
607 case ('MSH', 'msh', 'EXO','exo')
608 !Do nothing, already handled
609 case ('W')
610 if (neko_w_bc_label .eq. -1) then
611 neko_w_bc_label = current_internal_zone
612 current_internal_zone = current_internal_zone + 1
613 end if
614
615 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
616 re2v2_data_bc(i)%type, &
617 neko_w_bc_label, total_labeled_zone_offset, &
618 labeled_zone_offsets(neko_w_bc_label) .eq. 0)
619
620 labeled_zone_offsets(neko_w_bc_label) = 1
621
622 case ('v', 'V')
623 if (neko_v_bc_label .eq. -1) then
624 neko_v_bc_label = current_internal_zone
625 current_internal_zone = current_internal_zone + 1
626 end if
627
628 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
629 re2v2_data_bc(i)%type, &
630 neko_v_bc_label, total_labeled_zone_offset, &
631 labeled_zone_offsets(neko_v_bc_label) .eq. 0)
632
633 labeled_zone_offsets(neko_v_bc_label) = 1
634 case ('O', 'o')
635 if (neko_o_bc_label .eq. -1) then
636 neko_o_bc_label = current_internal_zone
637 current_internal_zone = current_internal_zone + 1
638 end if
639
640 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
641 re2v2_data_bc(i)%type, &
642 neko_o_bc_label, total_labeled_zone_offset, &
643 labeled_zone_offsets(neko_o_bc_label) .eq. 0)
644
645 labeled_zone_offsets(neko_o_bc_label) = 1
646 case ('SYM','sym')
647 if (neko_sym_bc_label .eq. -1) then
648 neko_sym_bc_label = current_internal_zone
649 current_internal_zone = current_internal_zone + 1
650 end if
651
652 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
653 re2v2_data_bc(i)%type, &
654 neko_sym_bc_label, total_labeled_zone_offset, &
655 labeled_zone_offsets(neko_sym_bc_label) .eq. 0)
656
657 labeled_zone_offsets(neko_sym_bc_label) = 1
658 case ('ON', 'on')
659 if (neko_on_bc_label .eq. -1) then
660 neko_on_bc_label = current_internal_zone
661 current_internal_zone = current_internal_zone + 1
662 end if
663
664 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
665 re2v2_data_bc(i)%type, &
666 neko_on_bc_label, total_labeled_zone_offset, &
667 labeled_zone_offsets(neko_on_bc_label) .eq. 0)
668
669 labeled_zone_offsets(neko_on_bc_label) = 1
670 case ('s', 'sl', 'sh', 'shl', 'S', 'SL', 'SH', 'SHL')
671 if (neko_shl_bc_label .eq. -1) then
672 neko_shl_bc_label = current_internal_zone
673 current_internal_zone = current_internal_zone + 1
674 end if
675
676 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
677 re2v2_data_bc(i)%type, &
678 neko_shl_bc_label, total_labeled_zone_offset, &
679 labeled_zone_offsets(neko_shl_bc_label) .eq. 0)
680
681 labeled_zone_offsets(neko_shl_bc_label) = 1
682 case ('P')
683 periodic = .true.
684 p_el_idx = int(re2v2_data_bc(i)%bc_data(1))
685 p_facet = facet_map(int(re2v2_data_bc(i)%bc_data(2)))
686 call msh%get_facet_ids(sym_facet, el_idx, pids)
687 call msh%mark_periodic_facet(sym_facet, el_idx, &
688 p_facet, p_el_idx, pids)
689 case default
690 write (*,*) trim(re2v2_data_bc(i)%type), 'bc type not supported yet'
691 write (*,*) re2v2_data_bc(i)%bc_data
692 end select
693 end do
694
695 !
696 ! Fix periodic condition for shared nodes
697 !
698 if (periodic) then
699 do j = 1, 3
700 do i = 1, nbcs
701 el_idx = re2v2_data_bc(i)%elem - dist%start_idx()
702 sym_facet = facet_map(int(re2v2_data_bc(i)%face))
703 select case(trim(re2v2_data_bc(i)%type))
704 case ('P')
705 p_el_idx = int(re2v2_data_bc(i)%bc_data(1))
706 p_facet = facet_map(int(re2v2_data_bc(i)%bc_data(2)))
707 call msh%create_periodic_ids(sym_facet, el_idx, &
708 p_facet, p_el_idx)
709 end select
710 end do
711 end do
712 end if
713 deallocate(re2v2_data_bc)
714
715 else ! V! format
716
717 !
718 ! Search for user labeled zones
719 !
720 do i = 1, nbcs
721
722 el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
723 sym_facet = facet_map(re2v1_data_bc(i)%face)
724
725 select case(trim(re2v1_data_bc(i)%type))
726 case ('MSH', 'msh', 'EXO','exo')
727
728 label = int(re2v1_data_bc(i)%bc_data(5))
729
730 if (label .lt. 1 .or. label .gt. neko_msh_max_zlbls) then
731 call neko_error('Invalid label id (valid range [1,...,20])')
732 end if
733
734 if (user_labeled_zones(label) .eq. 0) then
735 write (log_buf, "(A,I2,A,I3)") "Labeled zone ", label, &
736 " found."
737 call neko_log%message(log_buf)
738 user_labeled_zones(label) = 1
739 end if
740
741 ! Get the largest label possible in case we need to convert
742 ! re2 labels to labeled zones (see below).
743 total_labeled_zone_offset = max(total_labeled_zone_offset, label)
744
745 call msh%mark_labeled_facet(sym_facet, el_idx, label)
746
747 end select
748 end do
749
750 do i = 1, nbcs
751
752 el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
753 sym_facet = facet_map(re2v1_data_bc(i)%face)
754 select case(trim(re2v1_data_bc(i)%type))
755 case ('MSH', 'msh', 'EXO','exo')
756 !Do nothing, already handled
757 case ('W')
758 if (neko_w_bc_label .eq. -1) then
759 neko_w_bc_label = current_internal_zone
760 current_internal_zone = current_internal_zone + 1
761 end if
762
763 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
764 re2v1_data_bc(i)%type, &
765 neko_w_bc_label, total_labeled_zone_offset, &
766 labeled_zone_offsets(neko_w_bc_label) .eq. 0)
767
768 labeled_zone_offsets(neko_w_bc_label) = 1
769
770 case ('v', 'V')
771 if (neko_v_bc_label .eq. -1) then
772 neko_v_bc_label = current_internal_zone
773 current_internal_zone = current_internal_zone + 1
774 end if
775
776 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
777 re2v1_data_bc(i)%type, &
778 neko_v_bc_label, total_labeled_zone_offset, &
779 labeled_zone_offsets(neko_v_bc_label) .eq. 0)
780
781 labeled_zone_offsets(neko_v_bc_label) = 1
782 case ('O', 'o')
783 if (neko_o_bc_label .eq. -1) then
784 neko_o_bc_label = current_internal_zone
785 current_internal_zone = current_internal_zone + 1
786 end if
787
788 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
789 re2v1_data_bc(i)%type, &
790 neko_o_bc_label, total_labeled_zone_offset, &
791 labeled_zone_offsets(neko_o_bc_label) .eq. 0)
792
793 labeled_zone_offsets(neko_o_bc_label) = 1
794 case ('SYM','sym')
795 if (neko_sym_bc_label .eq. -1) then
796 neko_sym_bc_label = current_internal_zone
797 current_internal_zone = current_internal_zone + 1
798 end if
799
800 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
801 re2v1_data_bc(i)%type, &
802 neko_sym_bc_label, total_labeled_zone_offset, &
803 labeled_zone_offsets(neko_sym_bc_label) .eq. 0)
804
805 labeled_zone_offsets(neko_sym_bc_label) = 1
806 case ('ON', 'on')
807 if (neko_on_bc_label .eq. -1) then
808 neko_on_bc_label = current_internal_zone
809 current_internal_zone = current_internal_zone + 1
810 end if
811
812 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
813 re2v1_data_bc(i)%type, &
814 neko_on_bc_label, total_labeled_zone_offset, &
815 labeled_zone_offsets(neko_on_bc_label) .eq. 0)
816
817 labeled_zone_offsets(neko_on_bc_label) = 1
818 case ('s', 'sl', 'sh', 'shl', 'S', 'SL', 'SH', 'SHL')
819 if (neko_shl_bc_label .eq. -1) then
820 neko_shl_bc_label = current_internal_zone
821 current_internal_zone = current_internal_zone + 1
822 end if
823
824 call re2_file_mark_labeled_bc(msh, el_idx, sym_facet, &
825 re2v1_data_bc(i)%type, &
826 neko_shl_bc_label, total_labeled_zone_offset, &
827 labeled_zone_offsets(neko_shl_bc_label) .eq. 0)
828
829 labeled_zone_offsets(neko_shl_bc_label) = 1
830 case ('P')
831 periodic = .true.
832 p_el_idx = int(re2v1_data_bc(i)%bc_data(1))
833 p_facet = facet_map(int(re2v1_data_bc(i)%bc_data(2)))
834 call msh%get_facet_ids(sym_facet, el_idx, pids)
835 call msh%mark_periodic_facet(sym_facet, el_idx, &
836 p_facet, p_el_idx, pids)
837 case default
838 write (*,*) re2v1_data_bc(i)%type, 'bc type not supported yet'
839 write (*,*) re2v1_data_bc(i)%bc_data
840 end select
841 end do
842
843 !
844 ! Fix periodic condition for shared nodes
845 !
846 if (periodic) then
847 do j = 1, 3
848 do i = 1, nbcs
849 el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
850 sym_facet = facet_map(re2v1_data_bc(i)%face)
851 select case(trim(re2v1_data_bc(i)%type))
852 case ('P')
853 p_el_idx = int(re2v1_data_bc(i)%bc_data(1))
854 p_facet = facet_map(int(re2v1_data_bc(i)%bc_data(2)))
855 call msh%create_periodic_ids(sym_facet, el_idx, &
856 p_facet, p_el_idx)
857 end select
858 end do
859 end do
860 end if
861
862 deallocate(re2v1_data_bc)
863 end if
864
865 end subroutine re2_file_read_bcs
866
877 subroutine re2_file_mark_labeled_bc(msh, el_idx, facet, type, label, offset, print_info)
878 type(mesh_t), intent(inout) :: msh
879 integer, intent(in) :: el_idx
880 integer, intent(in) :: facet
881 character(len=*), intent(in) :: type
882 integer, intent(in) :: label
883 integer, intent(in) :: offset
884 logical, intent(in) :: print_info
885
886 integer :: mark_label
887 character(len=LOG_SIZE) :: log_buf
888
889 mark_label = offset + label
890
891 if (mark_label .lt. 1 .or. mark_label .gt. neko_msh_max_zlbls) then
892 call neko_error("You have reached the maximum amount of allowed labeled&
893 & zones (max allowed: 20). This happened when converting re2 internal labels&
894 & like e.g. 'w', 'V' or 'o' to labeled zones. Please reduce the number of&
895 & labeled zones that you have defined or make sure that they are labeled&
896 & from [1,...,20].")
897 end if
898
899 if (print_info) then
900 write (log_buf, "(A3,A,I2)") trim(type), " => Labeled index ", mark_label
901 call neko_log%message(log_buf)
902 end if
903 call msh%mark_labeled_facet(facet, el_idx, mark_label)
904
905 end subroutine re2_file_mark_labeled_bc
906
907 subroutine re2_file_add_point(htp, p, idx)
908 type(htable_pt_t), intent(inout) :: htp
909 type(point_t), intent(inout) :: p
910 integer, intent(inout) :: idx
911 integer :: tmp
912
913 if (htp%get(p, tmp) .gt. 0) then
914 idx = idx + 1
915 call htp%set(p, idx)
916 call p%set_id(idx)
917 else
918 call p%set_id(tmp)
919 end if
920
921 end subroutine re2_file_add_point
922
923end module re2_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 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
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
MPI derived types.
Definition mpi_types.f90:34
integer, public mpi_double_precision_size
Size of MPI type double precision.
Definition mpi_types.f90:63
type(mpi_datatype), public mpi_re2v2_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition mpi_types.f90:54
type(mpi_datatype), public mpi_re2v2_data_xy
MPI derived type for 2D NEKTON re2 data.
Definition mpi_types.f90:55
integer, public mpi_character_size
Size of MPI type character.
Definition mpi_types.f90:64
integer, public mpi_real_size
Size of MPI type real.
Definition mpi_types.f90:62
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:65
type(mpi_datatype), public mpi_re2v1_data_cv
MPI derived type for NEKTON re2 cv data.
Definition mpi_types.f90:51
type(mpi_datatype), public mpi_re2v2_data_cv
MPI derived type for NEKTON re2 cv data.
Definition mpi_types.f90:56
type(mpi_datatype), public mpi_re2v1_data_xy
MPI derived type for 2D NEKTON re2 data.
Definition mpi_types.f90:50
type(mpi_datatype), public mpi_re2v1_data_bc
MPI derived type for NEKTON re2 bc data.
Definition mpi_types.f90:52
type(mpi_datatype), public mpi_re2v2_data_bc
MPI derived type for NEKTON re2 bc data.
Definition mpi_types.f90:57
type(mpi_datatype), public mpi_re2v1_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition mpi_types.f90:49
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
Implements a point.
Definition point.f90:35
NEKTON mesh data in re2 format.
Definition re2_file.f90:35
subroutine re2_file_mark_labeled_bc(msh, el_idx, facet, type, label, offset, print_info)
Mark a labeled zone based on an offset.
Definition re2_file.f90:878
subroutine re2_file_read(this, data)
Load a binary NEKTON mesh from a re2 file.
Definition re2_file.f90:75
subroutine re2_file_write(this, data, t)
Definition re2_file.f90:223
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
subroutine re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
Definition re2_file.f90:526
subroutine re2_file_add_point(htp, p, idx)
Definition re2_file.f90:908
integer, public neko_v_bc_label
Definition re2_file.f90:58
integer, public neko_shl_bc_label
Definition re2_file.f90:62
subroutine re2_file_read_curve(msh, ncurve, dist, fh, mpi_offset, v2_format)
Definition re2_file.f90:424
integer, public neko_sym_bc_label
Definition re2_file.f90:60
subroutine re2_file_read_points(msh, ndim, nel, dist, fh, mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
Definition re2_file.f90:315
NEKTON re2 format.
Definition re2.f90:2
real(kind=sp), parameter re2_endian_test
NEKTION re2 endian test.
Definition re2.f90:10
integer, parameter re2_hdr_size
NEKTON re2 header size.
Definition re2.f90:7
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
NEKTON re2 bc data (version 1)
Definition re2.f90:39
NEKTON re2 curve data (version 1)
Definition re2.f90:31
NEKTON re2 element data (2d) (version 1)
Definition re2.f90:25
NEKTON re2 element data (3d) (version 1)
Definition re2.f90:18
NEKTON re2 bc data (version 2)
Definition re2.f90:73
NEKTON re2 curve data (version 2)
Definition re2.f90:65
NEKTON re2 element data (2d) (version 2)
Definition re2.f90:59
NEKTON re2 element data (3d) (version 2)
Definition re2.f90:52
Interface for NEKTON re2 files.
Definition re2_file.f90:65
#define max(a, b)
Definition tensor.cu:40