Neko  0.8.1
A portable framework for high-order spectral element flow simulations
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 !
35 module rea_file
36  use generic_file
37  use num_types
38  use utils
39  use mesh
40  use point
41  use map
42  use rea
43  use re2_file
44  use map_file
45  use comm
46  use datadist
47  use htable
48  use logger
49  implicit none
50  private
51 
53  type, public, extends(generic_file_t) :: rea_file_t
54  contains
55  procedure :: read => rea_file_read
56  procedure :: write => rea_file_write
57  end type rea_file_t
58 
59 contains
60 
62  subroutine rea_file_read(this, data)
63  class(rea_file_t) :: this
64  class(*), target, intent(inout) :: data
65  type(mesh_t), pointer :: msh
66  real(kind=dp), pointer :: params(:)
67  character(len=3), pointer :: cbc(:,:)
68  integer, allocatable :: curve_type(:,:)
69  logical, allocatable :: curve_element(:)
70  character(len=1) :: chtemp
71  integer :: ndim, nparam, nskip, nlogic, ncurve
72  integer :: nelgs, nelgv, i, j, ierr, l
73  integer :: el_idx, pt_idx
74  logical :: read_param, read_bcs, read_map
75  real(kind=dp) :: xc(8), yc(8), zc(8), curve(5)
76  real(kind=dp), allocatable :: bc_data(:,:,:)
77  real(kind=dp), allocatable :: curve_data(:,:,:)
78  type(point_t) :: p(8)
79  type(re2_file_t) :: re2_file
80  type(map_file_t) :: map_file
81  character(len=1024) :: re2_fname, map_fname
82  integer :: start_el, end_el, nel, edge
83  type(linear_dist_t) :: dist
84  type(map_t) :: nm
85  type(htable_pt_t) :: htp
86  integer :: sym_facet, pids(4), p_el_idx, p_facet
87  integer :: off
88  integer, parameter, dimension(6) :: facet_map = (/3, 2, 4, 1, 5, 6/)
89  logical :: curve_skip = .false.
90  character(len=LOG_SIZE) :: log_buf
91 
92  call this%check_exists()
93 
94  select type(data)
95  type is (rea_t)
96  call rea_free(data)
97  msh => data%msh
98  params => data%params
99  cbc => data%cbc
100  read_param = .true.
101  read_bcs = .true.
102  type is (mesh_t)
103  msh => data
104  read_param = .false.
105  read_bcs = .false.
106  class default
107  call neko_error('Invalid output data')
108  end select
109 
110  if (read_param .and. read_bcs .and. pe_size .gt. 1) then
111  call neko_error('Reading NEKTON session data only implemented in serial')
112  end if
113 
114 
115  open(unit=9,file=trim(this%fname), status='old', iostat=ierr)
116  call neko_log%message('Reading NEKTON file ' // this%fname)
117 
118  read(9, *)
119  read(9, *)
120  read(9, *) ndim
121  read(9, *) nparam
122 
123  if (.not. read_param) then
124  ! Skip parameters
125  do i = 1, nparam
126  read(9, *)
127  end do
128  else
129  allocate(params(nparam))
130  do i = 1, nparam
131  read(9, *) params(i)
132  end do
133  end if
134 
135  ! Skip passive scalars
136  read(9, *) nskip
137  do i = 1, nskip
138  read(9, *)
139  end do
140 
141  ! Skip logic switches
142  read(9, *) nlogic
143  do i = 1, nlogic
144  read(9, *)
145  end do
146 
147  ! Read mesh info
148  read(9, *)
149  read(9, *)
150  read(9, *) nelgs,ndim, nelgv
151  if (nelgs .lt. 0) then
152  re2_fname = trim(this%fname(1:scan(trim(this%fname), &
153  '.', back=.true.)))//'re2'
154  call re2_file%init(re2_fname)
155  call re2_file%read(msh)
156  else
157  write(log_buf,1) ndim, nelgv
158 1 format('gdim = ', i1, ', nelements =', i7)
159  call neko_log%message(log_buf)
160 
161  call filename_chsuffix(this%fname, map_fname, 'map')
162  inquire(file=map_fname, exist=read_map)
163  if (read_map) then
164  call map_init(nm, nelgv, 2**ndim)
165  call map_file%init(map_fname)
166  call map_file%read(nm)
167  else
168  call neko_log%warning('No NEKTON map file found')
169  end if
170 
171  ! Use a load-balanced linear distribution
172  dist = linear_dist_t(nelgv, pe_rank, pe_size, neko_comm)
173  nel = dist%num_local()
174  start_el = dist%start_idx() + 1
175  end_el = dist%end_idx() + 1
176 
177  call msh%init(ndim, dist)
178 
179  call htp%init(2**ndim * nel, ndim)
180 
181  el_idx = 1
182  pt_idx = 0
183  do i = 1, nelgv
184  read(9, *)
185  if (ndim .eq. 2) then
186  read(9, *) (xc(j),j=1,4)
187  read(9, *) (yc(j),j=1,4)
188  if (i .ge. start_el .and. i .le. end_el) then
189  do j = 1, 4
190  p(j) = point_t(real(xc(j),dp), real(yc(j),dp),real(0d0,dp))
191  call rea_file_add_point(htp, p(j), pt_idx)
192  end do
193  ! swap vertices to keep symmetric vertex numbering in neko
194  call msh%add_element(el_idx, p(1), p(2), p(4), p(3))
195  end if
196  else if (ndim .eq. 3) then
197  read(9, *) (xc(j),j=1,4)
198  read(9, *) (yc(j),j=1,4)
199  read(9, *) (zc(j),j=1,4)
200  read(9, *) (xc(j),j=5,8)
201  read(9, *) (yc(j),j=5,8)
202  read(9, *) (zc(j),j=5,8)
203  if (i .ge. start_el .and. i .le. end_el) then
204  do j = 1, 8
205  p(j) = point_t(real(xc(j),dp), real(yc(j),dp), real(zc(j),dp))
206  call rea_file_add_point(htp, p(j), pt_idx)
207  end do
208  ! swap vertices to keep symmetric vertex numbering in neko
209  call msh%add_element(el_idx, &
210  p(1), p(2), p(4), p(3), p(5), p(6), p(8), p(7))
211  end if
212  end if
213  if (i .ge. start_el .and. i .le. end_el) then
214  el_idx = el_idx + 1
215  end if
216  end do
217 
218  call htp%free()
219 
220  read(9, *)
221  read(9, *) ncurve
222  allocate(curve_data(5,8,nelgv))
223  allocate(curve_element(nelgv))
224  allocate(curve_type(8,nelgv))
225  do i = 1, nelgv
226  curve_element(i) = .false.
227  do j = 1, 8
228  curve_type(j,i) = 0
229  do l = 1, 5
230  curve_data(l,j,i) = 0d0
231  end do
232  end do
233  end do
234  do i = 1, ncurve
235  read(9, *) edge, el_idx, (curve(j),j=1,5), chtemp
236  do j = 1, 5
237  curve_data(j,edge,el_idx) = curve(j)
238  end do
239  curve_element(el_idx) = .true.
240  select case(trim(chtemp))
241  case ('s')
242  curve_type(edge,el_idx) = 1
243  curve_skip = .true.
244  case ('e')
245  curve_type(edge,el_idx) = 2
246  curve_skip = .true.
247  case ('C')
248  curve_type(edge,el_idx) = 3
249  case ('m')
250  curve_type(edge,el_idx) = 4
251  end select
252  end do
253  if (curve_skip) then
254  call neko_log%warning('Curve type: s, e are not supported, treating mesh as non-curved.')
255  else
256  do el_idx = 1, nelgv
257  if (curve_element(el_idx)) then
258  call msh%mark_curve_element(el_idx, &
259  curve_data(1,1,el_idx), curve_type(1,el_idx))
260  end if
261  end do
262  end if
263  deallocate(curve_data)
264  deallocate(curve_element)
265  deallocate(curve_type)
266 
267  ! Read fluid boundary conditions
268  read(9,*)
269  read(9,*)
270  if (.not. read_bcs) then ! Mark zones in the mesh
271  allocate(cbc(6,nelgv))
272  allocate(bc_data(6,2*ndim,nelgv))
273  off = 0
274  !Fix for different horrible .rea periodic bc formats.
275  if (nelgv .lt. 1000) off = 1
276  do i = 1, nelgv
277  if (i .ge. start_el .and. i .le. end_el) then
278  el_idx = i - start_el + 1
279  do j = 1, 2*ndim
280  read(9, *) cbc(j, i), (bc_data(l,j,i),l=1,6)
281  sym_facet = facet_map(j)
282  select case(trim(cbc(j,i)))
283  case ('W')
284  call msh%mark_wall_facet(sym_facet, el_idx)
285  case ('v', 'V')
286  call msh%mark_inlet_facet(sym_facet, el_idx)
287  case ('O', 'o')
288  call msh%mark_outlet_facet(sym_facet, el_idx)
289  case ('ON', 'on')
290  call msh%mark_outlet_normal_facet(sym_facet, el_idx)
291  case ('SYM')
292  call msh%mark_sympln_facet(sym_facet, el_idx)
293  case ('P')
294  p_el_idx = int(bc_data(2+off,j,i))
295  p_facet = facet_map(int(bc_data(3+off,j,i)))
296  call msh%get_facet_ids(sym_facet, el_idx, pids)
297  call msh%mark_periodic_facet(sym_facet, el_idx, &
298  p_facet, p_el_idx, pids)
299  end select
300  end do
301  end if
302  end do
303  do i = 1, nelgv
304  if (i .ge. start_el .and. i .le. end_el) then
305  el_idx = i - start_el + 1
306  do j = 1, 2*ndim
307  sym_facet = facet_map(j)
308  select case(trim(cbc(j,i)))
309  case ('P')
310  p_el_idx = int(bc_data(2+off,j,i))
311  p_facet = facet_map(int(bc_data(3+off,j,i)))
312  call msh%create_periodic_ids(sym_facet, el_idx, &
313  p_facet, p_el_idx)
314  end select
315  end do
316  end if
317  end do
318  do i = 1, nelgv
319  if (i .ge. start_el .and. i .le. end_el) then
320  el_idx = i - start_el + 1
321  do j = 1, 2*ndim
322  sym_facet = facet_map(j)
323  select case(trim(cbc(j,i)))
324  case ('P')
325  p_el_idx = int(bc_data(2+off,j,i))
326  p_facet = facet_map(int(bc_data(3+off,j,i)))
327  call msh%create_periodic_ids(sym_facet, el_idx, &
328  p_facet, p_el_idx)
329  end select
330  end do
331  end if
332  end do
333  do i = 1, nelgv
334  if (i .ge. start_el .and. i .le. end_el) then
335  el_idx = i - start_el + 1
336  do j = 1, 2*ndim
337  sym_facet = facet_map(j)
338  select case(trim(cbc(j,i)))
339  case ('P')
340  p_el_idx = int(bc_data(2+off,j,i))
341  p_facet = facet_map(int(bc_data(3+off,j,i)))
342  call msh%create_periodic_ids(sym_facet, el_idx, &
343  p_facet, p_el_idx)
344  end select
345  end do
346  end if
347  end do
348  deallocate(cbc)
349  deallocate(bc_data)
350  else ! Store bcs in a NEKTON session structure
351  allocate(cbc(6,nelgv))
352  do i = 1, nelgv
353  do j = 1, 2*ndim
354  read(9,'(a1, a3)') chtemp, cbc(j, i)
355  end do
356  end do
357  end if
358 
359  call msh%finalize()
360 
361  call neko_log%message('Done')
362  close(9)
363  endif
364 
365  end subroutine rea_file_read
366 
367  subroutine rea_file_write(this, data, t)
368  class(rea_file_t), intent(inout) :: this
369  class(*), target, intent(in) :: data
370  real(kind=rp), intent(in), optional :: t
371  end subroutine rea_file_write
372 
373  subroutine rea_file_add_point(htp, p, idx)
374  type(htable_pt_t), intent(inout) :: htp
375  type(point_t), intent(inout) :: p
376  integer, intent(inout) :: idx
377  integer :: tmp
378 
379  if (htp%get(p, tmp) .gt. 0) then
380  idx = idx + 1
381  call htp%set(p, idx)
382  call p%set_id(idx)
383  else
384  call p%set_id(tmp)
385  end if
386 
387  end subroutine rea_file_add_point
388 
389 end module rea_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 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
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 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
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:63
subroutine rea_file_add_point(htp, p, idx)
Definition: rea_file.f90:374
subroutine rea_file_write(this, data, t)
Definition: rea_file.f90:368
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 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
Interface for NEKTON re2 files.
Definition: re2_file.f90:55
NEKTON session data struct.
Definition: rea.f90:13
Interface for NEKTON ascii files.
Definition: rea_file.f90:53