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.
icbthm.F90 in NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/ICB – NEMO

source: NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/ICB/icbthm.F90 @ 14038

Last change on this file since 14038 was 14038, checked in by laurent, 4 years ago

Catch up with trunk at rev r14037

  • Property svn:keywords set to Id
File size: 15.0 KB
Line 
1MODULE icbthm
2   !!======================================================================
3   !!                       ***  MODULE  icbthm  ***
4   !! Icebergs:  thermodynamics routines for icebergs
5   !!======================================================================
6   !! History : 3.3.1 !  2010-01  (Martin&Adcroft) Original code
7   !!            -    !  2011-03  (Madec)          Part conversion to NEMO form
8   !!            -    !                            Removal of mapping from another grid
9   !!            -    !  2011-04  (Alderson)       Split into separate modules
10   !!            -    !  2011-05  (Alderson)       Use tmask instead of tmask_i
11   !!----------------------------------------------------------------------
12   !!----------------------------------------------------------------------
13   !!   icb_thm : initialise
14   !!             reference for equations - M = Martin + Adcroft, OM 34, 2010
15   !!----------------------------------------------------------------------
16   USE par_oce        ! NEMO parameters
17   USE dom_oce        ! NEMO domain
18   USE in_out_manager ! NEMO IO routines, numout in particular
19   USE lib_mpp        ! NEMO MPI routines, ctl_stop in particular
20   USE phycst         ! NEMO physical constants
21   USE sbc_oce
22   USE eosbn2         ! equation of state
23   USE lib_fortran, ONLY : DDPDD
24
25   USE icb_oce        ! define iceberg arrays
26   USE icbutl         ! iceberg utility routines
27   USE icbdia         ! iceberg budget routines
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   icb_thm ! routine called in icbstp.F90 module
33
34   !!----------------------------------------------------------------------
35   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
36   !! $Id$
37   !! Software governed by the CeCILL license (see ./LICENSE)
38   !!----------------------------------------------------------------------
39CONTAINS
40
41   SUBROUTINE icb_thm( kt )
42      !!----------------------------------------------------------------------
43      !!                  ***  ROUTINE icb_thm  ***
44      !!
45      !! ** Purpose :   compute the iceberg thermodynamics.
46      !!
47      !! ** Method  : - See Martin & Adcroft, Ocean Modelling 34, 2010
48      !!----------------------------------------------------------------------
49      INTEGER, INTENT(in) ::   kt   ! timestep number, just passed to icb_utl_print_berg
50      !
51      INTEGER  ::   ii, ij, jk, ikb
52      REAL(wp) ::   zM, zT, zW, zL, zSST, zVol, zLn, zWn, zTn, znVol, zIC, zDn, zD, zvb, zub, ztb
53      REAL(wp) ::   zMv, zMe, zMb, zmelt, zdvo, zdvob, zdva, zdM, zSs, zdMe, zdMb, zdMv
54      REAL(wp) ::   zSSS, zfzpt
55      REAL(wp) ::   zMnew, zMnew1, zMnew2, zheat_hcflux, zheat_latent, z1_12
56      REAL(wp) ::   zMbits, znMbits, zdMbitsE, zdMbitsM, zLbits, zAbits, zMbb
57      REAL(wp) ::   zxi, zyj, zff, z1_rday, z1_e1e2, zdt, z1_dt, z1_dt_e1e2, zdepw
58      REAL(wp), DIMENSION(jpk) :: ztoce, zuoce, zvoce, ze3t, zzMv
59      TYPE(iceberg), POINTER ::   this, next
60      TYPE(point)  , POINTER ::   pt
61      !
62      COMPLEX(dp), DIMENSION(jpi,jpj) :: cicb_melt, cicb_hflx
63      !!----------------------------------------------------------------------
64      !
65      !! initialiaze cicb_melt and cicb_heat
66      cicb_melt = CMPLX( 0.e0, 0.e0, dp ) 
67      cicb_hflx = CMPLX( 0.e0, 0.e0, dp ) 
68      !
69      z1_rday = 1._wp / rday
70      z1_12   = 1._wp / 12._wp
71      zdt     = berg_dt
72      z1_dt   = 1._wp / zdt
73      !
74      ! we're either going to ignore berg fresh water melt flux and associated heat
75      ! or we pass it into the ocean, so at this point we set them both to zero,
76      ! accumulate the contributions to them from each iceberg in the while loop following
77      ! and then pass them (or not) to the ocean
78      !
79      berg_grid%floating_melt(:,:) = 0._wp
80      ! calving_hflx re-used here as temporary workspace for the heat flux associated with melting
81      berg_grid%calving_hflx(:,:)  = 0._wp
82      !
83      this => first_berg
84      DO WHILE( ASSOCIATED(this) )
85         !
86         pt => this%current_point
87         nknberg = this%number(1)
88
89         CALL icb_utl_interp( pt%xi, pt%yj,            &   ! position
90             &                 pssu=pt%ssu, pua=pt%ua, &   ! oce/atm velocities
91             &                 pssv=pt%ssv, pva=pt%va, &   ! oce/atm velocities
92             &                 psst=pt%sst, pcn=pt%cn, &
93             &                 psss=pt%sss             )
94
95         IF ( nn_sample_rate > 0 .AND. MOD(kt-1,nn_sample_rate) == 0 ) THEN
96            CALL icb_utl_interp( pt%xi, pt%yj, pe1=pt%e1, pe2=pt%e2,                 &
97               &                 pui=pt%ui, pssh_i=pt%ssh_x, &
98               &                 pvi=pt%vi, pssh_j=pt%ssh_y, &
99               &                 phi=pt%hi,                  &
100               &                 plat=pt%lat, plon=pt%lon )
101         END IF
102         !
103         zSST = pt%sst
104         zSSS = pt%sss
105         CALL eos_fzp(zSSS,zfzpt)                       ! freezing point
106         zIC  = MIN( 1._wp, pt%cn + rn_sicn_shift )     ! Shift sea-ice concentration       !!gm ???
107         zM   = pt%mass
108         zT   = pt%thickness                               ! total thickness
109         zD   = rho_berg_1_oce * zT                        ! draught (keel depth)
110         zW   = pt%width
111         zL   = pt%length
112         zxi  = pt%xi                                      ! position in (i,j) referential
113         zyj  = pt%yj
114         ii  = INT( zxi + 0.5 )                            ! T-cell of the berg
115         ii  = mi1( ii )
116         ij  = INT( zyj + 0.5 )             
117         ij  = mj1( ij )
118         zVol = zT * zW * zL
119
120         ! Environment
121         ! default sst, ssu and ssv
122         ! ln_M2016: use temp, u and v profile
123         IF ( ln_M2016 ) THEN
124
125            ! load t, u, v and e3 profile at icb position
126            CALL icb_utl_interp( pt%xi, pt%yj, ptoce=ztoce, puoce=zuoce, pvoce=zvoce, pe3t=ze3t )
127           
128            !compute bottom level
129            CALL icb_utl_getkb( pt%kb, ze3t, zD )
130
131            ikb = MIN(pt%kb,mbkt(ii,ij))                             ! limit pt%kb by mbkt
132                                                                     ! => bottom temperature used to fill ztoce(mbkt:jpk)
133            ztb = ztoce(ikb)                                         ! basal temperature
134            zub = zuoce(ikb)
135            zvb = zvoce(ikb)
136         ELSE
137            ztb = pt%sst
138            zub = pt%ssu
139            zvb = pt%ssv
140         END IF
141
142         zdvob = SQRT( (pt%uvel-zub)**2 + (pt%vvel-zvb)**2 )        ! relative basal velocity
143         zdva  = SQRT( (pt%ua  -pt%ssu)**2 + (pt%va  -pt%ssv)**2 )  ! relative wind
144         zSs   = 1.5_wp * SQRT( zdva ) + 0.1_wp * zdva              ! Sea state      (eqn M.A9)
145         !
146         ! Melt rates in m/s (i.e. division by rday)
147         ! Buoyant convection at sides (eqn M.A10)
148         IF ( ln_M2016 ) THEN
149            ! averaging along all the iceberg draft
150            zzMv(:) = MAX( 7.62d-3*ztoce(:)+1.29d-3*(ztoce(:)**2), 0._wp ) * z1_rday
151            CALL icb_utl_zavg(zMv, zzMv, ze3t, zD, ikb )
152         ELSE
153            zMv = MAX( 7.62d-3*zSST+1.29d-3*(zSST**2), 0._wp ) * z1_rday
154         END IF
155         !
156         ! Basal turbulent melting     (eqn M.A7 )
157         IF ( zSST > zfzpt ) THEN                                                              ! Calculate basal melting only if SST above freezing point 
158            zMb = MAX( 0.58_wp*(zdvob**0.8_wp)*(ztb+4.0_wp)/(zL**0.2_wp) , 0._wp ) * z1_rday
159         ELSE
160            zMb = 0._wp                                                                        ! No basal melting if SST below freezing point     
161         ENDIF
162         !
163         ! Wave erosion                (eqn M.A8 )
164         zMe = MAX( z1_12*(zSST+2.)*zSs*(1._wp+COS(rpi*(zIC**3)))     , 0._wp ) * z1_rday
165
166         IF( ln_operator_splitting ) THEN      ! Operator split update of volume/mass
167            zTn    = MAX( zT - zMb*zdt , 0._wp )         ! new total thickness (m)
168            znVol  = zTn * zW * zL                       ! new volume (m^3)
169            zMnew1 = ( znVol / zVol ) * zM               ! new mass (kg)
170            zdMb   = zM - zMnew1                         ! mass lost to basal melting (>0) (kg)
171            !
172            zLn    = MAX( zL - zMv*zdt , 0._wp )         ! new length (m)
173            zWn    = MAX( zW - zMv*zdt , 0._wp )         ! new width (m)
174            znVol  = zTn * zWn * zLn                     ! new volume (m^3)
175            zMnew2 = ( znVol / zVol ) * zM               ! new mass (kg)
176            zdMv   = zMnew1 - zMnew2                     ! mass lost to buoyant convection (>0) (kg)
177            !
178            zLn    = MAX( zLn - zMe*zdt , 0._wp )        ! new length (m)
179            zWn    = MAX( zWn - zMe*zdt , 0._wp )        ! new width (m)
180            znVol  = zTn * zWn * zLn                     ! new volume (m^3)
181            zMnew  = ( znVol / zVol ) * zM               ! new mass (kg)
182            zdMe   = zMnew2 - zMnew                      ! mass lost to erosion (>0) (kg)
183            zdM    = zM - zMnew                          ! mass lost to all erosion and melting (>0) (kg)
184            !
185         ELSE                                         ! Update dimensions of berg
186            zLn = MAX( zL -(zMv+zMe)*zdt ,0._wp )        ! (m)
187            zWn = MAX( zW -(zMv+zMe)*zdt ,0._wp )        ! (m)
188            zTn = MAX( zT - zMb    *zdt ,0._wp )         ! (m)
189            ! Update volume and mass of berg
190            znVol = zTn*zWn*zLn                          ! (m^3)
191            zMnew = (znVol/zVol)*zM                      ! (kg)
192            zdM   = zM - zMnew                           ! (kg)
193            zdMb = (zM/zVol) * (zW*   zL ) *zMb*zdt      ! approx. mass loss to basal melting (kg)
194            zdMe = (zM/zVol) * (zT*(zW+zL)) *zMe*zdt     ! approx. mass lost to erosion (kg)
195            zdMv = (zM/zVol) * (zT*(zW+zL)) *zMv*zdt     ! approx. mass loss to buoyant convection (kg)
196         ENDIF
197
198         IF( rn_bits_erosion_fraction > 0._wp ) THEN     ! Bergy bits
199            !
200            zMbits   = pt%mass_of_bits                                               ! mass of bergy bits (kg)
201            zdMbitsE = rn_bits_erosion_fraction * zdMe                               ! change in mass of bits (kg)
202            znMbits  = zMbits + zdMbitsE                                             ! add new bergy bits to mass (kg)
203            zLbits   = MIN( zL, zW, zT, 40._wp )                                     ! assume bergy bits are smallest dimension or 40 meters
204            zAbits   = ( zMbits / rn_rho_bergs ) / zLbits                            ! Effective bottom area (assuming T=Lbits)
205            zMbb     = MAX( 0.58_wp*(zdvo**0.8_wp)*(zSST+2._wp) /   &
206               &                              ( zLbits**0.2_wp ) , 0._wp ) * z1_rday ! Basal turbulent melting (for bits)
207            zMbb     = rn_rho_bergs * zAbits * zMbb                                  ! in kg/s
208            zdMbitsM = MIN( zMbb*zdt , znMbits )                                     ! bergy bits mass lost to melting (kg)
209            znMbits  = znMbits-zdMbitsM                                              ! remove mass lost to bergy bits melt
210            IF( zMnew == 0._wp ) THEN                                                ! if parent berg has completely melted then
211               zdMbitsM = zdMbitsM + znMbits                                         ! instantly melt all the bergy bits
212               znMbits  = 0._wp
213            ENDIF
214         ELSE                                                     ! No bergy bits
215            zAbits   = 0._wp
216            zdMbitsE = 0._wp
217            zdMbitsM = 0._wp
218            znMbits  = pt%mass_of_bits                             ! retain previous value incase non-zero
219         ENDIF
220
221         ! use tmask rather than tmask_i when dealing with icebergs
222         IF( tmask(ii,ij,1) /= 0._wp ) THEN    ! Add melting to the grid and field diagnostics
223            z1_e1e2    = r1_e1e2t(ii,ij) * this%mass_scaling
224            z1_dt_e1e2 = z1_dt * z1_e1e2
225            !
226            ! iceberg melt
227            !! the use of DDPDD function for the cumulative sum is needed for reproducibility
228            zmelt    = ( zdM - ( zdMbitsE - zdMbitsM ) ) * z1_dt   ! kg/s
229            CALL DDPDD( CMPLX( zmelt * z1_e1e2, 0.e0, dp ), cicb_melt(ii,ij) )
230            !
231            ! iceberg heat flux
232            !! the use of DDPDD function for the cumulative sum is needed for reproducibility
233            !! NB. The src_calving_hflx field is currently hardwired to zero in icb_stp, which means that the
234            !!     heat density of the icebergs is zero and the heat content flux to the ocean from iceberg
235            !!     melting is always zero. Leaving the term in the code until such a time as this is fixed. DS.
236            zheat_hcflux = zmelt * pt%heat_density       ! heat content flux : kg/s x J/kg = J/s
237            zheat_latent = - zmelt * rLfus               ! latent heat flux:  kg/s x J/kg = J/s
238            CALL DDPDD( CMPLX( ( zheat_hcflux + zheat_latent ) * z1_e1e2, 0.e0, dp ), cicb_hflx(ii,ij) )
239            !
240            ! diagnostics
241            CALL icb_dia_melt( ii, ij, zMnew, zheat_hcflux, zheat_latent, this%mass_scaling,       &
242               &                       zdM, zdMbitsE, zdMbitsM, zdMb, zdMe,   &
243               &                       zdMv, z1_dt_e1e2 )
244         ELSE
245            WRITE(numout,*) 'icb_thm: berg ',this%number(:),' appears to have grounded  at ',narea,ii,ij
246            CALL icb_utl_print_berg( this, kt )
247            WRITE(numout,*) 'msk=',tmask(ii,ij,1), e1e2t(ii,ij)
248            CALL ctl_stop('icb_thm', 'berg appears to have grounded!')
249         ENDIF
250
251         ! Rolling
252         zDn = rho_berg_1_oce * zTn       ! draught (keel depth)
253         IF( zDn > 0._wp .AND. MAX(zWn,zLn) < SQRT( 0.92*(zDn**2) + 58.32*zDn ) ) THEN
254            zT  = zTn
255            zTn = zWn
256            zWn = zT
257         ENDIF
258
259         ! Store the new state of iceberg (with L>W)
260         pt%mass         = zMnew
261         pt%mass_of_bits = znMbits
262         pt%thickness    = zTn
263         pt%width        = MIN( zWn , zLn )
264         pt%length       = MAX( zWn , zLn )
265
266         next=>this%next
267
268!!gm  add a test to avoid over melting ?
269!!pm  I agree, over melting could break conservation (more melt than calving)
270
271         IF( zMnew <= 0._wp ) THEN       ! Delete the berg if completely melted
272            CALL icb_utl_delete( first_berg, this )
273            !
274         ELSE                            ! Diagnose mass distribution on grid
275            z1_e1e2 = r1_e1e2t(ii,ij) * this%mass_scaling
276            CALL icb_dia_size( ii, ij, zWn, zLn, zAbits,   &
277               &               this%mass_scaling, zMnew, znMbits, z1_e1e2 )
278         ENDIF
279         !
280         this=>next
281         !
282      END DO
283      !
284      berg_grid%floating_melt = REAL(cicb_melt,dp)    ! kg/m2/s
285      berg_grid%calving_hflx  = REAL(cicb_hflx,dp)
286      !
287      ! now use melt and associated heat flux in ocean (or not)
288      !
289      IF(.NOT. ln_passive_mode ) THEN
290         emp (:,:) = emp (:,:) - berg_grid%floating_melt(:,:)
291         qns (:,:) = qns (:,:) + berg_grid%calving_hflx (:,:) 
292      ENDIF
293      !
294   END SUBROUTINE icb_thm
295
296   !!======================================================================
297END MODULE icbthm
Note: See TracBrowser for help on using the repository browser.