source: LMDZ6/trunk/libf/phylmd/ocean_slab_mod.F90 @ 3409

Last change on this file since 3409 was 3102, checked in by oboucher, 7 years ago

Removing x permission from these files

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 39.6 KB
RevLine 
[3002]1!Completed
[781]2MODULE ocean_slab_mod
3!
4! This module is used for both surface ocean and sea-ice when using the slab ocean,
5! "ocean=slab".
6!
[2057]7
8  USE dimphy
9  USE indice_sol_mod
[2209]10  USE surface_data
[2656]11  USE mod_grid_phy_lmdz, ONLY: klon_glo
12  USE mod_phys_lmdz_mpi_data, ONLY: is_mpi_root
[2057]13
[781]14  IMPLICIT NONE
[996]15  PRIVATE
[2209]16  PUBLIC :: ocean_slab_init, ocean_slab_frac, ocean_slab_noice, ocean_slab_ice
[781]17
[2656]18!***********************************************************************************
[2209]19! Global saved variables
[2656]20!***********************************************************************************
21  ! number of slab vertical layers
22  INTEGER, PUBLIC, SAVE :: nslay
23  !$OMP THREADPRIVATE(nslay)
24  ! timestep for coupling (update slab temperature) in timesteps
[2057]25  INTEGER, PRIVATE, SAVE                           :: cpl_pas
26  !$OMP THREADPRIVATE(cpl_pas)
[2656]27  ! cyang = 1/heat capacity of top layer (rho.c.H)
[2057]28  REAL, PRIVATE, SAVE                              :: cyang
29  !$OMP THREADPRIVATE(cyang)
[2656]30  ! depth of slab layers (1 or 2)
[2057]31  REAL, ALLOCATABLE, DIMENSION(:), PRIVATE, SAVE   :: slabh
32  !$OMP THREADPRIVATE(slabh)
[2656]33  ! slab temperature
[2057]34  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC, SAVE   :: tslab
35  !$OMP THREADPRIVATE(tslab)
[2656]36  ! heat flux convergence due to Ekman
37  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC, SAVE   :: dt_ekman
38  !$OMP THREADPRIVATE(dt_ekman)
39  ! heat flux convergence due to horiz diffusion
40  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC, SAVE   :: dt_hdiff
41  !$OMP THREADPRIVATE(dt_hdiff)
[3002]42  ! heat flux convergence due to GM eddy advection
43  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC, SAVE   :: dt_gm
44  !$OMP THREADPRIVATE(dt_gm)
45  ! Heat Flux correction
46  REAL, ALLOCATABLE, DIMENSION(:,:), PUBLIC, SAVE   :: dt_qflux
47  !$OMP THREADPRIVATE(dt_qflux)
[2656]48  ! fraction of ocean covered by sea ice (sic / (oce+sic))
[2209]49  REAL, ALLOCATABLE, DIMENSION(:), PUBLIC, SAVE  :: fsic
50  !$OMP THREADPRIVATE(fsic)
[2656]51  ! temperature of the sea ice
[2209]52  REAL, ALLOCATABLE, DIMENSION(:), PUBLIC, SAVE  :: tice
53  !$OMP THREADPRIVATE(tice)
[2656]54  ! sea ice thickness, in kg/m2
[2209]55  REAL, ALLOCATABLE, DIMENSION(:), PUBLIC, SAVE  :: seaice
56  !$OMP THREADPRIVATE(seaice)
[2656]57  ! net surface heat flux, weighted by open ocean fraction
58  ! slab_bils accumulated over cpl_pas timesteps
[2057]59  REAL, ALLOCATABLE, DIMENSION(:), PRIVATE, SAVE  :: bils_cum
60  !$OMP THREADPRIVATE(bils_cum)
[2656]61  ! net heat flux into the ocean below the ice : conduction + solar radiation
[2209]62  REAL, ALLOCATABLE, DIMENSION(:), PUBLIC, SAVE  :: slab_bilg
63  !$OMP THREADPRIVATE(slab_bilg)
[2656]64  ! slab_bilg over cpl_pas timesteps
[2209]65  REAL, ALLOCATABLE, DIMENSION(:), PRIVATE, SAVE  :: bilg_cum
66  !$OMP THREADPRIVATE(bilg_cum)
[2656]67  ! wind stress saved over cpl_pas timesteps
68  REAL, ALLOCATABLE, DIMENSION(:), PRIVATE, SAVE  :: taux_cum
69  !$OMP THREADPRIVATE(taux_cum)
70  REAL, ALLOCATABLE, DIMENSION(:), PRIVATE, SAVE  :: tauy_cum
71  !$OMP THREADPRIVATE(tauy_cum)
[2057]72
[2656]73!***********************************************************************************
[2209]74! Parameters (could be read in def file: move to slab_init)
[2656]75!***********************************************************************************
[2209]76! snow and ice physical characteristics:
77    REAL, PARAMETER :: t_freeze=271.35 ! freezing sea water temp
78    REAL, PARAMETER :: t_melt=273.15   ! melting ice temp
79    REAL, PARAMETER :: sno_den=300. !mean snow density, kg/m3
[2656]80    REAL, PARAMETER :: ice_den=917. ! ice density
81    REAL, PARAMETER :: sea_den=1025. ! sea water density
82    REAL, PARAMETER :: ice_cond=2.17*ice_den !conductivity of ice
83    REAL, PARAMETER :: sno_cond=0.31*sno_den ! conductivity of snow
[2209]84    REAL, PARAMETER :: ice_cap=2067.   ! specific heat capacity, snow and ice
[2656]85    REAL, PARAMETER :: sea_cap=3995.   ! specific heat capacity, snow and ice
[2209]86    REAL, PARAMETER :: ice_lat=334000. ! freeze /melt latent heat snow and ice
87
88! control of snow and ice cover & freeze / melt (heights converted to kg/m2)
89    REAL, PARAMETER :: snow_min=0.05*sno_den !critical snow height 5 cm
90    REAL, PARAMETER :: snow_wfact=0.4 ! max fraction of falling snow blown into ocean
91    REAL, PARAMETER :: ice_frac_min=0.001
92    REAL, PARAMETER :: ice_frac_max=1. ! less than 1. if min leads fraction
93    REAL, PARAMETER :: h_ice_min=0.01*ice_den ! min ice thickness
94    REAL, PARAMETER :: h_ice_thin=0.15*ice_den ! thin ice thickness
95    ! below ice_thin, priority is melt lateral / grow height
96    ! ice_thin is also height of new ice
97    REAL, PARAMETER :: h_ice_thick=2.5*ice_den ! thin ice thickness
98    ! above ice_thick, priority is melt height / grow lateral
99    REAL, PARAMETER :: h_ice_new=1.*ice_den ! max height of new open ocean ice
100    REAL, PARAMETER :: h_ice_max=10.*ice_den ! max ice height
101
102! albedo  and radiation parameters
103    REAL, PARAMETER :: alb_sno_min=0.55 !min snow albedo
104    REAL, PARAMETER :: alb_sno_del=0.3  !max snow albedo = min + del
105    REAL, PARAMETER :: alb_ice_dry=0.75 !dry thick ice
106    REAL, PARAMETER :: alb_ice_wet=0.66 !melting thick ice
[2656]107    REAL, PARAMETER :: pen_frac=0.3 !fraction of shortwave penetrating into the
108    ! ice (no snow)
[2209]109    REAL, PARAMETER :: pen_ext=1.5 !extinction of penetrating shortwave (m-1)
110
[2656]111! horizontal transport
112   LOGICAL, PUBLIC, SAVE :: slab_hdiff
113   !$OMP THREADPRIVATE(slab_hdiff)
[3002]114   LOGICAL, PUBLIC, SAVE :: slab_gm
115   !$OMP THREADPRIVATE(slab_gm)
[2656]116   REAL, PRIVATE, SAVE    :: coef_hdiff ! coefficient for horizontal diffusion
117   !$OMP THREADPRIVATE(coef_hdiff)
118   INTEGER, PUBLIC, SAVE :: slab_ekman, slab_cadj ! Ekman, conv adjustment
119   !$OMP THREADPRIVATE(slab_ekman)
120   !$OMP THREADPRIVATE(slab_cadj)
[2209]121
[2656]122!***********************************************************************************
123
[781]124CONTAINS
125!
[2656]126!***********************************************************************************
[781]127!
[2057]128  SUBROUTINE ocean_slab_init(dtime, pctsrf_rst)
129  !, seaice_rst etc
[781]130
[2656]131    USE ioipsl_getin_p_mod, ONLY : getin_p
132    USE mod_phys_lmdz_transfert_para, ONLY : gather
133    USE slab_heat_transp_mod, ONLY : ini_slab_transp
[2057]134
135    ! Input variables
[2656]136!***********************************************************************************
[2057]137    REAL, INTENT(IN)                         :: dtime
138! Variables read from restart file
[2656]139    REAL, DIMENSION(klon, nbsrf), INTENT(IN) :: pctsrf_rst
140    ! surface fractions from start file
[2057]141
142! Local variables
[2656]143!************************************************************************************
[2057]144    INTEGER                :: error
[2656]145    REAL, DIMENSION(klon_glo) :: zmasq_glo
[2057]146    CHARACTER (len = 80)   :: abort_message
147    CHARACTER (len = 20)   :: modname = 'ocean_slab_intit'
148
[2656]149!***********************************************************************************
150! Define some parameters
151!***********************************************************************************
152! Number of slab layers
153    nslay=2
154    CALL getin_p('slab_layers',nslay)
155    print *,'number of slab layers : ',nslay
156! Layer thickness
157    ALLOCATE(slabh(nslay), stat = error)
158    IF (error /= 0) THEN
159       abort_message='Pb allocation slabh'
160       CALL abort_physic(modname,abort_message,1)
161    ENDIF
162    slabh(1)=50.
[3002]163    CALL getin_p('slab_depth',slabh(1))
[2656]164    IF (nslay.GT.1) THEN
165        slabh(2)=150.
166    END IF
167
168! cyang = 1/heat capacity of top layer (rho.c.H)
169    cyang=1/(slabh(1)*sea_den*sea_cap)
170
171! cpl_pas  coupling period (update of tslab and ice fraction)
172! pour un calcul a chaque pas de temps, cpl_pas=1
173    cpl_pas = NINT(86400./dtime * 1.0) ! une fois par jour
174    CALL getin_p('cpl_pas',cpl_pas)
175    print *,'cpl_pas',cpl_pas
176
177! Horizontal diffusion
178    slab_hdiff=.FALSE.
179    CALL getin_p('slab_hdiff',slab_hdiff)
180    coef_hdiff=25000.
181    CALL getin_p('coef_hdiff',coef_hdiff)
182! Ekman transport
183    slab_ekman=0
184    CALL getin_p('slab_ekman',slab_ekman)
[3002]185! GM eddy advection (2-layers only)
186    slab_gm=.FALSE.
187    CALL getin_p('slab_gm',slab_gm)
188    IF (slab_ekman.LT.2) THEN
189       slab_gm=.FALSE.
190    ENDIF
[2656]191! Convective adjustment
192    IF (nslay.EQ.1) THEN
193        slab_cadj=0
194    ELSE
195        slab_cadj=1
196    END IF
197    CALL getin_p('slab_cadj',slab_cadj)
198
199!************************************************************************************
[2057]200! Allocate surface fraction read from restart file
[2656]201!************************************************************************************
[2209]202    ALLOCATE(fsic(klon), stat = error)
[2057]203    IF (error /= 0) THEN
204       abort_message='Pb allocation tmp_pctsrf_slab'
[2311]205       CALL abort_physic(modname,abort_message,1)
[2057]206    ENDIF
[2209]207    fsic(:)=0.
[2656]208    !zmasq = continent fraction
[2209]209    WHERE (1.-zmasq(:)>EPSFRA)
210        fsic(:) = pctsrf_rst(:,is_sic)/(1.-zmasq(:))
211    END WHERE
[2057]212
[2656]213!************************************************************************************
214! Allocate saved fields
215!************************************************************************************
[2209]216    ALLOCATE(tslab(klon,nslay), stat=error)
[2311]217       IF (error /= 0) CALL abort_physic &
[2209]218         (modname,'pb allocation tslab', 1)
219
[2057]220    ALLOCATE(bils_cum(klon), stat = error)
221    IF (error /= 0) THEN
222       abort_message='Pb allocation slab_bils_cum'
[2311]223       CALL abort_physic(modname,abort_message,1)
[2057]224    ENDIF
225    bils_cum(:) = 0.0   
226
[2656]227    IF (version_ocean=='sicINT') THEN ! interactive sea ice
[2209]228        ALLOCATE(slab_bilg(klon), stat = error)
229        IF (error /= 0) THEN
230           abort_message='Pb allocation slab_bilg'
[2311]231           CALL abort_physic(modname,abort_message,1)
[2209]232        ENDIF
233        slab_bilg(:) = 0.0   
234        ALLOCATE(bilg_cum(klon), stat = error)
235        IF (error /= 0) THEN
236           abort_message='Pb allocation slab_bilg_cum'
[2311]237           CALL abort_physic(modname,abort_message,1)
[2209]238        ENDIF
239        bilg_cum(:) = 0.0   
240        ALLOCATE(tice(klon), stat = error)
241        IF (error /= 0) THEN
242           abort_message='Pb allocation slab_tice'
[2311]243           CALL abort_physic(modname,abort_message,1)
[2209]244        ENDIF
245        ALLOCATE(seaice(klon), stat = error)
246        IF (error /= 0) THEN
247           abort_message='Pb allocation slab_seaice'
[2311]248           CALL abort_physic(modname,abort_message,1)
[2209]249        ENDIF
250    END IF
251
[2656]252    IF (slab_hdiff) THEN !horizontal diffusion
253        ALLOCATE(dt_hdiff(klon,nslay), stat = error)
254        IF (error /= 0) THEN
255           abort_message='Pb allocation dt_hdiff'
256           CALL abort_physic(modname,abort_message,1)
257        ENDIF
258        dt_hdiff(:,:) = 0.0   
[2057]259    ENDIF
260
[3002]261    ALLOCATE(dt_qflux(klon,nslay), stat = error)
262    IF (error /= 0) THEN
263       abort_message='Pb allocation dt_qflux'
264       CALL abort_physic(modname,abort_message,1)
265    ENDIF
266    dt_qflux(:,:) = 0.0   
267
268    IF (slab_gm) THEN !GM advection
269        ALLOCATE(dt_gm(klon,nslay), stat = error)
270        IF (error /= 0) THEN
271           abort_message='Pb allocation dt_gm'
272           CALL abort_physic(modname,abort_message,1)
273        ENDIF
274        dt_gm(:,:) = 0.0   
275    ENDIF
276
[2656]277    IF (slab_ekman.GT.0) THEN ! ekman transport
278        ALLOCATE(dt_ekman(klon,nslay), stat = error)
279        IF (error /= 0) THEN
280           abort_message='Pb allocation dt_ekman'
281           CALL abort_physic(modname,abort_message,1)
282        ENDIF
283        dt_ekman(:,:) = 0.0   
284        ALLOCATE(taux_cum(klon), stat = error)
285        IF (error /= 0) THEN
286           abort_message='Pb allocation taux_cum'
287           CALL abort_physic(modname,abort_message,1)
288        ENDIF
289        taux_cum(:) = 0.0   
290        ALLOCATE(tauy_cum(klon), stat = error)
291        IF (error /= 0) THEN
292           abort_message='Pb allocation tauy_cum'
293           CALL abort_physic(modname,abort_message,1)
294        ENDIF
295        tauy_cum(:) = 0.0   
296    ENDIF
[2209]297
[2656]298! Initialize transport
299    IF (slab_hdiff.OR.(slab_ekman.GT.0)) THEN
300      CALL gather(zmasq,zmasq_glo)
[3002]301! Master thread/process only
302!$OMP MASTER 
[2656]303      IF (is_mpi_root) THEN
304          CALL ini_slab_transp(zmasq_glo)
305      END IF
306!$OMP END MASTER
307    END IF
308
[2057]309 END SUBROUTINE ocean_slab_init
310!
[2656]311!***********************************************************************************
[2057]312!
313  SUBROUTINE ocean_slab_frac(itime, dtime, jour, pctsrf_chg, is_modified)
314
[2656]315! this routine sends back the sea ice and ocean fraction to the main physics
316! routine. Called only with interactive sea ice
[1785]317
[996]318! Arguments
[2656]319!************************************************************************************
320    INTEGER, INTENT(IN)                        :: itime   ! current timestep
321    INTEGER, INTENT(IN)                        :: jour    !  day in year (not
322    REAL   , INTENT(IN)                        :: dtime   ! physics timestep (s)
[2057]323    REAL, DIMENSION(klon,nbsrf), INTENT(INOUT) :: pctsrf_chg  ! sub-surface fraction
[2656]324    LOGICAL, INTENT(OUT)                       :: is_modified ! true if pctsrf is
325                                                         ! modified at this time step
[781]326
[2209]327       pctsrf_chg(:,is_oce)=(1.-fsic(:))*(1.-zmasq(:))
328       pctsrf_chg(:,is_sic)=fsic(:)*(1.-zmasq(:))
[2057]329       is_modified=.TRUE.
[781]330
[996]331  END SUBROUTINE ocean_slab_frac
[781]332!
[2656]333!************************************************************************************
[781]334!
335  SUBROUTINE ocean_slab_noice( &
[996]336       itime, dtime, jour, knon, knindex, &
[2254]337       p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, temp_air, spechum, &
[1067]338       AcoefH, AcoefQ, BcoefH, BcoefQ, &
339       AcoefU, AcoefV, BcoefU, BcoefV, &
[2240]340       ps, u1, v1, gustiness, tsurf_in, &
[2209]341       radsol, snow, &
[1067]342       qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
[3002]343       tsurf_new, dflux_s, dflux_l, slab_bils)
[1067]344   
345    USE calcul_fluxs_mod
[3002]346    USE slab_heat_transp_mod, ONLY: divgrad_phy,slab_ekman1,slab_ekman2,slab_gmdiff
[2656]347    USE mod_phys_lmdz_para
[1785]348
[2254]349    INCLUDE "clesphys.h"
[781]350
[2656]351! This routine
352! (1) computes surface turbulent fluxes over points with some open ocean   
353! (2) reads additional Q-flux (everywhere)
354! (3) computes horizontal transport (diffusion & Ekman)
355! (4) updates slab temperature every cpl_pas ; creates new ice if needed.
356
357! Note :
358! klon total number of points
359! knon number of points with open ocean (varies with time)
360! knindex gives position of the knon points within klon.
361! In general, local saved variables have klon values
362! variables exchanged with PBL module have knon.
363
[781]364! Input arguments
[2656]365!***********************************************************************************
366    INTEGER, INTENT(IN)                  :: itime ! current timestep INTEGER,
367    INTEGER, INTENT(IN)                  :: jour  ! day in year (for Q-Flux)
368    INTEGER, INTENT(IN)                  :: knon  ! number of points
369    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex
370    REAL, INTENT(IN) :: dtime  ! timestep (s)
371    REAL, DIMENSION(klon), INTENT(IN)    :: p1lay
372    REAL, DIMENSION(klon), INTENT(IN)    :: cdragh, cdragq, cdragm
373    ! drag coefficients
374    REAL, DIMENSION(klon), INTENT(IN)    :: precip_rain, precip_snow
375    REAL, DIMENSION(klon), INTENT(IN)    :: temp_air, spechum ! near surface T, q
376    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefH, AcoefQ, BcoefH, BcoefQ
377    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefU, AcoefV, BcoefU, BcoefV
378    ! exchange coefficients for boundary layer scheme
379    REAL, DIMENSION(klon), INTENT(IN)    :: ps  ! surface pressure
380    REAL, DIMENSION(klon), INTENT(IN)    :: u1, v1, gustiness ! surface wind
381    REAL, DIMENSION(klon), INTENT(IN)    :: tsurf_in ! surface temperature
382    REAL, DIMENSION(klon), INTENT(INOUT) :: radsol ! net surface radiative flux
[781]383
384! In/Output arguments
[2656]385!************************************************************************************
386    REAL, DIMENSION(klon), INTENT(INOUT) :: snow ! in kg/m2
[781]387   
388! Output arguments
[2656]389!************************************************************************************
[781]390    REAL, DIMENSION(klon), INTENT(OUT)   :: qsurf
391    REAL, DIMENSION(klon), INTENT(OUT)   :: evap, fluxsens, fluxlat
[1067]392    REAL, DIMENSION(klon), INTENT(OUT)   :: flux_u1, flux_v1
[2656]393    REAL, DIMENSION(klon), INTENT(OUT)   :: tsurf_new ! new surface tempearture
[781]394    REAL, DIMENSION(klon), INTENT(OUT)   :: dflux_s, dflux_l     
[3002]395    REAL, DIMENSION(klon), INTENT(OUT)   :: slab_bils
[781]396
397! Local variables
[2656]398!************************************************************************************
399    INTEGER               :: i,ki,k
400    REAL                  :: t_cadj
401    !  for surface heat fluxes
[996]402    REAL, DIMENSION(klon) :: cal, beta, dif_grnd
[2656]403    ! for Q-Flux computation: d/dt SST, d/dt ice volume (kg/m2), surf fluxes
[3002]404    REAL, DIMENSION(klon) :: diff_sst, diff_siv
405    REAL, DIMENSION(klon,nslay) :: lmt_bils
[2656]406    ! for surface wind stress
[1067]407    REAL, DIMENSION(klon) :: u0, v0
408    REAL, DIMENSION(klon) :: u1_lay, v1_lay
[2656]409    ! for new ice creation
[2209]410    REAL                  :: e_freeze, h_new, dfsic
[2656]411    ! horizontal diffusion and Ekman local vars
412    ! dimension = global domain (klon_glo) instead of // subdomains
[3002]413    REAL, DIMENSION(klon_glo,nslay) :: dt_hdiff_glo,dt_ekman_glo,dt_gm_glo
[2656]414    ! dt_ekman_glo saved for diagnostic, dt_ekman_tmp used for time loop
415    REAL, DIMENSION(klon_glo,nslay) :: dt_hdiff_tmp, dt_ekman_tmp
416    REAL, DIMENSION(klon_glo,nslay) :: tslab_glo
417    REAL, DIMENSION(klon_glo) :: taux_glo,tauy_glo
[781]418
419!****************************************************************************************
[2656]420! 1) Surface fluxes calculation
421!   
[996]422!****************************************************************************************
[2656]423    cal(:)      = 0. ! infinite thermal inertia
424    beta(:)     = 1. ! wet surface
425    dif_grnd(:) = 0. ! no diffusion into ground
[781]426   
[1067]427! Suppose zero surface speed
428    u0(:)=0.0
429    v0(:)=0.0
430    u1_lay(:) = u1(:) - u0(:)
431    v1_lay(:) = v1(:) - v0(:)
432
[2656]433! Compute latent & sensible fluxes
[781]434    CALL calcul_fluxs(knon, is_oce, dtime, &
[2254]435         tsurf_in, p1lay, cal, beta, cdragh, cdragq, ps, &
[781]436         precip_rain, precip_snow, snow, qsurf,  &
[2240]437         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
[2254]438         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
[781]439         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
440
[2656]441! save total cumulated heat fluxes locally
442! radiative + turbulent + melt of falling snow
[2057]443    slab_bils(:)=0.
444    DO i=1,knon
445        ki=knindex(i)
[2209]446        slab_bils(ki)=(1.-fsic(ki))*(fluxlat(i)+fluxsens(i)+radsol(i) &
447                      -precip_snow(i)*ice_lat*(1.+snow_wfact*fsic(ki)))
[2057]448        bils_cum(ki)=bils_cum(ki)+slab_bils(ki)
449    END DO
450
[2656]451!  Compute surface wind stress
452    CALL calcul_flux_wind(knon, dtime, &
453         u0, v0, u1, v1, gustiness, cdragm, &
454         AcoefU, AcoefV, BcoefU, BcoefV, &
455         p1lay, temp_air, &
456         flux_u1, flux_v1) 
457
458! save cumulated wind stress
459    IF (slab_ekman.GT.0) THEN
460      DO i=1,knon
461          ki=knindex(i)
462          taux_cum(ki)=taux_cum(ki)+flux_u1(i)*(1.-fsic(ki))/cpl_pas
463          tauy_cum(ki)=tauy_cum(ki)+flux_v1(i)*(1.-fsic(ki))/cpl_pas
464      END DO
465    ENDIF
466
[781]467!****************************************************************************************
[2656]468! 2) Q-Flux : get global variables lmt_bils, diff_sst and diff_siv from file limit_slab.nc
[781]469!
470!****************************************************************************************
[2656]471    CALL limit_slab(itime, dtime, jour, lmt_bils, diff_sst, diff_siv)
[2209]472    ! lmt_bils and diff_sst,siv saved by limit_slab
[2057]473    ! qflux = total QFlux correction (in W/m2)
[3002]474    dt_qflux(:,1)=lmt_bils(:,1)+diff_sst(:)/cyang/86400.-diff_siv(:)*ice_den*ice_lat/86400.
475    IF (nslay.GT.1) THEN
476       dt_qflux(:,2:nslay)=lmt_bils(:,2:nslay)
477    END IF
[781]478
479!****************************************************************************************
[2656]480! 3) Recalculate new temperature (add Surf fluxes, Q-Flux, Ocean transport)
481!    Bring to freezing temp and make sea ice if necessary
482
[2057]483!***********************************************o*****************************************
484    tsurf_new=tsurf_in
485    IF (MOD(itime,cpl_pas).EQ.0) THEN ! time to update tslab & fraction
[2656]486! ***********************************
487!  Horizontal transport
488! ***********************************
489      IF (slab_ekman.GT.0) THEN
490          ! copy wind stress to global var
491          CALL gather(taux_cum,taux_glo)
492          CALL gather(tauy_cum,tauy_glo)
493      END IF
494
495      IF (slab_hdiff.OR.(slab_ekman.GT.0)) THEN
496        CALL gather(tslab,tslab_glo)
497      ! Compute horiz transport on one process only
[3002]498        IF (is_mpi_root .AND. is_omp_root) THEN ! Only master processus         
[2656]499          IF (slab_hdiff) THEN
500              dt_hdiff_glo(:,:)=0.
501          END IF
502          IF (slab_ekman.GT.0) THEN
503              dt_ekman_glo(:,:)=0.
504          END IF
[3002]505          IF (slab_gm) THEN
506              dt_gm_glo(:,:)=0.
507          END IF
[2656]508          DO i=1,cpl_pas ! time splitting for numerical stability
509            IF (slab_ekman.GT.0) THEN
510              SELECT CASE (slab_ekman)
511                CASE (1)
512                  CALL slab_ekman1(taux_glo,tauy_glo,tslab_glo,dt_ekman_tmp)
513                CASE (2)
[3002]514                  CALL slab_ekman2(taux_glo,tauy_glo,tslab_glo,dt_ekman_tmp,dt_hdiff_tmp,slab_gm)
[2656]515                CASE DEFAULT
516                  dt_ekman_tmp(:,:)=0.
517              END SELECT
518              dt_ekman_glo(:,:)=dt_ekman_glo(:,:)+dt_ekman_tmp(:,:)
519              ! convert dt_ekman from K.s-1.(kg.m-2) to K.s-1   
520              DO k=1,nslay
521                dt_ekman_tmp(:,k)=dt_ekman_tmp(:,k)/(slabh(k)*sea_den)
522              ENDDO
523              tslab_glo=tslab_glo+dt_ekman_tmp*dtime
[3002]524              IF (slab_gm) THEN ! Gent-McWilliams eddy advection
525                dt_gm_glo(:,:)=dt_gm_glo(:,:)+ dt_hdiff_tmp(:,:)
526                ! convert dt from K.s-1.(kg.m-2) to K.s-1   
527                DO k=1,nslay
528                  dt_hdiff_tmp(:,k)=dt_hdiff_tmp(:,k)/(slabh(k)*sea_den)
529                END DO
530                tslab_glo=tslab_glo+dt_hdiff_tmp*dtime
531              END IF
[2656]532            ENDIF
[3002]533! GM included in Ekman_2
534!            IF (slab_gm) THEN ! Gent-McWilliams eddy advection
535!              CALL slab_gmdiff(tslab_glo,dt_hdiff_tmp)
536!              ! convert dt_gm from K.m.s-1 to K.s-1   
537!              DO k=1,nslay
538!                dt_hdiff_tmp(:,k)=dt_hdiff_tmp(:,k)/slabh(k)
539!              END DO
540!              tslab_glo=tslab_glo+dt_hdiff_tmp*dtime
541!              dt_gm_glo(:,:)=dt_gm_glo(:,:)+ dt_hdiff_tmp(:,:)
542!            END IF
[2656]543            IF (slab_hdiff) THEN ! horizontal diffusion
544              ! laplacian of slab T
545              CALL divgrad_phy(nslay,tslab_glo,dt_hdiff_tmp)
546              ! multiply by diff coef and normalize to 50m slab equivalent
547              dt_hdiff_tmp=dt_hdiff_tmp*coef_hdiff*50./SUM(slabh)
548              dt_hdiff_glo(:,:)=dt_hdiff_glo(:,:)+ dt_hdiff_tmp(:,:)
549              tslab_glo=tslab_glo+dt_hdiff_tmp*dtime
550            END IF
551          END DO ! time splitting
552          IF (slab_hdiff) THEN
553            !dt_hdiff_glo saved in W/m2
554            DO k=1,nslay
555              dt_hdiff_glo(:,k)=dt_hdiff_glo(:,k)*slabh(k)*sea_den*sea_cap/cpl_pas
[2057]556            END DO
[2656]557          END IF
[3002]558          IF (slab_gm) THEN
559            !dt_hdiff_glo saved in W/m2
560            dt_gm_glo(:,:)=dt_gm_glo(:,:)*sea_cap/cpl_pas
561          END IF
[2656]562          IF (slab_ekman.GT.0) THEN
563            ! dt_ekman_glo saved in W/m2
564            dt_ekman_glo(:,:)=dt_ekman_glo(:,:)*sea_cap/cpl_pas
565          END IF
[3002]566        END IF ! master process
[2656]567!$OMP BARRIER
568        ! Send new fields back to all processes
569        CALL Scatter(tslab_glo,tslab)
570        IF (slab_hdiff) THEN
571          CALL Scatter(dt_hdiff_glo,dt_hdiff)
572        END IF
[3002]573        IF (slab_gm) THEN
574          CALL Scatter(dt_gm_glo,dt_gm)
575        END IF
[2656]576        IF (slab_ekman.GT.0) THEN
577          CALL Scatter(dt_ekman_glo,dt_ekman)
578          ! clear wind stress
579          taux_cum(:)=0.
580          tauy_cum(:)=0.
581        END IF
582      ENDIF ! transport
583
584! ***********************************
585! Other heat fluxes
586! ***********************************
587      ! Add read QFlux
[3002]588      DO k=1,nslay
589        tslab(:,k)=tslab(:,k)+dt_qflux(:,k)*cyang*dtime*cpl_pas &
590                   *slabh(1)/slabh(k)
591      END DO
[2656]592      ! Add cumulated surface fluxes
593      tslab(:,1)=tslab(:,1)+bils_cum(:)*cyang*dtime
594      ! Convective adjustment if 2 layers
595      IF ((nslay.GT.1).AND.(slab_cadj.GT.0)) THEN
596        DO i=1,klon
597          IF (tslab(i,2).GT.tslab(i,1)) THEN
598            ! mean (mass-weighted) temperature
599            t_cadj=SUM(tslab(i,:)*slabh(:))/SUM(slabh(:))
600            tslab(i,1)=t_cadj
601            tslab(i,2)=t_cadj
602          END IF
603        END DO
604      END IF
605! ***********************************
606! Update surface temperature and ice
607! ***********************************
608      SELECT CASE(version_ocean)
609      CASE('sicNO') ! no sea ice even below freezing !
610          DO i=1,knon
611              ki=knindex(i)
612              tsurf_new(i)=tslab(ki,1)
613          END DO
614      CASE('sicOBS') ! "realistic" case, for prescribed sea ice
615        ! tslab cannot be below freezing, or above it if there is sea ice
616          DO i=1,knon
617              ki=knindex(i)
618              IF ((tslab(ki,1).LT.t_freeze).OR.(fsic(ki).GT.epsfra)) THEN
619                  tslab(ki,1)=t_freeze
620              END IF
621              tsurf_new(i)=tslab(ki,1)
622          END DO
623      CASE('sicINT') ! interactive sea ice
624          DO i=1,knon
625              ki=knindex(i)
626              IF (fsic(ki).LT.epsfra) THEN ! Free of ice
627                  IF (tslab(ki,1).LT.t_freeze) THEN ! create new ice
628                      ! quantity of new ice formed
629                      e_freeze=(t_freeze-tslab(ki,1))/cyang/ice_lat
630                      ! new ice
631                      tice(ki)=t_freeze
632                      fsic(ki)=MIN(ice_frac_max,e_freeze/h_ice_thin)
633                      IF (fsic(ki).GT.ice_frac_min) THEN
634                          seaice(ki)=MIN(e_freeze/fsic(ki),h_ice_max)
635                          tslab(ki,1)=t_freeze
636                      ELSE
637                          fsic(ki)=0.
638                      END IF
639                      tsurf_new(i)=t_freeze
640                  ELSE
641                      tsurf_new(i)=tslab(ki,1)
642                  END IF
643              ELSE ! ice present
644                  tsurf_new(i)=t_freeze
645                  IF (tslab(ki,1).LT.t_freeze) THEN ! create new ice
646                      ! quantity of new ice formed over open ocean
647                      e_freeze=(t_freeze-tslab(ki,1))/cyang*(1.-fsic(ki)) &
648                               /(ice_lat+ice_cap/2.*(t_freeze-tice(ki)))
649                      ! new ice height and fraction
650                      h_new=MIN(h_ice_new,seaice(ki)) ! max new height ice_new
651                      dfsic=MIN(ice_frac_max-fsic(ki),e_freeze/h_new)
652                      h_new=MIN(e_freeze/dfsic,h_ice_max)
653                      ! update tslab to freezing over open ocean only
654                      tslab(ki,1)=tslab(ki,1)*fsic(ki)+t_freeze*(1.-fsic(ki))
655                      ! update sea ice
656                      seaice(ki)=(h_new*dfsic+seaice(ki)*fsic(ki)) &
657                                 /(dfsic+fsic(ki))
658                      fsic(ki)=fsic(ki)+dfsic
659                      ! update snow?
660                  END IF ! tslab below freezing
661              END IF ! sea ice present
662          END DO
663      END SELECT
664      bils_cum(:)=0.0! clear cumulated fluxes
[2057]665    END IF ! coupling time
666  END SUBROUTINE ocean_slab_noice
667!
[781]668!****************************************************************************************
[2209]669
670  SUBROUTINE ocean_slab_ice(   &
671       itime, dtime, jour, knon, knindex, &
672       tsurf_in, p1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, &
673       AcoefH, AcoefQ, BcoefH, BcoefQ, &
674       AcoefU, AcoefV, BcoefU, BcoefV, &
[2240]675       ps, u1, v1, gustiness, &
[2209]676       radsol, snow, qsurf, qsol, agesno, &
677       alb1_new, alb2_new, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
678       tsurf_new, dflux_s, dflux_l, swnet)
679
680   USE calcul_fluxs_mod
681
682   INCLUDE "YOMCST.h"
[2254]683   INCLUDE "clesphys.h"
[2209]684
685! Input arguments
[2057]686!****************************************************************************************
[2209]687    INTEGER, INTENT(IN)                  :: itime, jour, knon
688    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex
689    REAL, INTENT(IN)                     :: dtime
690    REAL, DIMENSION(klon), INTENT(IN)    :: tsurf_in
691    REAL, DIMENSION(klon), INTENT(IN)    :: p1lay
692    REAL, DIMENSION(klon), INTENT(IN)    :: cdragh, cdragm
693    REAL, DIMENSION(klon), INTENT(IN)    :: precip_rain, precip_snow
694    REAL, DIMENSION(klon), INTENT(IN)    :: temp_air, spechum
695    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefH, AcoefQ, BcoefH, BcoefQ
696    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefU, AcoefV, BcoefU, BcoefV
697    REAL, DIMENSION(klon), INTENT(IN)    :: ps
[2240]698    REAL, DIMENSION(klon), INTENT(IN)    :: u1, v1, gustiness
[2209]699    REAL, DIMENSION(klon), INTENT(IN)    :: swnet
700
701! In/Output arguments
702!****************************************************************************************
703    REAL, DIMENSION(klon), INTENT(INOUT)          :: snow, qsol
704    REAL, DIMENSION(klon), INTENT(INOUT)          :: agesno
705    REAL, DIMENSION(klon), INTENT(INOUT)          :: radsol
706
707! Output arguments
708!****************************************************************************************
709    REAL, DIMENSION(klon), INTENT(OUT)            :: qsurf
710    REAL, DIMENSION(klon), INTENT(OUT)            :: alb1_new  ! new albedo in visible SW interval
711    REAL, DIMENSION(klon), INTENT(OUT)            :: alb2_new  ! new albedo in near IR interval
712    REAL, DIMENSION(klon), INTENT(OUT)            :: evap, fluxsens, fluxlat
713    REAL, DIMENSION(klon), INTENT(OUT)            :: flux_u1, flux_v1
714    REAL, DIMENSION(klon), INTENT(OUT)            :: tsurf_new
715    REAL, DIMENSION(klon), INTENT(OUT)            :: dflux_s, dflux_l
716
717! Local variables
718!****************************************************************************************
719    INTEGER               :: i,ki
720    REAL, DIMENSION(klon) :: cal, beta, dif_grnd
721    REAL, DIMENSION(klon) :: u0, v0
722    REAL, DIMENSION(klon) :: u1_lay, v1_lay
723    ! intermediate heat fluxes:
724    REAL                  :: f_cond, f_swpen
725    ! for snow/ice albedo:
726    REAL                  :: alb_snow, alb_ice, alb_pond
727    REAL                  :: frac_snow, frac_ice, frac_pond
728    ! for ice melt / freeze
729    REAL                  :: e_melt, snow_evap, h_test
730    ! dhsic, dfsic change in ice mass, fraction.
731    REAL                  :: dhsic, dfsic, frac_mf
732
733!****************************************************************************************
[2057]734! 1) Flux calculation
735!****************************************************************************************
[2209]736! Suppose zero surface speed
737    u0(:)=0.0
738    v0(:)=0.0
739    u1_lay(:) = u1(:) - u0(:)
740    v1_lay(:) = v1(:) - v0(:)
741
742! set beta, cal, compute conduction fluxes inside ice/snow
743    slab_bilg(:)=0.
744    dif_grnd(:)=0.
745    beta(:) = 1.
746    DO i=1,knon
747    ki=knindex(i)
748        IF (snow(i).GT.snow_min) THEN
749            ! snow-layer heat capacity
750            cal(i)=2.*RCPD/(snow(i)*ice_cap)
751            ! snow conductive flux
752            f_cond=sno_cond*(tice(ki)-tsurf_in(i))/snow(i)
753            ! all shortwave flux absorbed
754            f_swpen=0.
755            ! bottom flux (ice conduction)
756            slab_bilg(ki)=ice_cond*(tice(ki)-t_freeze)/seaice(ki)
757            ! update ice temperature
758            tice(ki)=tice(ki)-2./ice_cap/(snow(i)+seaice(ki)) &
759                     *(slab_bilg(ki)+f_cond)*dtime
760       ELSE ! bare ice
761            ! ice-layer heat capacity
762            cal(i)=2.*RCPD/(seaice(ki)*ice_cap)
763            ! conductive flux
764            f_cond=ice_cond*(t_freeze-tice(ki))/seaice(ki)
765            ! penetrative shortwave flux...
766            f_swpen=swnet(i)*pen_frac*exp(-pen_ext*seaice(ki)/ice_den)
767            slab_bilg(ki)=f_swpen-f_cond
768        END IF
769        radsol(i)=radsol(i)+f_cond-f_swpen
770    END DO
771    ! weight fluxes to ocean by sea ice fraction
772    slab_bilg(:)=slab_bilg(:)*fsic(:)
773
[2057]774! calcul_fluxs (sens, lat etc)
[2209]775    CALL calcul_fluxs(knon, is_sic, dtime, &
[2254]776        tsurf_in, p1lay, cal, beta, cdragh, cdragh, ps, &
[2209]777        precip_rain, precip_snow, snow, qsurf,  &
[2240]778        radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
[2254]779        f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
[2209]780        tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
781    DO i=1,knon
782        IF (snow(i).LT.snow_min) tice(knindex(i))=tsurf_new(i)
783    END DO
784
[2057]785! calcul_flux_wind
[2209]786    CALL calcul_flux_wind(knon, dtime, &
[2240]787         u0, v0, u1, v1, gustiness, cdragm, &
[2209]788         AcoefU, AcoefV, BcoefU, BcoefV, &
789         p1lay, temp_air, &
790         flux_u1, flux_v1)
[781]791
[2057]792!****************************************************************************************
[2209]793! 2) Update snow and ice surface
[2057]794!****************************************************************************************
[2209]795! snow precip
796    DO i=1,knon
797        ki=knindex(i)
798        IF (precip_snow(i) > 0.) THEN
799            snow(i) = snow(i)+precip_snow(i)*dtime*(1.-snow_wfact*(1.-fsic(ki)))
800        END IF
801! snow and ice sublimation
802        IF (evap(i) > 0.) THEN
803           snow_evap = MIN (snow(i) / dtime, evap(i))
804           snow(i) = snow(i) - snow_evap * dtime
805           snow(i) = MAX(0.0, snow(i))
806           seaice(ki) = MAX(0.0,seaice(ki)-(evap(i)-snow_evap)*dtime)
807        ENDIF
[2656]808! Melt / Freeze snow from above if Tsurf>0
[2209]809        IF (tsurf_new(i).GT.t_melt) THEN
[2656]810            ! energy available for melting snow (in kg of melted snow /m2)
[2209]811            e_melt = MIN(MAX(snow(i)*(tsurf_new(i)-t_melt)*ice_cap/2. &
812               /(ice_lat+ice_cap/2.*(t_melt-tice(ki))),0.0),snow(i))
813            ! remove snow
814            IF (snow(i).GT.e_melt) THEN
815                snow(i)=snow(i)-e_melt
816                tsurf_new(i)=t_melt
817            ELSE ! all snow is melted
818                ! add remaining heat flux to ice
819                e_melt=e_melt-snow(i)
820                tice(ki)=tice(ki)+e_melt*ice_lat*2./(ice_cap*seaice(ki))
821                tsurf_new(i)=tice(ki)
822            END IF
823        END IF
824! melt ice from above if Tice>0
825        IF (tice(ki).GT.t_melt) THEN
826            ! quantity of ice melted (kg/m2)
827            e_melt=MAX(seaice(ki)*(tice(ki)-t_melt)*ice_cap/2. &
828             /(ice_lat+ice_cap/2.*(t_melt-t_freeze)),0.0)
829            ! melt from above, height only
830            dhsic=MIN(seaice(ki)-h_ice_min,e_melt)
831            e_melt=e_melt-dhsic
832            IF (e_melt.GT.0) THEN
833            ! lateral melt if ice too thin
834            dfsic=MAX(fsic(ki)-ice_frac_min,e_melt/h_ice_min*fsic(ki))
835            ! if all melted add remaining heat to ocean
836            e_melt=MAX(0.,e_melt*fsic(ki)-dfsic*h_ice_min)
837            slab_bilg(ki)=slab_bilg(ki)+ e_melt*ice_lat/dtime
838            ! update height and fraction
839            fsic(ki)=fsic(ki)-dfsic
840            END IF
841            seaice(ki)=seaice(ki)-dhsic
842            ! surface temperature at melting point
843            tice(ki)=t_melt
844            tsurf_new(i)=t_melt
845        END IF
846        ! convert snow to ice if below floating line
847        h_test=(seaice(ki)+snow(i))*ice_den-seaice(ki)*sea_den
848        IF (h_test.GT.0.) THEN !snow under water
849            ! extra snow converted to ice (with added frozen sea water)
850            dhsic=h_test/(sea_den-ice_den+sno_den)
851            seaice(ki)=seaice(ki)+dhsic
852            snow(i)=snow(i)-dhsic*sno_den/ice_den
853            ! available energy (freeze sea water + bring to tice)
854            e_melt=dhsic*(1.-sno_den/ice_den)*(ice_lat+ &
855                   ice_cap/2.*(t_freeze-tice(ki)))
856            ! update ice temperature
857            tice(ki)=tice(ki)+2.*e_melt/ice_cap/(snow(i)+seaice(ki))
858        END IF
859    END DO
860
[2057]861! New albedo
[2209]862    DO i=1,knon
863        ki=knindex(i)
864       ! snow albedo: update snow age
865        IF (snow(i).GT.0.0001) THEN
866             agesno(i)=(agesno(i) + (1.-agesno(i)/50.)*dtime/86400.)&
867                         * EXP(-1.*MAX(0.0,precip_snow(i))*dtime/5.)
868        ELSE
869            agesno(i)=0.0
870        END IF
871        ! snow albedo
872        alb_snow=alb_sno_min+alb_sno_del*EXP(-agesno(i)/50.)
873        ! ice albedo (varies with ice tkickness and temp)
874        alb_ice=MAX(0.0,0.13*LOG(100.*seaice(ki)/ice_den)+0.1)
875        IF (tice(ki).GT.t_freeze-0.01) THEN
876            alb_ice=MIN(alb_ice,alb_ice_wet)
877        ELSE
878            alb_ice=MIN(alb_ice,alb_ice_dry)
879        END IF
880        ! pond albedo
881        alb_pond=0.36-0.1*(2.0+MIN(0.0,MAX(tice(ki)-t_melt,-2.0)))
882        ! pond fraction
883        frac_pond=0.2*(2.0+MIN(0.0,MAX(tice(ki)-t_melt,-2.0)))
884        ! snow fraction
885        frac_snow=MAX(0.0,MIN(1.0-frac_pond,snow(i)/snow_min))
886        ! ice fraction
887        frac_ice=MAX(0.0,1.-frac_pond-frac_snow)
888        ! total albedo
889        alb1_new(i)=alb_snow*frac_snow+alb_ice*frac_ice+alb_pond*frac_pond
890    END DO
891    alb2_new(:) = alb1_new(:)
[2057]892
893!****************************************************************************************
[2209]894! 3) Recalculate new ocean temperature (add fluxes below ice)
[2057]895!    Melt / freeze from below
896!***********************************************o*****************************************
[2209]897    !cumul fluxes
898    bilg_cum(:)=bilg_cum(:)+slab_bilg(:)
899    IF (MOD(itime,cpl_pas).EQ.0) THEN ! time to update tslab & fraction
900        ! Add cumulated surface fluxes
901        tslab(:,1)=tslab(:,1)+bilg_cum(:)*cyang*dtime
902        DO i=1,knon
903            ki=knindex(i)
904            ! split lateral/top melt-freeze
905            frac_mf=MIN(1.,MAX(0.,(seaice(ki)-h_ice_thin)/(h_ice_thick-h_ice_thin)))
906            IF (tslab(ki,1).LE.t_freeze) THEN
907               ! ****** Form new ice from below *******
908               ! quantity of new ice
909                e_melt=(t_freeze-tslab(ki,1))/cyang &
910                       /(ice_lat+ice_cap/2.*(t_freeze-tice(ki)))
911               ! first increase height to h_thin
912               dhsic=MAX(0.,MIN(h_ice_thin-seaice(ki),e_melt/fsic(ki)))
913               seaice(ki)=dhsic+seaice(ki)
914               e_melt=e_melt-fsic(ki)*dhsic
915               IF (e_melt.GT.0.) THEN
916               ! frac_mf fraction used for lateral increase
917               dfsic=MIN(ice_frac_max-fsic(ki),e_melt*frac_mf/seaice(ki))
918               fsic(ki)=fsic(ki)+dfsic
919               e_melt=e_melt-dfsic*seaice(ki)
920               ! rest used to increase height
921               seaice(ki)=MIN(h_ice_max,seaice(ki)+e_melt/fsic(ki))
922               END IF
923               tslab(ki,1)=t_freeze
924           ELSE ! slab temperature above freezing
925               ! ****** melt ice from below *******
926               ! quantity of melted ice
927               e_melt=(tslab(ki,1)-t_freeze)/cyang &
928                       /(ice_lat+ice_cap/2.*(tice(ki)-t_freeze))
929               ! first decrease height to h_thick
930               dhsic=MAX(0.,MIN(seaice(ki)-h_ice_thick,e_melt/fsic(ki)))
931               seaice(ki)=seaice(ki)-dhsic
932               e_melt=e_melt-fsic(ki)*dhsic
933               IF (e_melt.GT.0) THEN
934               ! frac_mf fraction used for height decrease
935               dhsic=MAX(0.,MIN(seaice(ki)-h_ice_min,e_melt*frac_mf/fsic(ki)))
936               seaice(ki)=seaice(ki)-dhsic
937               e_melt=e_melt-fsic(ki)*dhsic
938               ! rest used to decrease fraction (up to 0!)
939               dfsic=MIN(fsic(ki),e_melt/seaice(ki))
940               ! keep remaining in ocean
941               e_melt=e_melt-dfsic*seaice(ki)
942               END IF
943               tslab(ki,1)=t_freeze+e_melt*ice_lat*cyang
944               fsic(ki)=fsic(ki)-dfsic
945           END IF
946        END DO
947        bilg_cum(:)=0.
948    END IF ! coupling time
949   
[2656]950    !tests ice fraction
[2209]951    WHERE (fsic.LT.ice_frac_min)
952        tslab(:,1)=tslab(:,1)-fsic*seaice*ice_lat*cyang
953        tice=t_melt
954        fsic=0.
955        seaice=0.
956    END WHERE
[2057]957
[2209]958  END SUBROUTINE ocean_slab_ice
[781]959!
960!****************************************************************************************
961!
[2057]962  SUBROUTINE ocean_slab_final
963
964!****************************************************************************************
965! Deallocate module variables
966!****************************************************************************************
967    IF (ALLOCATED(tslab)) DEALLOCATE(tslab)
[2209]968    IF (ALLOCATED(fsic)) DEALLOCATE(fsic)
[2656]969    IF (ALLOCATED(tice)) DEALLOCATE(tice)
970    IF (ALLOCATED(seaice)) DEALLOCATE(seaice)
[2209]971    IF (ALLOCATED(slab_bilg)) DEALLOCATE(slab_bilg)
972    IF (ALLOCATED(bilg_cum)) DEALLOCATE(bilg_cum)
973    IF (ALLOCATED(bils_cum)) DEALLOCATE(bils_cum)
[2656]974    IF (ALLOCATED(taux_cum)) DEALLOCATE(taux_cum)
975    IF (ALLOCATED(tauy_cum)) DEALLOCATE(tauy_cum)
976    IF (ALLOCATED(dt_ekman)) DEALLOCATE(dt_ekman)
977    IF (ALLOCATED(dt_hdiff)) DEALLOCATE(dt_hdiff)
[3002]978    IF (ALLOCATED(dt_gm)) DEALLOCATE(dt_gm)
979    IF (ALLOCATED(dt_qflux)) DEALLOCATE(dt_qflux)
[2057]980
981  END SUBROUTINE ocean_slab_final
982
[781]983END MODULE ocean_slab_mod
Note: See TracBrowser for help on using the repository browser.