updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_radiance / da_transform_xtoy_rttov_adj.inc
blob07f7a64973a65cf1d8ea882febe17f363b68328b
1 subroutine da_transform_xtoy_rttov_adj ( iv, jo_grad_y, jo_grad_x )
3    !---------------------------------------------------------------------------
4    ! Purpose: transform gradient from obs space to model grid space.
5    !
6    ! METHOD:  jo_grad_x = H^T jo_grad_y =  - H^T R^-1 ( d - H delta_x )
7    !           1. input gradient in obs space and reference state of RTTOV
8    !           2. call adjoint of RTM
9    !           3. adjoint of interpolation from model grid to obs loc
10    !---------------------------------------------------------------------------
12    implicit none
14    type (x_type), intent(inout)   :: jo_grad_x !
15    type (y_type),  intent(in)     :: jo_grad_y ! H' delta_x
16    type (iv_type), intent(in)     :: iv        ! O-B structure.
18 #ifdef RTTOV
20    integer                        :: i, j, k,c,l  ! Index dimension.
21    integer                        :: nlevels ! Number of obs levels.
22    integer                        :: num_rad  ! Number of radiance obs
23 #ifdef DM_PARALLEL
24    integer :: status(mpi_status_size) ! MPI status.
25 #endif
27    real, allocatable              :: model_mr(:)
29    integer            :: nchan, n
31    type(con_vars_type),  allocatable :: con_vars(:), con_vars_ad(:)
32    type(aux_vars_type),  allocatable :: aux_vars(:)
34    ! evenly distributed versions
35    type(con_vars_type),  allocatable  :: d_con_vars(:), d_con_vars_ad(:)
36    type(aux_vars_type),  allocatable  :: d_aux_vars(:)
37    real,    allocatable :: d_tb(:,:)
38    integer :: d_num_rad, l_num_rad,e,s,len,nl
39    real, allocatable :: temp_t(:,:),temp_q(:,:), temp_ps(:), rtemp(:,:)
40    real, allocatable :: temp_tb(:,:)
41    real, allocatable :: temp_t_jac(:,:,:), temp_q_jac(:,:,:), temp_ps_jac(:,:)
42    integer, allocatable :: itemp(:,:)
44    if (iv%num_inst < 1) return
46    if (trace_use) call da_trace_entry("da_transform_xtoy_rttov_adj")
48    do i = 1, iv%num_inst                 ! loop for sensor
50       if (iv%instid(i)%num_rad_glo < 1) cycle
52       num_rad  = iv%instid(i)%num_rad
53       nlevels  = iv%instid(i)%nlevels
54       nchan    = iv%instid(i)%nchan
56       if (num_rad > 0) then
58          allocate (model_mr(1:nlevels))
60          allocate (con_vars(num_rad))
61          allocate (con_vars_ad(num_rad))
62          allocate (aux_vars(num_rad))
64          !---------------------------------------------------------------
65          ! [1.0] assign tb = R^-1 Re :
66          !---------------------------------------------------------------
68          !---------------------------------------------
69          ! [2.0] get input of da_rttov_ad
70          !---------------------------------------------
72          do n = 1, num_rad
73             con_vars(n) % nlevels = nlevels 
75             allocate (con_vars(n)    % t(nlevels))
76             allocate (con_vars(n)    % q(nlevels))
77             allocate (con_vars_ad(n) % t(nlevels))
78             allocate (con_vars_ad(n) % q(nlevels))
80             con_vars_ad(n) % t(:) = 0.0
81             con_vars_ad(n) % q(:) = 0.0
82             con_vars_ad(n) % ps   = 0.0
84             con_vars(n) % t(:) = iv%instid(i)%t(:,n)
85             con_vars(n) % q(:) = iv%instid(i)%mr(:,n)
86             con_vars(n) % ps   = iv%instid(i)%ps(n)
88             if ( use_rttov_kmatrix ) then
89                allocate (con_vars(n) % t_jac(nchan,nlevels))
90                allocate (con_vars(n) % q_jac(nchan,nlevels))
91                allocate (con_vars(n) % ps_jac(nchan))
92                do k = 1, nlevels
93                   con_vars(n) % t_jac(:,k) = iv%instid(i)%t_jacobian(:,k,n)
94                   con_vars(n) % q_jac(:,k) = iv%instid(i)%q_jacobian(:,k,n)
95                end do
96                con_vars(n) % ps_jac(:)= iv%instid(i)%ps_jacobian(:,n)
97             end if
99             aux_vars(n) % t2m = iv%instid(i)%t2m(n)
100             aux_vars(n) % q2m = iv%instid(i)%mr2m(n)
101             aux_vars(n) % u10 = iv%instid(i)%u10(n)
102             aux_vars(n) % v10 = iv%instid(i)%v10(n)
104             aux_vars(n) % surftype = iv%instid(i)%surftype(n)
105             aux_vars(n) % surft    = iv%instid(i)%ts(n)
106             ! aux_vars(n) % fastem(:) = 0.0
108             aux_vars(n) % satzen  = iv%instid(i)%satzen(n)
109             aux_vars(n) % satazi  = iv%instid(i)%satazi(n)
110     
111             aux_vars(n) % solzen   = iv%instid(i)%solzen(n)
112             aux_vars(n) % solazi   = iv%instid(i)%solazi(n)
113             aux_vars(n) % elevation = iv%instid(i)%elevation(n)
114             aux_vars(n) % rlat     = iv%instid(i)%info%lat(1,n)
116          end do
117       end if
119       !-------------------------------------------------
120       ! [2.9] Distribute evenly across processors
121       !-------------------------------------------------
123       d_num_rad=num_tovs_after(i,myproc+1)
125       if (d_num_rad > 0) then
126          allocate (d_con_vars(d_num_rad))
127          allocate (d_con_vars_ad(d_num_rad))
128          allocate (d_aux_vars(d_num_rad))
129          allocate (d_tb(nchan,d_num_rad))
130          d_tb(:,:) = 0.0
131          d_con_vars(:) % nlevels = nlevels
132          do n = 1, d_num_rad
133             allocate (d_con_vars(n) % t(nlevels))
134             allocate (d_con_vars(n) % q(nlevels))
135             if ( use_rttov_kmatrix ) then
136                allocate (d_con_vars(n) % t_jac(nchan,nlevels))
137                allocate (d_con_vars(n) % q_jac(nchan,nlevels))
138                allocate (d_con_vars(n) % ps_jac(nchan))
139             end if
140             allocate (d_con_vars_ad(n) % t(nlevels))
141             allocate (d_con_vars_ad(n) % q(nlevels))
142          end do
144          ! Fill up with data that stays local
146          l_num_rad=min(num_rad,d_num_rad)
148          if (l_num_rad > 0) then
149             do n = 1, l_num_rad
150                d_con_vars(n) % t(:) = con_vars(n) % t(:) 
151                d_con_vars(n) % q(:) = con_vars(n) % q(:)
152                ! d_aux_vars(n) % fastem(:) = 0.0
153                if ( use_rttov_kmatrix ) then
154                   d_con_vars(n) % t_jac(:,:) = con_vars(n) % t_jac(:,:)
155                   d_con_vars(n) % q_jac(:,:) = con_vars(n) % q_jac(:,:)
156                   d_con_vars(n) % ps_jac(:)  = con_vars(n) % ps_jac(:)
157                end if
158             end do
159             d_con_vars(1:l_num_rad) % nlevels  = con_vars(1:l_num_rad) % nlevels
160             d_con_vars(1:l_num_rad) % ps       = con_vars(1:l_num_rad) % ps
161             d_aux_vars(1:l_num_rad) % t2m      = aux_vars(1:l_num_rad) % t2m
162             d_aux_vars(1:l_num_rad) % q2m      = aux_vars(1:l_num_rad) % q2m
163             d_aux_vars(1:l_num_rad) % u10      = aux_vars(1:l_num_rad) % u10
164             d_aux_vars(1:l_num_rad) % v10      = aux_vars(1:l_num_rad) % v10
165             d_aux_vars(1:l_num_rad) % surftype = aux_vars(1:l_num_rad) % surftype
166             d_aux_vars(1:l_num_rad) % surft    = aux_vars(1:l_num_rad) % surft
167             d_aux_vars(1:l_num_rad) % satzen   = aux_vars(1:l_num_rad) % satzen
168             d_aux_vars(1:l_num_rad) % satazi   = aux_vars(1:l_num_rad) % satazi
170             d_aux_vars(1:l_num_rad) % solzen   = aux_vars(1:l_num_rad) % solzen
171             d_aux_vars(1:l_num_rad) % solazi   = aux_vars(1:l_num_rad) % solazi
172             d_aux_vars(1:l_num_rad) % elevation   = aux_vars(1:l_num_rad) % elevation
173             d_aux_vars(1:l_num_rad) % rlat     = aux_vars(1:l_num_rad) % rlat
175             d_tb(:,1:l_num_rad) = jo_grad_y%instid(i)%tb(:,1:l_num_rad)
176         
177          end if
179          ! Get data from elsewhere
180 #ifdef DM_PARALLEL
181          do c=1,tovs_copy_count(i)
182             if (tovs_send_pe(i,c)==myproc) then
183                s=tovs_send_start(i,c)
184                len=tovs_send_count(i,c)
185                e=s+len-1
186                nl=len*nlevels
187                allocate(temp_t(nlevels,len))
188                do j=1,len
189                   temp_t(:,j)=con_vars(s+j-1) % t(:)
190                end do
191                call mpi_send( temp_t,nl, true_mpi_real, tovs_recv_pe(i,c), &
192                   c*11+1, comm, ierr)
194                allocate(temp_q(nlevels,len))
195                do j=1,len
196                   temp_q(:,j)=con_vars(s+j-1) % q(:)
197                end do
198                call mpi_send (temp_q,nl, true_mpi_real, tovs_recv_pe(i,c), &
199                   c*11+2, comm, ierr)
200     
201                allocate (temp_tb(nchan,len))
202                temp_tb(:,:)=jo_grad_y%instid(i)%tb(:,s:e)
203                call mpi_send (temp_tb,len*nchan,  &
204                   true_mpi_real, tovs_recv_pe(i,c), c*11+3, comm, ierr)
206                allocate (rtemp(len,12))
207                rtemp(:,1)= con_vars(s:e) % ps
208                rtemp(:,2)= aux_vars(s:e) % t2m
209                rtemp(:,3)= aux_vars(s:e) % q2m
210                rtemp(:,4)= aux_vars(s:e) % u10
211                rtemp(:,5)= aux_vars(s:e) % v10
212                rtemp(:,6)= aux_vars(s:e) % surft
213                rtemp(:,7)= aux_vars(s:e) % satzen
214                rtemp(:,8)= aux_vars(s:e) % satazi
215                rtemp(:,9)= aux_vars(s:e) % solzen
216                rtemp(:,10)= aux_vars(s:e) % solazi
217                rtemp(:,11)= aux_vars(s:e) % elevation
218                rtemp(:,12)= aux_vars(s:e) % rlat
219                call mpi_send(rtemp,len*12, true_mpi_real, tovs_recv_pe(i,c), &
220                   c*11+4, comm, ierr)
222                allocate (itemp(len,2))
223                itemp(:,1)= con_vars(s:e) % nlevels
224                itemp(:,2)= aux_vars(s:e) % surftype
225                call mpi_send(itemp,len*2, mpi_integer, tovs_recv_pe(i,c), &
226                   c*11+5, comm, ierr)
228                deallocate (temp_tb)
229                deallocate (temp_t)
230                deallocate (temp_q)
231                deallocate (itemp)
232                deallocate (rtemp)
234                if ( use_rttov_kmatrix ) then
235                   allocate(temp_t_jac(nchan,nlevels,len))
236                   do k = 1, nlevels
237                      do j=1,len
238                         temp_t_jac(:,k,j)=con_vars(s+j-1) % t_jac(:,k)
239                      end do
240                   end do
241                   call mpi_send( temp_t_jac, nl*nchan, true_mpi_real, tovs_recv_pe(i,c), &
242                      c*11+6, comm, ierr)
244                   allocate(temp_q_jac(nchan,nlevels,len))
245                   do k = 1, nlevels
246                      do j=1,len
247                         temp_q_jac(:,:,j)=con_vars(s+j-1) % q_jac(:,:)
248                      end do
249                   end do
250                   call mpi_send( temp_q_jac, nl*nchan, true_mpi_real, tovs_recv_pe(i,c), &
251                      c*11+7, comm, ierr)
253                   allocate(temp_ps_jac(nchan,len))
254                   do j=1,len
255                      temp_ps_jac(:,j)=con_vars(s+j-1) % ps_jac(:)
256                   end do
257                   call mpi_send( temp_ps_jac, len*nchan, true_mpi_real, tovs_recv_pe(i,c), &
258                      c*11+8, comm, ierr)
260                   deallocate (temp_t_jac)
261                   deallocate (temp_q_jac)
262                   deallocate (temp_ps_jac)
263                end if
265             end if
267             if (tovs_recv_pe(i,c)==myproc) then
268                s=tovs_recv_start(i,c)
269                len=tovs_send_count(i,c)
270                e=s+len-1
271                nl=len*nlevels
272                allocate(temp_t(nlevels,len))
273                call mpi_recv(temp_t,nl, true_mpi_real, tovs_send_pe(i,c), &
274                   c*11+1, comm, status, ierr)
275                do j=1,len
276                   d_con_vars(s+j-1) % t(:)=temp_t(:,j)
277                end do
279                allocate(temp_q(nlevels,len))
280                call mpi_recv(temp_q,nl, true_mpi_real, tovs_send_pe(i,c), &
281                   c*11+2, comm, status, ierr)
282                do j=1,len
283                   d_con_vars(s+j-1) % q(:)=temp_q(:,j)
284                end do
286                allocate (temp_tb(nchan,len))
287                call mpi_recv(temp_tb,len*nchan, true_mpi_real, &
288                   tovs_send_pe(i,c), c*11+3, comm, status, ierr)
289                d_tb(:,s:e)=temp_tb(:,:)
291                allocate (rtemp(len,12))
292                call mpi_recv(rtemp,len*12, true_mpi_real, &
293                   tovs_send_pe(i,c), c*11+4, comm, status, ierr)
294                d_con_vars(s:e) % ps       = rtemp(:,1)
295                d_aux_vars(s:e) % t2m      = rtemp(:,2)
296                d_aux_vars(s:e) % q2m      = rtemp(:,3)
297                d_aux_vars(s:e) % u10      = rtemp(:,4)
298                d_aux_vars(s:e) % v10      = rtemp(:,5)
299                d_aux_vars(s:e) % surft    = rtemp(:,6)
300                d_aux_vars(s:e) % satzen   = rtemp(:,7)
301                d_aux_vars(s:e) % satazi   = rtemp(:,8)
302                d_aux_vars(s:e) % solzen   = rtemp(:,9)
303                d_aux_vars(s:e) % solazi   = rtemp(:,10)
304                d_aux_vars(s:e) % elevation= rtemp(:,11)
305                d_aux_vars(s:e) % rlat     = rtemp(:,12)
307                allocate (itemp(len,2))
308                call mpi_recv(itemp,len*2, mpi_integer, &
309                   tovs_send_pe(i,c), c*11+5, comm, status, ierr)
310                d_con_vars(s:e) % nlevels  = itemp(:,1)
311                d_aux_vars(s:e) % surftype = itemp(:,2)
313                deallocate (temp_tb)
314                deallocate (temp_t)
315                deallocate (temp_q)
316                deallocate (itemp)
317                deallocate (rtemp)
319                if ( use_rttov_kmatrix ) then
320                   allocate(temp_t_jac(nchan,nlevels,len))
321                   call mpi_recv(temp_t_jac, nl*nchan, true_mpi_real, tovs_send_pe(i,c), &
322                      c*11+6, comm, status, ierr)
323                   do k = 1, nlevels
324                      do j=1,len
325                         d_con_vars(s+j-1) % t_jac(:,k)=temp_t_jac(:,k,j)
326                      end do
327                   end do
329                   allocate(temp_q_jac(nchan,nlevels,len))
330                   call mpi_recv(temp_q_jac, nl*nchan, true_mpi_real, tovs_send_pe(i,c), &
331                      c*11+7, comm, status, ierr)
332                   do k = 1, nlevels
333                      do j=1,len
334                         d_con_vars(s+j-1) % q_jac(:,k)=temp_q_jac(:,k,j)
335                      end do
336                   end do
338                   allocate (temp_ps_jac(nchan,len))
339                   call mpi_recv(temp_ps_jac, len*nchan, true_mpi_real, &
340                      tovs_send_pe(i,c), c*11+8, comm, status, ierr)
341                   do j=1,len
342                      d_con_vars(s+j-1) % ps_jac(:)=temp_ps_jac(:,j)
343                   end do
345                   deallocate (temp_t_jac)
346                   deallocate (temp_q_jac)
347                   deallocate (temp_ps_jac)
348                end if
349             end if
351          end do
352 #endif
354          if (d_num_rad > 0) then
355             if (tovs_batch) then
356                call da_rttov_ad (i, nchan, d_num_rad, d_con_vars, &
357                d_aux_vars, d_con_vars_ad, d_tb)
358             else
359                if (.not. use_rttov_kmatrix) then
360                   !$OMP PARALLEL DO  &
361                   !$OMP PRIVATE ( n )
362                    do n=1,d_num_rad
363                       call da_rttov_ad (i, nchan, 1, d_con_vars(n:n), &
364                          d_aux_vars(n:n), d_con_vars_ad(n:n), d_tb(:,n:n))
365                    end do
366                   !$OMP END PARALLEL DO
367                else
368                   do n=1,d_num_rad
369                      d_con_vars_ad(n) % t(:) = 0.0
370                      d_con_vars_ad(n) % q(:) = 0.0
371                      d_con_vars_ad(n) % ps   = 0.0
372                      do l = 1, nchan
373                         d_con_vars_ad(n) % ps = d_con_vars_ad(n) % ps + &
374                                          d_con_vars(n)%ps_jac(l) * d_tb(l,n)
375                         do k = 1,nlevels
376                            d_con_vars_ad(n) % t(k) = d_con_vars_ad(n) % t(k) + &
377                                          d_con_vars(n)%t_jac(l,k) * d_tb(l,n)
378                            d_con_vars_ad(n) % q(k) = d_con_vars_ad(n) % q(k) + &
379                                          d_con_vars(n)%q_jac(l,k) * d_tb(l,n)
380                         end do
381                      end do
382                   end do
383                 end if
384             end if
385          end if
387          ! Transfer data back
388        
389          if (l_num_rad > 0) then
390             do n = 1, l_num_rad
391                con_vars_ad(n) % t(:) = d_con_vars_ad(n) % t(:) 
392                con_vars_ad(n) % q(:) = d_con_vars_ad(n) % q(:)
393             end do
395             con_vars_ad(1:l_num_rad) % ps = d_con_vars_ad(1:l_num_rad) % ps
396          end if
397       end if ! d_num_rad > 0
399       ! Return the data to other processors. Note the meaning of send_pe and
400       ! recv_pe is swapped here
402 #ifdef DM_PARALLEL
403       do c=1,tovs_copy_count(i)
404          if (tovs_recv_pe(i,c)==myproc) then
405             s=tovs_recv_start(i,c)
406             len=tovs_send_count(i,c)
407             e=s+len-1
408             nl=len*nlevels
409             allocate(temp_t(nlevels,len))
410             do j=1,len
411                temp_t(:,j)=d_con_vars_ad(s+j-1) % t(:)
412             end do
413             call mpi_send (temp_t,nl, true_mpi_real, tovs_send_pe(i,c), &
414                c*11+9, comm, ierr)
416             allocate(temp_q(nlevels,len))
417             do j=1,len
418                temp_q(:,j)=d_con_vars_ad(s+j-1) % q(:)
419             end do
420             call mpi_send (temp_q,nl, true_mpi_real, tovs_send_pe(i,c), &
421                c*11+10, comm, ierr)
423             allocate(rtemp(len,1))
424             rtemp(:,1) = d_con_vars_ad(s:e) % ps
425             call mpi_send (rtemp,len, true_mpi_real, tovs_send_pe(i,c), &
426                c*11+11, comm, ierr)
428             deallocate (temp_t)
429             deallocate (temp_q)
430             deallocate (rtemp)
431          end if
433          if (tovs_send_pe(i,c)==myproc) then
434             s=tovs_send_start(i,c)
435             len=tovs_send_count(i,c)
436             e=s+len-1
437             nl=len*nlevels
438             allocate(temp_t(nlevels,len))
439             call mpi_recv (temp_t,nl, true_mpi_real, tovs_recv_pe(i,c), &
440                c*11+9, comm, status, ierr)
441             do j=1,len
442                con_vars_ad(s+j-1) % t(:)=temp_t(:,j)
443             end do
445             allocate(temp_q(nlevels,len))
446             call mpi_recv (temp_q,nl, true_mpi_real, tovs_recv_pe(i,c), &
447                c*11+10, comm, status, ierr)
448             do j=1,len
449                con_vars_ad(s+j-1) % q(:)=temp_q(:,j)
450             end do
452             allocate(rtemp(len,1))
453             call mpi_recv (rtemp,len, true_mpi_real, tovs_recv_pe(i,c), &
454                c*11+11, comm, status, ierr)
455             con_vars_ad(s:e) % ps=rtemp(:,1)
457             deallocate (temp_t)
458             deallocate (temp_q)
459             deallocate (rtemp)
460          end if
462       end do
463 #endif
464       if (d_num_rad > 0) then
465          do n=1,d_num_rad
466             deallocate (d_con_vars(n) % t)
467             deallocate (d_con_vars(n) % q)
468             deallocate (d_con_vars_ad(n) % t)
469             deallocate (d_con_vars_ad(n) % q)
470             if ( use_rttov_kmatrix ) then
471                deallocate (d_con_vars(n) % t_jac)
472                deallocate (d_con_vars(n) % q_jac)
473                deallocate (d_con_vars(n) % ps_jac)
474             end if
475          end do
477          deallocate (d_tb)
478          deallocate (d_con_vars)
479          deallocate (d_aux_vars)
480          deallocate (d_con_vars_ad)
481       end if
483       ! adjoint of convert to hPa
485       if (num_rad > 0) then
486          con_vars_ad(:)% ps = con_vars_ad(:)%ps * 0.01 
487       end if
489       do n=1,num_rad
490          ! 4.2 scale transform 
492          do k=1, nlevels
493             model_mr(k) = con_vars_ad(n) % q(k)
495             if (iv%instid(i)%info%zk(k,n) <= 0.0) then
496                con_vars_ad(n)%t(k)  = 0.0  !coefs(i) % ref_prfl_t(k,gas_id_watervapour)
497                con_vars_ad(n)%q(k)  = 0.0  !coefs(i) % ref_prfl_mr(k,gas_id_watervapour)
498             else
499                ! adjoint of q(kg/kg) to ppmv
501                con_vars_ad(n)%q(k) = model_mr(k) * q2ppmv
502             end if
503          end do
504       end do
506       if (num_rad > 0) then
507          allocate(temp_t(nlevels,num_rad))
508          do n=1,num_rad
509             temp_t(:,n) = con_vars_ad(n)% t(:)
510          end do
511          call da_interp_lin_3d_adj (jo_grad_x%t, iv%instid(i)%info, temp_t)
512          deallocate(temp_t)
514          allocate(temp_q(nlevels,num_rad))
515          do n=1,num_rad
516             temp_q(:,n) = con_vars_ad(n)% q(:)
517          end do
518          call da_interp_lin_3d_adj (jo_grad_x%q, iv%instid(i)%info, temp_q)
519          deallocate(temp_q)
521          allocate(temp_ps(num_rad))
522          do n=1,num_rad
523             temp_ps(n) = con_vars_ad(n)% ps
524          end do
525          call da_interp_lin_2d_adj (jo_grad_x% psfc, iv%instid(i)%info, 1, temp_ps)
526          deallocate(temp_ps)
527       end if
530       if (num_rad > 0) then
531          do n = 1, num_rad
532             deallocate (con_vars(n)    % t)
533             deallocate (con_vars(n)    % q)
534             deallocate (con_vars_ad(n) % t)
535             deallocate (con_vars_ad(n) % q)
536             if ( use_rttov_kmatrix ) then
537                deallocate (con_vars(n) % t_jac)
538                deallocate (con_vars(n) % q_jac)
539                deallocate (con_vars(n) % ps_jac)
540             end if
541          end do
543          deallocate (con_vars)
544          deallocate (aux_vars)
545          deallocate (con_vars_ad)
546          deallocate (model_mr)
547       end if
548    end do        ! end loop for sensor
550    if (trace_use) call da_trace_exit("da_transform_xtoy_rttov_adj")
551 #else
552     call da_error(__FILE__,__LINE__, &
553        (/"Must compile with $RTTOV option for radiances"/))
554 #endif
556 end subroutine da_transform_xtoy_rttov_adj