-
Notifications
You must be signed in to change notification settings - Fork 4
/
m_parameters.F90
268 lines (228 loc) · 7.79 KB
/
m_parameters.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
! File: m_parameters.F90
!
! Created: 6 August 2010
!
! Last modified: 6/8/2010
!
! Author: Pavel Sakov
! NERSC
!
! Purpose: Provide a simpl nml list-based parameter input into EnKF.
!
! Description: Provides code for reading parameters from a specified
! parameter file.
!
! Modifications: none
module m_parameters
#if defined(QMPI)
use qmpi
#else
use qmpi_fake
#endif
implicit none
integer, parameter, private :: STRLEN = 512
integer, parameter, private :: FID = 101
character(STRLEN), public :: PRMFNAME = "NONE"
integer, public :: ENSSIZE = 0
namelist /ensemble/ ENSSIZE
character(STRLEN), public :: METHODTAG = "NONE"
namelist /method/ METHODTAG
real, public :: LOCRAD = 0.0d0
character(STRLEN), public :: LOCFUNTAG = "GASPARI-COHN"
namelist /localisation/ LOCRAD, LOCFUNTAG
real, public :: INFL = 1.0d0
real, public :: RFACTOR1 = 1.0d0
real, public :: RFACTOR2 = 1.0d0
real, public :: KFACTOR = 2.0d0
namelist /moderation/ INFL, RFACTOR1, RFACTOR2, KFACTOR
character(STRLEN), public :: JMAPFNAME = "NONE"
character(STRLEN), public :: POINTFNAME = "NONE"
character(STRLEN), public :: MEANSSHFNAME = "NONE"
namelist /files/ JMAPFNAME, POINTFNAME, MEANSSHFNAME
integer, parameter, private :: NPRMESTMAX = 10
integer :: nprmest = 0
character(STRLEN), dimension(NPRMESTMAX), public :: PRMESTNAME
real, dimension(NPRMESTMAX), public :: PRMINFL
namelist /prmest/ PRMESTNAME, PRMINFL
public prm_read, prm_describe, prm_print, prm_getinfl, prm_prmestexists, ucase
contains
subroutine prm_read
integer :: ios, i
call getarg(1, PRMFNAME)
if (master) then
print *, 'EnKF: reading parameters from "', trim(PRMFNAME), '":'
end if
open(unit = FID, file = trim(PRMFNAME), form = "formatted",&
status = "old", iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: could not open "', trim(PRMFNAME), '", iostatus =', ios
stop
end if
end if
read(unit = FID, nml = method, iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: "', trim(PRMFNAME), '": could not read namelist "method"'
end if
stop
end if
rewind(FID)
read(unit = FID, nml = ensemble, iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: "', trim(PRMFNAME), '": could not read namelist "ensemble"'
end if
stop
end if
rewind(FID)
read(unit = FID, nml = localisation, iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: "', trim(PRMFNAME), '": could not read namelist "localisation"'
end if
stop
end if
rewind(FID)
read(unit = FID, nml = moderation, iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: "', trim(PRMFNAME), '": could not read namelist "moderation"'
end if
stop
end if
rewind(FID)
read(unit = FID, nml = files, iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: "', trim(PRMFNAME), '": could not read namelist "files"'
end if
stop
end if
rewind(FID)
do i = 1, NPRMESTMAX
PRMESTNAME(i) = ""
end do
read(unit = FID, nml = prmest, iostat = ios)
if (ios /= 0) then
if (master) then
print *, 'ERROR: "', trim(PRMFNAME), '": could not read namelist "prmest"'
end if
stop
end if
do i = 1, NPRMESTMAX
if (PRMESTNAME(i) == "") then
nprmest = i - 1
exit
end if
end do
rewind(FID)
close(FID)
call ucase(METHODTAG)
call ucase(LOCFUNTAG)
end subroutine prm_read
subroutine prm_describe
if (.not. master) then
return
end if
print '(a)', ' Example of EnKF parameter file:'
print *
print '(a)', '&method'
print '(a)', ' methodtag = "DEnKF"'
print '(a)', '/'
print '(a)', '&ensemble'
print '(a)', ' enssize = 0'
print '(a)', '/'
print '(a)', '&localisation'
print '(a)', ' locfuntag = "Gaspari-Cohn"'
print '(a)', ' locrad = 300.0'
print '(a)', '/'
print '(a)', '&moderation'
print '(a)', ' infl = 1.01 (<number>)'
print '(a)', ' rfactor1 = 1.0 (<number>)'
print '(a)', ' rfactor2 = 2.0 (<number>)'
print '(a)', ' kfactor = 2.0 (<number>)'
print '(a)', '/'
print '(a)', '&files'
print '(a)', ' jmapfname = "jmap.txt" (<file name>)'
print '(a)', ' pointfname = "point2nc.txt" (<file name>)'
print '(a)', ' meansshfname = "meanssh.uf" (<file name>)'
print *
print '(a)', 'Parameter options:'
print '(a)', ' method = "EnKF" | "DEnKF"*'
print '(a)', ' enssize = <number> (0* to use all available states)'
print '(a)', ' locfuntag = "Gaspari-Cohn"* | "Step" | "None"'
print '(a)', ' locrad = <support radius in km>'
print '(a)', ' infl = <multiple, for ensemble anomalies> (* 1.0)'
print '(a)', ' rfactor1 = <obs. error variance multiple> (* 1.0)'
print '(a)', ' rfactor2 = <additional multiple for updating ens. anomalies> (* 1.0)'
print '(a)', ' kfactor = <max. allowed increment in terms of ensemble spread> (* 2.0)'
print '(a)', ' jmapfname* = <file with j remapping> (* none)'
print '(a)', ' pointfname* = <file with point coordinates> (* none)'
print '(a)', ' meansshfname* = <file with mean SSH> (* none)'
end subroutine prm_describe
subroutine prm_print
integer :: i
if (.not. master) then
return
end if
print '(a)', ' EnKF parameters:'
print '(a)', ' method:'
print '(a, a, a)', ' methodtag = "', trim(METHODTAG), '"'
print '(a)', ' ensemble:'
print '(a, i0)', ' enssize = ', ENSSIZE
print '(a)', ' localisation:'
print '(a, f5.0)', ' locrad = ', LOCRAD
print '(a, a, a)', ' locfuntag = "', trim(LOCFUNTAG), '"'
print '(a)', ' moderation:'
print '(a, f5.3)', ' infl = ', INFL
print '(a, f3.1)', ' rfactor1 = ', RFACTOR1
print '(a, f3.1)', ' rfactor2 = ', RFACTOR2
print '(a, f3.1)', ' kfactor = ', KFACTOR
print '(a)', ' files:'
print '(a, a, a)', ' jmapfname = "', trim(JMAPFNAME), '"'
print '(a, a, a)', ' pointfname = "', trim(POINTFNAME), '"'
print '(a, a, a)', ' meansshfname = "', trim(MEANSSHFNAME), '"'
print '(a, i0, a)', ' prmest: ', nprmest, ' fields'
do i = 1, nprmest
print '(a, a, a, f5.3)', ' prmestname = "', trim(PRMESTNAME(i)), '", infl = ', PRMINFL(i)
end do
print *
end subroutine prm_print
function prm_getinfl(fldname)
real :: prm_getinfl
character(*), intent(in) :: fldname
integer :: i
prm_getinfl = INFL
do i = 1, nprmest
if (trim(fldname) == PRMESTNAME(i)) then
prm_getinfl = PRMINFL(i)
print '(a, a, a, f5.3)', ' "', trim(fldname), '": using inflation = ', prm_getinfl
return
end if
end do
end function prm_getinfl
function prm_prmestexists(varname)
logical :: prm_prmestexists
character(*), intent(in) :: varname
integer :: i
prm_prmestexists = .false.
do i = 1, nprmest
if (trim(varname) == PRMESTNAME(i)) then
prm_prmestexists = .true.
return
end if
end do
end function prm_prmestexists
! Shift a character string to upper case.
!
subroutine ucase(string)
character(*) :: string
integer :: i
do i = 1, len(string)
if (string(i:i) >= 'a' .and. string(i:i) <= 'z') then
string(i:i) = achar (iachar ( string(i:i) ) - 32)
end if
end do
end subroutine ucase
end module m_parameters