New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 15768 – NEMO

Changeset 15768


Ignore:
Timestamp:
2022-03-30T09:04:47+02:00 (2 years ago)
Author:
smasson
Message:

fix Incomplete initialisation in diahth.F90, see #2758

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/releases/r4.0/r4.0-HEAD/src/OCE/DIA/diahth.F90

    r15231 r15768  
    103103 
    104104      IF( kt == nit000 ) THEN 
    105          l_hth = .FALSE. 
    106          IF(   iom_use( 'mlddzt'   ) .OR. iom_use( 'mldr0_3'  ) .OR. iom_use( 'mldr0_1'  )    .OR.  &  
    107             &  iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' )    .OR.  &     
    108             &  iom_use( '20d'      ) .OR. iom_use( '26d'      ) .OR. iom_use( '28d'      )    .OR.  &     
    109             &  iom_use( 'hc300'    ) .OR. iom_use( 'hc700'    ) .OR. iom_use( 'hc2000'   )    .OR.  &     
    110             &  iom_use( 'pycndep'  ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    )  ) l_hth = .TRUE. 
     105         ! 
     106         l_hth = iom_use( 'mlddzt'   ) .OR. iom_use( 'mldr0_3'  ) .OR. iom_use( 'mldr0_1'  )    .OR.  &  
     107            &    iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' )    .OR.  &     
     108            &    iom_use( '20d'      ) .OR. iom_use( '26d'      ) .OR. iom_use( '28d'      )    .OR.  &     
     109            &    iom_use( 'hc300'    ) .OR. iom_use( 'hc700'    ) .OR. iom_use( 'hc2000'   )    .OR.  &     
     110            &    iom_use( 'pycndep'  ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    ) 
     111         ! 
    111112         !                                      ! allocate dia_hth array 
    112113         IF( l_hth ) THEN  
     
    121122      IF( l_hth ) THEN 
    122123         ! 
    123          IF( iom_use( 'mlddzt' ) .OR. iom_use( 'mldr0_3' ) .OR. iom_use( 'mldr0_1' ) ) THEN 
    124             ! initialization 
    125             ztinv  (:,:) = 0._wp   
    126             zdepinv(:,:) = 0._wp   
    127             zmaxdzT(:,:) = 0._wp   
     124         ! initialization 
     125         IF( iom_use( 'tinv'   ) )   ztinv  (:,:) = 0._wp   
     126         IF( iom_use( 'depti'  ) )   zdepinv(:,:) = 0._wp   
     127         IF( iom_use( 'mlddzt' ) )   zmaxdzT(:,:) = 0._wp   
     128         IF( iom_use( 'mlddzt' ) .OR. iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' )   & 
     129            &                    .OR. iom_use( 'mldr10_3' ) .OR. iom_use( 'pycndep'  ) ) THEN 
    128130            DO jj = 1, jpj 
    129131               DO ji = 1, jpi 
     
    134136                  zrho10_3(ji,jj) = zztmp 
    135137                  zpycn   (ji,jj) = zztmp 
    136                  END DO 
    137             END DO 
     138               END DO 
     139            END DO 
     140         ENDIF 
     141         IF( iom_use( 'mldr0_3' ) .OR. iom_use( 'mldr0_1' ) ) THEN 
    138142            IF( nla10 > 1 ) THEN  
    139143               DO jj = 1, jpj 
     
    145149               END DO 
    146150            ENDIF 
    147        
     151         ENDIF 
     152 
     153         IF( iom_use( 'mlddzt' ) .OR. iom_use( 'mldr0_3' ) .OR. iom_use( 'mldr0_1' ) ) THEN 
     154            ! ------------------------------------------------------------- ! 
     155            ! thermocline depth: strongest vertical gradient of temperature ! 
     156            ! turbocline depth (mixing layer depth): avt = zavt5            ! 
     157            ! MLD: rho = rho(1) + zrho3                                     ! 
     158            ! MLD: rho = rho(1) + zrho1                                     ! 
     159            ! ------------------------------------------------------------- ! 
     160            DO jk = jpkm1, 2, -1   ! loop from bottom to 2 
     161               DO jj = 1, jpj 
     162                  DO ji = 1, jpi 
     163                     ! 
     164                     zzdep = gdepw_n(ji,jj,jk) 
     165                     zztmp = ( tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) ) & 
     166                            & / zzdep * tmask(ji,jj,jk)   ! vertical gradient of temperature (dT/dz) 
     167                     zzdep = zzdep * tmask(ji,jj,1) 
     168 
     169                     IF( zztmp > zmaxdzT(ji,jj) ) THEN                         
     170                         zmaxdzT(ji,jj) = zztmp    
     171                         hth    (ji,jj) = zzdep                ! max and depth of dT/dz 
     172                     ENDIF 
     173                
     174                     IF( nla10 > 1 ) THEN  
     175                        zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1)                       ! delta rho(1) 
     176                        IF( zztmp > zrho3 )   zrho0_3(ji,jj) = zzdep                ! > 0.03 
     177                        IF( zztmp > zrho1 )   zrho0_1(ji,jj) = zzdep                ! > 0.01 
     178                     ENDIF 
     179                  END DO 
     180               END DO 
     181            END DO 
     182          
     183            CALL iom_put( 'mlddzt', hth )            ! depth of the thermocline 
     184            IF( nla10 > 1 ) THEN  
     185               CALL iom_put( 'mldr0_3', zrho0_3 )   ! MLD delta rho(surf) = 0.03 
     186               CALL iom_put( 'mldr0_1', zrho0_1 )   ! MLD delta rho(surf) = 0.01 
     187            ENDIF 
     188            ! 
     189         ENDIF 
     190         ! 
     191         IF(  iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' ) .OR.  &     
     192            &  iom_use( 'pycndep' ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    )  ) THEN 
     193      
    148194            ! Preliminary computation 
    149195            ! computation of zdelr = (dr/dT)(T,S,10m)*(-0.2 degC) 
     
    165211               END DO 
    166212            END DO 
    167  
    168             ! ------------------------------------------------------------- ! 
    169             ! thermocline depth: strongest vertical gradient of temperature ! 
    170             ! turbocline depth (mixing layer depth): avt = zavt5            ! 
    171             ! MLD: rho = rho(1) + zrho3                                     ! 
    172             ! MLD: rho = rho(1) + zrho1                                     ! 
    173             ! ------------------------------------------------------------- ! 
    174             DO jk = jpkm1, 2, -1   ! loop from bottom to 2 
    175                DO jj = 1, jpj 
    176                   DO ji = 1, jpi 
    177                      ! 
    178                      zzdep = gdepw_n(ji,jj,jk) 
    179                      zztmp = ( tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) ) & 
    180                             & / zzdep * tmask(ji,jj,jk)   ! vertical gradient of temperature (dT/dz) 
    181                      zzdep = zzdep * tmask(ji,jj,1) 
    182  
    183                      IF( zztmp > zmaxdzT(ji,jj) ) THEN                         
    184                          zmaxdzT(ji,jj) = zztmp    
    185                          hth    (ji,jj) = zzdep                ! max and depth of dT/dz 
    186                      ENDIF 
    187                 
    188                      IF( nla10 > 1 ) THEN  
    189                         zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1)                       ! delta rho(1) 
    190                         IF( zztmp > zrho3 )   zrho0_3(ji,jj) = zzdep                ! > 0.03 
    191                         IF( zztmp > zrho1 )   zrho0_1(ji,jj) = zzdep                ! > 0.01 
    192                      ENDIF 
    193                   END DO 
    194                END DO 
    195             END DO 
    196           
    197             CALL iom_put( 'mlddzt', hth )            ! depth of the thermocline 
    198             IF( nla10 > 1 ) THEN  
    199                CALL iom_put( 'mldr0_3', zrho0_3 )   ! MLD delta rho(surf) = 0.03 
    200                CALL iom_put( 'mldr0_1', zrho0_1 )   ! MLD delta rho(surf) = 0.01 
    201             ENDIF 
    202             ! 
    203          ENDIF 
    204          ! 
    205          IF(  iom_use( 'mld_dt02' ) .OR. iom_use( 'topthdep' ) .OR. iom_use( 'mldr10_3' ) .OR.  &     
    206             &  iom_use( 'pycndep' ) .OR. iom_use( 'tinv'     ) .OR. iom_use( 'depti'    )  ) THEN 
    207213            ! ------------------------------------------------------------- ! 
    208214            ! MLD: abs( tn - tn(10m) ) = ztem2                              ! 
Note: See TracChangeset for help on using the changeset viewer.