10 use fckit_configuration_module,
only: fckit_configuration
11 use fckit_log_module,
only : fckit_log
14 use missing_values_mod
32 character(len=MAXVARLEN),
public,
allocatable :: varin(:)
33 integer,
allocatable :: channels(:)
39 integer :: nchan_total
62 type(fckit_configuration),
intent(in) :: f_confOper
63 integer(c_int),
intent(in) :: channels(:)
65 type(fckit_configuration) :: f_confOpts
66 type(fckit_configuration) :: f_confLinOper
70 call f_confoper % get_or_die(
"obs options",f_confopts)
74 if ( f_confoper%has(
"linear obs operator") )
then
75 call f_confoper%get_or_die(
"linear obs operator",f_conflinoper)
84 write(
message,*)
'ufo_radiancerttov_setup error: H2O must be included in RTTOV Absorbers'
90 allocate(self%varin(nvars_in))
94 do jspec = 1, self%conf%ngas
95 self%varin(ind) = self%conf%Absorbers(jspec)
110 allocate(self%channels(
size(channels)))
111 self%channels(:) = channels(:)
120 integer(kind=jpim) :: errorstatus
122 self % ltraj = .false.
132 use fckit_mpi_module,
only: fckit_mpi_comm
138 type(c_ptr),
value,
intent(in) :: obss
141 real(c_double) :: missing
142 type(fckit_mpi_comm) :: f_comm
145 character(*),
parameter :: routine_name =
'ufo_radiancerttov_tlad_settraj'
146 type(rttov_chanprof),
allocatable :: chanprof(:)
149 integer(kind=jpim) :: errorstatus
151 integer :: i_inst, nchan_total, ichan, iprof, prof, iprof_rttov
152 integer :: nprof_sim, nprof_max_sim, ichan_sim
153 integer :: prof_start, prof_end
155 logical :: jacobian_needed
157 include
'rttov_k.interface'
160 missing = missing_value(missing)
163 call obsspace_get_comm(obss, f_comm)
172 self % nprofiles = geovals % nlocs
174 self % nlevels = geoval_temp % nval
178 write(
message,
'(A, A, I0, A, I0, A)') &
179 trim(routine_name),
': Allocating ', self % nprofiles,
' profiles with ', self % nlevels,
' levels'
181 call self % RTprof_K % alloc_profs(errorstatus, self % conf, self % nprofiles, self % nlevels, init=.true., asw=1)
184 write(
message,
'(A, A, I0, A, I0, A)') trim(routine_name),
': Creating RTTOV profiles from geovals'
186 call self % RTprof_K % setup(geovals,obss,self % conf)
198 trim(routine_name),
': Allocating Trajectory resources for RTTOV K: ', self % nprofiles *
nchan_inst,
' total channels'
199 call self % RTprof_K % alloc_profs_K(errorstatus, self % conf, self % nprofiles *
nchan_inst, self % nlevels, init=.true., asw=1)
202 allocate(self % RTprof_K % chanprof ( self % nprofiles *
nchan_inst ))
205 if(self % conf % prof_by_prof)
then
208 nprof_max_sim = max(1,self % conf % nchan_max_sim /
nchan_inst)
210 nprof_sim = min(nprof_max_sim, self % nprofiles)
213 nchan_sim = nprof_sim *
size(self%channels)
216 write(
message,
'(A,A,I0,A,I0,A)') &
217 trim(routine_name),
': Allocating resources for RTTOV direct (K): ', nprof_sim,
' and ',
nchan_sim,
' channels'
219 call self % RTprof_K % alloc_direct(errorstatus, self % conf, nprof_sim,
nchan_sim, self % nlevels, init=.true., asw=1)
221 write(
message,
'(A,A,I0,A,I0,A)') &
222 trim(routine_name),
': Allocating resources for RTTOV K code: ', nprof_sim,
' and ',
nchan_sim,
' channels'
224 call self % RTprof_K % alloc_k(errorstatus, self % conf, nprof_sim,
nchan_sim, self % nlevels, init=.true., asw=1)
226 prof_start = 1; prof_end = self % nprofiles
229 rttov_loop :
do while (prof_start <= prof_end)
232 nprof_sim = min(nprof_sim, prof_end - prof_start + 1)
233 nchan_sim = nprof_sim *
size(self%channels)
237 chanprof(:) % prof = 0
238 chanprof(:) % chan = 0
244 do iprof_rttov = 1, nprof_sim
245 errorstatus = errorstatus_success
248 iprof = prof_start + iprof_rttov - 1
251 if(any(self % conf % inspect == iprof))
call self % RTprof_K % print(self % conf, iprof, i_inst)
254 if(self % conf % RTTOV_profile_checkinput)
call self % RTprof_K % check(self % conf, iprof, i_inst, errorstatus)
256 if (errorstatus == errorstatus_success)
then
258 ichan_sim = ichan_sim + 1_jpim
259 chanprof(ichan_sim) % prof = iprof_rttov
260 chanprof(ichan_sim) % chan = self % channels(ichan)
261 self % RTprof_K % chanprof(nchan_total + ichan_sim) % prof = iprof
262 self % RTprof_K % chanprof(nchan_total + ichan_sim) % chan = self % channels(ichan)
270 call self % RTProf_K % init_emissivity(self % conf, prof_start)
279 self % conf % rttov_opts, &
280 self % RTprof_K % profiles(prof_start:prof_start + nprof_sim - 1), &
281 self % RTprof_K % profiles_k(nchan_total + 1 : nchan_total +
nchan_sim), &
282 self % conf % rttov_coef_array(i_inst), &
283 self % RTprof_K % transmission, &
284 self % RTprof_K % transmission_k, &
285 self % RTprof_K % radiance, &
286 self % RTprof_K % radiance_k, &
287 calcemis = self % RTprof_K % calcemis(1:
nchan_sim), &
288 emissivity = self % RTprof_K % emissivity(1:
nchan_sim), &
289 emissivity_k = self % RTprof_K % emissivity_k(1:
nchan_sim))
291 if ( errorstatus /= errorstatus_success )
then
292 write(
message,
'(A, A, 2I6)') trim(routine_name),
'after rttov_k: error ', errorstatus, i_inst, &
293 ' skipping profiles ', prof_start,
' -- ', prof_start + nprof_sim - 1
298 if(hofxdiags%nvar > 0)
call populate_hofxdiags(self % RTprof_K, chanprof, self % conf, prof_start, hofxdiags)
303 prof_start = prof_start + nprof_sim
305 self % nchan_total = nchan_total
306 deallocate (chanprof)
311 call self % RTprof_K % alloc_k(errorstatus, self % conf, -1, -1, -1, asw=0)
312 call self % RTprof_K % alloc_direct(errorstatus, self % conf, -1, -1, -1, asw=0)
313 call self % RTprof_K % alloc_profs(errorstatus, self % conf, -1, -1, asw=0)
318 self % ltraj = .true.
331 type(c_ptr),
value,
intent(in) :: obss
332 integer,
intent(in) :: nvars, nlocs
333 real(c_double),
intent(inout) :: hofx(nvars, nlocs)
335 character(len=*),
parameter :: myname_=
"ufo_radiancerttov_simobs_tl"
336 integer :: ichan, jchan, prof, jspec
338 type(
ufo_geoval),
pointer :: geoval_d, geoval_d2
344 if (.not. self % ltraj)
then
345 write(
message,*) myname_,
' trajectory wasnt set!'
350 if (geovals % nlocs /= self % nprofiles)
then
351 write(
message,*) myname_,
' error: nlocs inconsistent!'
364 if (geoval_d % nval /= self % nlevels)
then
365 write(
message,*) myname_,
' error: layers inconsistent!'
370 do ichan = 1, self % nchan_total,
size(self%channels)
371 prof = self % RTprof_K % chanprof(ichan) % prof
372 do jchan = 1,
size(self%channels)
373 hofx(jchan,prof) = hofx(jchan,prof) + &
374 sum(self % RTprof_K % profiles_k(ichan+jchan-1) % t(self % nlevels:1:-1) * geoval_d % vals(1:geoval_d % nval,prof))
378 do jspec = 1, self%conf%ngas
382 if (geoval_d % nval /= self % nlevels)
then
383 write(
message,*) myname_,
' error: layers inconsistent!'
390 do ichan = 1, self % nchan_total,
size(self%channels)
391 prof = self % RTprof_K % chanprof(ichan) % prof
392 do jchan = 1,
size(self%channels)
393 if(self%conf%Absorbers(jspec) ==
var_q)
then
394 hofx(jchan,prof) = hofx(jchan,prof) + &
395 sum(self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * &
396 geoval_d % vals(1:geoval_d % nval,prof)) * self%conf%scale_fac(gas_id_watervapour)
397 elseif(self%conf%Absorbers(jspec) ==
var_mixr)
then
398 hofx(jchan,prof) = hofx(jchan,prof) + &
399 sum(self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * &
400 geoval_d % vals(1:geoval_d % nval,prof)) * self%conf%scale_fac(gas_id_watervapour) /
g_to_kg
401 elseif(self%conf%Absorbers(jspec) ==
var_clw)
then
402 hofx(jchan,prof) = hofx(jchan,prof) + &
403 sum(self % RTprof_K % profiles_k(ichan+jchan-1) % clw(self % nlevels:1:-1) * &
404 geoval_d % vals(1:geoval_d % nval,prof))
424 do ichan = 1, self % nchan_total,
size(self%channels)
425 prof = self % RTprof_K % chanprof(ichan) % prof
426 do jchan = 1,
size(self%channels)
427 hofx(jchan,prof) = hofx(jchan,prof) + &
428 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % t * geoval_d % vals(1,prof)
435 do ichan = 1, self % nchan_total,
size(self%channels)
436 prof = self % RTprof_K % chanprof(ichan) % prof
437 do jchan = 1,
size(self%channels)
438 hofx(jchan,prof) = hofx(jchan,prof) + &
439 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % q * &
440 geoval_d % vals(1,prof) * self%conf%scale_fac(gas_id_watervapour)
448 do ichan = 1, self % nchan_total,
size(self%channels)
449 prof = self % RTprof_K % chanprof(ichan) % prof
450 do jchan = 1,
size(self%channels)
451 hofx(jchan,prof) = hofx(jchan,prof) + &
452 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % u * geoval_d % vals(1,prof) + &
453 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % v * geoval_d2 % vals(1,prof)
460 do ichan = 1, self % nchan_total,
size(self%channels)
461 prof = self % RTprof_K % chanprof(ichan) % prof
462 do jchan = 1,
size(self%channels)
463 hofx(jchan,prof) = hofx(jchan,prof) + &
464 self % RTprof_K % profiles_k(ichan+jchan-1) % skin % t * geoval_d % vals(1,prof)
479 type(c_ptr),
value,
intent(in) :: obss
480 integer,
intent(in) :: nvars, nlocs
481 real(c_double),
intent(in) :: hofx(nvars, nlocs)
483 type(
ufo_geoval),
pointer :: geoval_d, geoval_d2
485 real(c_double) :: missing
486 integer :: ichan, jchan, prof, jspec
488 character(len=*),
parameter :: myname_ =
"ufo_radiancerttov_simobs_ad"
491 missing = missing_value(missing)
497 if (.not. self % ltraj)
then
498 write(
message,*) myname_,
' trajectory wasnt set!'
503 if (geovals % nlocs /= self % nprofiles)
then
504 write(
message,*) myname_,
' error: nlocs inconsistent!'
512 do ichan = 1, self % nchan_total,
size(self%channels)
513 prof = self % RTprof_K % chanprof(ichan) % prof
514 do jchan = 1,
size(self%channels)
515 if (hofx(jchan, prof) /= missing)
then
516 geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
517 self % RTprof_K % profiles_k(ichan+jchan-1) % t(self % nlevels:1:-1) * hofx(jchan,prof)
526 do jspec = 1, self%conf%ngas
529 do ichan = 1, self % nchan_total,
size(self%channels)
530 prof = self % RTprof_K % chanprof(ichan) % prof
531 do jchan = 1,
size(self%channels)
532 if (hofx(jchan, prof) /= missing)
then
534 if(self%conf%Absorbers(jspec) ==
var_q)
then
535 geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
536 self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * hofx(jchan,prof) * &
537 self%conf%scale_fac(gas_id_watervapour)
538 elseif(self%conf%Absorbers(jspec) ==
var_mixr)
then
539 geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
540 (self % RTprof_K % profiles_k(ichan+jchan-1) % q(self % nlevels:1:-1) * hofx(jchan,prof)) * &
541 self%conf%scale_fac(gas_id_watervapour) /
g_to_kg
542 elseif(self%conf%Absorbers(jspec) ==
var_clw)
then
543 geoval_d % vals(:,prof) = geoval_d % vals(:,prof) + &
544 self % RTprof_K % profiles_k(ichan+jchan-1) % clw(self % nlevels:1:-1) * hofx(jchan,prof)
564 do ichan = 1, self % nchan_total,
size(self%channels)
565 prof = self % RTprof_K % chanprof(ichan) % prof
566 do jchan = 1,
size(self%channels)
567 if (hofx(jchan, prof) /= missing)
then
568 geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
569 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % t * hofx(jchan,prof)
577 do ichan = 1, self % nchan_total,
size(self%channels)
578 prof = self % RTprof_K % chanprof(ichan) % prof
579 do jchan = 1,
size(self%channels)
580 if (hofx(jchan, prof) /= missing)
then
581 geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
582 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % q * &
583 hofx(jchan,prof) * self%conf%scale_fac(gas_id_watervapour)
592 do ichan = 1, self % nchan_total,
size(self%channels)
593 prof = self % RTprof_K % chanprof(ichan) % prof
594 do jchan = 1,
size(self%channels)
595 if (hofx(jchan, prof) /= missing)
then
596 geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
597 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % u * hofx(jchan,prof)
599 geoval_d2 % vals(1,prof) = geoval_d2 % vals(1,prof) + &
600 self % RTprof_K % profiles_k(ichan+jchan-1) % s2m % v * hofx(jchan,prof)
608 do ichan = 1, self % nchan_total,
size(self%channels)
609 prof = self % RTprof_K % chanprof(ichan) % prof
610 do jchan = 1,
size(self%channels)
611 if (hofx(jchan, prof) /= missing)
then
612 geoval_d % vals(1,prof) = geoval_d % vals(1,prof) + &
613 self % RTprof_K % profiles_k(ichan+jchan-1) % skin % t * hofx(jchan,prof)
620 if (.not. geovals % linit ) geovals % linit=.true.
real(kind_real), parameter, public g_to_kg
real(kind_real), parameter, public zero
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module for radiancerttov tl/ad observation operator.
subroutine ufo_radiancerttov_tlad_delete(self)
subroutine ufo_radiancerttov_tlad_setup(self, f_confOper, channels)
subroutine ufo_radiancerttov_tlad_settraj(self, geovals, obss, hofxdiags)
subroutine ufo_radiancerttov_simobs_ad(self, geovals, obss, nvars, nlocs, hofx)
character(len=maxvarlen), dimension(1), parameter varin_default_tlad
subroutine ufo_radiancerttov_simobs_tl(self, geovals, obss, nvars, nlocs, hofx)
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
subroutine, public parse_hofxdiags(hofxdiags, jacobian_needed)
integer, public nchan_inst
subroutine, public rttov_conf_setup(conf, f_confOpts, f_confOper)
character(len=max_string), public message
subroutine, public rttov_conf_delete(conf)
subroutine, public populate_hofxdiags(RTProf, chanprof, conf, prof_start, hofxdiags)
integer, public nchan_sim
integer function, public ufo_vars_getindex(vars, varname)
character(len=maxvarlen), parameter, public var_clw
character(len=maxvarlen), parameter, public var_sfc_v10
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_sfc_u10
character(len=maxvarlen), parameter, public var_sfc_q2m
character(len=maxvarlen), parameter, public var_mixr
character(len=maxvarlen), parameter, public var_sfc_tskin
character(len=maxvarlen), parameter, public var_ts
character(len=maxvarlen), parameter, public var_sfc_t2m
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators
Fortran derived type for radiancerttov trajectory.