Opened 4 years ago

Closed 4 years ago

#848 closed bug (Fixed)

VERY IMPORTANT: Corrected bugs in tracer advection

Reported by: arango Owned by:
Priority: major Milestone: Release ROMS/TOMS 3.9
Component: Nonlinear Version: 3.9
Keywords: Cc:

Description

Corrected critical bug in some of the tracer advection options:

  • Corrected a bug in the HSIMT advection. Previously, I fixed the issue that three ghost points are needed for this scheme. At that time, I missed making a similar change in the computation of the inverse thickness (oHz). Therefore, in step3d_t.F we need to have
          logical :: LapplySrc, Lhsimt, Lmpdata
          ...
          
          Lhsimt =ANY(Hadvection(:,ng)%HSIMT).and.                          &
         &        ANY(Vadvection(:,ng)%HSIMT)
          Lmpdata=ANY(Hadvection(:,ng)%MPDATA).and.                         &
         &        ANY(Vadvection(:,ng)%MPDATA)
          ...
          
    !
    !  Compute reciprocal thickness, 1/Hz.
    !
          IF (Lmpdata.or.Lhsimt) THEN
            DO k=1,N(ng)
              DO j=Jstrm2,Jendp2
                DO i=Istrm2,Iendp2
                  oHz(i,j,k)=1.0_r8/Hz(i,j,k)
                END DO
              END DO
            END DO
          ELSE
            DO k=1,N(ng)
              DO j=Jstr,Jend
                DO i=Istr,Iend
                  oHz(i,j,k)=1.0_r8/Hz(i,j,k)
                END DO
              END DO
            END DO
          END IF
    
    I have been hunting on and off for a bug in the HSIMT scheme for awhile. Many thanks to Pierre St-Laurent for reporting this issue. It is always helpful when another set of eyes look into a code.
  • Corrected a weird bug in the tangent linear advection switches. In inp_par.F we need to have instead:
    #if defined TANGENT || defined TL_IOMS
    !
    !  Set tracer advection scheme switches for the tangent linear models
    !  (TLM and RPM) to the same values as the adjoint model.
    !
          DO ng=1,Ngrids
            DO i=1,NT(ng)
              tl_Hadvection(i,ng)%AKIMA4    = ad_Hadvection(i,ng)%AKIMA4
              tl_Hadvection(i,ng)%CENTERED2 = ad_Hadvection(i,ng)%CENTERED2
              tl_Hadvection(i,ng)%CENTERED4 = ad_Hadvection(i,ng)%CENTERED4
              tl_Hadvection(i,ng)%HSIMT     = ad_Hadvection(i,ng)%HSIMT
              tl_Hadvection(i,ng)%MPDATA    = ad_Hadvection(i,ng)%MPDATA
              tl_Hadvection(i,ng)%SPLINES   = ad_Hadvection(i,ng)%SPLINES
              tl_Hadvection(i,ng)%SPLIT_U3  = ad_Hadvection(i,ng)%SPLIT_U3
              tl_Hadvection(i,ng)%UPSTREAM3 = ad_Hadvection(i,ng)%UPSTREAM3
    !
              tl_Vadvection(i,ng)%AKIMA4    = ad_Vadvection(i,ng)%AKIMA4
              tl_Vadvection(i,ng)%CENTERED2 = ad_Vadvection(i,ng)%CENTERED2
              tl_Vadvection(i,ng)%CENTERED4 = ad_Vadvection(i,ng)%CENTERED4
              tl_Vadvection(i,ng)%HSIMT     = ad_Vadvection(i,ng)%HSIMT
              tl_Vadvection(i,ng)%MPDATA    = ad_Vadvection(i,ng)%MPDATA
              tl_Vadvection(i,ng)%SPLINES   = ad_Vadvection(i,ng)%SPLINES
              tl_Vadvection(i,ng)%SPLIT_U3  = ad_Vadvection(i,ng)%SPLIT_U3
              tl_Vadvection(i,ng)%UPSTREAM3 = ad_Vadvection(i,ng)%UPSTREAM3
            END DO
          END DO
    #endif
    
    Before, the assigment was to use the nonlinear values instead. As a consequence, the 4D-Var cost function was not decreasing correctly and the minimization failed because the broken symmetry between tangent linear (TLM) and adjoint (ADM) model operators. We need to have the same advection scheme in both TLM and ADM. Recall that in roms.in, the user specify just the values for the adjoint model.

Many thanks to Julia Levin for reporting this problem and Andy Moore for his help in debugging.

Change History (1)

comment:1 by arango, 4 years ago

Resolution: Fixed
Status: newclosed
Note: See TracTickets for help on using tickets.