Neko  0.8.1
A portable framework for high-order spectral element flow simulations
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 !
35 module re2_file
36  use generic_file
37  use num_types
38  use utils
39  use mesh
40  use point
41  use comm
42  use mpi_f08
43  use neko_mpi_types
44  use datadist
45  use re2
46  use map
47  use map_file
48  use htable
49  use logger
50  implicit none
51  private
52 
53 
55  type, public, extends(generic_file_t) :: re2_file_t
56  contains
57  procedure :: read => re2_file_read
58  procedure :: write => re2_file_write
59  end type re2_file_t
60 
61 contains
62 
64  subroutine re2_file_read(this, data)
65  class(re2_file_t) :: this
66  class(*), target, intent(inout) :: data
67  type(mesh_t), pointer :: msh
68  character(len=5) :: hdr_ver
69  character(len=54) :: hdr_str
70  character(len=80) :: hdr_full
71  integer :: nel, ndim, nelv, ierr, nBCre2
72  type(mpi_status) :: status
73  type(mpi_file) :: fh
74  integer (kind=MPI_OFFSET_KIND) :: mpi_offset
75  real(kind=sp) :: test
76  real(kind=dp) :: t2
77  integer :: ncurv, nbcs
78  type(linear_dist_t) :: dist
79  type(map_t) :: nm
80  type(map_file_t) :: map_file
81  character(len=1024) :: map_fname
82  logical :: read_map
83  integer :: re2_data_xy_size
84  integer :: re2_data_xyz_size
85  integer :: re2_data_cv_size
86  integer :: re2_data_bc_size
87  logical :: v2_format
88  character(len=LOG_SIZE) :: log_buf
89 
90  call this%check_exists()
91 
92  select type(data)
93  type is (mesh_t)
94  msh => data
95  class default
96  call neko_error('Invalid output data')
97  end select
98 
99  v2_format = .false.
100  open(unit=9,file=trim(this%fname), status='old', iostat=ierr)
101  call neko_log%message('Reading binary NEKTON file ' // this%fname)
102 
103  read(9,'(a80)') hdr_full
104  read(hdr_full, '(a5)') hdr_ver
105 
106  if (hdr_ver .eq. '#v004') then
107  read(hdr_full, '(a5,i16,i3,i16,i4,a36)') hdr_ver, nel, ndim, nelv, nbcre2, hdr_str
108  v2_format = .true.
109  else if (hdr_ver .eq. '#v002' .or. hdr_ver .eq. '#v003') then
110  read(hdr_full, '(a5,i9,i3,i9,a54)') hdr_ver, nel, ndim, nelv, hdr_str
111  v2_format = .true.
112  else if (hdr_ver .eq. '#v001') then
113  read(hdr_full, '(a5,i9,i3,i9,a54)') hdr_ver, nel, ndim, nelv, hdr_str
114  end if
115 
116  if (v2_format) then
117  call mpi_type_size(mpi_re2v2_data_xy, re2_data_xy_size, ierr)
118  call mpi_type_size(mpi_re2v2_data_xyz, re2_data_xyz_size, ierr)
119  call mpi_type_size(mpi_re2v2_data_cv, re2_data_cv_size, ierr)
120  call mpi_type_size(mpi_re2v2_data_bc, re2_data_bc_size, ierr)
121  else
122  call mpi_type_size(mpi_re2v1_data_xy, re2_data_xy_size, ierr)
123  call mpi_type_size(mpi_re2v1_data_xyz, re2_data_xyz_size, ierr)
124  call mpi_type_size(mpi_re2v1_data_cv, re2_data_cv_size, ierr)
125  call mpi_type_size(mpi_re2v1_data_bc, re2_data_bc_size, ierr)
126  end if
127 
128  write(log_buf,1) ndim, nelv
129 1 format('gdim = ', i1, ', nelements =', i7)
130  call neko_log%message(log_buf)
131  close(9)
132 
133  call filename_chsuffix(this%fname, map_fname,'map')
134 
135  inquire(file=map_fname, exist=read_map)
136  if (read_map) then
137  call map_init(nm, nelv, 2**ndim)
138  call map_file%init(map_fname)
139  call map_file%read(nm)
140  else
141  call neko_log%warning('No NEKTON map file found')
142  end if
143 
144  call mpi_file_open(neko_comm, trim(this%fname), &
145  mpi_mode_rdonly, mpi_info_null, fh, ierr)
146 
147  if (ierr .ne. 0) then
148  call neko_log%error("Can't open binary NEKTON file ")
149  end if
150  dist = linear_dist_t(nelv, pe_rank, pe_size, neko_comm)
151 
152 
153  call msh%init(ndim, dist)
154 
155  ! Set offset (header)
156  mpi_offset = re2_hdr_size * mpi_character_size
157 
158  call mpi_file_read_at_all(fh, mpi_offset, test, 1, mpi_real, status, ierr)
159  mpi_offset = mpi_offset + mpi_real_size
160 
161  if (abs(re2_endian_test - test) .gt. 1e-4) then
162  call neko_error('Invalid endian of re2 file, byte swap not implemented yet')
163  end if
164 
165  call re2_file_read_points(msh, ndim, nel, dist, fh, &
166  mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
167 
168 
169  ! Set offset to start of curved side data
171  if (ndim .eq. 2) then
172  mpi_offset = mpi_offset + int(dist%num_global(),i8) * int(re2_data_xy_size,i8)
173  else
174  mpi_offset = mpi_offset + int(dist%num_global(),i8) * int(re2_data_xyz_size,i8)
175  end if
176 
179  if (v2_format) then
180  call mpi_file_read_at_all(fh, mpi_offset, t2, 1, mpi_double_precision, status, ierr)
181  ncurv = int(t2)
182  mpi_offset = mpi_offset + mpi_double_precision_size
183  call re2_file_read_curve(msh, ncurv, dist, fh, mpi_offset, v2_format)
184  mpi_offset = mpi_offset + int(ncurv,i8) * int(re2_data_cv_size,i8)
185  call mpi_file_read_at_all(fh, mpi_offset, t2, 1, mpi_double_precision, status, ierr)
186  nbcs = int(t2)
187  mpi_offset = mpi_offset + mpi_double_precision_size
188 
189  call re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
190  else
191  call mpi_file_read_at_all(fh, mpi_offset, ncurv, 1, mpi_integer, status, ierr)
192  mpi_offset = mpi_offset + mpi_integer_size
193  call re2_file_read_curve(msh, ncurv, dist, fh, mpi_offset, v2_format)
194  mpi_offset = mpi_offset + int(ncurv,i8) * int(re2_data_cv_size,i8)
195  call mpi_file_read_at_all(fh, mpi_offset, nbcs, 1, mpi_integer, status, ierr)
196  mpi_offset = mpi_offset + mpi_integer_size
197 
198  call re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
199 
200  end if
201 
202  call mpi_file_close(fh, ierr)
203  call msh%finalize()
204 
205 
206  call neko_log%message('Done')
207 
208 
209  end subroutine re2_file_read
210 
211  subroutine re2_file_write(this, data, t)
212  class(re2_file_t), intent(inout) :: this
213  class(*), target, intent(in) :: data
214  real(kind=rp), intent(in), optional :: t
215  type(re2v1_xy_t), allocatable :: re2_data_xy(:)
216  type(re2v1_xyz_t), allocatable :: re2_data_xyz(:)
217  type(mesh_t), pointer :: msh
218  character(len=5), parameter :: RE2_HDR_VER = '#v001'
219  character(len=54), parameter :: RE2_HDR_STR = 'RE2 exported by NEKO'
220  integer :: i, j, ierr, nelgv
221  type(mpi_status) :: status
222  type(mpi_file) :: fh
223  integer (kind=MPI_OFFSET_KIND) :: mpi_offset
224  integer :: element_offset
225  integer :: re2_data_xy_size
226  integer :: re2_data_xyz_size
227  character(len=LOG_SIZE) :: log_buf
228 
229  select type(data)
230  type is (mesh_t)
231  msh => data
232  class default
233  call neko_error('Invalid output data')
234  end select
235 
236  call mpi_type_size(mpi_re2v1_data_xy, re2_data_xy_size, ierr)
237  call mpi_type_size(mpi_re2v1_data_xyz, re2_data_xyz_size, ierr)
238  call mpi_reduce(msh%nelv, nelgv, 1, &
239  mpi_integer, mpi_sum, 0, neko_comm, ierr)
240  element_offset = 0
241  call mpi_exscan(msh%nelv, element_offset, 1, &
242  mpi_integer, mpi_sum, neko_comm, ierr)
243 
244  call neko_log%message('Writing data as a binary NEKTON file ' // this%fname)
245 
246  if (pe_rank .eq. 0) then
247  open(unit=9,file=trim(this%fname), status='new', iostat=ierr)
248  write(9, '(a5,i9,i3,i9,a54)') re2_hdr_ver, nelgv, msh%gdim,&
249  nelgv, re2_hdr_str
250  close(9)
251  end if
252 
253  call mpi_barrier(neko_comm, ierr)
254  call mpi_file_open(neko_comm, trim(this%fname), &
255  mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
256  mpi_offset = re2_hdr_size * mpi_character_size
257 
258  call mpi_file_write_at(fh, mpi_offset, re2_endian_test, 1, &
259  mpi_real, status, ierr)
260  mpi_offset = mpi_offset + mpi_real_size
261 
262  if (msh%gdim .eq. 2) then
263  allocate(re2_data_xy(msh%nelv))
264  do i = 1, msh%nelv
265  re2_data_xy(i)%rgroup = 1.0 ! Not used
266  do j = 1, 4
267  re2_data_xy(i)%x(j) = real(msh%elements(i)%e%pts(j)%p%x(1))
268  re2_data_xy(i)%y(j) = real(msh%elements(i)%e%pts(j)%p%x(2))
269  end do
270  end do
271  mpi_offset = mpi_offset + element_offset * re2_data_xy_size
272  call mpi_file_write_at(fh, mpi_offset, &
273  re2_data_xy, msh%nelv, mpi_re2v1_data_xy, status, ierr)
274 
275  deallocate(re2_data_xy)
276 
277  else if (msh%gdim .eq. 3) then
278  allocate(re2_data_xyz(msh%nelv))
279  do i = 1, msh%nelv
280  re2_data_xyz(i)%rgroup = 1.0 ! Not used
281  do j = 1, 8
282  re2_data_xyz(i)%x(j) = real(msh%elements(i)%e%pts(j)%p%x(1))
283  re2_data_xyz(i)%y(j) = real(msh%elements(i)%e%pts(j)%p%x(2))
284  re2_data_xyz(i)%z(j) = real(msh%elements(i)%e%pts(j)%p%x(3))
285  end do
286  end do
287  mpi_offset = mpi_offset + element_offset * re2_data_xyz_size
288  call mpi_file_write_at(fh, mpi_offset, &
289  re2_data_xyz, msh%nelv, mpi_re2v1_data_xyz, status, ierr)
290 
291  deallocate(re2_data_xyz)
292  else
293  call neko_error("Invalid dimension of mesh")
294  end if
295 
296  call mpi_file_close(fh, ierr)
297  call neko_log%message('Done')
298 
300 
301  end subroutine re2_file_write
302 
303  subroutine re2_file_read_points(msh, ndim, nel, dist, fh, &
304  mpi_offset, re2_data_xy_size, re2_data_xyz_size, v2_format)
305  type(mesh_t), intent(inout) :: msh
306  integer (kind=MPI_OFFSET_KIND) :: mpi_offset
307  integer, intent(inout) :: ndim
308  integer, intent(inout) :: nel
309  type(mpi_file), intent(inout) :: fh
310  integer, intent(in) :: re2_data_xy_size
311  integer, intent(in) :: re2_data_xyz_size
312  logical, intent(in) :: v2_format
313  type(linear_dist_t) :: dist
314  integer :: element_offset
315  type(re2v1_xy_t), allocatable :: re2v1_data_xy(:)
316  type(re2v1_xyz_t), allocatable :: re2v1_data_xyz(:)
317  type(re2v2_xy_t), allocatable :: re2v2_data_xy(:)
318  type(re2v2_xyz_t), allocatable :: re2v2_data_xyz(:)
319  type(mpi_status) :: status
320  type(htable_pt_t) :: htp
321  type(point_t) :: p(8)
322  integer :: pt_idx, nelv
323  integer :: i, j, ierr
324 
325 
326  nelv = dist%num_local()
327  element_offset = dist%start_idx()
328 
329  call htp%init(2*nel, ndim)
330  pt_idx = 0
331  if (ndim .eq. 2) then
332  mpi_offset = mpi_offset + element_offset * re2_data_xy_size
333  if (.not. v2_format) then
334  allocate(re2v1_data_xy(nelv))
335  call mpi_file_read_at_all(fh, mpi_offset, &
336  re2v1_data_xy, nelv, mpi_re2v1_data_xy, status, ierr)
337  do i = 1, nelv
338  do j = 1, 4
339  p(j) = point_t(real(re2v1_data_xy(i)%x(j),dp), &
340  real(re2v1_data_xy(i)%y(j),dp), 0.0d0)
341  call re2_file_add_point(htp, p(j), pt_idx)
342  end do
343  if (nelv > 10) then
344  if(mod(i,nelv/10) .eq. 0) write(*,*) i, 'elements read'
345  end if
346  ! swap vertices to keep symmetric vertex numbering in neko
347  call msh%add_element(i, p(1), p(2), p(4), p(3))
348  end do
349  deallocate(re2v1_data_xy)
350  else
351  allocate(re2v2_data_xy(nelv))
352  call mpi_file_read_at_all(fh, mpi_offset, &
353  re2v2_data_xy, nelv, mpi_re2v2_data_xy, status, ierr)
354  do i = 1, nelv
355  do j = 1, 4
356  p(j) = point_t(re2v2_data_xy(i)%x(j), &
357  re2v2_data_xy(i)%y(j), 0.0d0)
358  call re2_file_add_point(htp, p(j), pt_idx)
359  end do
360  if (nelv > 10) then
361  if(mod(i,nelv/10) .eq. 0) write(*,*) i, 'elements read'
362  end if
363  ! swap vertices to keep symmetric vertex numbering in neko
364  call msh%add_element(i, p(1), p(2), p(4), p(3))
365  end do
366  deallocate(re2v2_data_xy)
367  end if
368  else if (ndim .eq. 3) then
369  mpi_offset = mpi_offset + element_offset * re2_data_xyz_size
370  if (.not. v2_format) then
371  allocate(re2v1_data_xyz(nelv))
372  call mpi_file_read_at_all(fh, mpi_offset, &
373  re2v1_data_xyz, nelv, mpi_re2v1_data_xyz, status, ierr)
374  do i = 1, nelv
375  do j = 1, 8
376  p(j) = point_t(real(re2v1_data_xyz(i)%x(j),dp), &
377  real(re2v1_data_xyz(i)%y(j),dp),&
378  real(re2v1_data_xyz(i)%z(j),dp))
379  call re2_file_add_point(htp, p(j), pt_idx)
380  end do
381  if (nelv > 100) then
382  if(mod(i,nelv/100) .eq. 0) write(*,*) i, 'elements read'
383  end if
384  ! swap vertices to keep symmetric vertex numbering in neko
385  call msh%add_element(i, &
386  p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
387  end do
388  deallocate(re2v1_data_xyz)
389  else
390  allocate(re2v2_data_xyz(nelv))
391  call mpi_file_read_at_all(fh, mpi_offset, &
392  re2v2_data_xyz, nelv, mpi_re2v2_data_xyz, status, ierr)
393  do i = 1, nelv
394  do j = 1, 8
395  p(j) = point_t(re2v2_data_xyz(i)%x(j), &
396  re2v2_data_xyz(i)%y(j),&
397  re2v2_data_xyz(i)%z(j))
398  call re2_file_add_point(htp, p(j), pt_idx)
399  end do
400  if (nelv > 100) then
401  if(mod(i,nelv/100) .eq. 0) write(*,*) i, 'elements read'
402  end if
403  ! swap vertices to keep symmetric vertex numbering in neko
404  call msh%add_element(i, &
405  p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
406  end do
407  deallocate(re2v2_data_xyz)
408  end if
409  end if
410 
411  call htp%free()
412  end subroutine re2_file_read_points
413 
414  subroutine re2_file_read_curve(msh, ncurve, dist, fh, mpi_offset, v2_format)
415  type(mesh_t), intent(inout) :: msh
416  integer (kind=MPI_OFFSET_KIND) :: mpi_offset
417  integer, intent(inout) :: ncurve
418  type(linear_dist_t) :: dist
419  type(mpi_file), intent(inout) :: fh
420  logical, intent(in) :: v2_format
421  type(mpi_status) :: status
422  integer :: i, j, l, ierr, el_idx, id
423  type(re2v1_curve_t), allocatable :: re2v1_data_curve(:)
424  type(re2v2_curve_t), allocatable :: re2v2_data_curve(:)
425  real(kind=dp), allocatable :: curve_data(:,:,:)
426  integer, allocatable :: curve_type(:,:)
427  logical, allocatable :: curve_element(:)
428  character(len=1) :: chtemp
429  logical :: curve_skip = .false.
430 
431  allocate(curve_data(5,12,msh%nelv))
432  allocate(curve_element(msh%nelv))
433  allocate(curve_type(12,msh%nelv))
434  do i = 1, msh%nelv
435  curve_element(i) = .false.
436  do j = 1, 12
437  curve_type(j,i) = 0
438  do l = 1, 5
439  curve_data(l,j,i) = 0d0
440  end do
441  end do
442  end do
443  write(*,*) 'reading curved data'
444  if (.not. v2_format) then
445  allocate(re2v1_data_curve(ncurve))
446  call mpi_file_read_at_all(fh, mpi_offset, re2v1_data_curve, ncurve, &
447  mpi_re2v1_data_cv, status, ierr)
448  else
449  allocate(re2v2_data_curve(ncurve))
450  call mpi_file_read_at_all(fh, mpi_offset, re2v2_data_curve, ncurve, &
451  mpi_re2v2_data_cv, status, ierr)
452  end if
453  !This can probably be made nicer...
454  do i = 1, ncurve
455  if(v2_format) then
456  el_idx = re2v2_data_curve(i)%elem - dist%start_idx()
457  id = re2v2_data_curve(i)%zone
458  chtemp = re2v2_data_curve(i)%type
459  do j = 1, 5
460  curve_data(j,id, el_idx) = re2v2_data_curve(i)%point(j)
461  enddo
462  else
463  el_idx = re2v1_data_curve(i)%elem - dist%start_idx()
464  id = re2v1_data_curve(i)%zone
465  chtemp = re2v1_data_curve(i)%type
466  do j = 1, 5
467  curve_data(j,id, el_idx) = real(re2v1_data_curve(i)%point(j),dp)
468  enddo
469  end if
470 
471  curve_element(el_idx) = .true.
472  !This might need to be extended
473  select case(trim(chtemp))
474  case ('s')
475  curve_type(id,el_idx) = 1
476  call neko_log%warning('curve type s not supported, treating mesh as non-curved')
477  curve_skip = .true.
478  exit
479  case ('e')
480  curve_type(id,el_idx) = 2
481  call neko_log%warning('curve type e not supported, treating mesh as non-curved')
482  curve_skip = .true.
483  exit
484  case ('C')
485  curve_type(id,el_idx) = 3
486  case ('m')
487  curve_type(id,el_idx) = 4
488  case default
489  write(*,*) chtemp, 'curve type not supported yet, treating mesh as non-curved',id, el_idx
490 
491  curve_skip = .true.
492  end select
493  end do
494 
495  if( v2_format) then
496  deallocate(re2v2_data_curve)
497  else
498  deallocate(re2v1_data_curve)
499  end if
500  if (.not. curve_skip) then
501  do el_idx = 1, msh%nelv
502  if (curve_element(el_idx)) then
503  call msh%mark_curve_element(el_idx, &
504  curve_data(1,1,el_idx), curve_type(1,el_idx))
505  end if
506  end do
507  end if
508 
509  deallocate(curve_data)
510  deallocate(curve_element)
511  deallocate(curve_type)
512  end subroutine re2_file_read_curve
513 
514 
515  subroutine re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
516  type(mesh_t), intent(inout) :: msh
517  integer (kind=MPI_OFFSET_KIND) :: mpi_offset
518  integer, intent(inout) :: nbcs
519  type(linear_dist_t) :: dist
520  type(mpi_file), intent(inout) :: fh
521  logical, intent(in) :: v2_format
522  type(mpi_status) :: status
523  integer :: pids(4)
524  integer :: sym_facet, label
525  integer :: p_el_idx, p_facet
526  integer :: i, j, ierr, el_idx
527  integer, parameter, dimension(6) :: facet_map = (/3, 2, 4, 1, 5, 6/)
528  logical :: periodic
529  type(re2v1_bc_t), allocatable :: re2v1_data_bc(:)
530  type(re2v2_bc_t), allocatable :: re2v2_data_bc(:)
531 
532  if (.not. v2_format) then
533  allocate(re2v1_data_bc(nbcs))
534  call mpi_file_read_at_all(fh, mpi_offset, re2v1_data_bc, nbcs, &
535  mpi_re2v1_data_bc, status, ierr)
536  else
537  allocate(re2v2_data_bc(nbcs))
538  call mpi_file_read_at_all(fh, mpi_offset, re2v2_data_bc, nbcs, &
539  mpi_re2v2_data_bc, status, ierr)
540  end if
541 
542  periodic = .false.
543 
545  if (v2_format) then ! V2 format
546  do i = 1, nbcs
547  el_idx = int(re2v2_data_bc(i)%elem) - dist%start_idx()
548  sym_facet = facet_map(int(re2v2_data_bc(i)%face))
549 
550  select case(trim(re2v2_data_bc(i)%type))
551  case ('W')
552  call msh%mark_wall_facet(sym_facet, el_idx)
553  case ('v', 'V')
554  call msh%mark_inlet_facet(sym_facet, el_idx)
555  case ('O', 'o')
556  call msh%mark_outlet_facet(sym_facet, el_idx)
557  case ('ON', 'on')
558  call msh%mark_outlet_normal_facet(sym_facet, el_idx)
559  case ('SYM')
560  call msh%mark_sympln_facet(sym_facet, el_idx)
561  case ('P')
562  periodic = .true.
563  p_el_idx = int(re2v2_data_bc(i)%bc_data(1))
564  p_facet = facet_map(int(re2v2_data_bc(i)%bc_data(2)))
565  call msh%get_facet_ids(sym_facet, el_idx, pids)
566  call msh%mark_periodic_facet(sym_facet, el_idx, &
567  p_facet, p_el_idx, pids)
568  case ('MSH', 'msh', 'EXO', 'exo')
569  label = int(re2v2_data_bc(i)%bc_data(5))
570  if (label .lt. 1 .or. label .gt. neko_msh_max_zlbls) then
571  call neko_error('Invalid label id (valid range [1,...,20])')
572  end if
573  call msh%mark_labeled_facet(sym_facet, el_idx, label)
574  case default
575  write (*,*) re2v2_data_bc(i)%type, 'bc type not supported yet'
576  write (*,*) re2v2_data_bc(i)%bc_data
577  end select
578  end do
579 
580  !
581  ! Fix periodic condition for shared nodes
582  !
583  if (periodic) then
584  do j = 1, 3
585  do i = 1, nbcs
586  el_idx = re2v2_data_bc(i)%elem - dist%start_idx()
587  sym_facet = facet_map(int(re2v2_data_bc(i)%face))
588  select case(trim(re2v2_data_bc(i)%type))
589  case ('P')
590  p_el_idx = int(re2v2_data_bc(i)%bc_data(1))
591  p_facet = facet_map(int(re2v2_data_bc(i)%bc_data(2)))
592  call msh%create_periodic_ids(sym_facet, el_idx, &
593  p_facet, p_el_idx)
594  end select
595  end do
596  end do
597  end if
598  deallocate(re2v2_data_bc)
599 
600  else ! V! format
601  do i = 1, nbcs
602  el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
603  sym_facet = facet_map(re2v1_data_bc(i)%face)
604  select case(trim(re2v1_data_bc(i)%type))
605  case ('W')
606  call msh%mark_wall_facet(sym_facet, el_idx)
607  case ('v', 'V')
608  call msh%mark_inlet_facet(sym_facet, el_idx)
609  case ('O', 'o')
610  call msh%mark_outlet_facet(sym_facet, el_idx)
611  case ('SYM')
612  call msh%mark_sympln_facet(sym_facet, el_idx)
613  case ('P')
614  periodic = .true.
615  p_el_idx = int(re2v1_data_bc(i)%bc_data(1))
616  p_facet = facet_map(int(re2v1_data_bc(i)%bc_data(2)))
617  call msh%get_facet_ids(sym_facet, el_idx, pids)
618  call msh%mark_periodic_facet(sym_facet, el_idx, &
619  p_facet, p_el_idx, pids)
620  case ('MSH', 'msh')
621  label = int(re2v1_data_bc(i)%bc_data(5))
622  call msh%mark_labeled_facet(sym_facet, el_idx, label)
623  case default
624  write (*,*) re2v1_data_bc(i)%type, 'bc type not supported yet'
625  write (*,*) re2v1_data_bc(i)%bc_data
626 
627  end select
628  end do
629 
630  !
631  ! Fix periodic condition for shared nodes
632  !
633  if (periodic) then
634  do j = 1, 3
635  do i = 1, nbcs
636  el_idx = re2v1_data_bc(i)%elem - dist%start_idx()
637  sym_facet = facet_map(re2v1_data_bc(i)%face)
638  select case(trim(re2v1_data_bc(i)%type))
639  case ('P')
640  p_el_idx = int(re2v1_data_bc(i)%bc_data(1))
641  p_facet = facet_map(int(re2v1_data_bc(i)%bc_data(2)))
642  call msh%create_periodic_ids(sym_facet, el_idx, &
643  p_facet, p_el_idx)
644  end select
645  end do
646  end do
647  end if
648 
649  deallocate(re2v1_data_bc)
650  end if
651 
652  end subroutine re2_file_read_bcs
653 
654 
655  subroutine re2_file_add_point(htp, p, idx)
656  type(htable_pt_t), intent(inout) :: htp
657  type(point_t), intent(inout) :: p
658  integer, intent(inout) :: idx
659  integer :: tmp
660 
661  if (htp%get(p, tmp) .gt. 0) then
662  idx = idx + 1
663  call htp%set(p, idx)
664  call p%set_id(idx)
665  else
666  call p%set_id(tmp)
667  end if
668 
669  end subroutine re2_file_add_point
670 
671 end module re2_file
double real
Definition: device_config.h:12
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:26
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
integer pe_size
MPI size of communicator.
Definition: comm.F90:29
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:61
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:55
MPI derived types.
Definition: mpi_types.f90:34
integer, public mpi_double_precision_size
Size of MPI type double precision.
Definition: mpi_types.f90:61
type(mpi_datatype), public mpi_re2v2_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition: mpi_types.f90:52
type(mpi_datatype), public mpi_re2v2_data_xy
MPI derived type for 2D NEKTON re2 data.
Definition: mpi_types.f90:53
integer, public mpi_character_size
Size of MPI type character.
Definition: mpi_types.f90:62
integer, public mpi_real_size
Size of MPI type real.
Definition: mpi_types.f90:60
integer, public mpi_integer_size
Size of MPI type integer.
Definition: mpi_types.f90:63
type(mpi_datatype), public mpi_re2v1_data_cv
MPI derived type for NEKTON re2 cv data.
Definition: mpi_types.f90:49
type(mpi_datatype), public mpi_re2v2_data_cv
MPI derived type for NEKTON re2 cv data.
Definition: mpi_types.f90:54
type(mpi_datatype), public mpi_re2v1_data_xy
MPI derived type for 2D NEKTON re2 data.
Definition: mpi_types.f90:48
type(mpi_datatype), public mpi_re2v1_data_bc
MPI derived type for NEKTON re2 bc data.
Definition: mpi_types.f90:50
type(mpi_datatype), public mpi_re2v2_data_bc
MPI derived type for NEKTON re2 bc data.
Definition: mpi_types.f90:55
type(mpi_datatype), public mpi_re2v1_data_xyz
MPI derived type for 3D NEKTON re2 data.
Definition: mpi_types.f90:47
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_read(this, data)
Load a binary NEKTON mesh from a re2 file.
Definition: re2_file.f90:65
subroutine re2_file_write(this, data, t)
Definition: re2_file.f90:212
subroutine re2_file_read_bcs(msh, nbcs, dist, fh, mpi_offset, v2_format)
Definition: re2_file.f90:516
subroutine re2_file_add_point(htp, p, idx)
Definition: re2_file.f90:656
subroutine re2_file_read_curve(msh, ncurve, dist, fh, mpi_offset, v2_format)
Definition: re2_file.f90:415
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:305
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 filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition: utils.f90:69
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:44
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:55