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 12680 for NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/ASM/asminc.F90 – NEMO

Ignore:
Timestamp:
2020-04-03T18:54:55+02:00 (4 years ago)
Author:
techene
Message:

dynatfQCO.F90, stepLF.F90 : fixed (remove pe3. from dyn_atf_qco input arguments), all : remove e3. tables and include gurvan's feedbacks

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/ASM/asminc.F90

    r12656 r12680  
    99   !!                 ! 2007-04  (A. Weaver)  Merge with OPAVAR/NEMOVAR 
    1010   !!   NEMO     3.3  ! 2010-05  (D. Lea)  Update to work with NEMO v3.2 
    11    !!             -   ! 2010-05  (D. Lea)  add calc_month_len routine based on day_init  
     11   !!             -   ! 2010-05  (D. Lea)  add calc_month_len routine based on day_init 
    1212   !!            3.4  ! 2012-10  (A. Weaver and K. Mogensen) Fix for direct initialization 
    1313   !!                 ! 2014-09  (D. Lea)  Local calc_date removed use routine from OBS 
     
    3131   USE zpshde          ! Partial step : Horizontal Derivative 
    3232   USE asmpar          ! Parameters for the assmilation interface 
    33    USE asmbkg          !  
     33   USE asmbkg          ! 
    3434   USE c1d             ! 1D initialization 
    3535   USE sbc_oce         ! Surface boundary condition variables. 
     
    4545   IMPLICIT NONE 
    4646   PRIVATE 
    47     
     47 
    4848   PUBLIC   asm_inc_init   !: Initialize the increment arrays and IAU weights 
    4949   PUBLIC   tra_asm_inc    !: Apply the tracer (T and S) increments 
     
    7272   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkg   , v_bkg      !: Background u- & v- velocity components 
    7373   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   t_bkginc, s_bkginc   !: Increment to the background T & S 
    74    REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkginc, v_bkginc   !: Increment to the u- & v-components  
     74   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkginc, v_bkginc   !: Increment to the u- & v-components 
    7575   REAL(wp), PUBLIC, DIMENSION(:)    , ALLOCATABLE ::   wgtiau               !: IAU weights for each time step 
    7676#if defined key_asminc 
     
    8080   INTEGER , PUBLIC ::   nitbkg      !: Time step of the background state used in the Jb term 
    8181   INTEGER , PUBLIC ::   nitdin      !: Time step of the background state for direct initialization 
    82    INTEGER , PUBLIC ::   nitiaustr   !: Time step of the start of the IAU interval  
     82   INTEGER , PUBLIC ::   nitiaustr   !: Time step of the start of the IAU interval 
    8383   INTEGER , PUBLIC ::   nitiaufin   !: Time step of the end of the IAU interval 
    84    !  
     84   ! 
    8585   INTEGER , PUBLIC ::   niaufn      !: Type of IAU weighing function: = 0   Constant weighting 
    86    !                                 !: = 1   Linear hat-like, centred in middle of IAU interval  
     86   !                                 !: = 1   Linear hat-like, centred in middle of IAU interval 
    8787   REAL(wp), PUBLIC ::   salfixmin   !: Ensure that the salinity is larger than this  value if (ln_salfix) 
    8888 
     
    106106      !!---------------------------------------------------------------------- 
    107107      !!                    ***  ROUTINE asm_inc_init  *** 
    108       !!           
     108      !! 
    109109      !! ** Purpose : Initialize the assimilation increment and IAU weights. 
    110110      !! 
    111111      !! ** Method  : Initialize the assimilation increment and IAU weights. 
    112112      !! 
    113       !! ** Action  :  
     113      !! ** Action  : 
    114114      !!---------------------------------------------------------------------- 
    115115      INTEGER, INTENT(in) ::  Kbb, Kmm, Krhs  ! time level indices 
     
    263263         ! 
    264264         !                                !--------------------------------------------------------- 
    265          IF( niaufn == 0 ) THEN           ! Constant IAU forcing  
     265         IF( niaufn == 0 ) THEN           ! Constant IAU forcing 
    266266            !                             !--------------------------------------------------------- 
    267267            DO jt = 1, iiauper 
     
    269269            END DO 
    270270            !                             !--------------------------------------------------------- 
    271          ELSEIF ( niaufn == 1 ) THEN      ! Linear hat-like, centred in middle of IAU interval  
     271         ELSEIF ( niaufn == 1 ) THEN      ! Linear hat-like, centred in middle of IAU interval 
    272272            !                             !--------------------------------------------------------- 
    273273            ! Compute the normalization factor 
    274274            znorm = 0._wp 
    275275            IF( MOD( iiauper, 2 ) == 0 ) THEN   ! Even number of time steps in IAU interval 
    276                imid = iiauper / 2  
     276               imid = iiauper / 2 
    277277               DO jt = 1, imid 
    278278                  znorm = znorm + REAL( jt ) 
     
    280280               znorm = 2.0 * znorm 
    281281            ELSE                                ! Odd number of time steps in IAU interval 
    282                imid = ( iiauper + 1 ) / 2         
     282               imid = ( iiauper + 1 ) / 2 
    283283               DO jt = 1, imid - 1 
    284284                  znorm = znorm + REAL( jt ) 
     
    307307             DO jt = 1, icycper 
    308308                ztotwgt = ztotwgt + wgtiau(jt) 
    309                 WRITE(numout,*) '         ', jt, '       ', wgtiau(jt)  
    310              END DO    
     309                WRITE(numout,*) '         ', jt, '       ', wgtiau(jt) 
     310             END DO 
    311311             WRITE(numout,*) '         ===================================' 
    312312             WRITE(numout,*) '         Time-integrated weight = ', ztotwgt 
    313313             WRITE(numout,*) '         ===================================' 
    314314          ENDIF 
    315           
     315 
    316316      ENDIF 
    317317 
     
    338338         CALL iom_open( c_asminc, inum ) 
    339339         ! 
    340          CALL iom_get( inum, 'time'       , zdate_inc   )  
     340         CALL iom_get( inum, 'time'       , zdate_inc   ) 
    341341         CALL iom_get( inum, 'z_inc_dateb', z_inc_dateb ) 
    342342         CALL iom_get( inum, 'z_inc_datef', z_inc_datef ) 
     
    345345         ! 
    346346         IF(lwp) THEN 
    347             WRITE(numout,*)  
     347            WRITE(numout,*) 
    348348            WRITE(numout,*) 'asm_inc_init : Assimilation increments valid between dates ', z_inc_dateb,' and ', z_inc_datef 
    349349            WRITE(numout,*) '~~~~~~~~~~~~' 
     
    359359            &                ' not agree with Direct Initialization time' ) 
    360360 
    361          IF ( ln_trainc ) THEN    
     361         IF ( ln_trainc ) THEN 
    362362            CALL iom_get( inum, jpdom_autoglo, 'bckint', t_bkginc, 1 ) 
    363363            CALL iom_get( inum, jpdom_autoglo, 'bckins', s_bkginc, 1 ) 
     
    371371         ENDIF 
    372372 
    373          IF ( ln_dyninc ) THEN    
    374             CALL iom_get( inum, jpdom_autoglo, 'bckinu', u_bkginc, 1 )               
    375             CALL iom_get( inum, jpdom_autoglo, 'bckinv', v_bkginc, 1 )               
     373         IF ( ln_dyninc ) THEN 
     374            CALL iom_get( inum, jpdom_autoglo, 'bckinu', u_bkginc, 1 ) 
     375            CALL iom_get( inum, jpdom_autoglo, 'bckinv', v_bkginc, 1 ) 
    376376            ! Apply the masks 
    377377            u_bkginc(:,:,:) = u_bkginc(:,:,:) * umask(:,:,:) 
     
    382382            WHERE( ABS( v_bkginc(:,:,:) ) > 1.0e+10 ) v_bkginc(:,:,:) = 0.0 
    383383         ENDIF 
    384          
     384 
    385385         IF ( ln_sshinc ) THEN 
    386386            CALL iom_get( inum, jpdom_autoglo, 'bckineta', ssh_bkginc, 1 ) 
     
    408408      IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN    ! Apply divergence damping filter 
    409409         !                                         !-------------------------------------- 
    410          ALLOCATE( zhdiv(jpi,jpj) )  
     410         ALLOCATE( zhdiv(jpi,jpj) ) 
    411411         ! 
    412412         DO jt = 1, nn_divdmp 
     
    427427                     &               + 0.2_wp * ( zhdiv(ji+1,jj) - zhdiv(ji  ,jj) ) * r1_e1u(ji,jj) * umask(ji,jj,jk) 
    428428                  v_bkginc(ji,jj,jk) = v_bkginc(ji,jj,jk)                         & 
    429                      &               + 0.2_wp * ( zhdiv(ji,jj+1) - zhdiv(ji,jj  ) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk)  
     429                     &               + 0.2_wp * ( zhdiv(ji,jj+1) - zhdiv(ji,jj  ) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk) 
    430430               END_2D 
    431431            END DO 
     
    433433         END DO 
    434434         ! 
    435          DEALLOCATE( zhdiv )  
     435         DEALLOCATE( zhdiv ) 
    436436         ! 
    437437      ENDIF 
     
    454454         CALL iom_open( c_asmdin, inum ) 
    455455         ! 
    456          CALL iom_get( inum, 'rdastp', zdate_bkg )  
     456         CALL iom_get( inum, 'rdastp', zdate_bkg ) 
    457457         ! 
    458458         IF(lwp) THEN 
    459             WRITE(numout,*)  
     459            WRITE(numout,*) 
    460460            WRITE(numout,*) '   ==>>>  Assimilation background state valid at : ', zdate_bkg 
    461461            WRITE(numout,*) 
     
    466466            &                ' not agree with Direct Initialization time' ) 
    467467         ! 
    468          IF ( ln_trainc ) THEN    
     468         IF ( ln_trainc ) THEN 
    469469            CALL iom_get( inum, jpdom_autoglo, 'tn', t_bkg ) 
    470470            CALL iom_get( inum, jpdom_autoglo, 'sn', s_bkg ) 
     
    473473         ENDIF 
    474474         ! 
    475          IF ( ln_dyninc ) THEN    
     475         IF ( ln_dyninc ) THEN 
    476476            CALL iom_get( inum, jpdom_autoglo, 'un', u_bkg ) 
    477477            CALL iom_get( inum, jpdom_autoglo, 'vn', v_bkg ) 
     
    501501      ! 
    502502   END SUBROUTINE asm_inc_init 
    503     
    504     
     503 
     504 
    505505   SUBROUTINE tra_asm_inc( kt, Kbb, Kmm, pts, Krhs ) 
    506506      !!---------------------------------------------------------------------- 
    507507      !!                    ***  ROUTINE tra_asm_inc  *** 
    508       !!           
     508      !! 
    509509      !! ** Purpose : Apply the tracer (T and S) assimilation increments 
    510510      !! 
    511511      !! ** Method  : Direct initialization or Incremental Analysis Updating 
    512512      !! 
    513       !! ** Action  :  
     513      !! ** Action  : 
    514514      !!---------------------------------------------------------------------- 
    515515      INTEGER                                  , INTENT(in   ) :: kt             ! Current time step 
     
    523523      !!---------------------------------------------------------------------- 
    524524      ! 
    525       ! freezing point calculation taken from oc_fz_pt (but calculated for all depths)  
    526       ! used to prevent the applied increments taking the temperature below the local freezing point  
     525      ! freezing point calculation taken from oc_fz_pt (but calculated for all depths) 
     526      ! used to prevent the applied increments taking the temperature below the local freezing point 
    527527      DO jk = 1, jpkm1 
    528528        CALL eos_fzp( pts(:,:,jk,jp_sal,Kmm), fzptnz(:,:,jk), gdept(:,:,jk,Kmm) ) 
     
    539539            ! 
    540540            IF(lwp) THEN 
    541                WRITE(numout,*)  
     541               WRITE(numout,*) 
    542542               WRITE(numout,*) 'tra_asm_inc : Tracer IAU at time step = ', kt,' with IAU weight = ', wgtiau(it) 
    543543               WRITE(numout,*) '~~~~~~~~~~~~' 
     
    549549                  ! Do not apply negative increments if the temperature will fall below freezing 
    550550                  WHERE(t_bkginc(:,:,jk) > 0.0_wp .OR. & 
    551                      &   pts(:,:,jk,jp_tem,Kmm) + pts(:,:,jk,jp_tem,Krhs) + t_bkginc(:,:,jk) * wgtiau(it) > fzptnz(:,:,jk) )  
    552                      pts(:,:,jk,jp_tem,Krhs) = pts(:,:,jk,jp_tem,Krhs) + t_bkginc(:,:,jk) * zincwgt   
     551                     &   pts(:,:,jk,jp_tem,Kmm) + pts(:,:,jk,jp_tem,Krhs) + t_bkginc(:,:,jk) * wgtiau(it) > fzptnz(:,:,jk) ) 
     552                     pts(:,:,jk,jp_tem,Krhs) = pts(:,:,jk,jp_tem,Krhs) + t_bkginc(:,:,jk) * zincwgt 
    553553                  END WHERE 
    554554               ELSE 
    555                   pts(:,:,jk,jp_tem,Krhs) = pts(:,:,jk,jp_tem,Krhs) + t_bkginc(:,:,jk) * zincwgt   
     555                  pts(:,:,jk,jp_tem,Krhs) = pts(:,:,jk,jp_tem,Krhs) + t_bkginc(:,:,jk) * zincwgt 
    556556               ENDIF 
    557557               IF (ln_salfix) THEN 
     
    559559                  ! minimum value salfixmin 
    560560                  WHERE(s_bkginc(:,:,jk) > 0.0_wp .OR. & 
    561                      &   pts(:,:,jk,jp_sal,Kmm) + pts(:,:,jk,jp_sal,Krhs) + s_bkginc(:,:,jk) * wgtiau(it) > salfixmin )  
     561                     &   pts(:,:,jk,jp_sal,Kmm) + pts(:,:,jk,jp_sal,Krhs) + s_bkginc(:,:,jk) * wgtiau(it) > salfixmin ) 
    562562                     pts(:,:,jk,jp_sal,Krhs) = pts(:,:,jk,jp_sal,Krhs) + s_bkginc(:,:,jk) * zincwgt 
    563563                  END WHERE 
     
    576576      ELSEIF ( ln_asmdin ) THEN        ! Direct Initialization 
    577577         !                             !-------------------------------------- 
    578          !             
     578         ! 
    579579         IF ( kt == nitdin_r ) THEN 
    580580            ! 
     
    584584            IF (ln_temnofreeze) THEN 
    585585               ! Do not apply negative increments if the temperature will fall below freezing 
    586                WHERE( t_bkginc(:,:,:) > 0.0_wp .OR. pts(:,:,:,jp_tem,Kmm) + t_bkginc(:,:,:) > fzptnz(:,:,:) )  
    587                   pts(:,:,:,jp_tem,Kmm) = t_bkg(:,:,:) + t_bkginc(:,:,:)    
     586               WHERE( t_bkginc(:,:,:) > 0.0_wp .OR. pts(:,:,:,jp_tem,Kmm) + t_bkginc(:,:,:) > fzptnz(:,:,:) ) 
     587                  pts(:,:,:,jp_tem,Kmm) = t_bkg(:,:,:) + t_bkginc(:,:,:) 
    588588               END WHERE 
    589589            ELSE 
    590                pts(:,:,:,jp_tem,Kmm) = t_bkg(:,:,:) + t_bkginc(:,:,:)    
     590               pts(:,:,:,jp_tem,Kmm) = t_bkg(:,:,:) + t_bkginc(:,:,:) 
    591591            ENDIF 
    592592            IF (ln_salfix) THEN 
    593593               ! Do not apply negative increments if the salinity will fall below a specified 
    594594               ! minimum value salfixmin 
    595                WHERE( s_bkginc(:,:,:) > 0.0_wp .OR. pts(:,:,:,jp_sal,Kmm) + s_bkginc(:,:,:) > salfixmin )  
    596                   pts(:,:,:,jp_sal,Kmm) = s_bkg(:,:,:) + s_bkginc(:,:,:)    
     595               WHERE( s_bkginc(:,:,:) > 0.0_wp .OR. pts(:,:,:,jp_sal,Kmm) + s_bkginc(:,:,:) > salfixmin ) 
     596                  pts(:,:,:,jp_sal,Kmm) = s_bkg(:,:,:) + s_bkginc(:,:,:) 
    597597               END WHERE 
    598598            ELSE 
    599                pts(:,:,:,jp_sal,Kmm) = s_bkg(:,:,:) + s_bkginc(:,:,:)    
     599               pts(:,:,:,jp_sal,Kmm) = s_bkg(:,:,:) + s_bkginc(:,:,:) 
    600600            ENDIF 
    601601 
     
    619619            DEALLOCATE( s_bkg    ) 
    620620         ENDIF 
    621          !   
     621         ! 
    622622      ENDIF 
    623623      ! Perhaps the following call should be in step 
     
    630630      !!---------------------------------------------------------------------- 
    631631      !!                    ***  ROUTINE dyn_asm_inc  *** 
    632       !!           
     632      !! 
    633633      !! ** Purpose : Apply the dynamics (u and v) assimilation increments. 
    634634      !! 
    635635      !! ** Method  : Direct initialization or Incremental Analysis Updating. 
    636636      !! 
    637       !! ** Action  :  
     637      !! ** Action  : 
    638638      !!---------------------------------------------------------------------- 
    639639      INTEGER                             , INTENT( in )  ::  kt             ! ocean time-step index 
     
    656656            ! 
    657657            IF(lwp) THEN 
    658                WRITE(numout,*)  
     658               WRITE(numout,*) 
    659659               WRITE(numout,*) 'dyn_asm_inc : Dynamics IAU at time step = ', kt,' with IAU weight = ', wgtiau(it) 
    660660               WRITE(numout,*) '~~~~~~~~~~~~' 
     
    676676      ELSEIF ( ln_asmdin ) THEN     ! Direct Initialization 
    677677         !                          !----------------------------------------- 
    678          !          
     678         ! 
    679679         IF ( kt == nitdin_r ) THEN 
    680680            ! 
     
    683683            ! Initialize the now fields with the background + increment 
    684684            puu(:,:,:,Kmm) = u_bkg(:,:,:) + u_bkginc(:,:,:) 
    685             pvv(:,:,:,Kmm) = v_bkg(:,:,:) + v_bkginc(:,:,:)   
     685            pvv(:,:,:,Kmm) = v_bkg(:,:,:) + v_bkginc(:,:,:) 
    686686            ! 
    687687            puu(:,:,:,Kbb) = puu(:,:,:,Kmm)         ! Update before fields 
     
    702702      !!---------------------------------------------------------------------- 
    703703      !!                    ***  ROUTINE ssh_asm_inc  *** 
    704       !!           
     704      !! 
    705705      !! ** Purpose : Apply the sea surface height assimilation increment. 
    706706      !! 
    707707      !! ** Method  : Direct initialization or Incremental Analysis Updating. 
    708708      !! 
    709       !! ** Action  :  
     709      !! ** Action  : 
    710710      !!---------------------------------------------------------------------- 
    711711      INTEGER, INTENT(IN) :: kt         ! Current time step 
     
    727727            ! 
    728728            IF(lwp) THEN 
    729                WRITE(numout,*)  
     729               WRITE(numout,*) 
    730730               WRITE(numout,*) 'ssh_asm_inc : SSH IAU at time step = ', & 
    731731                  &  kt,' with IAU weight = ', wgtiau(it) 
     
    779779      !!                  ***  ROUTINE ssh_asm_div  *** 
    780780      !! 
    781       !! ** Purpose :   ssh increment with z* is incorporated via a correction of the local divergence           
     781      !! ** Purpose :   ssh increment with z* is incorporated via a correction of the local divergence 
    782782      !!                across all the water column 
    783783      !! 
     
    795795      REAL(wp), DIMENSION(:,:)  , POINTER       ::   ztim     ! local array 
    796796      !!---------------------------------------------------------------------- 
    797       !  
     797      ! 
    798798#if defined key_asminc 
    799799      CALL ssh_asm_inc( kt, Kbb, Kmm ) !==   (calculate increments) 
    800800      ! 
    801       IF( ln_linssh ) THEN  
     801      IF( ln_linssh ) THEN 
    802802         phdivn(:,:,1) = phdivn(:,:,1) - ssh_iau(:,:) / e3t(:,:,1,Kmm) * tmask(:,:,1) 
    803       ELSE  
     803      ELSE 
    804804         ALLOCATE( ztim(jpi,jpj) ) 
    805805         ztim(:,:) = ssh_iau(:,:) / ( ht(:,:) + 1.0 - ssmask(:,:) ) 
    806          DO jk = 1, jpkm1                                  
    807             phdivn(:,:,jk) = phdivn(:,:,jk) - ztim(:,:) * tmask(:,:,jk)  
     806         DO jk = 1, jpkm1 
     807            phdivn(:,:,jk) = phdivn(:,:,jk) - ztim(:,:) * tmask(:,:,jk) 
    808808         END DO 
    809809         ! 
     
    818818      !!---------------------------------------------------------------------- 
    819819      !!                    ***  ROUTINE seaice_asm_inc  *** 
    820       !!           
     820      !! 
    821821      !! ** Purpose : Apply the sea ice assimilation increment. 
    822822      !! 
    823823      !! ** Method  : Direct initialization or Incremental Analysis Updating. 
    824824      !! 
    825       !! ** Action  :  
     825      !! ** Action  : 
    826826      !! 
    827827      !!---------------------------------------------------------------------- 
     
    844844            ! 
    845845            it = kt - nit000 + 1 
    846             zincwgt = wgtiau(it)      ! IAU weight for the current time step  
     846            zincwgt = wgtiau(it)      ! IAU weight for the current time step 
    847847            ! note this is not a tendency so should not be divided by rdt (as with the tracer and other increments) 
    848848            ! 
    849849            IF(lwp) THEN 
    850                WRITE(numout,*)  
     850               WRITE(numout,*) 
    851851               WRITE(numout,*) 'seaice_asm_inc : sea ice conc IAU at time step = ', kt,' with IAU weight = ', wgtiau(it) 
    852852               WRITE(numout,*) '~~~~~~~~~~~~' 
     
    866866            ! 
    867867            ! Nudge sea ice depth to bring it up to a required minimum depth 
    868             WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin )  
    869                zhicifinc(:,:) = (zhicifmin - hm_i(:,:)) * zincwgt     
     868            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin ) 
     869               zhicifinc(:,:) = (zhicifmin - hm_i(:,:)) * zincwgt 
    870870            ELSEWHERE 
    871871               zhicifinc(:,:) = 0.0_wp 
     
    907907            zofrld (:,:) = 1._wp - at_i(:,:) 
    908908            zohicif(:,:) = hm_i(:,:) 
    909             !  
     909            ! 
    910910            ! Initialize the now fields the background + increment 
    911911            at_i(:,:) = 1. - MIN( MAX( 1.-at_i(:,:) - seaice_bkginc(:,:), 0.0_wp), 1.0_wp) 
    912             at_i_b(:,:) = at_i(:,:)  
     912            at_i_b(:,:) = at_i(:,:) 
    913913            fr_i(:,:) = at_i(:,:)        ! adjust ice fraction 
    914914            ! 
     
    916916            ! 
    917917            ! Nudge sea ice depth to bring it up to a required minimum depth 
    918             WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin )  
     918            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin ) 
    919919               zhicifinc(:,:) = zhicifmin - hm_i(:,:) 
    920920            ELSEWHERE 
     
    946946!#if defined defined key_si3 || defined key_cice 
    947947! 
    948 !            IF (ln_seaicebal ) THEN        
     948!            IF (ln_seaicebal ) THEN 
    949949!             !! balancing salinity increments 
    950950!             !! simple case from limflx.F90 (doesn't include a mass flux) 
     
    958958! 
    959959!             DO jj = 1, jpj 
    960 !               DO ji = 1, jpi  
     960!               DO ji = 1, jpi 
    961961!           ! calculate change in ice and snow mass per unit area 
    962962!           ! positive values imply adding salt to the ocean (results from ice formation) 
     
    969969! 
    970970!           ! prevent small mld 
    971 !           ! less than 10m can cause salinity instability  
     971!           ! less than 10m can cause salinity instability 
    972972!                 IF (mld < 10) mld=10 
    973973! 
    974 !           ! set to bottom of a level  
     974!           ! set to bottom of a level 
    975975!                 DO jk = jpk-1, 2, -1 
    976 !                   IF ((mld > gdepw(ji,jj,jk)) .and. (mld < gdepw(ji,jj,jk+1))) THEN  
     976!                   IF ((mld > gdepw(ji,jj,jk)) .and. (mld < gdepw(ji,jj,jk+1))) THEN 
    977977!                     mld=gdepw(ji,jj,jk+1) 
    978978!                     jkmax=jk 
     
    981981! 
    982982!            ! avoid applying salinity balancing in shallow water or on land 
    983 !            !  
     983!            ! 
    984984! 
    985985!            ! dsal_ocn (psu kg m^-2) / (kg m^-3 * m) 
     
    992992! 
    993993!           ! put increments in for levels in the mixed layer 
    994 !           ! but prevent salinity below a threshold value  
    995 ! 
    996 !                   DO jk = 1, jkmax               
    997 ! 
    998 !                     IF (dsal_ocn > 0.0_wp .or. sb(ji,jj,jk)+dsal_ocn > sal_thresh) THEN  
     994!           ! but prevent salinity below a threshold value 
     995! 
     996!                   DO jk = 1, jkmax 
     997! 
     998!                     IF (dsal_ocn > 0.0_wp .or. sb(ji,jj,jk)+dsal_ocn > sal_thresh) THEN 
    999999!                           sb(ji,jj,jk) = sb(ji,jj,jk) + dsal_ocn 
    10001000!                           sn(ji,jj,jk) = sn(ji,jj,jk) + dsal_ocn 
     
    10071007!      ! 
    10081008!      !! Adjust fsalt. A +ve fsalt means adding salt to ocean 
    1009 !      !!           fsalt(ji,jj) =  fsalt(ji,jj) + zpmess     ! adjust fsalt   
    1010 !      !!                
    1011 !      !!           emps(ji,jj) = emps(ji,jj) + zpmess        ! or adjust emps (see icestp1d)  
     1009!      !!           fsalt(ji,jj) =  fsalt(ji,jj) + zpmess     ! adjust fsalt 
     1010!      !! 
     1011!      !!           emps(ji,jj) = emps(ji,jj) + zpmess        ! or adjust emps (see icestp1d) 
    10121012!      !!                                                     ! E-P (kg m-2 s-2) 
    10131013!      !            emp(ji,jj) = emp(ji,jj) + zpmess          ! E-P (kg m-2 s-2) 
     
    10221022      ! 
    10231023   END SUBROUTINE seaice_asm_inc 
    1024     
     1024 
    10251025   !!====================================================================== 
    10261026END MODULE asminc 
Note: See TracChangeset for help on using the changeset viewer.