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 6761 – NEMO

Changeset 6761


Ignore:
Timestamp:
2016-06-30T15:21:49+02:00 (8 years ago)
Author:
kingr
Message:

Merged branches/UKMO/dev_r5518_v3.4_asm_nemovar_community@6754

Location:
branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/ASM/asmbkg.F90

    r6757 r6761  
    120120#endif 
    121121            CALL iom_rstput( kt, nitbkg_r, inum, 'gcx'    , gcx               ) 
     122            CALL iom_rstput( kt, nitbkg_r, inum, 'avt'    , avt               ) 
    122123            ! 
    123124            CALL iom_close( inum ) 
  • branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90

    r6757 r6761  
    4040#endif 
    4141   USE sbc_oce          ! Surface boundary condition variables. 
     42   USE zdfmxl, ONLY :  &   
     43   &  hmld_tref,       &    
     44#if defined key_karaml 
     45   &  hmld_kara,       & 
     46   &  ln_kara,         & 
     47#endif    
     48   &  hmld,            &  
     49   &  hmlp,            & 
     50   &  hmlpt 
     51#if defined key_bdy  
     52   USE bdy_oce, ONLY: bdytmask   
     53#endif   
    4254 
    4355   IMPLICIT NONE 
     
    88100   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   seaice_bkginc         ! Increment to the background sea ice conc 
    89101 
     102   INTEGER :: mld_choice        = 4   !: choice of mld criteria to use for physics assimilation 
     103                                      !: 1) hmld      - Turbocline/mixing depth                           [W points] 
     104                                      !: 2) hmlp      - Density criterion (0.01 kg/m^3 change from 10m)   [W points] 
     105                                      !: 3) hmld_kara - Kara MLD                                          [Interpolated] 
     106                                      !: 4) hmld_tref - Temperature criterion (0.2 K change from surface) [T points] 
     107 
     108 
    90109   !! * Substitutions 
    91110#  include "domzgr_substitute.h90" 
     
    119138      INTEGER :: iitiaustr_date  ! Date YYYYMMDD of IAU interval start time step 
    120139      INTEGER :: iitiaufin_date  ! Date YYYYMMDD of IAU interval final time step 
     140      INTEGER :: isurfstat       ! Local integer for status of reading surft variable 
    121141      ! 
    122142      REAL(wp) :: znorm        ! Normalization factor for IAU weights 
     
    127147      REAL(wp) :: zdate_inc    ! Time axis in increments file 
    128148      ! 
     149      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: &  
     150          &       t_bkginc_2d  ! file for reading in 2D   
     151      !                        ! temperature increments  
     152      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: &  
     153          &       z_mld     ! Mixed layer depth  
     154           
    129155      REAL(wp), POINTER, DIMENSION(:,:) ::   hdiv   ! 2D workspace 
     156      ! 
     157      LOGICAL :: lk_surft      ! Logical: T => Increments file contains surft variable  
     158                               !               so only apply surft increments. 
    130159      !! 
    131160      NAMELIST/nam_asminc/ ln_bkgwri,                                      & 
     
    133162         &                 ln_asmdin, ln_asmiau,                           & 
    134163         &                 nitbkg, nitdin, nitiaustr, nitiaufin, niaufn,   & 
    135          &                 ln_salfix, salfixmin, nn_divdmp 
     164         &                 ln_salfix, salfixmin, nn_divdmp, mld_choice 
    136165      !!---------------------------------------------------------------------- 
    137166 
     
    139168      ! Read Namelist nam_asminc : assimilation increment interface 
    140169      !----------------------------------------------------------------------- 
     170 
     171      ! Set default values 
     172      ln_bkgwri = .FALSE. 
     173      ln_trainc = .FALSE. 
     174      ln_dyninc = .FALSE. 
     175      ln_sshinc = .FALSE. 
    141176      ln_seaiceinc = .FALSE. 
     177      ln_asmdin = .FALSE. 
     178      ln_asmiau = .TRUE. 
     179      ln_salfix = .FALSE. 
    142180      ln_temnofreeze = .FALSE. 
     181      salfixmin = -9999 
     182      nitbkg    = 0 
     183      nitdin    = 0       
     184      nitiaustr = 1 
     185      nitiaufin = 150 
     186      niaufn    = 0 
    143187 
    144188      REWIND( numnam_ref )              ! Namelist nam_asminc in reference namelist : Assimilation increment 
     
    171215         WRITE(numout,*) '      Logical switch for ensuring that the sa > salfixmin      ln_salfix = ', ln_salfix 
    172216         WRITE(numout,*) '      Minimum salinity after applying the increments           salfixmin = ', salfixmin 
     217         WRITE(numout,*) '      Choice of MLD for physics assimilation                  mld_choice = ', mld_choice 
    173218      ENDIF 
    174219 
     
    327372      !-------------------------------------------------------------------- 
    328373 
    329       ALLOCATE( t_bkginc(jpi,jpj,jpk) ) 
    330       ALLOCATE( s_bkginc(jpi,jpj,jpk) ) 
    331       ALLOCATE( u_bkginc(jpi,jpj,jpk) ) 
    332       ALLOCATE( v_bkginc(jpi,jpj,jpk) ) 
    333       ALLOCATE( ssh_bkginc(jpi,jpj)   ) 
    334       ALLOCATE( seaice_bkginc(jpi,jpj)) 
     374      IF ( ln_trainc ) THEN 
     375         ALLOCATE( t_bkginc(jpi,jpj,jpk) ) 
     376         ALLOCATE( s_bkginc(jpi,jpj,jpk) ) 
     377         t_bkginc(:,:,:) = 0.0 
     378         s_bkginc(:,:,:) = 0.0 
     379      ENDIF 
     380      IF ( ln_dyninc ) THEN  
     381         ALLOCATE( u_bkginc(jpi,jpj,jpk) ) 
     382         ALLOCATE( v_bkginc(jpi,jpj,jpk) ) 
     383         u_bkginc(:,:,:) = 0.0 
     384         v_bkginc(:,:,:) = 0.0 
     385      ENDIF 
     386      IF ( ln_sshinc ) THEN 
     387         ALLOCATE( ssh_bkginc(jpi,jpj)   ) 
     388         ssh_bkginc(:,:) = 0.0 
     389      ENDIF 
     390      IF ( ln_seaiceinc ) THEN  
     391         ALLOCATE( seaice_bkginc(jpi,jpj)) 
     392         seaice_bkginc(:,:) = 0.0 
     393      ENDIF 
    335394#if defined key_asminc 
    336395      ALLOCATE( ssh_iau(jpi,jpj)      ) 
    337 #endif 
    338       t_bkginc(:,:,:) = 0.0 
    339       s_bkginc(:,:,:) = 0.0 
    340       u_bkginc(:,:,:) = 0.0 
    341       v_bkginc(:,:,:) = 0.0 
    342       ssh_bkginc(:,:) = 0.0 
    343       seaice_bkginc(:,:) = 0.0 
    344 #if defined key_asminc 
    345396      ssh_iau(:,:)    = 0.0 
    346397#endif 
     
    378429 
    379430         IF ( ln_trainc ) THEN    
    380             CALL iom_get( inum, jpdom_autoglo, 'bckint', t_bkginc, 1 ) 
    381             CALL iom_get( inum, jpdom_autoglo, 'bckins', s_bkginc, 1 ) 
    382             ! Apply the masks 
    383             t_bkginc(:,:,:) = t_bkginc(:,:,:) * tmask(:,:,:) 
    384             s_bkginc(:,:,:) = s_bkginc(:,:,:) * tmask(:,:,:) 
    385             ! Set missing increments to 0.0 rather than 1e+20 
    386             ! to allow for differences in masks 
    387             WHERE( ABS( t_bkginc(:,:,:) ) > 1.0e+10 ) t_bkginc(:,:,:) = 0.0 
    388             WHERE( ABS( s_bkginc(:,:,:) ) > 1.0e+10 ) s_bkginc(:,:,:) = 0.0 
     431             
     432            !Test if the increments file contains the surft variable. 
     433            isurfstat = iom_varid( inum, 'bckinsurft', ldstop = .FALSE. ) 
     434            IF ( isurfstat == -1 ) THEN 
     435               lk_surft = .FALSE. 
     436            ELSE 
     437               lk_surft = .TRUE. 
     438               CALL ctl_warn( ' Applying 2D temperature increment to bottom of ML: ', & 
     439            &                 ' bckinsurft found in increments file.' ) 
     440            ENDIF              
     441 
     442            IF (lk_surft) THEN  
     443                 
     444               ALLOCATE(z_mld(jpi,jpj))  
     445               SELECT CASE(mld_choice)  
     446               CASE(1)  
     447                  z_mld = hmld  
     448               CASE(2)  
     449                  z_mld = hmlp  
     450               CASE(3)  
     451#if defined key_karaml 
     452                  IF ( ln_kara ) THEN 
     453                     z_mld = hmld_kara 
     454                  ELSE 
     455                     CALL ctl_stop("Kara mixed layer not calculated as ln_kara=.false.") 
     456                  ENDIF 
     457#else 
     458                  CALL ctl_stop("Kara mixed layer not defined in current version of NEMO")  ! JW: Safety feature, should be removed 
     459                                                                                            ! once the Kara mixed layer is available  
     460#endif 
     461               CASE(4)  
     462                  z_mld = hmld_tref  
     463               END SELECT        
     464                       
     465               ALLOCATE( t_bkginc_2d(jpi,jpj) )  
     466               CALL iom_get( inum, jpdom_autoglo, 'bckinsurft', t_bkginc_2d, 1)  
     467#if defined key_bdy                 
     468               DO jk = 1,jpkm1  
     469                  WHERE( z_mld(:,:) > fsdepw(:,:,jk) )  
     470                     t_bkginc(:,:,jk) = t_bkginc_2d(:,:) * 0.5 * & 
     471                     &       ( 1 + cos( (fsdept(:,:,jk)/z_mld(:,:) ) * rpi ) ) 
     472                      
     473                     t_bkginc(:,:,jk) = t_bkginc(:,:,jk) * bdytmask(:,:)  
     474                  ELSEWHERE  
     475                     t_bkginc(:,:,jk) = 0.  
     476                  ENDWHERE  
     477               ENDDO  
     478#else  
     479               t_bkginc(:,:,:) = 0.  
     480#endif                 
     481               s_bkginc(:,:,:) = 0.  
     482                 
     483               DEALLOCATE(z_mld, t_bkginc_2d)  
     484             
     485            ELSE  
     486                
     487               CALL iom_get( inum, jpdom_autoglo, 'bckint', t_bkginc, 1 ) 
     488               CALL iom_get( inum, jpdom_autoglo, 'bckins', s_bkginc, 1 ) 
     489               ! Apply the masks 
     490               t_bkginc(:,:,:) = t_bkginc(:,:,:) * tmask(:,:,:) 
     491               s_bkginc(:,:,:) = s_bkginc(:,:,:) * tmask(:,:,:) 
     492               ! Set missing increments to 0.0 rather than 1e+20 
     493               ! to allow for differences in masks 
     494               WHERE( ABS( t_bkginc(:,:,:) ) > 1.0e+10 ) t_bkginc(:,:,:) = 0.0 
     495               WHERE( ABS( s_bkginc(:,:,:) ) > 1.0e+10 ) s_bkginc(:,:,:) = 0.0 
     496          
     497            ENDIF 
     498          
    389499         ENDIF 
    390500 
  • branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90

    r6758 r6761  
    7474      INTEGER, INTENT(in) ::   kt                      ! time step 
    7575      !  
    76       INTEGER             ::   jk                      ! dummy loop indice 
     76      INTEGER             ::   ji, jj, jk              ! dummy loop indices 
    7777      REAL(wp)            ::   z2dt, z1_rau0           ! local scalars 
    7878      !!---------------------------------------------------------------------- 
     
    9494      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog) 
    9595      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt 
     96 
     97#if defined key_asminc 
     98      !                                                ! Include the IAU weighted SSH increment 
     99      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN 
     100         CALL ssh_asm_inc( kt ) 
     101#if defined key_vvl 
     102! Don't directly adjust ssh but change hdivn at all levels instead 
     103! In trasbc also add in the heat and salt content associated with these changes at each level   
     104        DO jk = 1, jpkm1                                  
     105                 hdivn(:,:,jk) = hdivn(:,:,jk) - ( ssh_iau(:,:) / ( ht_0(:,:) + 1.0 - ssmask(:,:) ) ) * ( e3t_0(:,:,jk) / fse3t_n(:,:,jk) ) * tmask(:,:,jk)  
     106        END DO 
     107        CALL lbc_lnk( hdivn, 'T', 1. ) ! Not sure that's necessary 
     108      ENDIF 
     109#endif 
     110#endif 
    96111 
    97112      !                                           !------------------------------! 
     
    123138#endif 
    124139 
    125 #if defined key_asminc 
    126       !                                                ! Include the IAU weighted SSH increment 
    127       IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN 
    128          CALL ssh_asm_inc( kt ) 
    129          ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:) 
    130       ENDIF 
    131 #endif 
    132  
    133140      !                                           !------------------------------! 
    134141      !                                           !           outputs            ! 
  • branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/TRA/trasbc.F90

    r6758 r6761  
    3333   USE timing          ! Timing 
    3434   USE eosbn2 
     35#if defined key_asminc    
     36   USE asminc          ! Assimilation increment 
     37#endif 
    3538 
    3639   IMPLICIT NONE 
     
    278281         END DO   
    279282      ENDIF 
     283 
     284#if defined key_asminc 
     285! WARNING: THIS MAY WELL NOT BE REQUIRED - WE DON'T WANT TO CHANGE T&S BUT THIS MAY COMPENSATE ANOTHER TERM... 
     286! Rate of change in e3t for each level is ssh_iau*e3t_0/ht_0 
     287! Contribution to tsa should be rate of change in level / per m of ocean? (hence the division by fse3t_n) 
     288      IF( ln_sshinc ) THEN         ! input of heat and salt due to assimilation 
     289         DO jj = 2, jpj  
     290            DO ji = fs_2, fs_jpim1 
     291               zdep = ssh_iau(ji,jj) / ( ht_0(ji,jj) + 1.0 - ssmask(ji, jj) ) 
     292               DO jk = 1, jpkm1 
     293                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)   & 
     294                                        &            + tsn(ji,jj,jk,jp_tem) * zdep * ( e3t_0(ji,jj,jk) / fse3t_n(ji,jj,jk) ) 
     295                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal)   & 
     296                                        &            + tsn(ji,jj,jk,jp_sal) * zdep * ( e3t_0(ji,jj,jk) / fse3t_n(ji,jj,jk) ) 
     297               END DO 
     298            END DO   
     299         END DO   
     300      ENDIF 
     301#endif 
    280302  
    281303      IF( l_trdtra )   THEN                      ! send trends for further diagnostics 
  • branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90

    r6758 r6761  
    2727 
    2828   PUBLIC   zdf_mxl       ! called by step.F90 
     29   PUBLIC   zdf_mxl_tref  ! called by asminc.F90 
    2930   PUBLIC   zdf_mxl_alloc ! Used in zdf_tke_init 
    3031 
     
    3334   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m] 
    3435   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlpt   !: mixed layer depth at t-points        [m] 
     36   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmld_tref  !: mixed layer depth at t-points - temperature criterion [m] 
    3537 
    3638   REAL(wp), PUBLIC ::   rho_c = 0.01_wp    !: density criterion for mixed layer depth 
     
    5254      zdf_mxl_alloc = 0      ! set to zero if no array to be allocated 
    5355      IF( .NOT. ALLOCATED( nmln ) ) THEN 
    54          ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc ) 
     56         ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), & 
     57         &                           hmld_tref(jpi,jpj), STAT= zdf_mxl_alloc ) 
    5558         ! 
    5659         IF( lk_mpp             )   CALL mpp_sum ( zdf_mxl_alloc ) 
     
    146149   END SUBROUTINE zdf_mxl 
    147150 
     151 
     152   SUBROUTINE zdf_mxl_tref() 
     153      !!---------------------------------------------------------------------- 
     154      !!                  ***  ROUTINE zdf_mxl_tref  *** 
     155      !!                    
     156      !! ** Purpose :   Compute the mixed layer depth with temperature criteria. 
     157      !! 
     158      !! ** Method  :   The temperature-defined mixed layer depth is required 
     159      !!                   when assimilating SST in a 2D analysis.  
     160      !! 
     161      !! ** Action  :   hmld_tref 
     162      !!---------------------------------------------------------------------- 
     163      ! 
     164      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
     165      REAL(wp) ::   t_ref               ! Reference temperature   
     166      REAL(wp) ::   temp_c = 0.2        ! temperature criterion for mixed layer depth   
     167      !!---------------------------------------------------------------------- 
     168      ! 
     169      ! Initialise array 
     170      IF( zdf_mxl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl_tref : unable to allocate arrays' ) 
     171       
     172      !For the AMM model assimiation uses a temperature based mixed layer depth   
     173      !This is defined here   
     174      DO jj = 1, jpj   
     175         DO ji = 1, jpi   
     176           hmld_tref(ji,jj)=fsdept(ji,jj,1  )    
     177           IF(ssmask(ji,jj) > 0.)THEN   
     178             t_ref=tsn(ji,jj,1,jp_tem)  
     179             DO jk=2,jpk   
     180               IF(ssmask(ji,jj)==0.)THEN   
     181                  hmld_tref(ji,jj)=fsdept(ji,jj,jk )   
     182                  EXIT   
     183               ELSEIF( ABS(tsn(ji,jj,jk,jp_tem)-t_ref) < temp_c)THEN   
     184                  hmld_tref(ji,jj)=fsdept(ji,jj,jk )   
     185               ELSE   
     186                  EXIT   
     187               ENDIF   
     188             ENDDO   
     189           ENDIF   
     190         ENDDO   
     191      ENDDO 
     192 
     193   END SUBROUTINE zdf_mxl_tref 
     194 
    148195   !!====================================================================== 
    149196END MODULE zdfmxl 
  • branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/nemogcm.F90

    r6758 r6761  
    473473 
    474474      !                                     ! Assimilation increments 
    475       IF( lk_asminc     )   CALL asm_inc_init   ! Initialize assimilation increments 
     475      IF( lk_asminc ) THEN  
     476#if defined key_shelf  
     477         CALL  zdf_mxl_tref()     ! Initialization of hmld_tref 
     478#endif  
     479         CALL asm_inc_init     ! Initialize assimilation increments  
     480      ENDIF 
     481            
    476482      IF(lwp) WRITE(numout,*) 'Euler time step switch is ', neuler 
    477483      ! 
Note: See TracChangeset for help on using the changeset viewer.