MODULE ad_htobs_mod ! !git $Id$ !svn $Id: ad_htobs.F 1151 2023-02-09 03:08:53Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group Andrew M. Moore ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This routine computes the adjoint observation operator, ! ! ! ! transpose(H) * X ! ! ! ! which loads the observation vector into the adjoint forcing arrays. ! ! ! ! The observation screening and quality control variable "ObsScale" ! ! is not modified in this routine. Their values are the ones set in ! ! obs_write.F during the running of the nonlinear model. ! ! ! !======================================================================= ! implicit none CONTAINS ! !*********************************************************************** SUBROUTINE ad_htobs (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_ocean ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! integer :: IminS, ImaxS, JminS, JmaxS integer :: LBi, UBi, LBj, UBj, LBij, UBij ! ! Set horizontal starting and ending indices for automatic private ! storage arrays. ! IminS=BOUNDS(ng)%Istr(tile)-3 ImaxS=BOUNDS(ng)%Iend(tile)+3 JminS=BOUNDS(ng)%Jstr(tile)-3 JmaxS=BOUNDS(ng)%Jend(tile)+3 ! ! Determine array lower and upper bounds in the I- and J-directions. ! LBi=BOUNDS(ng)%LBi(tile) UBi=BOUNDS(ng)%UBi(tile) LBj=BOUNDS(ng)%LBj(tile) UBj=BOUNDS(ng)%UBj(tile) ! ! Set array lower and upper bounds for MIN(I,J) directions and ! MAX(I,J) directions. ! LBij=BOUNDS(ng)%LBij UBij=BOUNDS(ng)%UBij ! CALL ad_htobs_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & GRID(ng) % rmask, & & GRID(ng) % umask, & & GRID(ng) % vmask, & & GRID(ng) % z_r, & & GRID(ng) % z_v, & & OCEAN(ng) % f_u, & & OCEAN(ng) % f_v, & & OCEAN(ng) % f_t, & & OCEAN(ng) % f_ubar, & & OCEAN(ng) % f_vbar, & & OCEAN(ng) % f_zeta) RETURN END SUBROUTINE ad_htobs ! !*********************************************************************** SUBROUTINE ad_htobs_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & rmask, umask, vmask, & & z_r, z_v, & & f_u, f_v, f_t, & & f_ubar, f_vbar, f_zeta) !*********************************************************************** ! USE mod_param USE mod_parallel USE mod_fourdvar USE mod_iounits USE mod_ncparam USE mod_scalars ! USE distribute_mod, ONLY : mp_collect USE mp_exchange_mod, ONLY : ad_mp_exchange2d USE mp_exchange_mod, ONLY : ad_mp_exchange3d USE ad_extract_obs_mod, ONLY : ad_extract_obs2d USE ad_extract_obs_mod, ONLY : ad_extract_obs3d ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS ! real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: z_r(LBi:,LBj:,:) real(r8), intent(inout) :: z_v(LBi:,LBj:,:) real(r8), intent(inout) :: f_u(LBi:,LBj:,:) real(r8), intent(inout) :: f_v(LBi:,LBj:,:) real(r8), intent(inout) :: f_t(LBi:,LBj:,:,:) real(r8), intent(inout) :: f_ubar(LBi:,LBj:) real(r8), intent(inout) :: f_vbar(LBi:,LBj:) real(r8), intent(inout) :: f_zeta(LBi:,LBj:) ! ! Local variable declarations. ! integer :: Mstr, Mend, ObsSum, ObsVoid integer :: Ncollect integer :: i, ie, iobs, is, j integer :: itrc, k real(r8) :: angle real(r8), parameter :: IniVal = 0.0_r8 real(r8) :: ad_uradial(Mobs), ad_vradial(Mobs) ! !----------------------------------------------------------------------- ! Set lower and upper tile bounds and staggered variables bounds for ! this horizontal domain partition. Notice that if tile=-1, it will ! set the values for the global grid. !----------------------------------------------------------------------- ! integer :: Istr, IstrB, IstrP, IstrR, IstrT, IstrM, IstrU integer :: Iend, IendB, IendP, IendR, IendT integer :: Jstr, JstrB, JstrP, JstrR, JstrT, JstrM, JstrV integer :: Jend, JendB, JendP, JendR, JendT integer :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1 integer :: Iendp1, Iendp2, Iendp2i, Iendp3 integer :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1 integer :: Jendp1, Jendp2, Jendp2i, Jendp3 ! Istr =BOUNDS(ng) % Istr (tile) IstrB =BOUNDS(ng) % IstrB (tile) IstrM =BOUNDS(ng) % IstrM (tile) IstrP =BOUNDS(ng) % IstrP (tile) IstrR =BOUNDS(ng) % IstrR (tile) IstrT =BOUNDS(ng) % IstrT (tile) IstrU =BOUNDS(ng) % IstrU (tile) Iend =BOUNDS(ng) % Iend (tile) IendB =BOUNDS(ng) % IendB (tile) IendP =BOUNDS(ng) % IendP (tile) IendR =BOUNDS(ng) % IendR (tile) IendT =BOUNDS(ng) % IendT (tile) Jstr =BOUNDS(ng) % Jstr (tile) JstrB =BOUNDS(ng) % JstrB (tile) JstrM =BOUNDS(ng) % JstrM (tile) JstrP =BOUNDS(ng) % JstrP (tile) JstrR =BOUNDS(ng) % JstrR (tile) JstrT =BOUNDS(ng) % JstrT (tile) JstrV =BOUNDS(ng) % JstrV (tile) Jend =BOUNDS(ng) % Jend (tile) JendB =BOUNDS(ng) % JendB (tile) JendP =BOUNDS(ng) % JendP (tile) JendR =BOUNDS(ng) % JendR (tile) JendT =BOUNDS(ng) % JendT (tile) ! Istrm3 =BOUNDS(ng) % Istrm3 (tile) ! Istr-3 Istrm2 =BOUNDS(ng) % Istrm2 (tile) ! Istr-2 Istrm1 =BOUNDS(ng) % Istrm1 (tile) ! Istr-1 IstrUm2=BOUNDS(ng) % IstrUm2(tile) ! IstrU-2 IstrUm1=BOUNDS(ng) % IstrUm1(tile) ! IstrU-1 Iendp1 =BOUNDS(ng) % Iendp1 (tile) ! Iend+1 Iendp2 =BOUNDS(ng) % Iendp2 (tile) ! Iend+2 Iendp2i=BOUNDS(ng) % Iendp2i(tile) ! Iend+2 interior Iendp3 =BOUNDS(ng) % Iendp3 (tile) ! Iend+3 Jstrm3 =BOUNDS(ng) % Jstrm3 (tile) ! Jstr-3 Jstrm2 =BOUNDS(ng) % Jstrm2 (tile) ! Jstr-2 Jstrm1 =BOUNDS(ng) % Jstrm1 (tile) ! Jstr-1 JstrVm2=BOUNDS(ng) % JstrVm2(tile) ! JstrV-2 JstrVm1=BOUNDS(ng) % JstrVm1(tile) ! JstrV-1 Jendp1 =BOUNDS(ng) % Jendp1 (tile) ! Jend+1 Jendp2 =BOUNDS(ng) % Jendp2 (tile) ! Jend+2 Jendp2i=BOUNDS(ng) % Jendp2i(tile) ! Jend+2 interior Jendp3 =BOUNDS(ng) % Jendp3 (tile) ! Jend+3 ! !----------------------------------------------------------------------- ! Compute model minus observations adjoint misfit forcing. The ! representer coefficients (or its approximation PSI) has been ! already loaded into vector ADmodVal in the conjugate gradient ! or read in. ! ! The adjoint of this operator is tricky in parallel (tile partitions) ! because we need to avoid adding observation contributions to ghost ! points. In the case that an observation is located between neighbor ! tiles, both tiles need to process it and the contribution to f_var ! (an adjoint variable) is only done in the tile where (i,j) or ! (i,j,k) is not a ghost point. ! ! Alternatively, only one tile process such observation and the ! ad_mp_exchange*d routine is used to add the contribution to the ! correct (i,j) or (i,j,k) point. This is the strategy used here. ! ! The processing flag used to reject (ObsVetting=0) or accept ! (ObsVetting=1) observations is computed here but it is never ! used. The observation screening and quality control variable ! (ObsScale) is only computed in routine obs_write. !----------------------------------------------------------------------- ! IF (ProcessObs(ng)) THEN ! ! Set starting and ending indices of representer coefficient vector to ! proccess. The adjoint forcing is only computed for current time ! survey observations. ! Mstr=NstrObs(ng) Mend=NendObs(ng) ! ! Initialize observation reject/accept processing flag. ! DO iobs=Mstr,Mend ObsVetting(iobs)=IniVal END DO ! ! Free-surface. ! DO i=LBi,UBi DO j=LBj,UBj f_zeta(i,j)=0.0_r8 END DO END DO IF (FOURDVAR(ng)%ObsCount(isFsur).gt.0) THEN CALL ad_extract_obs2d (ng, 0, Lm(ng)+1, 0, Mm(ng)+1, & & LBi, UBi, LBj, UBj, & & ObsState2Type(isFsur), & & Mobs, Mstr, Mend, & & rXmin(ng), rXmax(ng), & & rYmin(ng), rYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, & & f_zeta, & & rmask, & & ADmodVal) CALL ad_mp_exchange2d (ng, tile, iADM, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & f_zeta) END IF ! ! 2D u-momentum component. ! DO i=LBi,UBi DO j=LBj,UBj f_ubar(i,j)=0.0_r8 END DO END DO IF (FOURDVAR(ng)%ObsCount(isUbar).gt.0) THEN CALL ad_extract_obs2d (ng, 1, Lm(ng)+1, 0, Mm(ng)+1, & & LBi, UBi, LBj, UBj, & & ObsState2Type(isUbar), & & Mobs, Mstr, Mend, & & uXmin(ng), uXmax(ng), & & uYmin(ng), uYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, & & f_ubar, & & umask, & & ADmodVal) CALL ad_mp_exchange2d (ng, tile, iADM, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & f_ubar) END IF ! ! 2D v-momentum component. ! DO i=LBi,UBi DO j=LBj,UBj f_vbar(i,j)=0.0_r8 END DO END DO IF (FOURDVAR(ng)%ObsCount(isVbar).gt.0) THEN CALL ad_extract_obs2d (ng, 0, Lm(ng)+1, 1, Mm(ng)+1, & & LBi, UBi, LBj, UBj, & & ObsState2Type(isVbar), & & Mobs, Mstr, Mend, & & vXmin(ng), vXmax(ng), & & vYmin(ng), vYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, & & f_vbar, & & vmask, & & ADmodVal) CALL ad_mp_exchange2d (ng, tile, iADM, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & f_vbar) END IF ! ! 3D u-momentum component. ! DO k=1,N(ng) DO i=LBi,UBi DO j=LBj,UBj f_u(i,j,k)=0.0_r8 END DO END DO END DO IF (FOURDVAR(ng)%ObsCount(isUvel).gt.0) THEN DO k=1,N(ng) DO j=Jstr-1,Jend+1 DO i=IstrU-1,Iend+1 z_v(i,j,k)=0.5_r8*(z_r(i-1,j,k)+ & & z_r(i ,j,k)) END DO END DO END DO CALL ad_extract_obs3d (ng, 1, Lm(ng)+1, 0, Mm(ng)+1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ObsState2Type(isUvel), & & Mobs, Mstr, Mend, & & uXmin(ng), uXmax(ng), & & uYmin(ng), uYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, Zobs, & & f_u, z_v, & & umask, & & ADmodVal) CALL ad_mp_exchange3d (ng, tile, iADM, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & f_u) END IF ! ! 3D v-momentum component. ! DO k=1,N(ng) DO i=LBi,UBi DO j=LBj,UBj f_v(i,j,k)=0.0_r8 END DO END DO END DO IF (FOURDVAR(ng)%ObsCount(isVvel).gt.0) THEN DO k=1,N(ng) DO j=JstrV-1,Jend+1 DO i=Istr-1,Iend+1 z_v(i,j,k)=0.5_r8*(z_r(i,j-1,k)+ & & z_r(i,j ,k)) END DO END DO END DO CALL ad_extract_obs3d (ng, 0, Lm(ng)+1, 1, Mm(ng)+1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ObsState2Type(isVvel), & & Mobs, Mstr, Mend, & & vXmin(ng), vXmax(ng), & & vYmin(ng), vYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, Zobs, & & f_v, z_v, & & vmask, & & ADmodVal) CALL ad_mp_exchange3d (ng, tile, iADM, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & f_v) END IF ! ! Radial Velocity. The observations are in terms of radial speed and ! angle (stored in obs_meta). The observation angle converts the ! velocity components to geographical EAST and North components. ! By default, the radial velocity observations are processed as ! magnitude and heading angle (obs_meta; radians) in the navigation ! convention: an azimuth that is clockwise from TRUE North. ! ! In curvilinear coordinates, the radial forward problem is: ! ! radial = u * SIN(obs_meta + angler) + v * COS(obs_meta + angler) ! ! In the adjoint, u and v are given by: ! ! f_v = f_v + ADmodVal * COS(obs_meta + angler) ! f_u = f_u + ADmodVal * SIN(obs_meta + angler) ! IF (FOURDVAR(ng)%ObsCount(isRadial).gt.0) THEN DO iobs=Mstr,Mend ad_uradial(iobs)=IniVal ad_vradial(iobs)=IniVal END DO DO iobs=Mstr,Mend IF (ObsType(iobs).eq.ObsState2Type(isRadial)) THEN angle=ObsMeta(iobs)+ObsAngler(iobs) ad_uradial(iobs)=ad_uradial(iobs)+ & & ADmodVal(iobs)*SIN(angle) ad_vradial(iobs)=ad_vradial(iobs)+ & & ADmodVal(iobs)*COS(angle) ADmodVal(iobs)=0.0_r8 END IF END DO DO k=1,N(ng) DO j=JstrV-1,Jend+1 DO i=Istr-1,Iend+1 z_v(i,j,k)=0.5_r8*(z_r(i,j-1,k)+ & & z_r(i,j ,k)) END DO END DO END DO CALL ad_extract_obs3d (ng, 0, Lm(ng)+1, 1, Mm(ng)+1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ObsState2Type(isRadial), & & Mobs, Mstr, Mend, & & vXmin(ng), vXmax(ng), & & vYmin(ng), vYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, Zobs, & & f_v, z_v, & & vmask, & & ad_vradial) DO k=1,N(ng) DO j=Jstr-1,Jend+1 DO i=IstrU-1,Iend+1 z_v(i,j,k)=0.5_r8*(z_r(i-1,j,k)+ & & z_r(i ,j,k)) END DO END DO END DO CALL ad_extract_obs3d (ng, 1, Lm(ng)+1, 0, Mm(ng)+1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ObsState2Type(isRadial), & & Mobs, Mstr, Mend, & & uXmin(ng), uXmax(ng), & & uYmin(ng), uYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, Zobs, & & f_u, z_v, & & umask, & & ad_uradial) END IF ! ! Exchange adjoint velocites forcing terms, after all the 3D velocity ! observations are processed. ! IF ((FOURDVAR(ng)%ObsCount(isUvel).gt.0).or. & & (FOURDVAR(ng)%ObsCount(isVvel).gt.0).or. & & (FOURDVAR(ng)%ObsCount(isRadial).gt.0)) THEN CALL ad_mp_exchange3d (ng, tile, iADM, 2, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, EWperiodic(ng), & & NSperiodic(ng), f_u, f_v) END IF ! ! Tracer type variables. ! DO itrc=1,NT(ng) DO k=1,N(ng) DO i=LBi,UBi DO j=LBj,UBj f_t(i,j,k,itrc)=0.0_r8 END DO END DO END DO IF (FOURDVAR(ng)%ObsCount(isTvar(itrc)).gt.0) THEN CALL ad_extract_obs3d (ng, 0, Lm(ng)+1, 0, Mm(ng)+1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ObsState2Type(isTvar(itrc)), & & Mobs, Mstr, Mend, & & rXmin(ng), rXmax(ng), & & rYmin(ng), rYmax(ng), & & time(ng), dt(ng), & & ObsType, ObsVetting, & & Tobs, Xobs, Yobs, Zobs, & & f_t(:,:,:,itrc), z_r, & & rmask, & & ADmodVal) CALL ad_mp_exchange3d (ng, tile, iADM, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & f_t(:,:,:,itrc)) END IF END DO ! !----------------------------------------------------------------------- ! For debugging purposes, collect all observations reject/accept ! processing flag. !----------------------------------------------------------------------- ! Ncollect=Mend-Mstr+1 CALL mp_collect (ng, model, Ncollect, IniVal, & & ObsVetting(Mstr:)) ! !----------------------------------------------------------------------- ! Set counters for the number of rejected observations for each state ! variable. Although unnecessary, the counters are recomputed here to ! check if "ObsScale" changed from its initial values. !----------------------------------------------------------------------- ! DO iobs=Mstr,Mend IF (ObsScale(iobs).lt.1.0) THEN IF (ObsType(iobs).eq.ObsState2Type(isFsur)) THEN FOURDVAR(ng)%ObsReject(isFsur)= & & FOURDVAR(ng)%ObsReject(isFsur)+1 ELSE IF (ObsType(iobs).eq.ObsState2Type(isUbar)) THEN FOURDVAR(ng)%ObsReject(isUbar)= & & FOURDVAR(ng)%ObsReject(isUbar)+1 ELSE IF (ObsType(iobs).eq.ObsState2Type(isVbar)) THEN FOURDVAR(ng)%ObsReject(isVbar)= & & FOURDVAR(ng)%ObsReject(isVbar)+1 ELSE IF (ObsType(iobs).eq.ObsState2Type(isUvel)) THEN FOURDVAR(ng)%ObsReject(isUvel)= & & FOURDVAR(ng)%ObsReject(isUvel)+1 ELSE IF (ObsType(iobs).eq.ObsState2Type(isVvel)) THEN FOURDVAR(ng)%ObsReject(isVvel)= & & FOURDVAR(ng)%ObsReject(isVvel)+1 ELSE IF (ObsType(iobs).eq.ObsState2Type(isRadial)) THEN FOURDVAR(ng)%ObsReject(isRadial)= & & FOURDVAR(ng)%ObsReject(isRadial)+1 ELSE DO itrc=1,NT(ng) IF (ObsType(iobs).eq.ObsState2Type(isTvar(itrc))) THEN i=isTvar(itrc) FOURDVAR(ng)%ObsReject(i)=FOURDVAR(ng)%ObsReject(i)+1 END IF END DO END IF END IF END DO ! ! Load total available and rejected observations into structure ! array. ! DO i=1,NobsVar(ng) FOURDVAR(ng)%ObsCount(0)=FOURDVAR(ng)%ObsCount(0)+ & & FOURDVAR(ng)%ObsCount(i) FOURDVAR(ng)%ObsReject(0)=FOURDVAR(ng)%ObsReject(0)+ & & FOURDVAR(ng)%ObsReject(i) END DO ! ! Report. ! IF (DOMAIN(ng)%SouthWest_Test(tile)) THEN IF (Master) THEN ObsSum=0 ObsVoid=0 is=NstrObs(ng) DO i=1,NobsVar(ng) IF (FOURDVAR(ng)%ObsCount(i).gt.0) THEN ie=is+FOURDVAR(ng)%ObsCount(i)-1 WRITE (stdout,10) TRIM(ObsName(i)), is, ie, & & ie-is+1, FOURDVAR(ng)%ObsReject(i) is=ie+1 ObsSum=ObsSum+FOURDVAR(ng)%ObsCount(i) ObsVoid=ObsVoid+FOURDVAR(ng)%ObsReject(i) END IF END DO WRITE (stdout,20) ObsSum, ObsVoid, & & FOURDVAR(ng)%ObsCount(0), & & FOURDVAR(ng)%ObsReject(0) WRITE (stdout,30) time_code(ng), NstrObs(ng), NendObs(ng), & & iic(ng) 10 FORMAT (10x,a,t25,4(1x,i10)) 20 FORMAT (/,10x,'Total',t47,2(1x,i10), & & /,10x,'Obs Tally',t47,2(1x,i10),/) 30 FORMAT (3x,' AD_HTOBS - Computed adjoint observations ', & & 'forcing,',t68,a,/,19x,'(Observation ', & & 'records = ',i7.7,' - ',i7.7,', iic = ',i7.7,')') END IF END IF END IF RETURN END SUBROUTINE ad_htobs_tile END MODULE ad_htobs_mod