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 8698 for trunk/NEMOGCM/NEMO/OPA_SRC – NEMO

Ignore:
Timestamp:
2017-11-13T15:24:50+01:00 (7 years ago)
Author:
davestorkey
Message:

Bug fixes for tracer trends diagnostics in trunk. See ticket #1877.
These fixes correspond to the changes to the 3.6 branch at revisions 8102 and 8305.

Location:
trunk/NEMOGCM/NEMO/OPA_SRC
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/NEMO/OPA_SRC/IOM/iom.F90

    r8573 r8698  
    15971597      ! frequency of the call of iom_put (attribut: freq_op) 
    15981598      f_op%timestep = 1        ;  f_of%timestep = 0  ; CALL iom_set_field_attr('field_definition', freq_op=f_op, freq_offset=f_of) 
     1599      f_op%timestep = 2        ;  f_of%timestep = 0  ; CALL iom_set_field_attr('trendT_even'     , freq_op=f_op, freq_offset=f_of) 
     1600      f_op%timestep = 2        ;  f_of%timestep = -1 ; CALL iom_set_field_attr('trendT_odd'      , freq_op=f_op, freq_offset=f_of) 
    15991601      f_op%timestep = nn_fsbc  ;  f_of%timestep = 0  ; CALL iom_set_field_attr('SBC'             , freq_op=f_op, freq_offset=f_of) 
    16001602      f_op%timestep = nn_fsbc  ;  f_of%timestep = 0  ; CALL iom_set_field_attr('SBC_scalar'      , freq_op=f_op, freq_offset=f_of) 
  • trunk/NEMOGCM/NEMO/OPA_SRC/TRA/tranxt.F90

    r7753 r8698  
    121121      IF( l_trdtra )   THEN                     
    122122         CALL wrk_alloc( jpi, jpj, jpk, ztrdt, ztrds ) 
    123          ztrdt(:,:,jk) = 0._wp 
    124          ztrds(:,:,jk) = 0._wp 
     123         ztrdt(:,:,jpk) = 0._wp 
     124         ztrds(:,:,jpk) = 0._wp 
    125125         IF( ln_traldf_iso ) THEN              ! diagnose the "pure" Kz diffusive trend  
    126126            CALL trd_tra( kt, 'TRA', jp_tem, jptra_zdfp, ztrdt ) 
     
    128128         ENDIF 
    129129         ! total trend for the non-time-filtered variables.  
    130             zfact = 1.0 / rdt 
     130         zfact = 1.0 / rdt 
     131         ! G Nurser 23 Mar 2017. Recalculate trend as Delta(e3t*T)/e3tn; e3tn cancel from tsn terms 
    131132         DO jk = 1, jpkm1 
    132             ztrdt(:,:,jk) = ( tsa(:,:,jk,jp_tem) - tsn(:,:,jk,jp_tem) ) * zfact  
    133             ztrds(:,:,jk) = ( tsa(:,:,jk,jp_sal) - tsn(:,:,jk,jp_sal) ) * zfact  
     133            ztrdt(:,:,jk) = ( tsa(:,:,jk,jp_tem)*e3t_a(:,:,jk) / e3t_n(:,:,jk) - tsn(:,:,jk,jp_tem)) * zfact 
     134            ztrds(:,:,jk) = ( tsa(:,:,jk,jp_sal)*e3t_a(:,:,jk) / e3t_n(:,:,jk) - tsn(:,:,jk,jp_sal)) * zfact 
    134135         END DO 
    135136         CALL trd_tra( kt, 'TRA', jp_tem, jptra_tot, ztrdt ) 
    136137         CALL trd_tra( kt, 'TRA', jp_sal, jptra_tot, ztrds ) 
    137          ! Store now fields before applying the Asselin filter  
    138          ! in order to calculate Asselin filter trend later. 
    139          ztrdt(:,:,:) = tsn(:,:,:,jp_tem)  
    140          ztrds(:,:,:) = tsn(:,:,:,jp_sal) 
     138         IF( ln_linssh ) THEN  
     139            ! Store now fields before applying the Asselin filter  
     140            ! in order to calculate Asselin filter trend later. 
     141            ztrdt(:,:,:) = tsn(:,:,:,jp_tem)  
     142            ztrds(:,:,:) = tsn(:,:,:,jp_sal) 
     143         ENDIF 
    141144      ENDIF 
    142145 
     
    147150            END DO 
    148151         END DO 
     152         IF (l_trdtra .AND. .NOT. ln_linssh) THEN  ! Zero Asselin filter contribution must be explicitly written out since for vvl 
     153                                                   ! Asselin filter is output by tra_nxt_vvl that is not called on this time step 
     154            ztrdt(:,:,:) = 0._wp 
     155            ztrds(:,:,:) = 0._wp 
     156            CALL trd_tra( kt, 'TRA', jp_tem, jptra_atf, ztrdt ) 
     157            CALL trd_tra( kt, 'TRA', jp_sal, jptra_atf, ztrds ) 
     158         END IF 
    149159         ! 
    150160      ELSE                                            ! Leap-Frog + Asselin filter time stepping 
     
    162172      ENDIF      
    163173      ! 
    164       IF( l_trdtra ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt      
     174      IF( l_trdtra .AND. ln_linssh ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt      
     175         zfact = 1._wp / r2dt              
    165176         DO jk = 1, jpkm1 
    166             zfact = 1._wp / r2dt              
    167177            ztrdt(:,:,jk) = ( tsb(:,:,jk,jp_tem) - ztrdt(:,:,jk) ) * zfact 
    168178            ztrds(:,:,jk) = ( tsb(:,:,jk,jp_sal) - ztrds(:,:,jk) ) * zfact 
     
    170180         CALL trd_tra( kt, 'TRA', jp_tem, jptra_atf, ztrdt ) 
    171181         CALL trd_tra( kt, 'TRA', jp_sal, jptra_atf, ztrds ) 
    172          CALL wrk_dealloc( jpi, jpj, jpk, ztrdt, ztrds ) 
    173182      END IF 
     183      IF( l_trdtra ) CALL wrk_dealloc( jpi, jpj, jpk, ztrdt, ztrds ) 
    174184      ! 
    175185      !                        ! control print 
     
    259269      LOGICAL  ::   ll_traqsr, ll_rnf, ll_isf   ! local logical 
    260270      INTEGER  ::   ji, jj, jk, jn              ! dummy loop indices 
    261       REAL(wp) ::   zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar 
     271      REAL(wp) ::   zfact, zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar 
    262272      REAL(wp) ::   zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d   !   -      - 
     273      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztrd_atf 
    263274      !!---------------------------------------------------------------------- 
    264275      ! 
     
    279290      ENDIF 
    280291      ! 
     292      IF( ( l_trdtra .and. cdtype == 'TRA' ) .OR. ( l_trdtrc .and. cdtype == 'TRC' ) )   THEN 
     293         CALL wrk_alloc( jpi, jpj, jpk, kjpt, ztrd_atf ) 
     294         ztrd_atf(:,:,:,:) = 0.0_wp 
     295      ENDIF 
     296      zfact = 1._wp / r2dt 
    281297      DO jn = 1, kjpt       
    282298         DO jk = 1, jpkm1 
     
    331347                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta 
    332348                  ! 
     349                  IF( ( l_trdtra .and. cdtype == 'TRA' ) .OR. ( l_trdtrc .and. cdtype == 'TRC' ) ) THEN 
     350                     ztrd_atf(ji,jj,jk,jn) = (ztc_f - ztc_n) * zfact/ze3t_n 
     351                  ENDIF 
     352                  ! 
    333353               END DO 
    334354            END DO 
     
    337357      END DO 
    338358      ! 
     359      IF( l_trdtra .and. cdtype == 'TRA' ) THEN  
     360         CALL trd_tra( kt, cdtype, jp_tem, jptra_atf, ztrd_atf(:,:,:,jp_tem) ) 
     361         CALL trd_tra( kt, cdtype, jp_sal, jptra_atf, ztrd_atf(:,:,:,jp_sal) ) 
     362         CALL wrk_dealloc( jpi, jpj, jpk, kjpt, ztrd_atf ) 
     363      ENDIF 
     364      IF( l_trdtrc .and. cdtype == 'TRC' ) THEN 
     365         DO jn = 1, kjpt 
     366            CALL trd_tra( kt, cdtype, jn, jptra_atf, ztrd_atf(:,:,:,jn) ) 
     367         END DO 
     368         CALL wrk_dealloc( jpi, jpj, jpk, kjpt, ztrd_atf ) 
     369      ENDIF 
     370      ! 
    339371   END SUBROUTINE tra_nxt_vvl 
    340372 
  • trunk/NEMOGCM/NEMO/OPA_SRC/TRA/trazdf.F90

    r7753 r8698  
    8989      IF( l_trdtra )   THEN                      ! save the vertical diffusive trends for further diagnostics 
    9090         DO jk = 1, jpkm1 
    91             ztrdt(:,:,jk) = ( ( tsa(:,:,jk,jp_tem) - tsb(:,:,jk,jp_tem) ) / r2dt ) - ztrdt(:,:,jk) 
    92             ztrds(:,:,jk) = ( ( tsa(:,:,jk,jp_sal) - tsb(:,:,jk,jp_sal) ) / r2dt ) - ztrds(:,:,jk) 
     91            ztrdt(:,:,jk) = ( ( tsa(:,:,jk,jp_tem)*e3t_a(:,:,jk) - tsb(:,:,jk,jp_tem)*e3t_b(:,:,jk) ) & 
     92                 & / (e3t_n(:,:,jk)*r2dt) ) - ztrdt(:,:,jk) 
     93            ztrds(:,:,jk) = ( ( tsa(:,:,jk,jp_sal)*e3t_a(:,:,jk) - tsb(:,:,jk,jp_sal)*e3t_b(:,:,jk) ) & 
     94                 & / (e3t_n(:,:,jk)*r2dt) ) - ztrds(:,:,jk) 
    9395         END DO 
    9496!!gm this should be moved in trdtra.F90 and done on all trends 
  • trunk/NEMOGCM/NEMO/OPA_SRC/TRD/trdtra.F90

    r7646 r8698  
    104104                                 ztrds(:,:,:) = 0._wp 
    105105                                 CALL trd_tra_mng( trdt, ztrds, ktrd, kt ) 
     106         CASE( jptra_evd )   ;   avt_evd(:,:,:) = ptrd(:,:,:) * tmask(:,:,:) 
    106107         CASE DEFAULT                 ! other trends: masked trends 
    107108            trdt(:,:,:) = ptrd(:,:,:) * tmask(:,:,:)              ! mask & store 
     
    311312!!gm Rq: mask the trends already masked in trd_tra, but lbc_lnk should probably be added 
    312313      ! 
     314      ! Trends evaluated every time step that could go to the standard T file and can be output every ts into a 1ts file if 1ts output is selected 
    313315      SELECT CASE( ktrd ) 
    314       CASE( jptra_xad  )   ;   CALL iom_put( "ttrd_xad" , ptrdx )        ! x- horizontal advection 
    315                                CALL iom_put( "strd_xad" , ptrdy ) 
    316       CASE( jptra_yad  )   ;   CALL iom_put( "ttrd_yad" , ptrdx )        ! y- horizontal advection 
    317                                CALL iom_put( "strd_yad" , ptrdy ) 
    318       CASE( jptra_zad  )   ;   CALL iom_put( "ttrd_zad" , ptrdx )        ! z- vertical   advection 
    319                                CALL iom_put( "strd_zad" , ptrdy ) 
    320                                IF( ln_linssh ) THEN                   ! cst volume : adv flux through z=0 surface 
    321                                   CALL wrk_alloc( jpi, jpj, z2dx, z2dy ) 
    322                                   z2dx(:,:) = wn(:,:,1) * tsn(:,:,1,jp_tem) / e3t_n(:,:,1) 
    323                                   z2dy(:,:) = wn(:,:,1) * tsn(:,:,1,jp_sal) / e3t_n(:,:,1) 
    324                                   CALL iom_put( "ttrd_sad", z2dx ) 
    325                                   CALL iom_put( "strd_sad", z2dy ) 
    326                                   CALL wrk_dealloc( jpi, jpj, z2dx, z2dy ) 
    327                                ENDIF 
    328       CASE( jptra_totad  ) ;   CALL iom_put( "ttrd_totad" , ptrdx )      ! total   advection 
    329                                CALL iom_put( "strd_totad" , ptrdy ) 
    330       CASE( jptra_ldf  )   ;   CALL iom_put( "ttrd_ldf" , ptrdx )        ! lateral diffusion 
    331                                CALL iom_put( "strd_ldf" , ptrdy ) 
    332       CASE( jptra_zdf  )   ;   CALL iom_put( "ttrd_zdf" , ptrdx )        ! vertical diffusion (including Kz contribution) 
    333                                CALL iom_put( "strd_zdf" , ptrdy ) 
    334       CASE( jptra_zdfp )   ;   CALL iom_put( "ttrd_zdfp", ptrdx )        ! PURE vertical diffusion (no isoneutral contribution) 
    335                                CALL iom_put( "strd_zdfp", ptrdy ) 
    336       CASE( jptra_evd )    ;   CALL iom_put( "ttrd_evd", ptrdx )         ! EVD trend (convection) 
    337                                CALL iom_put( "strd_evd", ptrdy ) 
    338       CASE( jptra_dmp  )   ;   CALL iom_put( "ttrd_dmp" , ptrdx )        ! internal restoring (damping) 
    339                                CALL iom_put( "strd_dmp" , ptrdy ) 
    340       CASE( jptra_bbl  )   ;   CALL iom_put( "ttrd_bbl" , ptrdx )        ! bottom boundary layer 
    341                                CALL iom_put( "strd_bbl" , ptrdy ) 
    342       CASE( jptra_npc  )   ;   CALL iom_put( "ttrd_npc" , ptrdx )        ! static instability mixing 
    343                                CALL iom_put( "strd_npc" , ptrdy ) 
    344       CASE( jptra_nsr  )   ;   CALL iom_put( "ttrd_qns" , ptrdx(:,:,1) )        ! surface forcing + runoff (ln_rnf=T) 
    345                                CALL iom_put( "strd_cdt" , ptrdy(:,:,1) )        ! output as 2D surface fields 
    346       CASE( jptra_qsr  )   ;   CALL iom_put( "ttrd_qsr" , ptrdx )        ! penetrative solar radiat. (only on temperature) 
    347       CASE( jptra_bbc  )   ;   CALL iom_put( "ttrd_bbc" , ptrdx )        ! geothermal heating   (only on temperature) 
    348       CASE( jptra_atf  )   ;   CALL iom_put( "ttrd_atf" , ptrdx )        ! asselin time Filter 
    349                                CALL iom_put( "strd_atf" , ptrdy ) 
    350       CASE( jptra_tot  )   ;   CALL iom_put( "ttrd_tot" , ptrdx )        ! model total trend 
    351                                CALL iom_put( "strd_tot" , ptrdy ) 
     316      ! This total trend is done every time step 
     317      CASE( jptra_tot  )   ;   CALL iom_put( "ttrd_tot" , ptrdx )           ! model total trend 
     318         CALL iom_put( "strd_tot" , ptrdy ) 
    352319      END SELECT 
     320 
     321      ! These trends are done every second time step. When 1ts output is selected must go different (2ts) file from standard T-file 
     322      IF( MOD( kt, 2 ) == 0 ) THEN 
     323         SELECT CASE( ktrd ) 
     324         CASE( jptra_xad  )   ;   CALL iom_put( "ttrd_xad" , ptrdx )        ! x- horizontal advection 
     325            CALL iom_put( "strd_xad" , ptrdy ) 
     326         CASE( jptra_yad  )   ;   CALL iom_put( "ttrd_yad" , ptrdx )        ! y- horizontal advection 
     327            CALL iom_put( "strd_yad" , ptrdy ) 
     328         CASE( jptra_zad  )   ;   CALL iom_put( "ttrd_zad" , ptrdx )        ! z- vertical   advection 
     329            CALL iom_put( "strd_zad" , ptrdy ) 
     330            IF( ln_linssh ) THEN                   ! cst volume : adv flux through z=0 surface 
     331               CALL wrk_alloc( jpi, jpj, z2dx, z2dy ) 
     332               z2dx(:,:) = wn(:,:,1) * tsn(:,:,1,jp_tem) / e3t_n(:,:,1) 
     333               z2dy(:,:) = wn(:,:,1) * tsn(:,:,1,jp_sal) / e3t_n(:,:,1) 
     334               CALL iom_put( "ttrd_sad", z2dx ) 
     335               CALL iom_put( "strd_sad", z2dy ) 
     336               CALL wrk_dealloc( jpi, jpj, z2dx, z2dy ) 
     337            ENDIF 
     338         CASE( jptra_totad  ) ;   CALL iom_put( "ttrd_totad" , ptrdx )      ! total   advection 
     339            CALL iom_put( "strd_totad" , ptrdy ) 
     340         CASE( jptra_ldf  )   ;   CALL iom_put( "ttrd_ldf" , ptrdx )        ! lateral diffusion 
     341            CALL iom_put( "strd_ldf" , ptrdy ) 
     342         CASE( jptra_zdf  )   ;   CALL iom_put( "ttrd_zdf" , ptrdx )        ! vertical diffusion (including Kz contribution) 
     343            CALL iom_put( "strd_zdf" , ptrdy ) 
     344         CASE( jptra_zdfp )   ;   CALL iom_put( "ttrd_zdfp", ptrdx )        ! PURE vertical diffusion (no isoneutral contribution) 
     345            CALL iom_put( "strd_zdfp", ptrdy ) 
     346         CASE( jptra_evd )    ;   CALL iom_put( "ttrd_evd", ptrdx )         ! EVD trend (convection) 
     347            CALL iom_put( "strd_evd", ptrdy ) 
     348         CASE( jptra_dmp  )   ;   CALL iom_put( "ttrd_dmp" , ptrdx )        ! internal restoring (damping) 
     349            CALL iom_put( "strd_dmp" , ptrdy ) 
     350         CASE( jptra_bbl  )   ;   CALL iom_put( "ttrd_bbl" , ptrdx )        ! bottom boundary layer 
     351            CALL iom_put( "strd_bbl" , ptrdy ) 
     352         CASE( jptra_npc  )   ;   CALL iom_put( "ttrd_npc" , ptrdx )        ! static instability mixing 
     353            CALL iom_put( "strd_npc" , ptrdy ) 
     354         CASE( jptra_bbc  )   ;   CALL iom_put( "ttrd_bbc" , ptrdx )        ! geothermal heating   (only on temperature) 
     355         CASE( jptra_nsr  )   ;   CALL iom_put( "ttrd_qns" , ptrdx(:,:,1) ) ! surface forcing + runoff (ln_rnf=T) 
     356            CALL iom_put( "strd_cdt" , ptrdy(:,:,1) )        ! output as 2D surface fields 
     357         CASE( jptra_qsr  )   ;   CALL iom_put( "ttrd_qsr" , ptrdx )        ! penetrative solar radiat. (only on temperature) 
     358         END SELECT 
     359         ! the Asselin filter trend  is also every other time step but needs to be lagged one time step 
     360         ! Even when 1ts output is selected can go to the same (2ts) file as the trends plotted every even time step. 
     361      ELSE IF( MOD( kt, 2 ) == 1 ) THEN 
     362         SELECT CASE( ktrd ) 
     363         CASE( jptra_atf  )   ;   CALL iom_put( "ttrd_atf" , ptrdx )        ! asselin time Filter 
     364            CALL iom_put( "strd_atf" , ptrdy ) 
     365         END SELECT 
     366      END IF 
    353367      ! 
    354368   END SUBROUTINE trd_tra_iom 
Note: See TracChangeset for help on using the changeset viewer.