Neko  0.9.99
A portable framework for high-order spectral element flow simulations
utils.f90
Go to the documentation of this file.
1 ! Copyright (c) 2019-2021, 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 utils
36  use, intrinsic :: iso_fortran_env, only: error_unit, output_unit
37  implicit none
38  private
39 
40  integer, parameter :: neko_fname_len = 1024
41 
42  interface neko_error
43  module procedure neko_error_plain, neko_error_msg
44  end interface neko_error
45 
50 
51 
52 contains
53 
55  pure function filename_suffix_pos(fname) result(suffix_pos)
56  character(len=*), intent(in) :: fname
57  integer :: suffix_pos
58  suffix_pos = scan(trim(fname), '.', back = .true.)
59  end function filename_suffix_pos
60 
62  pure function filename_tslash_pos(fname) result(tslash_pos)
63  character(len=*), intent(in) :: fname
64  integer :: tslash_pos
65  tslash_pos = scan(trim(fname), '/', back = .true.)
66  end function filename_tslash_pos
67 
69  subroutine filename_suffix(fname, suffix)
70  character(len=*) :: fname
71  character(len=*) :: suffix
72  suffix = trim(fname(filename_suffix_pos(fname) + 1:len_trim(fname)))
73  end subroutine filename_suffix
74 
76  subroutine filename_chsuffix(fname, new_fname, new_suffix)
77  character(len=*) :: fname
78  character(len=*) :: new_fname
79  character(len=*) :: new_suffix
80  integer :: suffix_pos
81 
82  suffix_pos = filename_suffix_pos(fname)
83  new_fname = trim(fname(1:suffix_pos))//new_suffix
84 
85  end subroutine filename_chsuffix
86 
93  function extract_fld_file_index(fld_filename, default_index) result(index)
94  character(len=*), intent(in) :: fld_filename
95  integer, intent(in) :: default_index
96 
97  character(len=80) :: suffix
98  integer :: index, fpos, i
99  logical :: valid
100 
101  call filename_suffix(fld_filename, suffix)
102 
103  valid = .true.
104 
105  ! This value will be modified when reading the file name extension
106  ! e.g. "field0.f00035" will set sample_idx = 35
107  index = default_index
108 
109  !
110  ! Try to extract the index of the field file from the suffix "fxxxxx"
111  !
112  fpos = scan(trim(suffix), 'f')
113  if (fpos .eq. 1) then
114  ! Make sure that the suffix only contains integers from 0 to 9
115  do i = 2, len(trim(suffix))
116  if (.not. (iachar(suffix(i:i)) >= iachar('0') &
117  .and. iachar(suffix(i:i)) <= iachar('9'))) then
118  valid = .false.
119  end if
120  end do
121  else
122  valid = .false.
123  end if
124 
125  ! Must be exactly 6 characters long, i.e. an 'f' with 5 integers after
126  if (len(trim(suffix)) .ne. 6) valid = .false.
127 
128  if (valid) read (suffix(2:), "(I5.5)") index
129 
130  end function extract_fld_file_index
131 
135  function split_string(string, delimiter) result(split_str)
136  character(len=*) :: string
137  character(len=*) :: delimiter
138  character(len=100), allocatable :: split_str(:)
139  integer :: length, i, i2, offset, j
140  i = 0
141  offset = 1
142  length = 1
143  if (len(trim(string)) .eq. 0) then
144  allocate(split_str(1))
145  split_str(1) = trim(string)
146  return
147  end if
148  do while (.true.)
149  i = scan(string(offset:), delimiter, back = .false.)
150  if (i .eq. 0) exit
151  length = length + 1
152  offset = offset + i
153  end do
154 
155  allocate(split_str(length))
156  i = 0
157  j = 1
158  offset = 1
159  do while (.true.)
160  i2 = scan(trim(string(offset:)), delimiter, back = .false.)
161  if (i2 .eq. 0) then
162  split_str(j) = trim(string(offset:))
163  exit
164  end if
165  split_str(j) = trim(string(offset:offset+i2-2))
166  offset = offset+i2
167  j = j + 1
168  end do
169  end function split_string
170 
173  pure function linear_index(i, j, k, l, lx, ly, lz) result(index)
174  integer, intent(in) :: i, j, k, l, lx, ly, lz
175  integer :: index
176 
177  index = (i + lx * ((j - 1) + ly * ((k - 1) + lz * ((l - 1)))))
178  end function linear_index
179 
182  pure function nonlinear_index(linear_index, lx, ly, lz) result(index)
183  integer, intent(in) :: linear_index, lx, ly, lz
184  integer :: index(4)
185  integer :: lin_idx
186  lin_idx = linear_index -1
187  index(4) = lin_idx/(lx*ly*lz)
188  index(3) = (lin_idx-(lx*ly*lz)*index(4))/(lx*ly)
189  index(2) = (lin_idx-(lx*ly*lz)*index(4)-(lx*ly)*index(3))/lx
190  index(1) = (lin_idx-(lx*ly*lz)*index(4)-(lx*ly)*index(3)-lx*index(2))
191  index(1) = index(1) + 1
192  index(2) = index(2) + 1
193  index(3) = index(3) + 1
194  index(4) = index(4) + 1
195 
196  end function nonlinear_index
197 
198  pure function index_is_on_facet(i, j, k, lx, ly, lz, facet) result(is_on)
199  integer, intent(in) :: i, j, k, lx, ly, lz, facet
200  logical :: is_on
201 
202  is_on = .false.
203  select case (facet)
204  case (1)
205  if (i .eq. 1) is_on = .true.
206  case (2)
207  if (i .eq. lx) is_on = .true.
208  case (3)
209  if (j .eq. 1) is_on = .true.
210  case (4)
211  if (j .eq. ly) is_on = .true.
212  case (5)
213  if (k .eq. 1) is_on = .true.
214  case (6)
215  if (k .eq. lz) is_on = .true.
216  end select
217 
218  end function index_is_on_facet
219 
222  subroutine neko_error_plain(error_code)
223  integer, optional :: error_code
224 
225  if (present(error_code)) then
226  write(error_unit, *) '*** ERROR ***', error_code
227  error stop
228  else
229  write(error_unit, *) '*** ERROR ***'
230  error stop
231  end if
232 
233  end subroutine neko_error_plain
234 
237  subroutine neko_error_msg(error_msg)
238  character(len=*) :: error_msg
239  write(error_unit, *) '*** ERROR: ', error_msg, ' ***'
240  error stop
241  end subroutine neko_error_msg
242 
244  subroutine neko_warning(warning_msg)
245  character(len=*) :: warning_msg
246  write(output_unit, *) '*** WARNING: ', warning_msg, ' ***'
247  end subroutine neko_warning
248 
254  function concat_string_array(array, sep, prepend) result(result)
255  character(len=*), intent(in) :: array(:)
256  character(len=*), intent(in) :: sep
257  logical, intent(in) :: prepend
258  character(:), allocatable :: result
259  integer :: i
260 
261  result = trim(array(1))
262  do i = 2, size(array)
263  result = result // sep // trim(array(i))
264  end do
265 
266  if (prepend .eqv. .true.) then
267  result = sep // result
268  end if
269 
270  end function concat_string_array
271 
272 end module utils
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Utilities.
Definition: utils.f90:35
integer function, public extract_fld_file_index(fld_filename, default_index)
Extracts the index of a field file. For example, "myfield.f00045" will return 45. If the suffix of th...
Definition: utils.f90:94
character(:) function, allocatable, public concat_string_array(array, sep, prepend)
Concatenate an array of strings into one string with array items separated by spaces.
Definition: utils.f90:255
character(len=100) function, dimension(:), allocatable, public split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
Definition: utils.f90:136
pure logical function, public index_is_on_facet(i, j, k, lx, ly, lz, facet)
Definition: utils.f90:199
subroutine neko_error_msg(error_msg)
Reports an error and stops execution.
Definition: utils.f90:238
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
Definition: utils.f90:174
integer, parameter, public neko_fname_len
Definition: utils.f90:40
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition: utils.f90:245
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition: utils.f90:77
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition: utils.f90:70
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Definition: utils.f90:63
subroutine neko_error_plain(error_code)
Reports an error and stops execution.
Definition: utils.f90:223
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition: utils.f90:56