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 8586 for branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

Ignore:
Timestamp:
2017-10-04T09:19:23+02:00 (7 years ago)
Author:
gm
Message:

#1911 (ENHANCE-09): PART I.3 - phasing with branch dev_r8183_ICEMODEL revision 8575

Location:
branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ASM
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ASM/asmbkg.F90

    r8215 r8586  
    3636   USE asmpar             ! Parameters for the assmilation interface 
    3737   USE zdfmxl             ! mixed layer depth 
    38 #if defined key_lim2 
    39    USE ice_2 
    40 #endif 
    4138#if defined key_lim3 
    4239   USE ice 
     
    142139            CALL iom_rstput( kt, nitdin_r, inum, 'sn'     , tsn(:,:,:,jp_sal) ) 
    143140            CALL iom_rstput( kt, nitdin_r, inum, 'sshn'   , sshn              ) 
    144 #if defined key_lim2 || defined key_lim3 
    145             IF( nn_ice == 2  .OR.  nn_ice == 3 ) THEN 
    146                IF( ALLOCATED(frld) ) THEN 
    147                   CALL iom_rstput( kt, nitdin_r, inum, 'iceconc', 1._wp - frld(:,:)   ) 
     141#if defined key_lim3 
     142            IF( nn_ice == 2 ) THEN 
     143               IF( ALLOCATED(at_i) ) THEN 
     144                  CALL iom_rstput( kt, nitdin_r, inum, 'iceconc', at_i(:,:)   ) 
    148145               ELSE 
    149                   CALL ctl_warn('Ice concentration not written to background as ice variable frld not allocated on this timestep') 
     146                  CALL ctl_warn('asm_bkg_wri: Ice concentration not written to background ',   & 
     147                     &          'as ice variable at_i not allocated on this timestep') 
    150148               ENDIF 
    151149            ENDIF 
  • branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90

    r6140 r8586  
    2222   !!   seaice_asm_inc : Apply the seaice increment 
    2323   !!---------------------------------------------------------------------- 
    24    USE wrk_nemo         ! Memory Allocation 
    25    USE par_oce          ! Ocean space and time domain variables 
    26    USE dom_oce          ! Ocean space and time domain 
    27    USE domvvl           ! domain: variable volume level 
    28    USE oce              ! Dynamics and active tracers defined in memory 
    29    USE ldfdyn           ! lateral diffusion: eddy viscosity coefficients 
    30    USE eosbn2           ! Equation of state - in situ and potential density 
    31    USE zpshde           ! Partial step : Horizontal Derivative 
    32    USE iom              ! Library to read input files 
    33    USE asmpar           ! Parameters for the assmilation interface 
    34    USE c1d              ! 1D initialization 
    35    USE in_out_manager   ! I/O manager 
    36    USE lib_mpp          ! MPP library 
    37 #if defined key_lim2 
    38    USE ice_2            ! LIM2 
    39 #endif 
    40    USE sbc_oce          ! Surface boundary condition variables. 
    41    USE diaobs, ONLY: calc_date     ! Compute the calendar date on a given step 
     24   USE oce             ! Dynamics and active tracers defined in memory 
     25   USE par_oce         ! Ocean space and time domain variables 
     26   USE dom_oce         ! Ocean space and time domain 
     27   USE domvvl          ! domain: variable volume level 
     28   USE ldfdyn          ! lateral diffusion: eddy viscosity coefficients 
     29   USE eosbn2          ! Equation of state - in situ and potential density 
     30   USE zpshde          ! Partial step : Horizontal Derivative 
     31   USE asmpar          ! Parameters for the assmilation interface 
     32   USE c1d             ! 1D initialization 
     33   USE sbc_oce         ! Surface boundary condition variables. 
     34   USE diaobs   , ONLY : calc_date     ! Compute the calendar date on a given step 
     35#if defined key_lim3 
     36   USE ice      , ONLY : hm_i, at_i, at_i_b 
     37#endif 
     38   ! 
     39   USE in_out_manager  ! I/O manager 
     40   USE iom             ! Library to read input files 
     41   USE lib_mpp         ! MPP library 
    4242 
    4343   IMPLICIT NONE 
     
    8686   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   ssh_bkg, ssh_bkginc   ! Background sea surface height and its increment 
    8787   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   seaice_bkginc         ! Increment to the background sea ice conc 
     88#if defined key_cice && defined key_asminc 
     89   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   ndaice_da             ! ice increment tendency into CICE 
     90#endif 
    8891 
    8992   !! * Substitutions 
     
    124127      REAL(wp) :: zdate_inc    ! Time axis in increments file 
    125128      ! 
    126       REAL(wp), POINTER, DIMENSION(:,:) ::   hdiv   ! 2D workspace 
     129      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zhdiv   ! 2D workspace 
    127130      !! 
    128131      NAMELIST/nam_asminc/ ln_bkgwri,                                      & 
     
    170173      ENDIF 
    171174 
    172       nitbkg_r    = nitbkg    + nit000 - 1  ! Background time referenced to nit000 
    173       nitdin_r    = nitdin    + nit000 - 1  ! Background time for DI referenced to nit000 
    174       nitiaustr_r = nitiaustr + nit000 - 1  ! Start of IAU interval referenced to nit000 
    175       nitiaufin_r = nitiaufin + nit000 - 1  ! End of IAU interval referenced to nit000 
    176  
    177       iiauper = nitiaufin_r - nitiaustr_r + 1  ! IAU interval length 
    178       icycper = nitend      - nit000      + 1  ! Cycle interval length 
    179  
    180       ! Date of final time step 
    181       CALL calc_date( nitend, ditend_date ) 
    182  
    183       ! Background time for Jb referenced to ndate0 
    184       CALL calc_date( nitbkg_r, ditbkg_date ) 
    185  
    186       ! Background time for DI referenced to ndate0 
    187       CALL calc_date( nitdin_r, ditdin_date ) 
    188  
    189       ! IAU start time referenced to ndate0 
    190       CALL calc_date( nitiaustr_r, ditiaustr_date ) 
    191  
    192       ! IAU end time referenced to ndate0 
    193       CALL calc_date( nitiaufin_r, ditiaufin_date ) 
     175      nitbkg_r    = nitbkg    + nit000 - 1            ! Background time referenced to nit000 
     176      nitdin_r    = nitdin    + nit000 - 1            ! Background time for DI referenced to nit000 
     177      nitiaustr_r = nitiaustr + nit000 - 1            ! Start of IAU interval referenced to nit000 
     178      nitiaufin_r = nitiaufin + nit000 - 1            ! End of IAU interval referenced to nit000 
     179 
     180      iiauper     = nitiaufin_r - nitiaustr_r + 1     ! IAU interval length 
     181      icycper     = nitend      - nit000      + 1     ! Cycle interval length 
     182 
     183      CALL calc_date( nitend     , ditend_date    )   ! Date of final time step 
     184      CALL calc_date( nitbkg_r   , ditbkg_date    )   ! Background time for Jb referenced to ndate0 
     185      CALL calc_date( nitdin_r   , ditdin_date    )   ! Background time for DI referenced to ndate0 
     186      CALL calc_date( nitiaustr_r, ditiaustr_date )   ! IAU start time referenced to ndate0 
     187      CALL calc_date( nitiaufin_r, ditiaufin_date )   ! IAU end time referenced to ndate0 
    194188 
    195189      IF(lwp) THEN 
     
    263257         ALLOCATE( wgtiau( icycper ) ) 
    264258 
    265          wgtiau(:) = 0.0 
     259         wgtiau(:) = 0._wp 
    266260 
    267261         IF ( niaufn == 0 ) THEN 
     
    339333      ALLOCATE( ssh_iau(jpi,jpj)      ) 
    340334#endif 
    341       t_bkginc(:,:,:) = 0.0 
    342       s_bkginc(:,:,:) = 0.0 
    343       u_bkginc(:,:,:) = 0.0 
    344       v_bkginc(:,:,:) = 0.0 
    345       ssh_bkginc(:,:) = 0.0 
    346       seaice_bkginc(:,:) = 0.0 
     335#if defined key_cice && defined key_asminc 
     336      ALLOCATE( ndaice_da(jpi,jpj)      ) 
     337#endif 
     338      t_bkginc     (:,:,:) = 0._wp 
     339      s_bkginc     (:,:,:) = 0._wp 
     340      u_bkginc     (:,:,:) = 0._wp 
     341      v_bkginc     (:,:,:) = 0._wp 
     342      ssh_bkginc   (:,:)   = 0._wp 
     343      seaice_bkginc(:,:)   = 0._wp 
    347344#if defined key_asminc 
    348       ssh_iau(:,:)    = 0.0 
     345      ssh_iau      (:,:)   = 0._wp 
     346#endif 
     347#if defined key_cice && defined key_asminc 
     348      ndaice_da    (:,:)   = 0._wp 
    349349#endif 
    350350      IF ( ( ln_trainc ).OR.( ln_dyninc ).OR.( ln_sshinc ).OR.( ln_seaiceinc ) ) THEN 
     
    432432      IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN 
    433433         ! 
    434          CALL wrk_alloc( jpi,jpj,   hdiv )  
     434         ALLOCATE( zhdiv(jpi,jpj) )  
    435435         ! 
    436436         DO jt = 1, nn_divdmp 
    437437            ! 
    438             DO jk = 1, jpkm1           ! hdiv = e1e1 * div 
    439                hdiv(:,:) = 0._wp 
     438            DO jk = 1, jpkm1           ! zhdiv = e1e1 * div 
     439               zhdiv(:,:) = 0._wp 
    440440               DO jj = 2, jpjm1 
    441441                  DO ji = fs_2, fs_jpim1   ! vector opt. 
    442                      hdiv(ji,jj) = (  e2u(ji  ,jj) * e3u_n(ji  ,jj,jk) * u_bkginc(ji  ,jj,jk)    & 
    443                         &           - e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) * u_bkginc(ji-1,jj,jk)    & 
    444                         &           + e1v(ji,jj  ) * e3v_n(ji,jj  ,jk) * v_bkginc(ji,jj  ,jk)    & 
    445                         &           - e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) * v_bkginc(ji,jj-1,jk)  ) / e3t_n(ji,jj,jk) 
     442                     zhdiv(ji,jj) = (  e2u(ji  ,jj) * e3u_n(ji  ,jj,jk) * u_bkginc(ji  ,jj,jk)    & 
     443                        &            - e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) * u_bkginc(ji-1,jj,jk)    & 
     444                        &            + e1v(ji,jj  ) * e3v_n(ji,jj  ,jk) * v_bkginc(ji,jj  ,jk)    & 
     445                        &            - e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) * v_bkginc(ji,jj-1,jk)  ) / e3t_n(ji,jj,jk) 
    446446                  END DO 
    447447               END DO 
    448                CALL lbc_lnk( hdiv, 'T', 1. )   ! lateral boundary cond. (no sign change) 
     448               CALL lbc_lnk( zhdiv, 'T', 1. )   ! lateral boundary cond. (no sign change) 
    449449               ! 
    450450               DO jj = 2, jpjm1 
    451451                  DO ji = fs_2, fs_jpim1   ! vector opt. 
    452452                     u_bkginc(ji,jj,jk) = u_bkginc(ji,jj,jk)                         & 
    453                         &               + 0.2_wp * ( hdiv(ji+1,jj) - hdiv(ji  ,jj) ) * r1_e1u(ji,jj) * umask(ji,jj,jk) 
     453                        &               + 0.2_wp * ( zhdiv(ji+1,jj) - zhdiv(ji  ,jj) ) * r1_e1u(ji,jj) * umask(ji,jj,jk) 
    454454                     v_bkginc(ji,jj,jk) = v_bkginc(ji,jj,jk)                         & 
    455                         &               + 0.2_wp * ( hdiv(ji,jj+1) - hdiv(ji,jj  ) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk)  
     455                        &               + 0.2_wp * ( zhdiv(ji,jj+1) - zhdiv(ji,jj  ) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk)  
    456456                  END DO 
    457457               END DO 
     
    460460         END DO 
    461461         ! 
    462          CALL wrk_dealloc( jpi,jpj,   hdiv )  
     462         DEALLOCATE( zhdiv )  
    463463         ! 
    464464      ENDIF 
     
    800800      INTEGER  ::   it 
    801801      REAL(wp) ::   zincwgt   ! IAU weight for current time step 
    802 #if defined key_lim2 
     802#if defined key_lim3 
    803803      REAL(wp), DIMENSION(jpi,jpj) ::   zofrld, zohicif, zseaicendg, zhicifinc  ! LIM 
    804804      REAL(wp) ::   zhicifmin = 0.5_wp      ! ice minimum depth in metres 
     
    822822            ENDIF 
    823823            ! 
    824             ! Sea-ice : LIM-3 case (to add) 
    825             ! 
    826 #if defined key_lim2 
    827             ! Sea-ice : LIM-2 case 
    828             zofrld (:,:) = frld(:,:) 
    829             zohicif(:,:) = hicif(:,:) 
    830             ! 
    831             frld  = MIN( MAX( frld (:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp) 
    832             pfrld = MIN( MAX( pfrld(:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp) 
    833             fr_i(:,:) = 1.0_wp - frld(:,:)        ! adjust ice fraction 
    834             ! 
    835             zseaicendg(:,:) = zofrld(:,:) - frld(:,:)   ! find out actual sea ice nudge applied 
     824            ! Sea-ice : LIM-3 case 
     825            ! 
     826#if defined key_lim3 
     827            zofrld (:,:) = 1._wp - at_i(:,:) 
     828            zohicif(:,:) = hm_i(:,:) 
     829            ! 
     830            at_i  (:,:) = 1. - MIN( MAX( 1.-at_i  (:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp) 
     831            at_i_b(:,:) = 1. - MIN( MAX( 1.-at_i_b(:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp) 
     832            fr_i(:,:) = at_i(:,:)        ! adjust ice fraction 
     833            ! 
     834            zseaicendg(:,:) = zofrld(:,:) - (1. - at_i(:,:))   ! find out actual sea ice nudge applied 
    836835            ! 
    837836            ! Nudge sea ice depth to bring it up to a required minimum depth 
    838             WHERE( zseaicendg(:,:) > 0.0_wp .AND. hicif(:,:) < zhicifmin )  
    839                zhicifinc(:,:) = (zhicifmin - hicif(:,:)) * zincwgt     
     837            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin )  
     838               zhicifinc(:,:) = (zhicifmin - hm_i(:,:)) * zincwgt     
    840839            ELSEWHERE 
    841840               zhicifinc(:,:) = 0.0_wp 
     
    843842            ! 
    844843            ! nudge ice depth 
    845             hicif (:,:) = hicif (:,:) + zhicifinc(:,:) 
    846             phicif(:,:) = phicif(:,:) + zhicifinc(:,:)        
     844            hm_i (:,:) = hm_i (:,:) + zhicifinc(:,:) 
    847845            ! 
    848846            ! seaice salinity balancing (to add) 
     
    873871            neuler = 0                    ! Force Euler forward step 
    874872            ! 
    875             ! Sea-ice : LIM-3 case (to add) 
    876             ! 
    877 #if defined key_lim2 
    878             ! Sea-ice : LIM-2 case. 
    879             zofrld(:,:)=frld(:,:) 
    880             zohicif(:,:)=hicif(:,:) 
     873            ! Sea-ice : LIM-3 case 
     874            ! 
     875#if defined key_lim3 
     876            zofrld (:,:) = 1._wp - at_i(:,:) 
     877            zohicif(:,:) = hm_i(:,:) 
    881878            !  
    882879            ! Initialize the now fields the background + increment 
    883             frld (:,:) = MIN( MAX( frld(:,:) - seaice_bkginc(:,:), 0.0_wp), 1.0_wp) 
    884             pfrld(:,:) = frld(:,:)  
    885             fr_i (:,:) = 1.0_wp - frld(:,:)                ! adjust ice fraction 
    886             zseaicendg(:,:) = zofrld(:,:) - frld(:,:)      ! find out actual sea ice nudge applied 
     880            at_i(:,:) = 1. - MIN( MAX( 1.-at_i(:,:) - seaice_bkginc(:,:), 0.0_wp), 1.0_wp) 
     881            at_i_b(:,:) = at_i(:,:)  
     882            fr_i(:,:) = at_i(:,:)        ! adjust ice fraction 
     883            ! 
     884            zseaicendg(:,:) = zofrld(:,:) - (1. - at_i(:,:))   ! find out actual sea ice nudge applied 
    887885            ! 
    888886            ! Nudge sea ice depth to bring it up to a required minimum depth 
    889             WHERE( zseaicendg(:,:) > 0.0_wp .AND. hicif(:,:) < zhicifmin )  
    890                zhicifinc(:,:) = (zhicifmin - hicif(:,:)) * zincwgt     
     887            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin )  
     888               zhicifinc(:,:) = (zhicifmin - hm_i(:,:)) * zincwgt     
    891889            ELSEWHERE 
    892                zhicifinc(:,:) = 0._wp 
     890               zhicifinc(:,:) = 0.0_wp 
    893891            END WHERE 
    894892            ! 
    895893            ! nudge ice depth 
    896             hicif (:,:) = hicif (:,:) + zhicifinc(:,:) 
    897             phicif(:,:) = phicif(:,:)        
     894            hm_i (:,:) = hm_i (:,:) + zhicifinc(:,:) 
    898895            ! 
    899896            ! seaice salinity balancing (to add) 
     
    917914         ENDIF 
    918915 
    919 !#if defined defined key_lim2 || defined key_cice 
     916!#if defined defined key_lim3 || defined key_cice 
    920917! 
    921918!            IF (ln_seaicebal ) THEN        
Note: See TracChangeset for help on using the changeset viewer.