12 use fckit_configuration_module,
only: fckit_configuration
15 use missing_values_mod
31 character(len=MAXVARLEN),
public,
allocatable :: varin(:)
32 integer,
allocatable :: channels(:)
38 type(crtm_atmosphere_type),
allocatable :: atm_k(:,:)
39 type(crtm_surface_type),
allocatable :: sfc_k(:,:)
41 logical,
allocatable :: skip_profiles(:)
61 type(fckit_configuration),
intent(in) :: f_confOper
62 integer(c_int),
intent(in) :: channels(:)
66 type(fckit_configuration) :: f_confOpts,f_confLinOper
67 character(max_string) :: err_msg
69 call f_confoper%get_or_die(
"obs options",f_confopts)
72 if ( f_confoper%has(
"linear obs operator") )
then
73 call f_confoper%get_or_die(
"linear obs operator",f_conflinoper)
82 nvars_in =
size(
varin_default) + self%conf%n_Absorbers + self%conf%n_Clouds + self%conf%n_Surfaces
83 allocate(self%varin(nvars_in))
88 do jspec = 1, self%conf%n_Absorbers
89 self%varin(ind) = self%conf%Absorbers(jspec)
92 do jspec = 1, self%conf%n_Clouds
93 self%varin(ind) = self%conf%Clouds(jspec,1)
96 do jspec = 1, self%conf%n_Surfaces
97 self%varin(ind) = self%conf%Surfaces(jspec)
102 trim(self%conf_traj%sfc_wind_geovars) /=
"vector")
then
103 write(err_msg,*)
'ufo_radiancecrtm_tlad_setup error: sfc_wind_geovars not supported in tlad --> ', self%conf_traj%sfc_wind_geovars
104 call abor1_ftn(err_msg)
108 allocate(self%channels(
size(channels)))
109 self%channels(:) = channels(:)
125 if (
allocated(self%atm_k))
then
126 call crtm_atmosphere_destroy(self%atm_K)
127 deallocate(self%atm_k)
130 if (
allocated(self%sfc_k))
then
131 call crtm_surface_destroy(self%sfc_K)
132 deallocate(self%sfc_k)
135 if (
allocated(self%Skip_Profiles))
deallocate(self%Skip_Profiles)
142 use fckit_mpi_module,
only: fckit_mpi_comm
143 use fckit_log_module,
only: fckit_log
144 use ieee_arithmetic,
only: ieee_is_nan
151 type(c_ptr),
value,
intent(in) :: obss
155 character(*),
parameter :: PROGRAM_NAME =
'ufo_radiancecrtm_tlad_settraj'
156 character(255) :: message, version
157 character(max_string) :: err_msg
158 integer :: err_stat, alloc_stat
161 integer :: jvar, ichannel, jchannel, jprofile, jlevel, jspec
162 real(c_double) :: missing
163 type(fckit_mpi_comm) :: f_comm
167 type(crtm_channelinfo_type) :: chinfo(self%conf_traj%n_Sensors)
168 type(crtm_geometry_type),
allocatable :: geo(:)
171 type(crtm_atmosphere_type),
allocatable :: atm(:)
172 type(crtm_surface_type),
allocatable :: sfc(:)
173 type(crtm_rtsolution_type),
allocatable :: rts(:,:)
174 type(crtm_options_type),
allocatable :: Options(:)
177 type(crtm_rtsolution_type),
allocatable :: rts_K(:,:)
180 type(crtm_geometry_type),
allocatable :: geo_hf(:)
181 type(crtm_atmosphere_type),
allocatable :: atm_Ka(:,:)
182 type(crtm_surface_type),
allocatable :: sfc_Ka(:,:)
183 type(crtm_rtsolution_type),
allocatable :: rtsa(:,:)
184 type(crtm_rtsolution_type),
allocatable :: rts_Ka(:,:)
188 character(len=MAXVARLEN) :: varstr
189 character(len=MAXVARLEN),
dimension(hofxdiags%nvar) :: &
190 ystr_diags, xstr_diags
191 character(10),
parameter :: jacobianstr =
"_jacobian_"
192 integer :: str_pos(4), ch_diags(hofxdiags%nvar),numNaN
194 real(kind_real) :: total_od, secant_term, wfunc_max
195 real(kind_real),
allocatable :: TmpVar(:)
196 real(kind_real),
allocatable :: Tao(:)
197 real(kind_real),
allocatable :: Wfunc(:)
199 call obsspace_get_comm(obss, f_comm)
203 self%n_Profiles = geovals%nlocs
205 self%n_Layers = temp%nval
223 err_stat = crtm_init( self%conf_traj%SENSOR_ID, chinfo, &
224 file_path=trim(self%conf_traj%COEFFICIENT_PATH), &
225 irwatercoeff_file=trim(self%conf_traj%IRwaterCoeff_File), &
226 irlandcoeff_file=trim(self%conf_traj%IRlandCoeff_File), &
227 irsnowcoeff_file=trim(self%conf_traj%IRsnowCoeff_File), &
228 iricecoeff_file=trim(self%conf_traj%IRiceCoeff_File), &
229 viswatercoeff_file=trim(self%conf_traj%VISwaterCoeff_File), &
230 vislandcoeff_file=trim(self%conf_traj%VISlandCoeff_File), &
231 vissnowcoeff_file=trim(self%conf_traj%VISsnowCoeff_File), &
232 visicecoeff_file=trim(self%conf_traj%VISiceCoeff_File), &
233 mwwatercoeff_file=trim(self%conf_traj%MWwaterCoeff_File), &
235 message =
'Error initializing CRTM (setTraj)'
240 sensor_loop:
do n = 1, self%conf_traj%n_Sensors
245 err_stat = crtm_channelinfo_subset(chinfo(n), self%channels, reset=.false.)
246 message =
'Error subsetting channels'
251 self%n_Channels = crtm_channelinfo_n_channels(chinfo(n))
256 allocate( geo( self%n_Profiles ) , &
257 atm( self%n_Profiles ) , &
258 sfc( self%n_Profiles ) , &
259 rts( self%n_Channels, self%n_Profiles ) , &
260 self%atm_K( self%n_Channels, self%n_Profiles ) , &
261 self%sfc_K( self%n_Channels, self%n_Profiles ) , &
262 rts_k( self%n_Channels, self%n_Profiles ) , &
263 options( self%n_Profiles ) , &
265 message =
'Error allocating structure arrays (setTraj)'
270 call crtm_atmosphere_create( atm, self%n_Layers, self%conf_traj%n_Absorbers, self%conf_traj%n_Clouds, self%conf_traj%n_Aerosols )
271 if ( any(.NOT. crtm_atmosphere_associated(atm)) )
THEN
272 message =
'Error allocating CRTM Forward Atmosphere structure (setTraj)'
273 CALL display_message( program_name, message, failure )
277 if (self%n_Layers > 0)
CALL crtm_rtsolution_create(rts, self%n_Layers )
281 call crtm_surface_create(sfc, self%n_Channels)
282 IF ( any(.NOT. crtm_surface_associated(sfc)) )
THEN
283 message =
'Error allocating CRTM Surface structure (setTraj)'
284 CALL display_message( program_name, message, failure )
291 call crtm_atmosphere_create( self%atm_K, self%n_Layers, self%conf_traj%n_Absorbers, &
292 self%conf_traj%n_Clouds, self%conf_traj%n_Aerosols )
293 if ( any(.NOT. crtm_atmosphere_associated(self%atm_K)) )
THEN
294 message =
'Error allocating CRTM K-matrix Atmosphere structure (setTraj)'
295 CALL display_message( program_name, message, failure )
302 call crtm_surface_create(self%sfc_K, self%n_Channels)
303 IF ( any(.NOT. crtm_surface_associated(self%sfc_K)) )
THEN
304 message =
'Error allocating CRTM K-matrix Surface structure (setTraj)'
305 CALL display_message( program_name, message, failure )
311 call load_atm_data(self%N_PROFILES,self%N_LAYERS,geovals,atm,self%conf_traj)
312 call load_sfc_data(self%N_PROFILES,self%n_Channels,self%channels,geovals,sfc,chinfo,obss,self%conf_traj)
313 if (
cmp_strings(self%conf%SENSOR_ID(n),
'gmi_gpm'))
then
314 allocate( geo_hf( self%n_Profiles ))
322 call crtm_atmosphere_zero( self%atm_K )
323 call crtm_surface_zero( self%sfc_K )
328 rts_k%Radiance = zero
329 rts_k%Brightness_Temperature = one
331 if (
allocated(self%Skip_Profiles))
deallocate(self%Skip_Profiles)
332 allocate(self%Skip_Profiles(self%n_Profiles))
334 profile_loop:
do jprofile = 1, self%n_Profiles
335 options(jprofile)%Skip_Profile = self%Skip_Profiles(jprofile)
337 do jlevel = atm(jprofile)%n_layers, 1, -1
338 if ( atm(jprofile)%level_pressure(jlevel) <= atm(jprofile)%level_pressure(jlevel-1) )
then
339 options(jprofile)%Skip_Profile = .true.
347 err_stat = crtm_k_matrix( atm , &
356 message =
'Error calling CRTM (setTraj) K-Matrix Model for '//trim(self%conf_traj%SENSOR_ID(n))
358 if (
cmp_strings(self%conf%SENSOR_ID(n),
'gmi_gpm'))
then
359 allocate( atm_ka( self%n_Channels, self%n_Profiles ), &
360 sfc_ka( self%n_Channels, self%n_Profiles ), &
361 rts_ka( self%n_Channels, self%n_Profiles ), &
362 rtsa( self%n_Channels, self%n_Profiles ), &
364 message =
'Error allocating K structure arrays rtsa, atm_Ka ......'
373 call crtm_atmosphere_zero( self%atm_K )
374 call crtm_surface_zero( self%sfc_K )
377 rts_k%Radiance = zero
378 rts_k%Brightness_Temperature = one
381 err_stat = crtm_k_matrix( atm , &
390 message =
'Error calling CRTM (setTraj, geo_hf) K-Matrix Model for '&
391 //trim(self%conf_traj%SENSOR_ID(n))
394 do lch = 1,
size(self%channels)
395 if ( self%channels(lch) <= 9 )
then
396 self%atm_K(lch,:) = atm_ka(lch,:)
397 self%sfc_K(lch,:) = sfc_ka(lch,:)
398 rts_k(lch,:) = rts_ka(lch,:)
399 rts(lch,:) = rtsa(lch,:)
402 deallocate(atm_ka,sfc_ka,rts_ka,rtsa)
409 do jprofile = 1, self%n_Profiles
410 do jchannel = 1,
size(self%channels)
411 do jlevel = 1, self%atm_K(jchannel,jprofile)%n_layers
412 if (ieee_is_nan(self%atm_K(jchannel,jprofile)%Temperature(jlevel)))
then
413 self%Skip_Profiles(jprofile) = .true.
415 write(message,*) numnan,
'th NaN in Jacobian Profiles'
416 call fckit_log%info(message)
429 do jvar = 1, hofxdiags%nvar
430 varstr = hofxdiags%variables(jvar)
431 str_pos(4) = len_trim(varstr)
432 if (str_pos(4) < 1) cycle
433 str_pos(3) = index(varstr,
"_",back=.true.)
434 read(varstr(str_pos(3)+1:str_pos(4)),*, err=999) ch_diags(jvar)
435 999 str_pos(1) = index(varstr,jacobianstr) - 1
436 if (str_pos(1) == 0)
then
437 write(err_msg,*)
'ufo_radiancecrtm_simobs: _jacobian_ must be // &
438 & preceded by dependent variable in config: ', &
439 & hofxdiags%variables(jvar)
440 call abor1_ftn(err_msg)
441 else if (str_pos(1) > 0)
then
443 ystr_diags(jvar) = varstr(1:str_pos(1))
444 str_pos(2) = str_pos(1) + len(jacobianstr) + 1
445 str_pos(4) = str_pos(3) - str_pos(2)
446 xstr_diags(jvar)(1:str_pos(4)) = varstr(str_pos(2):str_pos(3)-1)
447 xstr_diags(jvar)(str_pos(4)+1:) =
""
450 xstr_diags(jvar) =
""
451 ystr_diags(jvar)(1:str_pos(3)-1) = varstr(1:str_pos(3)-1)
452 ystr_diags(jvar)(str_pos(3):) =
""
453 if (ch_diags(jvar) < 0) ystr_diags(jvar) = varstr
458 missing = missing_value(missing)
462 do jvar = 1, hofxdiags%nvar
464 if (len(trim(hofxdiags%variables(jvar))) < 1) cycle
466 if (ch_diags(jvar) > 0)
then
467 if (
size(pack(self%channels,self%channels==ch_diags(jvar))) /= 1)
then
468 write(err_msg,*)
'ufo_radiancecrtm_simobs: mismatch between// &
469 & h(x) channels(', self%channels,
') and// &
470 & ch_diags(jvar) = ', ch_diags(jvar)
471 call abor1_ftn(err_msg)
475 do ichannel = 1,
size(self%channels)
476 if (ch_diags(jvar) == self%channels(ichannel))
then
482 if (
allocated(hofxdiags%geovals(jvar)%vals)) &
483 deallocate(hofxdiags%geovals(jvar)%vals)
490 select case(ystr_diags(jvar))
493 hofxdiags%geovals(jvar)%nval = self%n_Layers
494 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
495 hofxdiags%geovals(jvar)%vals = missing
496 do jprofile = 1, self%n_Profiles
497 if (.not.self%Skip_Profiles(jprofile))
then
498 do jlevel = 1, hofxdiags%geovals(jvar)%nval
499 hofxdiags%geovals(jvar)%vals(jlevel,jprofile) = &
500 rts(jchannel,jprofile) % layer_optical_depth(jlevel)
507 hofxdiags%geovals(jvar)%nval = 1
508 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
509 hofxdiags%geovals(jvar)%vals = missing
510 do jprofile = 1, self%n_Profiles
511 if (.not.self%Skip_Profiles(jprofile))
then
512 hofxdiags%geovals(jvar)%vals(1,jprofile) = &
513 rts(jchannel,jprofile) % Brightness_Temperature
519 hofxdiags%geovals(jvar)%nval = self%n_Layers
520 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
521 hofxdiags%geovals(jvar)%vals = missing
522 allocate(tmpvar(self%n_Profiles))
523 call obsspace_get_db(obss,
"MetaData",
"sensor_zenith_angle", tmpvar)
524 do jprofile = 1, self%n_Profiles
525 if (.not.self%Skip_Profiles(jprofile))
then
526 secant_term = one/cos(tmpvar(jprofile)*
deg2rad)
528 do jlevel = 1, self%n_Layers
529 total_od = total_od + rts(jchannel,jprofile) % layer_optical_depth(jlevel)
530 hofxdiags%geovals(jvar)%vals(jlevel,jprofile) = &
531 exp(-min(limit_exp,total_od*secant_term))
539 hofxdiags%geovals(jvar)%nval = self%n_Layers
540 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
541 hofxdiags%geovals(jvar)%vals = missing
542 allocate(tmpvar(self%n_Profiles))
543 allocate(tao(self%n_Layers))
544 call obsspace_get_db(obss,
"MetaData",
"sensor_zenith_angle", tmpvar)
545 do jprofile = 1, self%n_Profiles
546 if (.not.self%Skip_Profiles(jprofile))
then
548 secant_term = one/cos(tmpvar(jprofile)*
deg2rad)
550 do jlevel = 1, self%n_Layers
551 total_od = total_od + rts(jchannel,jprofile) % layer_optical_depth(jlevel)
552 tao(jlevel) = exp(-min(limit_exp,total_od*secant_term))
555 do jlevel = self%n_Layers-1, 1, -1
556 hofxdiags%geovals(jvar)%vals(jlevel,jprofile) = &
557 abs( (tao(jlevel+1)-tao(jlevel))/ &
558 (log(atm(jprofile)%pressure(jlevel+1))- &
559 log(atm(jprofile)%pressure(jlevel))) )
561 hofxdiags%geovals(jvar)%vals(self%n_Layers,jprofile) = &
562 hofxdiags%geovals(jvar)%vals(self%n_Layers-1,jprofile)
570 hofxdiags%geovals(jvar)%nval = 1
571 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
572 hofxdiags%geovals(jvar)%vals = missing
573 allocate(tmpvar(self%n_Profiles))
574 allocate(tao(self%n_Layers))
575 allocate(wfunc(self%n_Layers))
576 call obsspace_get_db(obss,
"MetaData",
"sensor_zenith_angle", tmpvar)
577 do jprofile = 1, self%n_Profiles
578 if (.not.self%Skip_Profiles(jprofile))
then
580 secant_term = one/cos(tmpvar(jprofile)*
deg2rad)
582 do jlevel = 1, self%n_Layers
583 total_od = total_od + rts(jchannel,jprofile) % layer_optical_depth(jlevel)
584 tao(jlevel) = exp(-min(limit_exp,total_od*secant_term))
587 do jlevel = self%n_Layers-1, 1, -1
589 abs( (tao(jlevel+1)-tao(jlevel))/ &
590 (log(atm(jprofile)%pressure(jlevel+1))- &
591 log(atm(jprofile)%pressure(jlevel))) )
593 wfunc(self%n_Layers) = wfunc(self%n_Layers-1)
596 do jlevel = self%n_Layers-1, 1, -1
597 if (wfunc(jlevel) > wfunc_max)
then
598 wfunc_max = wfunc(jlevel)
599 hofxdiags%geovals(jvar)%vals(1,jprofile) = jlevel
609 write(err_msg,*)
'ufo_radiancecrtm_simobs: //&
610 & ObsDiagnostic is unsupported, ', &
611 & hofxdiags%variables(jvar)
613 hofxdiags%geovals(jvar)%nval = 1
614 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
615 hofxdiags%geovals(jvar)%vals = missing
617 else if (ystr_diags(jvar) ==
var_tb)
then
619 select case (xstr_diags(jvar))
622 hofxdiags%geovals(jvar)%nval = 1
623 allocate(hofxdiags%geovals(jvar)%vals(hofxdiags%geovals(jvar)%nval,self%n_Profiles))
624 hofxdiags%geovals(jvar)%vals = missing
625 do jprofile = 1, self%n_Profiles
626 if (.not.self%Skip_Profiles(jprofile))
then
627 hofxdiags%geovals(jvar)%vals(1,jprofile) = &
628 rts_k(jchannel,jprofile) % surface_emissivity
632 write(err_msg,*)
'ufo_radiancecrtm_simobs: //&
633 & ObsDiagnostic is unsupported, ', &
634 & hofxdiags%variables(jvar)
635 call abor1_ftn(err_msg)
638 write(err_msg,*)
'ufo_radiancecrtm_simobs: //&
639 & ObsDiagnostic is unsupported, ', &
640 & hofxdiags%variables(jvar)
641 call abor1_ftn(err_msg)
647 call crtm_geometry_destroy(geo)
648 call crtm_atmosphere_destroy(atm)
649 call crtm_rtsolution_destroy(rts_k)
650 call crtm_rtsolution_destroy(rts)
651 call crtm_surface_destroy(sfc)
656 deallocate(geo, atm, sfc, rts, rts_k, options, stat = alloc_stat)
657 if(
allocated(geo_hf))
deallocate(geo_hf)
658 message =
'Error deallocating structure arrays (setTraj)'
667 err_stat = crtm_destroy( chinfo )
668 message =
'Error destroying CRTM (setTraj)'
684 type(c_ptr),
value,
intent(in) :: obss
685 integer,
intent(in) :: nvars, nlocs
686 real(c_double),
intent(inout) :: hofx(nvars, nlocs)
688 character(len=*),
parameter :: myname_=
"ufo_radiancecrtm_simobs_tl"
689 character(max_string) :: err_msg
690 integer :: jprofile, jchannel, jlevel, jspec, ispec
697 if (.not. self%ltraj)
then
698 write(err_msg,*) myname_,
' trajectory wasnt set!'
699 call abor1_ftn(err_msg)
703 if (geovals%nlocs /= self%n_Profiles)
then
704 write(err_msg,*) myname_,
' error: nlocs inconsistent!'
705 call abor1_ftn(err_msg)
710 hofx(:,:) = 0.0_kind_real
719 if (geoval_d%nval /= self%n_Layers)
then
720 write(err_msg,*) myname_,
' error: layers inconsistent!'
721 call abor1_ftn(err_msg)
725 do jprofile = 1, self%n_Profiles
726 if (.not.self%Skip_Profiles(jprofile))
then
727 do jchannel = 1,
size(self%channels)
728 do jlevel = 1, geoval_d%nval
729 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
730 self%atm_K(jchannel,jprofile)%Temperature(jlevel) * &
731 geoval_d%vals(jlevel,jprofile)
740 do jspec = 1, self%conf%n_Absorbers
747 do jprofile = 1, self%n_Profiles
748 if (.not.self%Skip_Profiles(jprofile))
then
749 do jchannel = 1,
size(self%channels)
750 do jlevel = 1, geoval_d%nval
751 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
752 self%atm_K(jchannel,jprofile)%Absorber(jlevel,ispec) * &
753 geoval_d%vals(jlevel,jprofile)
763 do jspec = 1, self%conf%n_Clouds
770 do jprofile = 1, self%n_Profiles
771 if (.not.self%Skip_Profiles(jprofile))
then
772 do jchannel = 1,
size(self%channels)
773 do jlevel = 1, geoval_d%nval
774 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
775 self%atm_K(jchannel,jprofile)%Cloud(ispec)%Water_Content(jlevel) * &
776 geoval_d%vals(jlevel,jprofile)
786 do jspec = 1, self%conf%n_Surfaces
790 select case(self%conf%Surfaces(jspec))
795 do jprofile = 1, self%n_Profiles
796 if (.not.self%Skip_Profiles(jprofile))
then
797 do jchannel = 1,
size(self%channels)
799 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
800 self%sfc_K(jchannel,jprofile)%water_temperature * &
801 geoval_d%vals(jlevel,jprofile)
807 do jprofile = 1, self%n_Profiles
808 if (.not.self%Skip_Profiles(jprofile))
then
809 do jchannel = 1,
size(self%channels)
811 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
812 self%sfc_K(jchannel,jprofile)%wind_speed * &
813 geoval_d%vals(jlevel,jprofile)
819 do jprofile = 1, self%n_Profiles
820 if (.not.self%Skip_Profiles(jprofile))
then
821 do jchannel = 1,
size(self%channels)
823 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
824 self%sfc_K(jchannel,jprofile)%wind_direction * &
825 geoval_d%vals(jlevel,jprofile)
831 do jprofile = 1, self%n_Profiles
832 if (.not.self%Skip_Profiles(jprofile))
then
833 do jchannel = 1,
size(self%channels)
835 hofx(jchannel, jprofile) = hofx(jchannel, jprofile) + &
836 self%sfc_K(jchannel,jprofile)%salinity * &
837 geoval_d%vals(jlevel,jprofile)
856 type(c_ptr),
value,
intent(in) :: obss
857 integer,
intent(in) :: nvars, nlocs
858 real(c_double),
intent(in) :: hofx(nvars, nlocs)
860 character(len=*),
parameter :: myname_=
"ufo_radiancecrtm_simobs_ad"
861 character(max_string) :: err_msg
862 integer :: jprofile, jchannel, jlevel, jspec, ispec
864 real(c_double) :: missing
871 if (.not. self%ltraj)
then
872 write(err_msg,*) myname_,
' trajectory wasnt set!'
873 call abor1_ftn(err_msg)
877 if (geovals%nlocs /= self%n_Profiles)
then
878 write(err_msg,*) myname_,
' error: nlocs inconsistent!'
879 call abor1_ftn(err_msg)
883 missing = missing_value(missing)
892 do jprofile = 1, self%n_Profiles
893 if (.not.self%Skip_Profiles(jprofile))
then
894 do jchannel = 1,
size(self%channels)
895 if (hofx(jchannel, jprofile) /= missing)
then
896 do jlevel = 1, geoval_d%nval
897 geoval_d%vals(jlevel,jprofile) = geoval_d%vals(jlevel,jprofile) + &
898 self%atm_K(jchannel,jprofile)%Temperature(jlevel) * &
899 hofx(jchannel, jprofile)
909 do jspec = 1, self%conf%n_Absorbers
916 do jprofile = 1, self%n_Profiles
917 if (.not.self%Skip_Profiles(jprofile))
then
918 do jchannel = 1,
size(self%channels)
919 if (hofx(jchannel, jprofile) /= missing)
then
920 do jlevel = 1, geoval_d%nval
921 geoval_d%vals(jlevel,jprofile) = geoval_d%vals(jlevel,jprofile) + &
922 self%atm_K(jchannel,jprofile)%Absorber(jlevel,ispec) * &
923 hofx(jchannel, jprofile)
934 do jspec = 1, self%conf%n_Clouds
941 do jprofile = 1, self%n_Profiles
942 if (.not.self%Skip_Profiles(jprofile))
then
943 do jchannel = 1,
size(self%channels)
944 if (hofx(jchannel, jprofile) /= missing)
then
945 do jlevel = 1, geoval_d%nval
946 geoval_d%vals(jlevel,jprofile) = geoval_d%vals(jlevel,jprofile) + &
947 self%atm_K(jchannel,jprofile)%Cloud(ispec)%Water_Content(jlevel) * &
948 hofx(jchannel, jprofile)
959 do jspec = 1, self%conf%n_Surfaces
963 select case(self%conf%Surfaces(jspec))
967 do jprofile = 1, self%n_Profiles
968 if (.not.self%Skip_Profiles(jprofile))
then
969 do jchannel = 1,
size(self%channels)
970 if (hofx(jchannel, jprofile) /= missing)
then
972 geoval_d%vals(jlevel, jprofile) = geoval_d%vals(jlevel,jprofile) + &
973 self%sfc_K(jchannel,jprofile)%water_temperature * &
974 hofx(jchannel,jprofile)
981 do jprofile = 1, self%n_Profiles
982 if (.not.self%Skip_Profiles(jprofile))
then
983 do jchannel = 1,
size(self%channels)
984 if (hofx(jchannel, jprofile) /= missing)
then
986 geoval_d%vals(jlevel, jprofile) = geoval_d%vals(jlevel,jprofile) + &
987 self%sfc_K(jchannel,jprofile)%wind_speed * &
988 hofx(jchannel,jprofile)
995 do jprofile = 1, self%n_Profiles
996 if (.not.self%Skip_Profiles(jprofile))
then
997 do jchannel = 1,
size(self%channels)
998 if (hofx(jchannel, jprofile) /= missing)
then
1000 geoval_d%vals(jlevel, jprofile) = geoval_d%vals(jlevel,jprofile) + &
1001 self%sfc_K(jchannel,jprofile)%wind_direction * &
1002 hofx(jchannel,jprofile)
1009 do jprofile = 1, self%n_Profiles
1010 if (.not.self%Skip_Profiles(jprofile))
then
1011 do jchannel = 1,
size(self%channels)
1012 if (hofx(jchannel, jprofile) /= missing)
then
1014 geoval_d%vals(jlevel, jprofile) = geoval_d%vals(jlevel,jprofile) + &
1015 self%sfc_K(jchannel,jprofile)%salinity * &
1016 hofx(jchannel,jprofile)
1028 if (.not. geovals%linit ) geovals%linit=.true.
real(kind_real), parameter, public deg2rad
Fortran module to provide code shared between nonlinear and tlm/adm radiance calculations.
subroutine, public crtm_conf_setup(conf, f_confOpts, f_confOper)
subroutine, public load_geom_data(obss, geo, geo_hf)
subroutine, public crtm_comm_stat_check(stat, PROGRAM_NAME, message, f_comm)
subroutine, public load_sfc_data(n_Profiles, n_Channels, channels, geovals, sfc, chinfo, obss, conf)
subroutine, public ufo_crtm_skip_profiles(n_Profiles, n_Channels, channels, obss, Skip_Profiles)
subroutine, public crtm_conf_delete(conf)
subroutine, public load_atm_data(n_Profiles, n_Layers, geovals, atm, conf)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module to handle tl/ad for radiancecrtm observations.
subroutine ufo_radiancecrtm_simobs_ad(self, geovals, obss, nvars, nlocs, hofx)
subroutine ufo_radiancecrtm_tlad_setup(self, f_confOper, channels)
character(len=maxvarlen), dimension(1), parameter varin_default
subroutine ufo_radiancecrtm_tlad_delete(self)
subroutine ufo_radiancecrtm_simobs_tl(self, geovals, obss, nvars, nlocs, hofx)
subroutine ufo_radiancecrtm_tlad_settraj(self, geovals, obss, hofxdiags)
Fortran module with various useful routines.
logical function, public cmp_strings(str1, str2)
character(len=maxvarlen), parameter, public var_sfc_emiss
integer function, public ufo_vars_getindex(vars, varname)
character(len=maxvarlen), parameter, public var_sfc_wtmp
character(len=maxvarlen), parameter, public var_sfc_wdir
character(len=maxvarlen), parameter, public var_lvl_transmit
character(len=maxvarlen), parameter, public var_tb
character(len=maxvarlen), parameter, public var_lvl_weightfunc
character(len=maxvarlen), parameter, public var_sfc_sss
character(len=maxvarlen), parameter, public var_ts
character(len=maxvarlen), parameter, public var_pmaxlev_weightfunc
character(len=maxvarlen), parameter, public var_opt_depth
character(len=maxvarlen), parameter, public var_sfc_wspeed
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators
Fortran derived type for radiancecrtm trajectory.