!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer,
! version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer
!version 3 or (at your option) any later version that should be found at:
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software;
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec),
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
!
subroutine tovs_rttov_ad 1,18
#if defined (doc)
!
!**s/r tovs_rttov_ad - Adjoint of computation of radiance with rttov_ad
! (adapted from part of code of avtov)
!
!
!author : j. halle *cmda/aes april 19, 2005
!
!revision 001 : a. beaulne *cmda/smc june 2006
! -addition of ozone and IR surface emissivities
!revision 002 : r. sarrazin cmda april 2008
! -adapt to CSR
!revision 003 : s. heilliette
! -adapt to IASI
! -------------------
! purpose:
!
!arguments
!
#endif
use mod_tovs
implicit none
!implicits
#include "pardim.cdk"
#include "comdimo.cdk"
#include "comlun.cdk"
#include "comcst.cdk"
#include "comoabdy.cdk"
#include "comoahdr.cdk"
#include "comoba.cdk"
#include "commvohr.cdk"
#include "partov.cdk"
#include "comtov.cdk"
#include "comtovst.cdk"
#include "cvcord.cdk"
real*8, parameter :: q_mixratio_to_ppmv = 1.60771704e+6
integer :: jpol, joff1, joff2, joff3, ilev
integer :: nbtout, no_id, isurface
integer :: nlevels
integer :: nfrequencies
integer :: nchannels
integer :: len_nchannels1
integer :: len_nfrequencies1
integer :: nfrequencies_max, nchannels_max, nbtout_max
integer :: alloc_status(60)
integer :: errorstatus(jppf) ! rttov error return code
integer :: omp_get_num_threads
integer :: nthreads, max_nthreads
integer :: j, i, krtid, io, jf, iobs, iobs1
integer :: ibegin, ibeginob, ilast, ilastob, jj
integer :: ival, iplatform, isat, knpf
integer :: jo, jdata, idata, idatend, idatyp
integer :: jk, jn, jl
integer :: instrum, ilen, istride, ith
integer :: sensor_type ! sensor type (1=infrared; 2=microwave; 3=high resolution)
integer :: knpf_tot, ioffset, nrank, ichn
integer, allocatable :: ncan (:) ! number of channels for each profile
integer, allocatable :: knpf_th (:)
integer, allocatable :: iptobs (:)
integer, allocatable :: nfrequencies_th (:)
integer, allocatable :: nchannels_th (:)
integer, allocatable :: nbtout_th (:)
integer, allocatable :: lprofiles_th (:,:)
integer, allocatable :: channels_th (:,:)
integer, allocatable :: polarisations_th(:,:)
integer, allocatable :: errorstatus_th (:,:)
integer, allocatable :: iptobs_th (:,:)
real*8, allocatable :: surfem1_th (:,:)
real*8, allocatable :: emissivity_th (:,:)
real*8, allocatable :: emissivity_ad_th (:,:)
type(radiance_type) , allocatable :: radiancedata_d_th (:) ! radiances full structure buffer used in rttov calls
type(radiance_type) , allocatable :: radiancedata_ad_th (:) ! ad radiances full structure buffer used in rttov calls
type(transmission_type), allocatable :: transmission_th (:) ! transmission
type(transmission_type), allocatable :: transmission_ad_th (:) ! transmission ad
type( profile_type ) , allocatable :: profilesdata_th (:,:) ! profiles buffer used in rttov calls
type( profile_type ) , allocatable :: profilesdata_ad_th (:,:) ! ad profiles buffer used in rttov calls
logical :: assim, lcloud, switchrad
logical, allocatable :: calcemis_th (:,:)
real*8, allocatable :: surfem1 (:)
if(nobtov.eq.0) return ! exit if there are not tovs data
! 1. Get number of threads available and allocate memory for some variables
! . ----------------------------------------------------------------------
!
100 continue
no_id = 1
lcloud =.false.
switchrad =.true.
!$omp parallel
max_nthreads = omp_get_num_threads()
!$omp end parallel
alloc_status(:) = 0
allocate ( knpf_th (max_nthreads) ,stat=alloc_status(1))
allocate ( iptobs_th (jppf,max_nthreads) ,stat=alloc_status(2))
allocate ( errorstatus_th (jppf,max_nthreads) ,stat=alloc_status(3))
allocate ( nfrequencies_th(max_nthreads) ,stat=alloc_status(4))
allocate ( nchannels_th (max_nthreads) ,stat=alloc_status(5))
allocate ( nbtout_th (max_nthreads) ,stat=alloc_status(6))
allocate ( iptobs (jppf*max_nthreads) ,stat=alloc_status(7))
if( any(alloc_status /= 0) ) then
write(nulout,*) ' tovs_rttov_ad: alloc_status = ', alloc_status(:)
write(nulout,'(" tovs_rttov_ad: arrays #1 memory allocation error")')
call abort3d
(nulout,'tovs_rttov_ad ')
end if
!
! 2. Computation of adjoint hx for tovs data only
! . --------------------------------------------
200 continue
! Loop over all sensors specified by user
do 295 krtid = 1, nsensors
alloc_status(:) = 0
sensor_type = coef(krtid) % id_sensor
iplatform = coef(krtid) % id_platform
isat = coef(krtid) % id_sat
instrum = coef(krtid) % id_inst
! loop over all obs.
knpf = 0
do 290 iobs = 1, nobtov
jo = lobsno(iobs)
! Currently processed sensor?
if ( lsensor(iobs) .eq. krtid ) then
knpf = knpf + 1
iptobs(knpf) = iobs
endif
if ( knpf .le. 0 ) go to 290
if ( knpf .ne. jppf*max_nthreads .and. iobs .ne. nobtov ) go to 290
! . 2.1 Separate profiles according to the number of threads and
! . calculate the actual number of threads which will be used.
! . ----------------------------------------------------------
ith = 0
knpf_th (:) = 0
iptobs_th(:,:) = 0
istride = min(knpf,jppf)
do jn = 1, knpf,istride
ilen = min (knpf-jn+1,istride)
ith = ith + 1
knpf_th(ith) = ilen
do i = 1, ilen
iptobs_th(i,ith) = iptobs(jn+i-1)
enddo
enddo
! set nthreads to actual number of threads which will be used.
nthreads = min(max_nthreads,ith)
! . 2.2 Prepare all input variables required by rttov_ad.
! for the purposes of using openmp; an additional dimension
! is added to most input variables, corresponding
! to the number of threads.
! . ---------------------------------------------------------
allocate (ncan(knpf),stat= alloc_status(1))
ncan(:) = nchan(krtid)
! compute max possible values for nfrequencies, nchannels and nbtout, using knpf, i.e. maximum number of profiles.
call rttov_setupchan(knpf,ncan,coef(krtid),nfrequencies_max,nchannels_max,nbtout_max)
allocate ( surfem1_th (nchannels_max,nthreads) ,stat=alloc_status(2))
allocate ( lprofiles_th (nfrequencies_max,nthreads) ,stat=alloc_status(3))
allocate ( emissivity_th (nchannels_max,nthreads) ,stat=alloc_status(4))
allocate ( emissivity_ad_th(nchannels_max,nthreads) ,stat=alloc_status(5))
allocate ( channels_th (nfrequencies_max,nthreads) ,stat=alloc_status(6))
allocate ( polarisations_th(nchannels_max*3,nthreads) ,stat=alloc_status(7))
allocate ( calcemis_th (nchannels_max,nthreads) ,stat=alloc_status(8))
allocate ( surfem1 (nchannels_max) ,stat=alloc_status(9))
! get AIRS & IASI ir emissivities
if (instrum.eq.11 .or. instrum.eq.16) then
surfem1(:) = 0.
do jn = 1, knpf
iobs1 = iptobs(jn)
jo = lobsno(iobs1)
idata = mobhdr(ncmrln,jo)
idatend = mobhdr(ncmnlv,jo) + idata - 1
do jdata = idata, idatend
if(mobdata(ncmass,jdata).eq.1) then
ichn = nint(robdata8(ncmppp,jdata))
ichn = max(0,min(ichn,jpch+1))
do nrank = 1, nchan(krtid)
if ( ichn == ichan(nrank,krtid) ) exit
end do
surfem1 ( nrank + (jn-1)*nchan(krtid) ) = robdata8(ncmprl,jdata)
end if
end do
end do
end if
! loop on threads
knpf_tot = 0
do ith = 1, nthreads
call rttov_setupchan(knpf_th(ith),ncan,coef(krtid),nfrequencies_th(ith), &
& nchannels_th(ith),nbtout_th(ith))
knpf_tot = knpf_tot + knpf_th(ith)
surfem1_th (:,ith) = 0.
do j = 1 , knpf_th(ith)
len_nchannels1 = nchannels_th(ith)/knpf_th(ith)
len_nfrequencies1 = nfrequencies_th(ith)/knpf_th(ith)
joff1=len_nchannels1*(j-1)
joff2=len_nfrequencies1*(j-1)
ioffset=len_nfrequencies1*(knpf_tot-knpf_th(ith))
isurface = profiles(iptobs_th(j,ith)) % skin % surftype
if (sensor_type .eq. 2) then
if ( isurface .eq. 0 .or. &
isurface .eq. 2 ) then
calcemis_th(joff1+1:joff1+len_nchannels1 ,ith) = .false.
surfem1_th (joff2+1:joff2+len_nfrequencies1,ith) = 0.75
else
calcemis_th(joff1+1:joff1+len_nchannels1 ,ith) = .true.
surfem1_th (joff2+1:joff2+len_nfrequencies1,ith) = 0.
endif
elseif ((instrum .eq. 11) .or. (instrum .eq. 16) ) then
calcemis_th(joff1+1:joff1+len_nchannels1 ,ith) = .false.
surfem1_th (joff2+1:joff2+len_nfrequencies1,ith) = &
& surfem1(joff2+1+ioffset:joff2+len_nfrequencies1+ioffset)
elseif ((instrum .eq. 20) .or. (instrum .eq. 21) .or. &
& (instrum .eq. 22) .or. (instrum .eq. 24) ) then
calcemis_th(joff1+1:joff1+len_nchannels1 ,ith) = .true.
surfem1_th (joff2+1:joff2+len_nfrequencies1,ith) = 0.
else
call abort3d
(nulout,'tovs_rttov_ad. invalid sensor type')
endif
enddo
call rttov_setupindex (ncan,knpf_th(ith),nfrequencies_th(ith),nchannels_th(ith),&
& nbtout_th(ith),coef(krtid),surfem1_th(:,ith),lprofiles_th (:,ith),&
& channels_th(:,ith),polarisations_th(:,ith),emissivity_th(:,ith))
enddo
allocate ( radiancedata_d_th(nthreads) ,stat=alloc_status(10) ) ! radiances full structure buffer used in rttov calls
allocate ( radiancedata_ad_th(nthreads) ,stat=alloc_status(11)) ! ad radiances full structure buffer used in rttov calls
allocate ( transmission_th (nthreads) ,stat=alloc_status(12)) ! transmission
allocate ( transmission_ad_th (nthreads) ,stat=alloc_status(13)) ! transmission ad
allocate ( profilesdata_th(istride,nthreads) ,stat=alloc_status(14)) ! profilesdata
allocate ( profilesdata_ad_th(istride,nthreads),stat=alloc_status(15)) ! profilesdata ad
do ith = 1, nthreads
! allocate transmittance structures
call tovs_allocate_transmission
(transmission_th(ith) ,nchannels_th(ith),coef(krtid) % nlevels)
call tovs_allocate_transmission
(transmission_ad_th(ith),nchannels_th(ith),coef(krtid) % nlevels)
! allocate radiance structures
call tovs_allocate_radiance
(radiancedata_d_th(ith) ,nchannels_th(ith),coef(krtid) % nlevels,nbtout_th(ith))
call tovs_allocate_radiance
(radiancedata_ad_th(ith),nchannels_th(ith),coef(krtid) % nlevels,nbtout_th(ith))
! allocate profile structures
do j = 1, knpf_th(ith)
call tovs_allocate_profile
(profilesdata_th(j,ith) ,coef(no_id) % nlevels)
call tovs_allocate_profile
(profilesdata_ad_th(j,ith),coef(no_id) % nlevels)
end do
if( any(alloc_status /= 0) ) then
write(nulout,*) ' tovs_rttov_ad: alloc_status = ', alloc_status(:)
write(nulout,'(" tovs_rttov_ad: arrays #2 memory allocation error")')
call abort3d
(nulout,'tovs_rttov_ad ')
end if
!.. fill profilesdata arrays
do j = 1 , knpf_th(ith)
profilesdata_th(j,ith) % nlevels = profiles(iptobs_th(j,ith)) % nlevels
profilesdata_th(j,ith) % ozone_data = profiles(iptobs_th(j,ith)) % ozone_data
profilesdata_th(j,ith) % co2_data = profiles(iptobs_th(j,ith)) % co2_data
profilesdata_th(j,ith) % clw_data = profiles(iptobs_th(j,ith)) % clw_data
profilesdata_th(j,ith) % p(:) = profiles(iptobs_th(j,ith)) % p(:)
profilesdata_th(j,ith) % t(:) = profiles(iptobs_th(j,ith)) % t(:)
profilesdata_th(j,ith) % q(:) = profiles(iptobs_th(j,ith)) % q(:)
profilesdata_th(j,ith) % o3(:) = profiles(iptobs_th(j,ith)) % o3(:)
profilesdata_th(j,ith) % ctp = profiles(iptobs_th(j,ith)) % ctp
profilesdata_th(j,ith) % cfraction = profiles(iptobs_th(j,ith)) % cfraction
profilesdata_th(j,ith) % ozone_data = profiles(iptobs_th(j,ith)) % ozone_data
profilesdata_th(j,ith) % co2_data = profiles(iptobs_th(j,ith)) % co2_data
profilesdata_th(j,ith) % clw_data = profiles(iptobs_th(j,ith)) % clw_data
profilesdata_th(j,ith) % zenangle = profiles(iptobs_th(j,ith)) % zenangle
profilesdata_th(j,ith) % azangle = profiles(iptobs_th(j,ith)) % azangle
profilesdata_th(j,ith) % skin % surftype = profiles(iptobs_th(j,ith)) % skin % surftype
profilesdata_th(j,ith) % skin % t = profiles(iptobs_th(j,ith)) % skin % t
profilesdata_th(j,ith) % skin % fastem(:) = profiles(iptobs_th(j,ith)) % skin % fastem(:)
profilesdata_th(j,ith) % s2m % t = profiles(iptobs_th(j,ith)) % s2m % t
profilesdata_th(j,ith) % s2m % q = profiles(iptobs_th(j,ith)) % s2m % q
profilesdata_th(j,ith) % s2m % p = profiles(iptobs_th(j,ith)) % s2m % p
profilesdata_th(j,ith) % s2m % u = profiles(iptobs_th(j,ith)) % s2m % u
profilesdata_th(j,ith) % s2m % v = profiles(iptobs_th(j,ith)) % s2m % v
end do
!.. climatological moisture clip for profilesdata and profilesdata_tl arrays
do j = 1 , knpf_th(ith)
do jl = 1 , coef(no_id) % nlevels
if ( profiles(iptobs_th(j,ith)) % q(jl) &
.le. oqmin(jl)*q_mixratio_to_ppmv ) then
profilesdata_th(j,ith) % q(jl) = oqmin(jl)*q_mixratio_to_ppmv
elseif ( profiles(iptobs_th(j,ith)) % q(jl) &
.ge. oqmax(jl)*q_mixratio_to_ppmv ) then
profilesdata_th(j,ith) % q(jl) = oqmax(jl)*q_mixratio_to_ppmv
endif
enddo
end do
!.. fill profilesdata_ad arrays
do j = 1 , knpf_th(ith)
profilesdata_ad_th(j,ith) % nlevels = profiles(iptobs_th(j,ith)) % nlevels
profilesdata_ad_th(j,ith) % ozone_data = profiles(iptobs_th(j,ith)) % ozone_data
profilesdata_ad_th(j,ith) % co2_data = profiles(iptobs_th(j,ith)) % co2_data
profilesdata_ad_th(j,ith) % clw_data = profiles(iptobs_th(j,ith)) % clw_data
profilesdata_ad_th(j,ith) % p(:) = 0.
profilesdata_ad_th(j,ith) % t(:) = 0.
profilesdata_ad_th(j,ith) % q(:) = 0.
profilesdata_ad_th(j,ith) % o3(:) = 0.
profilesdata_ad_th(j,ith) % ctp = 0.
profilesdata_ad_th(j,ith) % cfraction = 0.
profilesdata_ad_th(j,ith) % zenangle = profiles(iptobs_th(j,ith)) % zenangle
profilesdata_ad_th(j,ith) % azangle = profiles(iptobs_th(j,ith)) % azangle
profilesdata_ad_th(j,ith) % skin % surftype = profiles(iptobs_th(j,ith)) % skin % surftype
profilesdata_ad_th(j,ith) % skin % t = 0.
profilesdata_ad_th(j,ith) % skin % fastem(:) = 0.
profilesdata_ad_th(j,ith) % s2m % t = 0.
profilesdata_ad_th(j,ith) % s2m % q = 0.
profilesdata_ad_th(j,ith) % s2m % p = 0.
profilesdata_ad_th(j,ith) % s2m % u = 0.
profilesdata_ad_th(j,ith) % s2m % v = 0.
end do
radiancedata_d_th(ith) % clear (:) = 0.
radiancedata_d_th(ith) % cloudy (:) = 0.
radiancedata_d_th(ith) % total (:) = 0.
radiancedata_d_th(ith) % bt_clear (:) = 0.
radiancedata_d_th(ith) % upclear (:) = 0.
radiancedata_d_th(ith) % dnclear (:) = 0.
radiancedata_d_th(ith) % reflclear(:) = 0.
radiancedata_d_th(ith) % overcast (:,:) = 0.
radiancedata_d_th(ith) % downcld (:,:) = 0.
radiancedata_d_th(ith) % out (:) = 0.
radiancedata_d_th(ith) % bt (:) = 0.
radiancedata_d_th(ith) % out_clear(:) = 0.
radiancedata_d_th(ith) % total_out(:) = 0.
radiancedata_d_th(ith) % clear_out(:) = 0.
!.. fill radiancedata_ad_th arrays
do jn = 1, knpf_th(ith)
io = iptobs_th(jn,ith)
joff1=nbtout_th(ith)/knpf_th(ith)*(jn-1)
radiancedata_ad_th(ith) % out(joff1+1:joff1+nbtout_th(ith)/knpf_th(ith)) = radiance_ad(io) % out(:)
enddo
radiancedata_ad_th(ith) % clear (:) = 0.
radiancedata_ad_th(ith) % cloudy (:) = 0.
radiancedata_ad_th(ith) % total (:) = 0.
radiancedata_ad_th(ith) % bt_clear (:) = 0.
radiancedata_ad_th(ith) % upclear (:) = 0.
radiancedata_ad_th(ith) % dnclear (:) = 0.
radiancedata_ad_th(ith) % reflclear(:) = 0.
radiancedata_ad_th(ith) % overcast (:,:) = 0.
radiancedata_ad_th(ith) % downcld (:,:) = 0.
radiancedata_ad_th(ith) % bt (:) = 0.
radiancedata_ad_th(ith) % out_clear(:) = 0.
radiancedata_ad_th(ith) % total_out(:) = 0.
radiancedata_ad_th(ith) % clear_out(:) = 0.
!.. fill transmission_ad_th arrays
transmission_ad_th(ith) % tau_surf (:) = 0.
transmission_ad_th(ith) % tau_layer (:,:) = 0.
transmission_ad_th(ith) % od_singlelayer(:,:) = 0.
!.. fill transmission_th arrays
transmission_th(ith) % tau_surf (:) = 0.
transmission_th(ith) % tau_layer (:,:) = 0.
transmission_th(ith) % od_singlelayer(:,:) = 0.
enddo
! . 2.3 Compute ad radiance with rttov_ad
! . ---------------------------------
errorstatus_th (:,:) = 0
emissivity_ad_th(:,:) = 0.0
!$omp parallel
!$omp do private(ith)
do ith = 1, nthreads
call rttov_ad( &
errorstatus_th(:,ith), & ! out
nfrequencies_th(ith), & ! in
nchannels_th(ith), & ! in
nbtout_th(ith), & ! in
knpf_th(ith), & ! in
channels_th(:,ith), & ! in
polarisations_th(:,ith), & ! in
lprofiles_th(:,ith), & ! in
profilesdata_th(:,ith), & ! in
coef(krtid), & ! in
lcloud, & ! in
switchrad, & ! in
calcemis_th(:,ith), & ! in
emissivity_th(:,ith), & ! inout
profilesdata_ad_th(:,ith), & ! in
emissivity_ad_th(:,ith), & ! inout
transmission_th(ith), & ! inout
transmission_ad_th(ith), & ! inout
radiancedata_d_th(ith), & ! inout
radiancedata_ad_th(ith) ) ! inout
enddo
!$omp end do
!$omp end parallel
!.. store results from rttov_ad into profiles_ad
do ith = 1, nthreads
do j = 1 , knpf_th(ith)
io = iptobs_th(j,ith)
profiles_ad(io) % t(:) = profilesdata_ad_th(j,ith) % t(:)
profiles_ad(io) % q(:) = profilesdata_ad_th(j,ith) % q(:)
profiles_ad(io) % skin % t = profilesdata_ad_th(j,ith) % skin % t
profiles_ad(io) % s2m % t = profilesdata_ad_th(j,ith) % s2m % t
profiles_ad(io) % s2m % q = profilesdata_ad_th(j,ith) % s2m % q
profiles_ad(io) % s2m % p = profilesdata_ad_th(j,ith) % s2m % p
profiles_ad(io) % s2m % u = profilesdata_ad_th(j,ith) % s2m % u
profiles_ad(io) % s2m % v = profilesdata_ad_th(j,ith) % s2m % v
end do
enddo
!.. adjoint of climatological moisture clip for profilesdata arrays
do ith = 1, nthreads
do j = 1 , knpf_th(ith)
io = iptobs_th(j,ith)
do jl = 1 , coef(no_id) % nlevels
if ( profiles(iptobs_th(j,ith)) % q(jl) &
.le. oqmin(jl)*q_mixratio_to_ppmv ) then
profiles_ad(io) % q(jl) = 0.0
elseif ( profiles(iptobs_th(j,ith)) % q(jl) &
.ge. oqmax(jl)*q_mixratio_to_ppmv ) then
profiles_ad(io) % q(jl) = 0.0
endif
enddo
end do
enddo
! de-allocate memory
alloc_status(:) = 0
do ith = 1, nthreads
call tovs_deallocate_transmission
(transmission_th(ith))
call tovs_deallocate_transmission
(transmission_ad_th(ith))
call tovs_deallocate_radiance
(radiancedata_d_th(ith))
call tovs_deallocate_radiance
(radiancedata_ad_th(ith))
do jn = 1, knpf_th(ith)
call tovs_deallocate_profile
(profilesdata_th(jn,ith))
call tovs_deallocate_profile
(profilesdata_ad_th(jn,ith))
enddo
enddo
deallocate ( profilesdata_th ,stat=alloc_status(1) )
deallocate ( profilesdata_ad_th,stat=alloc_status(2) )
deallocate ( radiancedata_d_th ,stat=alloc_status(3) ) ! radiances full structure buffer used in rttov calls
deallocate ( radiancedata_ad_th,stat=alloc_status(4) ) ! ad radiances full structure buffer used in rttov calls
deallocate ( transmission_th ,stat=alloc_status(5) ) ! transmission
deallocate ( transmission_ad_th,stat=alloc_status(6) ) ! transmission ad
deallocate ( ncan ,stat=alloc_status(7) )
deallocate ( surfem1_th ,stat=alloc_status(8) )
deallocate ( lprofiles_th ,stat=alloc_status(9) )
deallocate ( emissivity_th ,stat=alloc_status(10))
deallocate ( emissivity_ad_th ,stat=alloc_status(11))
deallocate ( channels_th ,stat=alloc_status(12))
deallocate ( polarisations_th ,stat=alloc_status(13))
deallocate ( calcemis_th ,stat=alloc_status(14))
deallocate ( surfem1 ,stat=alloc_status(15))
if( any(alloc_status /= 0) ) then
write(nulout,*) ' tovs_rttov_ad: alloc_status = ', alloc_status(:)
write(nulout,'(" tovs_rttov_ad: arrays #2 memory deallocation error")')
call abort3d
(nulout,'tovs_rttov_ad ')
end if
! next bunch !
knpf = 0
290 continue
295 continue
! 3. Close up
! . --------
! deallocate memory
alloc_status(:) = 0
deallocate ( knpf_th ,stat=alloc_status(1))
deallocate ( iptobs_th ,stat=alloc_status(2))
deallocate ( errorstatus_th ,stat=alloc_status(3))
deallocate ( nfrequencies_th,stat=alloc_status(4))
deallocate ( nchannels_th ,stat=alloc_status(5))
deallocate ( nbtout_th ,stat=alloc_status(6))
deallocate ( iptobs ,stat=alloc_status(7))
if( any(alloc_status /= 0) ) then
write(nulout,*) ' tovs_rttov_ad: alloc_status = ', alloc_status(:)
write(nulout,'(" tovs_rttov_ad: arrays #1 memory deallocation error")')
call abort3d
(nulout,'tovs_rttov_ad ')
end if
return
end subroutine tovs_rttov_ad