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.
iceupdate.F90 in NEMO/branches/UKMO/NEMO_4.0.4_remade_heat_balance/src/ICE – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.4_remade_heat_balance/src/ICE/iceupdate.F90 @ 14615

Last change on this file since 14615 was 14615, checked in by dancopsey, 3 years ago

Added missing ampersand

File size: 28.4 KB
Line 
1MODULE iceupdate
2   !!======================================================================
3   !!                       ***  MODULE iceupdate   ***
4   !!  Sea-ice :   computation of the flux at the sea ice/ocean interface
5   !!======================================================================
6   !! History :  4.0  !  2018     (many people)      SI3 [aka Sea Ice cube]
7   !!----------------------------------------------------------------------
8#if defined key_si3
9   !!----------------------------------------------------------------------
10   !!   'key_si3'                                       SI3 sea-ice model
11   !!----------------------------------------------------------------------
12   !!   ice_update_alloc : allocate the iceupdate arrays
13   !!   ice_update_init  : initialisation
14   !!   ice_update_flx   : updates mass, heat and salt fluxes at the ocean surface
15   !!   ice_update_tau   : update i- and j-stresses, and its modulus at the ocean surface
16   !!----------------------------------------------------------------------
17   USE oce     , ONLY : sshn, sshb
18   USE phycst         ! physical constants
19   USE dom_oce        ! ocean domain
20   USE ice            ! sea-ice: variables
21   USE sbc_ice        ! Surface boundary condition: ice   fields
22   USE sbc_oce        ! Surface boundary condition: ocean fields
23   USE sbccpl         ! Surface boundary condition: coupled interface
24   USE icealb         ! sea-ice: albedo parameters
25   USE traqsr         ! add penetration of solar flux in the calculation of heat budget
26   USE icectl         ! sea-ice: control prints
27   USE zdfdrg  , ONLY : ln_drgice_imp
28   !
29   USE in_out_manager ! I/O manager
30   USE iom            ! I/O manager library
31   USE lib_mpp        ! MPP library
32   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
33   USE lbclnk         ! lateral boundary conditions (or mpp links)
34   USE timing         ! Timing
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   ice_update_init   ! called by ice_init
40   PUBLIC   ice_update_flx    ! called by ice_stp
41   PUBLIC   ice_update_tau    ! called by ice_stp
42
43   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   utau_oce, vtau_oce   ! air-ocean surface i- & j-stress     [N/m2]
44   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   tmod_io              ! modulus of the ice-ocean velocity   [m/s]
45
46   !! * Substitutions
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
49   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
50   !! $Id$
51   !! Software governed by the CeCILL license (see ./LICENSE)
52   !!----------------------------------------------------------------------
53CONTAINS
54
55   INTEGER FUNCTION ice_update_alloc()
56      !!-------------------------------------------------------------------
57      !!             ***  ROUTINE ice_update_alloc ***
58      !!-------------------------------------------------------------------
59      ALLOCATE( utau_oce(jpi,jpj), vtau_oce(jpi,jpj), tmod_io(jpi,jpj), STAT=ice_update_alloc )
60      !
61      CALL mpp_sum( 'iceupdate', ice_update_alloc )
62      IF( ice_update_alloc /= 0 )   CALL ctl_stop( 'STOP', 'ice_update_alloc: failed to allocate arrays' )
63      !
64   END FUNCTION ice_update_alloc
65
66
67   SUBROUTINE ice_update_flx( kt )
68      !!-------------------------------------------------------------------
69      !!                ***  ROUTINE ice_update_flx ***
70      !! 
71      !! ** Purpose :   Update the surface ocean boundary condition for heat
72      !!                salt and mass over areas where sea-ice is non-zero
73      !!         
74      !! ** Action  : - computes the heat and freshwater/salt fluxes
75      !!                at the ice-ocean interface.
76      !!              - Update the ocean sbc
77      !!     
78      !! ** Outputs : - qsr     : sea heat flux:     solar
79      !!              - qns     : sea heat flux: non solar
80      !!              - emp     : freshwater budget: volume flux
81      !!              - sfx     : salt flux
82      !!              - fr_i    : ice fraction
83      !!              - tn_ice  : sea-ice surface temperature
84      !!              - alb_ice : sea-ice albedo (recomputed only for coupled mode)
85      !!
86      !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90.
87      !!              Tartinville et al. 2001 Ocean Modelling, 3, 95-108.
88      !!              These refs are now obsolete since everything has been revised
89      !!              The ref should be Rousset et al., 2015
90      !!---------------------------------------------------------------------
91      INTEGER, INTENT(in) ::   kt   ! number of iteration
92      !
93      INTEGER  ::   ji, jj, jl, jk   ! dummy loop indices
94      REAL(wp) ::   zqsr             ! New solar flux received by the ocean
95      REAL(wp), DIMENSION(jpi,jpj) ::   z2d                  ! 2D workspace
96      !!---------------------------------------------------------------------
97      IF( ln_timing )   CALL timing_start('ice_update')
98
99      IF( kt == nit000 .AND. lwp ) THEN
100         WRITE(numout,*)
101         WRITE(numout,*)'ice_update_flx: update fluxes (mass, salt and heat) at the ice-ocean interface'
102         WRITE(numout,*)'~~~~~~~~~~~~~~'
103      ENDIF
104
105      ! Net heat flux on top of the ice-ocean (W.m-2)
106      !----------------------------------------------
107      qt_atm_oi(:,:) = qns_tot(:,:) + qsr_tot(:,:) 
108
109      ! --- case we bypass ice thermodynamics --- !
110      IF( .NOT. ln_icethd ) THEN   ! we suppose ice is impermeable => ocean is isolated from atmosphere
111         qt_atm_oi  (:,:)   = ( 1._wp - at_i_b(:,:) ) * ( qns_oce(:,:) + qsr_oce(:,:) ) + qemp_oce(:,:)
112         qt_oce_ai  (:,:)   = ( 1._wp - at_i_b(:,:) ) *   qns_oce(:,:)                  + qemp_oce(:,:)
113         emp_ice    (:,:)   = 0._wp
114         qemp_ice   (:,:)   = 0._wp
115         qevap_ice  (:,:,:) = 0._wp
116      ENDIF
117     
118      DO jj = 1, jpj
119         DO ji = 1, jpi
120
121            ! Solar heat flux reaching the ocean (max) = zqsr (W.m-2)
122            !---------------------------------------------------
123            zqsr = qsr_tot(ji,jj) - SUM( a_i_b(ji,jj,:) * ( qsr_ice(ji,jj,:) - qtr_ice_bot(ji,jj,:) ) )
124
125            ! Total heat flux reaching the ocean = qt_oce_ai (W.m-2)
126            ! It is calculated using the total heat into the top of the ocean and ice (qt_atm_oi) minus energy taken up heating the ice interior (hfx_dif),
127            ! plus energy released by ice dynamics (hfx_dyn),
128            ! minus energy used up in top melt of snow and ice (hfx_snw & hfx_sum),
129            ! minus energy needed for melt/growth along the ocean/ice interface (hfx_bom, hfx_bog, hfx_lam, hfx_opw),
130            ! plus energy stored in snow precipitation (hfx_spr). Even though this isn't a real heat flux it is included in the surface
131            ! energy balance as zqemp_ice in sbccpl.F90 so needs to be cancelled out here.
132            ! plus the residual (hfx_res)
133            !---------------------------------------------------
134            qt_oce_ai(ji,jj) = qt_atm_oi(ji,jj) - hfx_dif(ji,jj) + hfx_dyn(ji,jj)                                   &
135               &                                - hfx_snw(ji,jj) - hfx_sum(ji,jj)                                   &
136               &                                - hfx_bom(ji,jj) - hfx_bog(ji,jj) - hfx_lam(ji,jj) - hfx_opw(ji,jj) &
137               &                                + hfx_spr(ji,jj)                                                    &
138               &                                + hfx_res(ji,jj)             
139           
140            ! New qsr and qns used to compute the oceanic heat flux at the next time step
141            !----------------------------------------------------------------------------
142            ! if warming and some ice remains, then we suppose that the whole solar flux has been consumed to melt the ice
143            ! else ( cooling or no ice left ), then we suppose that     no    solar flux has been consumed
144            !
145            IF( fhld(ji,jj) > 0._wp .AND. at_i(ji,jj) > 0._wp ) THEN   !-- warming and some ice remains
146               !                                        solar flux transmitted thru the 1st level of the ocean (i.e. not used by sea-ice)
147               qsr(ji,jj) = ( 1._wp - at_i_b(ji,jj) ) * qsr_oce(ji,jj) * ( 1._wp - frq_m(ji,jj) ) &
148                  !                                   + solar flux transmitted thru ice and the 1st ocean level (also not used by sea-ice)
149                  &             + SUM( a_i_b(ji,jj,:) * qtr_ice_bot(ji,jj,:) ) * ( 1._wp - frq_m(ji,jj) )
150               !
151            ELSE                                                       !-- cooling or no ice left
152               qsr(ji,jj) = zqsr
153            ENDIF
154            !
155            ! the non-solar is simply derived from the solar flux
156            qns(ji,jj) = qt_oce_ai(ji,jj) - zqsr             
157
158            ! Mass flux at the atm. surface       
159            !-----------------------------------
160            wfx_sub(ji,jj) = wfx_snw_sub(ji,jj) + wfx_ice_sub(ji,jj)
161
162            ! Mass flux at the ocean surface     
163            !------------------------------------
164            ! ice-ocean  mass flux
165            wfx_ice(ji,jj) = wfx_bog(ji,jj) + wfx_bom(ji,jj) + wfx_sum(ji,jj) + wfx_sni(ji,jj)   &
166               &           + wfx_opw(ji,jj) + wfx_dyn(ji,jj) + wfx_res(ji,jj) + wfx_lam(ji,jj) + wfx_pnd(ji,jj)
167
168            ! snw-ocean mass flux
169            wfx_snw(ji,jj) = wfx_snw_sni(ji,jj) + wfx_snw_dyn(ji,jj) + wfx_snw_sum(ji,jj)
170
171            ! total mass flux at the ocean/ice interface
172            fmmflx(ji,jj) =                - wfx_ice(ji,jj) - wfx_snw(ji,jj) - wfx_err_sub(ji,jj)   ! ice-ocean mass flux saved at least for biogeochemical model
173            emp   (ji,jj) = emp_oce(ji,jj) - wfx_ice(ji,jj) - wfx_snw(ji,jj) - wfx_err_sub(ji,jj)   ! atm-ocean + ice-ocean mass flux
174
175            ! Salt flux at the ocean surface     
176            !------------------------------------------
177            sfx(ji,jj) = sfx_bog(ji,jj) + sfx_bom(ji,jj) + sfx_sum(ji,jj) + sfx_sni(ji,jj) + sfx_opw(ji,jj)   &
178               &       + sfx_res(ji,jj) + sfx_dyn(ji,jj) + sfx_bri(ji,jj) + sfx_sub(ji,jj) + sfx_lam(ji,jj)
179           
180            ! Mass of snow and ice per unit area   
181            !----------------------------------------
182            snwice_mass_b(ji,jj) = snwice_mass(ji,jj)       ! save mass from the previous ice time step
183            !                                               ! new mass per unit area
184            snwice_mass  (ji,jj) = tmask(ji,jj,1) * ( rhos * vt_s(ji,jj) + rhoi * vt_i(ji,jj)  ) 
185            !                                               ! time evolution of snow+ice mass
186            snwice_fmass (ji,jj) = ( snwice_mass(ji,jj) - snwice_mass_b(ji,jj) ) * r1_rdtice
187           
188         END DO
189      END DO
190
191      ! Storing the transmitted variables
192      !----------------------------------
193      fr_i  (:,:)   = at_i(:,:)             ! Sea-ice fraction           
194      tn_ice(:,:,:) = t_su(:,:,:)           ! Ice surface temperature                     
195
196      ! Snow/ice albedo (only if sent to coupler, useless in forced mode)
197      !------------------------------------------------------------------
198      CALL ice_alb( t_su, h_i, h_s, ln_pnd_alb, a_ip_eff, h_ip, cloud_fra, alb_ice ) ! ice albedo
199
200      !
201      IF( lrst_ice ) THEN                       !* write snwice_mass fields in the restart file
202         CALL update_rst( 'WRITE', kt )
203      ENDIF
204      !
205      ! output all fluxes
206      !------------------
207      !
208      ! --- salt fluxes [kg/m2/s] --- !
209      !                           ! sfxice =  sfxbog + sfxbom + sfxsum + sfxsni + sfxopw + sfxres + sfxdyn + sfxbri + sfxsub + sfxlam
210      IF( iom_use('sfxice'  ) )   CALL iom_put( 'sfxice', sfx     * 1.e-03 )   ! salt flux from total ice growth/melt
211      IF( iom_use('sfxbog'  ) )   CALL iom_put( 'sfxbog', sfx_bog * 1.e-03 )   ! salt flux from bottom growth
212      IF( iom_use('sfxbom'  ) )   CALL iom_put( 'sfxbom', sfx_bom * 1.e-03 )   ! salt flux from bottom melting
213      IF( iom_use('sfxsum'  ) )   CALL iom_put( 'sfxsum', sfx_sum * 1.e-03 )   ! salt flux from surface melting
214      IF( iom_use('sfxlam'  ) )   CALL iom_put( 'sfxlam', sfx_lam * 1.e-03 )   ! salt flux from lateral melting
215      IF( iom_use('sfxsni'  ) )   CALL iom_put( 'sfxsni', sfx_sni * 1.e-03 )   ! salt flux from snow ice formation
216      IF( iom_use('sfxopw'  ) )   CALL iom_put( 'sfxopw', sfx_opw * 1.e-03 )   ! salt flux from open water formation
217      IF( iom_use('sfxdyn'  ) )   CALL iom_put( 'sfxdyn', sfx_dyn * 1.e-03 )   ! salt flux from ridging rafting
218      IF( iom_use('sfxbri'  ) )   CALL iom_put( 'sfxbri', sfx_bri * 1.e-03 )   ! salt flux from brines
219      IF( iom_use('sfxres'  ) )   CALL iom_put( 'sfxres', sfx_res * 1.e-03 )   ! salt flux from undiagnosed processes
220      IF( iom_use('sfxsub'  ) )   CALL iom_put( 'sfxsub', sfx_sub * 1.e-03 )   ! salt flux from sublimation
221
222      ! --- mass fluxes [kg/m2/s] --- !
223      CALL iom_put( 'emp_oce', emp_oce )   ! emp over ocean (taking into account the snow blown away from the ice)
224      CALL iom_put( 'emp_ice', emp_ice )   ! emp over ice   (taking into account the snow blown away from the ice)
225
226      !                           ! vfxice = vfxbog + vfxbom + vfxsum + vfxsni + vfxopw + vfxdyn + vfxres + vfxlam + vfxpnd
227      CALL iom_put( 'vfxice'    , wfx_ice     )   ! mass flux from total ice growth/melt
228      CALL iom_put( 'vfxbog'    , wfx_bog     )   ! mass flux from bottom growth
229      CALL iom_put( 'vfxbom'    , wfx_bom     )   ! mass flux from bottom melt
230      CALL iom_put( 'vfxsum'    , wfx_sum     )   ! mass flux from surface melt
231      CALL iom_put( 'vfxlam'    , wfx_lam     )   ! mass flux from lateral melt
232      CALL iom_put( 'vfxsni'    , wfx_sni     )   ! mass flux from snow-ice formation
233      CALL iom_put( 'vfxopw'    , wfx_opw     )   ! mass flux from growth in open water
234      CALL iom_put( 'vfxdyn'    , wfx_dyn     )   ! mass flux from dynamics (ridging)
235      CALL iom_put( 'vfxres'    , wfx_res     )   ! mass flux from undiagnosed processes
236      CALL iom_put( 'vfxpnd'    , wfx_pnd     )   ! mass flux from melt ponds
237      CALL iom_put( 'vfxsub'    , wfx_ice_sub )   ! mass flux from ice sublimation (ice-atm.)
238      CALL iom_put( 'vfxsub_err', wfx_err_sub )   ! "excess" of sublimation sent to ocean     
239
240      IF ( iom_use( 'vfxthin' ) ) THEN   ! mass flux from ice growth in open water + thin ice (<20cm) => comparable to observations 
241         WHERE( hm_i(:,:) < 0.2 .AND. hm_i(:,:) > 0. ) ; z2d = wfx_bog
242         ELSEWHERE                                     ; z2d = 0._wp
243         END WHERE
244         CALL iom_put( 'vfxthin', wfx_opw + z2d )
245      ENDIF
246
247      !                            ! vfxsnw = vfxsnw_sni + vfxsnw_dyn + vfxsnw_sum
248      CALL iom_put( 'vfxsnw'     , wfx_snw     )   ! mass flux from total snow growth/melt
249      CALL iom_put( 'vfxsnw_sum' , wfx_snw_sum )   ! mass flux from snow melt at the surface
250      CALL iom_put( 'vfxsnw_sni' , wfx_snw_sni )   ! mass flux from snow melt during snow-ice formation
251      CALL iom_put( 'vfxsnw_dyn' , wfx_snw_dyn )   ! mass flux from dynamics (ridging)
252      CALL iom_put( 'vfxsnw_sub' , wfx_snw_sub )   ! mass flux from snow sublimation (ice-atm.)
253      CALL iom_put( 'vfxsnw_pre' , wfx_spr     )   ! snow precip
254
255      ! --- heat fluxes [W/m2] --- !
256      !                              ! qt_atm_oi - qt_oce_ai = hfxdhc - ( dihctrp + dshctrp )
257      IF( iom_use('qsr_oce'    ) )   CALL iom_put( 'qsr_oce'    , qsr_oce * ( 1._wp - at_i_b )                               )   !     solar flux at ocean surface
258      IF( iom_use('qns_oce'    ) )   CALL iom_put( 'qns_oce'    , qns_oce * ( 1._wp - at_i_b ) + qemp_oce                    )   ! non-solar flux at ocean surface
259      IF( iom_use('qsr_ice'    ) )   CALL iom_put( 'qsr_ice'    , SUM( qsr_ice * a_i_b, dim=3 )                              )   !     solar flux at ice surface
260      IF( iom_use('qns_ice'    ) )   CALL iom_put( 'qns_ice'    , SUM( qns_ice * a_i_b, dim=3 ) + qemp_ice                   )   ! non-solar flux at ice surface
261      IF( iom_use('qtr_ice_bot') )   CALL iom_put( 'qtr_ice_bot', SUM( qtr_ice_bot * a_i_b, dim=3 )                          )   !     solar flux transmitted thru ice
262      IF( iom_use('qtr_ice_top') )   CALL iom_put( 'qtr_ice_top', SUM( qtr_ice_top * a_i_b, dim=3 )                          )   !     solar flux transmitted thru ice surface
263      IF( iom_use('qt_oce'     ) )   CALL iom_put( 'qt_oce'     ,      ( qsr_oce + qns_oce ) * ( 1._wp - at_i_b ) + qemp_oce )
264      IF( iom_use('qt_ice'     ) )   CALL iom_put( 'qt_ice'     , SUM( ( qns_ice + qsr_ice ) * a_i_b, dim=3 )     + qemp_ice )
265      IF( iom_use('qt_oce_ai'  ) )   CALL iom_put( 'qt_oce_ai'  , qt_oce_ai * tmask(:,:,1)                                   )   ! total heat flux at the ocean   surface: interface oce-(ice+atm)
266      IF( iom_use('qt_atm_oi'  ) )   CALL iom_put( 'qt_atm_oi'  , qt_atm_oi * tmask(:,:,1)                                   )   ! total heat flux at the oce-ice surface: interface atm-(ice+oce)
267      IF( iom_use('qemp_oce'   ) )   CALL iom_put( 'qemp_oce'   , qemp_oce                                                   )   ! Downward Heat Flux from E-P over ocean
268      IF( iom_use('qemp_ice'   ) )   CALL iom_put( 'qemp_ice'   , qemp_ice                                                   )   ! Downward Heat Flux from E-P over ice
269
270      ! heat fluxes from ice transformations
271      !                            ! hfxdhc = hfxbog + hfxbom + hfxsum + hfxopw + hfxdif + hfxsnw - ( hfxthd + hfxdyn + hfxres + hfxsub + hfxspr )
272      CALL iom_put ('hfxbog'     , hfx_bog     )   ! heat flux used for ice bottom growth
273      CALL iom_put ('hfxbom'     , hfx_bom     )   ! heat flux used for ice bottom melt
274      CALL iom_put ('hfxlam'     , hfx_lam     )   ! heat flux used for ice lateral melt
275      CALL iom_put ('hfxsum'     , hfx_sum     )   ! heat flux used for ice surface melt
276      CALL iom_put ('hfxopw'     , hfx_opw     )   ! heat flux used for ice formation in open water
277      CALL iom_put ('hfxdif'     , hfx_dif     )   ! heat flux used for ice temperature change
278      CALL iom_put ('hfxsnw'     , hfx_snw     )   ! heat flux used for snow melt
279      CALL iom_put ('hfxerr'     , hfx_err_dif )   ! heat flux error after heat diffusion
280
281      ! heat fluxes associated with mass exchange (freeze/melt/precip...)
282      CALL iom_put ('hfxthd'     , hfx_thd     )   
283      CALL iom_put ('hfxdyn'     , hfx_dyn     )   
284      CALL iom_put ('hfxres'     , hfx_res     )   
285      CALL iom_put ('hfxsub'     , hfx_sub     )   
286      CALL iom_put ('hfxspr'     , hfx_spr     )   ! Heat flux from snow precip heat content
287
288      ! other heat fluxes
289      IF( iom_use('hfxsensib'  ) )   CALL iom_put( 'hfxsensib'  ,     -qsb_ice_bot * at_i_b         )   ! Sensible oceanic heat flux
290      IF( iom_use('hfxcndbot'  ) )   CALL iom_put( 'hfxcndbot'  , SUM( qcn_ice_bot * a_i_b, dim=3 ) )   ! Bottom conduction flux
291      IF( iom_use('hfxcndtop'  ) )   CALL iom_put( 'hfxcndtop'  , SUM( qcn_ice_top * a_i_b, dim=3 ) )   ! Surface conduction flux
292
293      ! controls
294      !---------
295#if ! defined key_agrif
296      IF( ln_icediachk )   CALL ice_cons_final('iceupdate')                                       ! conservation
297#endif
298      IF( ln_icectl    )   CALL ice_prt       (kt, iiceprt, jiceprt, 3, 'Final state ice_update') ! prints
299      IF( ln_ctl       )   CALL ice_prt3D     ('iceupdate')                                       ! prints
300      IF( ln_timing    )   CALL timing_stop   ('ice_update')                                      ! timing
301      !
302   END SUBROUTINE ice_update_flx
303
304
305   SUBROUTINE ice_update_tau( kt, pu_oce, pv_oce )
306      !!-------------------------------------------------------------------
307      !!                ***  ROUTINE ice_update_tau ***
308      !! 
309      !! ** Purpose : Update the ocean surface stresses due to the ice
310      !!         
311      !! ** Action  : * at each ice time step (every nn_fsbc time step):
312      !!                - compute the modulus of ice-ocean relative velocity
313      !!                  (*rho*Cd) at T-point (C-grid) or I-point (B-grid)
314      !!                      tmod_io = rhoco * | U_ice-U_oce |
315      !!                - update the modulus of stress at ocean surface
316      !!                      taum = (1-a) * taum + a * tmod_io * | U_ice-U_oce |
317      !!              * at each ocean time step (every kt):
318      !!                  compute linearized ice-ocean stresses as
319      !!                      Utau = tmod_io * | U_ice - pU_oce |
320      !!                using instantaneous current ocean velocity (usually before)
321      !!
322      !!    NB: - ice-ocean rotation angle no more allowed
323      !!        - here we make an approximation: taum is only computed every ice time step
324      !!          This avoids mutiple average to pass from T -> U,V grids and next from U,V grids
325      !!          to T grid. taum is used in TKE and GLS, which should not be too sensitive to this approximaton...
326      !!
327      !! ** Outputs : - utau, vtau   : surface ocean i- and j-stress (u- & v-pts) updated with ice-ocean fluxes
328      !!              - taum         : modulus of the surface ocean stress (T-point) updated with ice-ocean fluxes
329      !!---------------------------------------------------------------------
330      INTEGER ,                     INTENT(in) ::   kt               ! ocean time-step index
331      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pu_oce, pv_oce   ! surface ocean currents
332      !
333      INTEGER  ::   ji, jj   ! dummy loop indices
334      REAL(wp) ::   zat_u, zutau_ice, zu_t, zmodt   ! local scalar
335      REAL(wp) ::   zat_v, zvtau_ice, zv_t, zrhoco  !   -      -
336      REAL(wp) ::   zflagi                          !   -      -
337      !!---------------------------------------------------------------------
338      IF( ln_timing )   CALL timing_start('ice_update_tau')
339
340      IF( kt == nit000 .AND. lwp ) THEN
341         WRITE(numout,*)
342         WRITE(numout,*)'ice_update_tau: update stress at the ice-ocean interface'
343         WRITE(numout,*)'~~~~~~~~~~~~~~'
344      ENDIF
345
346      zrhoco = rau0 * rn_cio
347      !
348      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !==  Ice time-step only  ==!   (i.e. surface module time-step)
349         DO jj = 2, jpjm1                             !* update the modulus of stress at ocean surface (T-point)
350            DO ji = fs_2, fs_jpim1
351               !                                               ! 2*(U_ice-U_oce) at T-point
352               zu_t = u_ice(ji,jj) + u_ice(ji-1,jj) - u_oce(ji,jj) - u_oce(ji-1,jj)   
353               zv_t = v_ice(ji,jj) + v_ice(ji,jj-1) - v_oce(ji,jj) - v_oce(ji,jj-1) 
354               !                                              ! |U_ice-U_oce|^2
355               zmodt =  0.25_wp * (  zu_t * zu_t + zv_t * zv_t  )
356               !                                               ! update the ocean stress modulus
357               taum(ji,jj) = ( 1._wp - at_i(ji,jj) ) * taum(ji,jj) + at_i(ji,jj) * zrhoco * zmodt
358               tmod_io(ji,jj) = zrhoco * SQRT( zmodt )          ! rhoco * |U_ice-U_oce| at T-point
359            END DO
360         END DO
361         CALL lbc_lnk_multi( 'iceupdate', taum, 'T', 1., tmod_io, 'T', 1. )
362         !
363         utau_oce(:,:) = utau(:,:)                    !* save the air-ocean stresses at ice time-step
364         vtau_oce(:,:) = vtau(:,:)
365         !
366      ENDIF
367      !
368      !                                      !==  every ocean time-step  ==!
369      IF ( ln_drgice_imp ) THEN
370         ! Save drag with right sign to update top drag in the ocean implicit friction
371         rCdU_ice(:,:) = -r1_rau0 * tmod_io(:,:) * at_i(:,:) * tmask(:,:,1) 
372         zflagi = 0._wp
373      ELSE
374         zflagi = 1._wp
375      ENDIF
376      !
377      DO jj = 2, jpjm1                                !* update the stress WITHOUT an ice-ocean rotation angle
378         DO ji = fs_2, fs_jpim1   ! Vect. Opt.   
379            ! ice area at u and v-points
380            zat_u  = ( at_i(ji,jj) * tmask(ji,jj,1) + at_i (ji+1,jj    ) * tmask(ji+1,jj  ,1) )  &
381               &     / MAX( 1.0_wp , tmask(ji,jj,1) + tmask(ji+1,jj  ,1) )
382            zat_v  = ( at_i(ji,jj) * tmask(ji,jj,1) + at_i (ji  ,jj+1  ) * tmask(ji  ,jj+1,1) )  &
383               &     / MAX( 1.0_wp , tmask(ji,jj,1) + tmask(ji  ,jj+1,1) )
384            !                                                   ! linearized quadratic drag formulation
385            zutau_ice   = 0.5_wp * ( tmod_io(ji,jj) + tmod_io(ji+1,jj) ) * ( u_ice(ji,jj) - zflagi * pu_oce(ji,jj) )
386            zvtau_ice   = 0.5_wp * ( tmod_io(ji,jj) + tmod_io(ji,jj+1) ) * ( v_ice(ji,jj) - zflagi * pv_oce(ji,jj) )
387            !                                                   ! stresses at the ocean surface
388            utau(ji,jj) = ( 1._wp - zat_u ) * utau_oce(ji,jj) + zat_u * zutau_ice
389            vtau(ji,jj) = ( 1._wp - zat_v ) * vtau_oce(ji,jj) + zat_v * zvtau_ice
390         END DO
391      END DO
392      CALL lbc_lnk_multi( 'iceupdate', utau, 'U', -1., vtau, 'V', -1. )   ! lateral boundary condition
393      !
394      IF( ln_timing )   CALL timing_stop('ice_update_tau')
395     
396   END SUBROUTINE ice_update_tau
397
398
399   SUBROUTINE ice_update_init
400      !!-------------------------------------------------------------------
401      !!                  ***  ROUTINE ice_update_init  ***
402      !!             
403      !! ** Purpose :   allocate ice-ocean stress fields and read restarts
404      !!                containing the snow & ice mass
405      !!
406      !!-------------------------------------------------------------------
407      INTEGER  ::   ji, jj, jk               ! dummy loop indices
408      REAL(wp) ::   zcoefu, zcoefv, zcoeff   ! local scalar
409      !!-------------------------------------------------------------------
410      !
411      IF(lwp) WRITE(numout,*)
412      IF(lwp) WRITE(numout,*) 'ice_update_init: ice-ocean stress init'
413      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~'
414      !
415      !                                      ! allocate ice_update array
416      IF( ice_update_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'ice_update_init : unable to allocate standard arrays' )
417      !
418      CALL update_rst( 'READ' )  !* read or initialize all required files
419      !
420   END SUBROUTINE ice_update_init
421
422
423   SUBROUTINE update_rst( cdrw, kt )
424      !!---------------------------------------------------------------------
425      !!                   ***  ROUTINE rhg_evp_rst  ***
426      !!                     
427      !! ** Purpose :   Read or write RHG file in restart file
428      !!
429      !! ** Method  :   use of IOM library
430      !!----------------------------------------------------------------------
431      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! 'READ'/'WRITE' flag
432      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
433      !
434      INTEGER  ::   iter   ! local integer
435      INTEGER  ::   id1    ! local integer
436      !!----------------------------------------------------------------------
437      !
438      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
439         !                                   ! ---------------
440         IF( ln_rstart ) THEN                   !* Read the restart file
441            !
442            id1 = iom_varid( numrir, 'snwice_mass' , ldstop = .FALSE. )
443            !
444            IF( id1 > 0 ) THEN                       ! fields exist
445               CALL iom_get( numrir, jpdom_autoglo, 'snwice_mass'  , snwice_mass   )
446               CALL iom_get( numrir, jpdom_autoglo, 'snwice_mass_b', snwice_mass_b )
447            ELSE                                     ! start from rest
448               IF(lwp) WRITE(numout,*) '   ==>>   previous run without snow-ice mass output then set it'
449               snwice_mass  (:,:) = tmask(:,:,1) * ( rhos * vt_s(:,:) + rhoi * vt_i(:,:) )
450               snwice_mass_b(:,:) = snwice_mass(:,:)
451            ENDIF
452         ELSE                                   !* Start from rest
453            IF(lwp) WRITE(numout,*) '   ==>>   start from rest: set the snow-ice mass'
454            snwice_mass  (:,:) = tmask(:,:,1) * ( rhos * vt_s(:,:) + rhoi * vt_i(:,:) )
455            snwice_mass_b(:,:) = snwice_mass(:,:)
456         ENDIF
457         !
458      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
459         !                                   ! -------------------
460         IF(lwp) WRITE(numout,*) '---- update-rst ----'
461         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
462         !
463         CALL iom_rstput( iter, nitrst, numriw, 'snwice_mass'  , snwice_mass   )
464         CALL iom_rstput( iter, nitrst, numriw, 'snwice_mass_b', snwice_mass_b )
465         !
466      ENDIF
467      !
468   END SUBROUTINE update_rst
469
470#else
471   !!----------------------------------------------------------------------
472   !!   Default option         Dummy module           NO SI3 sea-ice model
473   !!----------------------------------------------------------------------
474#endif 
475
476   !!======================================================================
477END MODULE iceupdate
Note: See TracBrowser for help on using the repository browser.