0.6.2-dev0
FORCES
FORtran lib for Comp. Env. Sys.
Loading...
Searching...
No Matches
mo_sentinel.f90
Go to the documentation of this file.
1!> \file mo_sentinel.f90
2!> \brief \copybrief mo_sentinel
3!> \details \copydetails mo_sentinel
4
5!> \brief Module to handle sentinels.
6!> \details This module provides standard sentinels for all used types as well as set, get and check routines.
7!! Sentinels are used to mark variables in order to check if they were accessed by an algorithm.
8!! A common use case is to check if variables were set by reading a fortran namelist.
9!> \version 0.1
10!> \authors Sebastian Mueller
11!> \date Sep 2022
12!> \copyright Copyright 2005-\today, the CHS Developers, Sabine Attinger: All rights reserved.
13!! FORCES is released under the LGPLv3+ license \license_note
15
16 use mo_kind, only : i1, i2, i4, i8, sp, dp, spc, dpc
17 use, intrinsic :: ieee_arithmetic, only : ieee_value, ieee_quiet_nan, ieee_is_nan
18
19 implicit none
20 private
21
22 ! sentinel values
23 public :: sentinel_i1
24 public :: sentinel_i2
25 public :: sentinel_i4
26 public :: sentinel_i8
27 public :: sentinel_sp
28 public :: sentinel_dp
29 public :: sentinel_spc
30 public :: sentinel_dpc
31 public :: sentinel_char
32
33 ! sentinel setters
34 public set_sentinel
35 ! sentinel getters
36 public get_sentinel
37 ! sentinel checker
38 public check_sentinel
39
40 !> \brief Set variable to sentinel value.
41 !> \details This routine will set sentinal values to the given variable.
42 !> \param[in] "integer/real/complex/character :: value" value to be set to sentinel
43 interface set_sentinel
44 module procedure :: set_sentinel_i1
45 module procedure :: set_sentinel_i2
46 module procedure :: set_sentinel_i4
47 module procedure :: set_sentinel_i8
48 module procedure :: set_sentinel_sp
49 module procedure :: set_sentinel_dp
50 module procedure :: set_sentinel_spc
51 module procedure :: set_sentinel_dpc
52 module procedure :: set_sentinel_char
53 end interface
54
55 !> \brief Get sentinel values of the given kind.
56 !> \details This routine will return sentinal values of the same kind and shape as the given variable.
57 !> \param[in] "integer/real/complex/character :: value" value to determine sentinel
58 !> \retval "integer/real/complex/character :: sentinel_value" sentinal values of the same kind and shape
59 interface get_sentinel
60 module procedure :: get_sentinel_i1
61 module procedure :: get_sentinel_i2
62 module procedure :: get_sentinel_i4
63 module procedure :: get_sentinel_i8
64 module procedure :: get_sentinel_sp
65 module procedure :: get_sentinel_dp
66 module procedure :: get_sentinel_spc
67 module procedure :: get_sentinel_dpc
68 module procedure :: get_sentinel_char
69 end interface
70
71 !> \brief Check given variable to be equal to the sentinel value.
72 !> \details This routine will check the given variable to be equal to the specified sentinel value.
73 !> \param[in] "integer/real/complex/character :: value" value to be checked
74 !> \retval "logical :: is_sentinel" result
76 module procedure :: check_sentinel_i1
77 module procedure :: check_sentinel_i2
78 module procedure :: check_sentinel_i4
79 module procedure :: check_sentinel_i8
80 module procedure :: check_sentinel_sp
81 module procedure :: check_sentinel_dp
82 module procedure :: check_sentinel_spc
83 module procedure :: check_sentinel_dpc
84 module procedure :: check_sentinel_char
85 end interface
86
87 ! use maxval(empty) intrinsic to get lowest possible number for each integer kind
88 integer(i1), parameter :: empty_i1(0) = 0_i1
89 integer(i2), parameter :: empty_i2(0) = 0_i2
90 integer(i4), parameter :: empty_i4(0) = 0_i4
91 integer(i8), parameter :: empty_i8(0) = 0_i8
92 character, parameter :: sent_char = achar(0) ! NULL character
93
94contains
95
96 !> \brief Default sentinel value (Null character) for characters.
97 !> \return Default sentinel value (Null character) for characters.
98 pure character function sentinel_char()
99 sentinel_char = sent_char
100 end function sentinel_char
101
102 !> \brief Default sentinel value (-huge()-1) for i1.
103 !> \return Default sentinel value (-huge()-1) for i1.
104 pure integer(i1) function sentinel_i1()
105 sentinel_i1 = maxval(empty_i1)
106 end function sentinel_i1
107
108 !> \brief Default sentinel value (-huge()-1) for i2.
109 !> \return Default sentinel value (-huge()-1) for i2.
110 pure integer(i2) function sentinel_i2()
111 sentinel_i2 = maxval(empty_i2)
112 end function sentinel_i2
113
114 !> \brief Default sentinel value (-huge()-1) for i4.
115 !> \return Default sentinel value (-huge()-1) for i4.
116 pure integer(i4) function sentinel_i4()
117 sentinel_i4 = maxval(empty_i4)
118 end function sentinel_i4
119
120 !> \brief Default sentinel value (-huge()-1) for i8.
121 !> \return Default sentinel value (-huge()-1) for i8.
122 pure integer(i8) function sentinel_i8()
123 sentinel_i8 = maxval(empty_i8)
124 end function sentinel_i8
125
126 !> \brief Default sentinel value (NaN) for sp.
127 !> \return Default sentinel value (NaN) for sp.
128 pure real(sp) function sentinel_sp()
129 sentinel_sp = ieee_value(0._sp, ieee_quiet_nan)
130 end function sentinel_sp
131
132 !> \brief Default sentinel value (NaN) for dp.
133 !> \return Default sentinel value (NaN) for dp.
134 pure real(dp) function sentinel_dp()
135 sentinel_dp = ieee_value(0._dp, ieee_quiet_nan)
136 end function sentinel_dp
137
138 !> \brief Default sentinel value (NaN) for spc.
139 !> \return Default sentinel value (NaN) for spc.
140 pure complex(spc) function sentinel_spc()
141 real(sp) :: sentinel
142 sentinel = sentinel_sp()
143 sentinel_spc = cmplx(sentinel, sentinel, spc)
144 end function sentinel_spc
145
146 !> \brief Default sentinel value (NaN) for dpc.
147 !> \return Default sentinel value (NaN) for dpc.
148 pure complex(dpc) function sentinel_dpc()
149 real(dp) :: sentinel
150 sentinel = sentinel_dp()
151 sentinel_dpc = cmplx(sentinel, sentinel, dpc)
152 end function sentinel_dpc
153
154
155 !> \brief Set sentinel value for i1.
156 elemental subroutine set_sentinel_i1(value)
157 integer(i1), intent(inout) :: value !< value to be set to sentinel
158 value = sentinel_i1()
159 end subroutine set_sentinel_i1
160
161 !> \brief Set sentinel value for i2.
162 elemental subroutine set_sentinel_i2(value)
163 integer(i2), intent(inout) :: value !< value to be set to sentinel
164 value = sentinel_i2()
165 end subroutine set_sentinel_i2
166
167 !> \brief Set sentinel value for i4.
168 elemental subroutine set_sentinel_i4(value)
169 integer(i4), intent(inout) :: value !< value to be set to sentinel
170 value = sentinel_i4()
171 end subroutine set_sentinel_i4
172
173 !> \brief Set sentinel value for i8.
174 elemental subroutine set_sentinel_i8(value)
175 integer(i8), intent(inout) :: value !< value to be set to sentinel
176 value = sentinel_i8()
177 end subroutine set_sentinel_i8
178
179 !> \brief Set sentinel value for sp.
180 elemental subroutine set_sentinel_sp(value)
181 real(sp), intent(inout) :: value !< value to be set to sentinel
182 value = sentinel_sp()
183 end subroutine set_sentinel_sp
184
185 !> \brief Set sentinel value for dp.
186 elemental subroutine set_sentinel_dp(value)
187 real(dp), intent(inout) :: value !< value to be set to sentinel
188 value = sentinel_dp()
189 end subroutine set_sentinel_dp
190
191 !> \brief Set sentinel value for spc.
192 elemental subroutine set_sentinel_spc(value)
193 complex(spc), intent(inout) :: value !< value to be set to sentinel
194 value = sentinel_spc()
195 end subroutine set_sentinel_spc
196
197 !> \brief Set sentinel value for dpc.
198 elemental subroutine set_sentinel_dpc(value)
199 complex(dpc), intent(inout) :: value !< value to be set to sentinel
200 value = sentinel_dpc()
201 end subroutine set_sentinel_dpc
202
203 !> \brief Set sentinel value for char.
204 elemental subroutine set_sentinel_char(value)
205 character(len=*), intent(inout) :: value !< value to be set to sentinel
206 value = sentinel_char()
207 end subroutine set_sentinel_char
208
209
210 !> \brief Get sentinel value for char.
211 elemental character function get_sentinel_char(value)
212 character(*), intent(in) :: value !< value to determine sentinel
214 end function get_sentinel_char
215
216 !> \brief Get sentinel value for i1.
217 elemental integer(i1) function get_sentinel_i1(value)
218 integer(i1), intent(in) :: value !< value to determine sentinel
220 end function get_sentinel_i1
221
222 !> \brief Get sentinel value for i2.
223 elemental integer(i2) function get_sentinel_i2(value)
224 integer(i2), intent(in) :: value !< value to determine sentinel
226 end function get_sentinel_i2
227
228 !> \brief Get sentinel value for i4.
229 elemental integer(i4) function get_sentinel_i4(value)
230 integer(i4), intent(in) :: value !< value to determine sentinel
232 end function get_sentinel_i4
233
234 !> \brief Get sentinel value for i8.
235 elemental integer(i8) function get_sentinel_i8(value)
236 integer(i8), intent(in) :: value !< value to determine sentinel
238 end function get_sentinel_i8
239
240 !> \brief Get sentinel value for sp.
241 elemental real(sp) function get_sentinel_sp(value)
242 real(sp), intent(in) :: value !< value to determine sentinel
244 end function get_sentinel_sp
245
246 !> \brief Get sentinel value for dp.
247 elemental real(dp) function get_sentinel_dp(value)
248 real(dp), intent(in) :: value !< value to determine sentinel
250 end function get_sentinel_dp
251
252 !> \brief Get sentinel value for spc.
253 elemental complex(spc) function get_sentinel_spc(value)
254 complex(spc), intent(in) :: value !< value to determine sentinel
256 end function get_sentinel_spc
257
258 !> \brief Get sentinel value for dpc.
259 elemental complex(dpc) function get_sentinel_dpc(value)
260 complex(dpc), intent(in) :: value !< value to determine sentinel
262 end function get_sentinel_dpc
263
264
265 !> \brief Check for sentinel value for string.
266 elemental logical function check_sentinel_char(value)
267 character(*), intent(in) :: value !< value to determine sentinel
268 ! if trimmed string has not exactly one char, return false
269 if ( len_trim(value) /= 1 ) then
270 check_sentinel_char = .false.
271 else
272 check_sentinel_char = value(1:1) == sent_char(1:1)
273 end if
274 end function check_sentinel_char
275
276 !> \brief Check for sentinel value for i1.
277 elemental logical function check_sentinel_i1(value)
278 integer(i1), intent(in) :: value !< value to determine sentinel
279 integer(i1) :: sentinel
280 sentinel = sentinel_i1()
281 check_sentinel_i1 = sentinel == value
282 end function check_sentinel_i1
283
284 !> \brief Check for sentinel value for i2.
285 elemental logical function check_sentinel_i2(value)
286 integer(i2), intent(in) :: value !< value to determine sentinel
287 integer(i2) :: sentinel
288 sentinel = sentinel_i2()
289 check_sentinel_i2 = sentinel == value
290 end function check_sentinel_i2
291
292 !> \brief Check for sentinel value for i4.
293 elemental logical function check_sentinel_i4(value)
294 integer(i4), intent(in) :: value !< value to determine sentinel
295 integer(i4) :: sentinel
296 sentinel = sentinel_i4()
297 check_sentinel_i4 = sentinel == value
298 end function check_sentinel_i4
299
300 !> \brief Check for sentinel value for i8.
301 elemental logical function check_sentinel_i8(value)
302 integer(i8), intent(in) :: value !< value to determine sentinel
303 integer(i8) :: sentinel
304 sentinel = sentinel_i8()
305 check_sentinel_i8 = sentinel == value
306 end function check_sentinel_i8
307
308 !> \brief Check for sentinel value for sp.
309 elemental logical function check_sentinel_sp(value)
310 real(sp), intent(in) :: value !< value to determine sentinel
311 check_sentinel_sp = ieee_is_nan(value)
312 end function check_sentinel_sp
313
314 !> \brief Check for sentinel value for dp.
315 elemental logical function check_sentinel_dp(value)
316 real(dp), intent(in) :: value !< value to determine sentinel
317 check_sentinel_dp = ieee_is_nan(value)
318 end function check_sentinel_dp
319
320 !> \brief Check for sentinel value for spc.
321 elemental logical function check_sentinel_spc(value)
322 complex(spc), intent(in) :: value !< value to determine sentinel
323 check_sentinel_spc = ieee_is_nan(real(value)) .or. ieee_is_nan(aimag(value))
324 end function check_sentinel_spc
325
326 !> \brief Check for sentinel value for dpc.
327 elemental logical function check_sentinel_dpc(value)
328 complex(dpc), intent(in) :: value !< value to determine sentinel
329 check_sentinel_dpc = ieee_is_nan(real(value)) .or. ieee_is_nan(aimag(value))
330 end function check_sentinel_dpc
331
332end module mo_sentinel
Check given variable to be equal to the sentinel value.
Get sentinel values of the given kind.
Set variable to sentinel value.
Define number representations.
Definition mo_kind.F90:17
integer, parameter sp
Single Precision Real Kind.
Definition mo_kind.F90:44
integer, parameter i4
4 Byte Integer Kind
Definition mo_kind.F90:40
integer, parameter i8
8 Byte Integer Kind
Definition mo_kind.F90:42
integer, parameter i2
2 Byte Integer Kind
Definition mo_kind.F90:38
integer, parameter i1
1 Byte Integer Kind
Definition mo_kind.F90:36
integer, parameter dpc
Double Precision Complex Kind.
Definition mo_kind.F90:52
integer, parameter dp
Double Precision Real Kind.
Definition mo_kind.F90:46
integer, parameter spc
Single Precision Complex Kind.
Definition mo_kind.F90:50
Module to handle sentinels.
elemental logical function check_sentinel_dp(value)
Check for sentinel value for dp.
elemental integer(i4) function get_sentinel_i4(value)
Get sentinel value for i4.
pure integer(i1) function, public sentinel_i1()
Default sentinel value (-huge()-1) for i1.
pure character function, public sentinel_char()
Default sentinel value (Null character) for characters.
pure complex(dpc) function, public sentinel_dpc()
Default sentinel value (NaN) for dpc.
elemental integer(i8) function get_sentinel_i8(value)
Get sentinel value for i8.
elemental subroutine set_sentinel_i4(value)
Set sentinel value for i4.
elemental character function get_sentinel_char(value)
Get sentinel value for char.
elemental logical function check_sentinel_i2(value)
Check for sentinel value for i2.
elemental logical function check_sentinel_i8(value)
Check for sentinel value for i8.
pure complex(spc) function, public sentinel_spc()
Default sentinel value (NaN) for spc.
elemental subroutine set_sentinel_i1(value)
Set sentinel value for i1.
elemental real(dp) function get_sentinel_dp(value)
Get sentinel value for dp.
elemental complex(dpc) function get_sentinel_dpc(value)
Get sentinel value for dpc.
elemental integer(i1) function get_sentinel_i1(value)
Get sentinel value for i1.
pure integer(i2) function, public sentinel_i2()
Default sentinel value (-huge()-1) for i2.
pure real(sp) function, public sentinel_sp()
Default sentinel value (NaN) for sp.
elemental logical function check_sentinel_sp(value)
Check for sentinel value for sp.
elemental subroutine set_sentinel_char(value)
Set sentinel value for char.
pure integer(i8) function, public sentinel_i8()
Default sentinel value (-huge()-1) for i8.
elemental logical function check_sentinel_spc(value)
Check for sentinel value for spc.
elemental complex(spc) function get_sentinel_spc(value)
Get sentinel value for spc.
elemental integer(i2) function get_sentinel_i2(value)
Get sentinel value for i2.
elemental subroutine set_sentinel_i2(value)
Set sentinel value for i2.
elemental subroutine set_sentinel_i8(value)
Set sentinel value for i8.
pure integer(i4) function, public sentinel_i4()
Default sentinel value (-huge()-1) for i4.
elemental logical function check_sentinel_dpc(value)
Check for sentinel value for dpc.
elemental subroutine set_sentinel_dpc(value)
Set sentinel value for dpc.
elemental logical function check_sentinel_i4(value)
Check for sentinel value for i4.
elemental logical function check_sentinel_i1(value)
Check for sentinel value for i1.
elemental logical function check_sentinel_char(value)
Check for sentinel value for string.
elemental subroutine set_sentinel_dp(value)
Set sentinel value for dp.
elemental subroutine set_sentinel_spc(value)
Set sentinel value for spc.
elemental subroutine set_sentinel_sp(value)
Set sentinel value for sp.
pure real(dp) function, public sentinel_dp()
Default sentinel value (NaN) for dp.
elemental real(sp) function get_sentinel_sp(value)
Get sentinel value for sp.