8 use,
intrinsic:: iso_c_binding
15 use missing_values_mod
20 use fckit_log_module,
only : fckit_log
39 use fckit_configuration_module,
only: fckit_configuration
42 type(fckit_configuration),
intent(in) :: f_conf
52 real(kind_real),
intent(inout) :: hofx(:)
53 type(c_ptr),
value,
intent(in) :: obss
54 character(len=*),
parameter :: myname =
"ufo_gnssro_bndnbam_simobs"
55 character(max_string) :: err_msg
56 integer :: nrecs, nlocs
57 integer,
parameter :: nlevAdd = 13
58 integer,
parameter :: ngrd = 80
59 integer :: iobs, k, igrd, irec, icount, kk
60 integer :: nlev, nlev1, nlevExt, nlevCheck
61 type(
ufo_geoval),
pointer :: t, q, gph, prs, zs
62 real(kind_real),
allocatable :: gesT(:,:), gesZ(:,:), gesP(:,:), gesQ(:,:), gesTv(:,:), gesZs(:)
63 real(kind_real),
allocatable :: obsLat(:), obsImpP(:),obsLocR(:), obsGeoid(:), obsValue(:)
64 integer(c_size_t),
allocatable :: obsRecnum(:)
65 real(kind_real),
allocatable :: temperature(:)
68 real(kind_real) :: grids(ngrd)
69 real(kind_real),
allocatable :: refIndex(:), refXrad(:), geomz(:)
70 real(kind_real),
allocatable :: ref(:), radius(:)
71 real(kind_real) :: sIndx
74 integer,
allocatable :: nlocs_begin(:)
75 integer,
allocatable :: nlocs_end(:)
76 real(c_double) :: missing
77 integer,
allocatable :: super_refraction_flag(:), super(:), obs_max(:)
78 real(kind_real),
allocatable :: toss_max(:)
80 real(kind_real) :: gradRef, obsImpH
81 integer,
allocatable :: LayerIdx(:)
83 write(err_msg,*) myname,
": begin"
84 call fckit_log%info(err_msg)
86 nlocs = obsspace_get_nlocs(obss)
87 nrecs = obsspace_get_nrecs(obss)
88 write(err_msg,*) myname,
': nlocs from gelvals and hofx, nrecs', geovals%nlocs, nlocs, nrecs
89 call fckit_log%info(err_msg)
90 missing = missing_value(missing)
92 allocate(temperature(nlocs))
94 allocate(super_refraction_flag(nlocs))
95 super_refraction_flag = 0
96 allocate(layeridx(nlocs))
102 if (geovals%nlocs /=
size(hofx))
then
103 write(err_msg,*) myname,
': nlocs inconsistent!', geovals%nlocs,
size(hofx)
104 call abor1_ftn(err_msg)
112 if (self%roconf%vertlayer .eq.
"mass")
then
115 else if (self%roconf%vertlayer .eq.
"full")
then
121 write(err_msg,*) myname,
': vertlayer has to be mass of full, '//new_line(
'a')// &
122 ' will use full layer anyway'
123 call fckit_log%info(err_msg)
129 allocate(gesp(nlev1,nlocs))
130 allocate(gesz(nlev1,nlocs))
131 allocate(gest(nlev,nlocs))
132 allocate(gestv(nlev,nlocs))
133 allocate(gesq(nlev,nlocs))
134 allocate(geszs(nlocs))
138 if (prs%vals(1,1) .lt. prs%vals(prs%nval,1) )
then
140 write(err_msg,
'(a)')
' ufo_gnssro_bndnbam_simobs:'//new_line(
'a')// &
141 ' Model vertical height profile is in descending order,'//new_line(
'a')// &
142 ' but bndNBAM requires it to be ascending order, need flip'
143 call fckit_log%info(err_msg)
145 gest(k,:) = t%vals(nlev-k+1,:)
146 gesq(k,:) = q%vals(nlev-k+1,:)
147 gestv(k,:)= gest(k,:)*(1+gesq(k,:)*(rv_over_rd-1))
150 gesp(k,:) = prs%vals(nlev1-k+1,:)
151 gesz(k,:) = gph%vals(nlev1-k+1,:)
155 gest(k,:) = t%vals(k,:)
156 gesq(k,:) = q%vals(k,:)
157 gestv(k,:) = gest(k,:)*(1+gesq(k,:)*(rv_over_rd-1))
161 gesp(k,:) = prs%vals(k,:)
162 gesz(k,:) = gph%vals(k,:)
165 geszs(:) = zs%vals(1,:)
170 if ( nlev1 /= nlev )
then
172 gesq(k,:) = half* (gesq(k,:) + gesq(k-1,:))
173 gestv(k,:) = half* (gestv(k,:) + gestv(k-1,:))
177 gest(k,:) = gestv(k,:)/(1+ gesq(k,:)*(rv_over_rd-1))
182 allocate(obslat(nlocs))
183 allocate(obsimpp(nlocs))
184 allocate(obslocr(nlocs))
185 allocate(obsgeoid(nlocs))
186 allocate(obsvalue(nlocs))
187 allocate(obsrecnum(nlocs))
188 allocate(nlocs_begin(nrecs))
189 allocate(nlocs_end(nrecs))
191 call obsspace_get_db(obss,
"MetaData",
"latitude", obslat)
192 call obsspace_get_db(obss,
"MetaData",
"impact_parameter", obsimpp)
193 call obsspace_get_db(obss,
"MetaData",
"earth_radius_of_curvature", obslocr)
194 call obsspace_get_db(obss,
"MetaData",
"geoid_height_above_reference_ellipsoid", obsgeoid)
195 call obsspace_get_db(obss,
"ObsValue",
"bending_angle", obsvalue)
196 call obsspace_get_recnum(obss, obsrecnum)
202 if (obsrecnum(iobs+1) /= obsrecnum(iobs))
then
204 nlocs_end(icount-1)= iobs
205 nlocs_begin(icount) = iobs+1
208 nlocs_end(nrecs)= nlocs
209 if (icount /= nrecs)
then
210 write(err_msg,*) myname,
': record number is not consistent :', icount, nrecs
211 call abor1_ftn(err_msg)
214 nlevext = nlev + nlevadd
215 nlevcheck = min(23, nlev)
219 grids(igrd+1) = igrd *
ds
223 allocate(geomz(nlev))
224 allocate(radius(nlev))
225 allocate(ref(nlevext))
226 allocate(refindex(nlev))
227 allocate(refxrad(0:nlevext+1))
229 allocate(super(nlocs))
230 allocate(toss_max(nrecs))
231 allocate(obs_max(nrecs))
239 rec_loop:
do irec = 1, nrecs
241 obs_loop:
do icount = nlocs_begin(irec), nlocs_end(irec)
246 call geop2geometric(obslat(iobs), gesz(k,iobs)-geszs(iobs), geomz(k))
247 radius(k) = geomz(k) + geszs(iobs) + obsgeoid(iobs) + obslocr(iobs)
250 ref(k), self%roconf%use_compress)
251 refindex(k) = one + (
r1em6*ref(k))
252 refxrad(k) = refindex(k) * radius(k)
258 if (sindx < one .or. sindx > float(nlev)) cycle obs_loop
261 layeridx(iobs) = min(max(1, int(sindx)), nlev)
265 wi=min(max(1,indx),nlev)
266 wi2=max(1,min(indx+1,nlev))
268 wf=max(zero,min(wf,one))
269 temperature(iobs)=gestv(wi,iobs)*(one-wf)+gestv(wi2,iobs)*wf
273 if(trim(self%roconf%super_ref_qc) ==
"NBAM")
then
275 obsimph = (obsimpp(iobs) - obslocr(iobs)) *
r1em3
277 if (obsimph <=
six)
then
278 do k = nlevcheck, 1, -1
281 gradref = 1000.0 * (ref(k+1)-ref(k))/(radius(k+1)-radius(k))
283 if (abs(gradref) >= 0.75*
crit_gradrefr .and. obsimpp(iobs) <= refxrad(k+2))
then
284 super_refraction_flag(iobs) = 1
291 if (self%roconf%sr_steps > 1 &
292 .and. obsvalue(iobs) >= 0.03 )
then
294 do k = nlevcheck, 1, -1
295 gradref = 1000.0 * (ref(k+1)-ref(k))/(radius(k+1)-radius(k))
297 .and. super(iobs) == 0 &
298 .and. toss_max(irec) <= obsvalue(iobs))
then
300 toss_max(irec)= max(toss_max(irec), obsvalue(iobs))
309 else if(trim(self%roconf%super_ref_qc) ==
"ECMWF")
then
313 if (refxrad(k) - refxrad(k-1) < 10.0)
THEN
319 if (obsimpp(iobs) < refxrad(sr_hgt_idx))
then
320 super_refraction_flag(iobs) = 1
325 write(err_msg,*) myname,
': super refraction method has to be NBAM or ECMWF!'
326 call abor1_ftn(err_msg)
329 if (super_refraction_flag(iobs) .eq. 0)
then
331 obslat(iobs), obsgeoid(iobs), obslocr(iobs), obsimpp(iobs), &
333 nlev, nlevext, nlevadd, nlevcheck, &
334 radius(1:nlev),ref(1:nlevext),refindex(1:nlev),refxrad(0:nlevext), &
341 if (trim(self%roconf%super_ref_qc) ==
"NBAM" .and. self%roconf%sr_steps > 1 )
then
342 rec_loop2:
do irec = 1, nrecs
344 if (obs_max(irec) > 0 )
then
346 obs_loop2:
do k = nlocs_begin(irec), nlocs_end(irec)
347 obsimph = (obsimpp(k) - obslocr(k)) *
r1em3
348 if (obsimph<=
six .and. obsimpp(k)<=obsimpp(obs_max(irec)) .and. &
349 hofx(k)/=missing .and. super_refraction_flag(k)==0)
then
350 super_refraction_flag(k)=2
374 deallocate(obsrecnum)
375 deallocate(nlocs_begin)
376 deallocate(nlocs_end)
381 write(err_msg,*) myname,
": complete"
382 call fckit_log%info(err_msg)
386 call obsspace_put_db(obss,
"MetaData",
"temperature", temperature)
388 call obsspace_put_db(obss,
"SRflag",
"bending_angle", super_refraction_flag)
390 call obsspace_put_db(obss,
"LayerIdx",
"bending_angle", layeridx)
391 deallocate(super_refraction_flag)
392 deallocate(temperature)