source: LMDZ6/branches/LMDZ-INCA-Dyn/libf/phylmd/surf_ocean_mod.F90 @ 5213

Last change on this file since 5213 was 3395, checked in by lguez, 6 years ago

Indent file.

  • 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: 11.5 KB
Line 
1!
2! $Id: surf_ocean_mod.F90 3395 2018-09-25 15:22:13Z evignon $
3!
4MODULE surf_ocean_mod
5
6  IMPLICIT NONE
7
8CONTAINS
9  !
10  !******************************************************************************
11  !
12  SUBROUTINE surf_ocean(rlon, rlat, swnet, lwnet, alb1, &
13       windsp, rmu0, fder, tsurf_in, &
14       itime, dtime, jour, knon, knindex, &
15       p1lay, z1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, &
16       AcoefH, AcoefQ, BcoefH, BcoefQ, &
17       AcoefU, AcoefV, BcoefU, BcoefV, &
18       ps, u1, v1, gustiness, rugoro, pctsrf, &
19       snow, qsurf, agesno, &
20       z0m, z0h, SFRWL, alb_dir_new, alb_dif_new, evap, fluxsens, fluxlat, &
21       tsurf_new, dflux_s, dflux_l, lmt_bils, &
22       flux_u1, flux_v1)
23
24    use albedo, only: alboc, alboc_cd
25    USE dimphy, ONLY: klon, zmasq
26    USE surface_data, ONLY     : type_ocean
27    USE ocean_forced_mod, ONLY : ocean_forced_noice
28    USE ocean_slab_mod, ONLY   : ocean_slab_noice
29    USE ocean_cpl_mod, ONLY    : ocean_cpl_noice
30    USE indice_sol_mod, ONLY : nbsrf, is_oce
31    USE limit_read_mod
32    !
33    ! This subroutine will make a call to ocean_XXX_noice according to the ocean mode (force,
34    ! slab or couple). The calculations of albedo and rugosity for the ocean surface are
35    ! done in here because they are identical for the different modes of ocean.
36
37
38    INCLUDE "YOMCST.h"
39
40    include "clesphys.h"
41    ! for cycle_diurne and for iflag_z0_oce==-1 (prescribed z0)
42
43    ! Input variables
44    !******************************************************************************
45    INTEGER, INTENT(IN)                      :: itime, jour, knon
46    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
47    REAL, INTENT(IN)                         :: dtime
48    REAL, DIMENSION(klon), INTENT(IN)        :: rlon, rlat
49    REAL, DIMENSION(klon), INTENT(IN)        :: swnet  ! net shortwave radiation at surface 
50    REAL, DIMENSION(klon), INTENT(IN)        :: lwnet  ! net longwave radiation at surface 
51    REAL, DIMENSION(klon), INTENT(IN)        :: alb1   ! albedo in visible SW interval
52    REAL, DIMENSION(klon), INTENT(IN)        :: windsp
53    REAL, DIMENSION(klon), INTENT(IN)        :: rmu0 
54    REAL, DIMENSION(klon), INTENT(IN)        :: fder
55    REAL, DIMENSION(klon), INTENT(IN)        :: tsurf_in
56    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay,z1lay ! pression (Pa) et altitude (m) du premier niveau
57    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh
58    REAL, DIMENSION(klon), INTENT(IN)        :: cdragm
59    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow
60    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
61    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
62    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
63    REAL, DIMENSION(klon), INTENT(IN)        :: ps
64    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
65    REAL, DIMENSION(klon), INTENT(IN)        :: rugoro
66    REAL, DIMENSION(klon,nbsrf), INTENT(IN)  :: pctsrf
67
68    ! In/Output variables
69    !******************************************************************************
70    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
71    REAL, DIMENSION(klon), INTENT(INOUT)     :: qsurf
72    REAL, DIMENSION(klon), INTENT(INOUT)     :: agesno
73    REAL, DIMENSION(klon), INTENT(inOUT):: z0h
74
75    ! Output variables
76    !******************************************************************************
77    REAL, DIMENSION(klon), INTENT(OUT)       :: z0m
78    !albedo SB >>>
79    !    REAL, DIMENSION(klon), INTENT(OUT)       :: alb1_new  ! new albedo in visible SW interval
80    !    REAL, DIMENSION(klon), INTENT(OUT)       :: alb2_new  ! new albedo in near IR interval
81    REAL, DIMENSION(6), INTENT(IN)          :: SFRWL
82    REAL, DIMENSION(klon,nsw), INTENT(OUT)       :: alb_dir_new,alb_dif_new
83    !albedo SB <<<     
84    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
85    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new
86    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l     
87    REAL, DIMENSION(klon), INTENT(OUT)       :: lmt_bils
88    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
89
90    ! Local variables
91    !******************************************************************************
92    INTEGER               :: i, k
93    REAL                  :: tmp
94    REAL, PARAMETER       :: cepdu2=(0.1)**2
95    REAL, DIMENSION(klon) :: alb_eau, z0_lim
96    REAL, DIMENSION(klon) :: radsol
97    REAL, DIMENSION(klon) :: cdragq ! Cdrag pour l'evaporation
98    CHARACTER(len=20),PARAMETER :: modname="surf_ocean"
99
100    ! End definition
101    !******************************************************************************
102
103
104    !******************************************************************************
105    ! Calculate total net radiance at surface
106    !
107    !******************************************************************************
108    radsol(1:klon) = 0.0 ! initialisation a priori inutile
109    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
110
111    !******************************************************************************
112    ! Cdragq computed from cdrag
113    ! The difference comes only from a factor (f_z0qh_oce) on z0, so that
114    ! it can be computed inside surf_ocean
115    ! More complicated appraches may require the propagation through
116    ! pbl_surface of an independant cdragq variable.
117    !******************************************************************************
118
119    IF ( f_z0qh_oce .ne. 1.) THEN
120       ! Si on suit les formulations par exemple de Tessel, on
121       ! a z0h=0.4*nu/u*, z0q=0.62*nu/u*, d'ou f_z0qh_oce=0.62/0.4=1.55
122       cdragq(1:knon)=cdragh(1:knon)*                                      &
123            log(z1lay(1:knon)/z0h(1:knon))/log(z1lay(1:knon)/(f_z0qh_oce*z0h(1:knon)))
124    ELSE
125       cdragq(1:knon)=cdragh(1:knon)
126    ENDIF
127
128
129    !******************************************************************************
130    ! Switch according to type of ocean (couple, slab or forced)
131    !******************************************************************************
132    SELECT CASE(type_ocean)
133    CASE('couple')
134       CALL ocean_cpl_noice( &
135            swnet, lwnet, alb1, &
136            windsp, fder, &
137            itime, dtime, knon, knindex, &
138            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow,temp_air,spechum,&
139            AcoefH, AcoefQ, BcoefH, BcoefQ, &
140            AcoefU, AcoefV, BcoefU, BcoefV, &
141            ps, u1, v1, gustiness, &
142            radsol, snow, agesno, &
143            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
144            tsurf_new, dflux_s, dflux_l)
145
146    CASE('slab')
147       CALL ocean_slab_noice( &
148            itime, dtime, jour, knon, knindex, &
149            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, temp_air, spechum,&
150            AcoefH, AcoefQ, BcoefH, BcoefQ, &
151            AcoefU, AcoefV, BcoefU, BcoefV, &
152            ps, u1, v1, gustiness, tsurf_in, &
153            radsol, snow, &
154            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
155            tsurf_new, dflux_s, dflux_l, lmt_bils)
156
157    CASE('force')
158       CALL ocean_forced_noice( &
159            itime, dtime, jour, knon, knindex, &
160            p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, &
161            temp_air, spechum, &
162            AcoefH, AcoefQ, BcoefH, BcoefQ, &
163            AcoefU, AcoefV, BcoefU, BcoefV, &
164            ps, u1, v1, gustiness, &
165            radsol, snow, agesno, &
166            qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
167            tsurf_new, dflux_s, dflux_l)
168    END SELECT
169
170    !******************************************************************************
171    ! fcodron: compute lmt_bils  forced case (same as wfbils_oce / 1.-contfracatm)
172    !******************************************************************************
173    IF (type_ocean.NE.'slab') THEN
174       lmt_bils(1:klon)=0.
175       DO i=1,knon
176          lmt_bils(knindex(i))=(swnet(i)+lwnet(i)+fluxsens(i)+fluxlat(i)) &
177               *pctsrf(knindex(i),is_oce)/(1.-zmasq(knindex(i)))
178       END DO
179    END IF
180
181    !******************************************************************************
182    ! Calculate ocean surface albedo
183    !******************************************************************************
184    !albedo SB >>>
185    IF (iflag_albedo==0) THEN
186       !--old parametrizations of ocean surface albedo
187       !
188       IF (iflag_cycle_diurne.GE.1) THEN
189          !
190          CALL alboc_cd(rmu0,alb_eau)
191          !
192          !--ad-hoc correction for model radiative balance tuning
193          !--now outside alboc_cd routine
194          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
195          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.0),1.0)
196          !
197       ELSE
198          !
199          CALL alboc(REAL(jour),rlat,alb_eau)
200          !--ad-hoc correction for model radiative balance tuning
201          !--now outside alboc routine
202          alb_eau(1:klon) = fmagic*alb_eau(1:klon) + pmagic
203          alb_eau(1:klon)=MIN(MAX(alb_eau(1:klon),0.04),0.60)
204          !
205       ENDIF
206       !
207       DO i =1, knon
208          DO  k=1,nsw
209             alb_dir_new(i,k) = alb_eau(knindex(i))
210          ENDDO
211       ENDDO
212       !IM 09122015 next line corresponds to the old way of doing in LMDZ5A/IPSLCM5A versions
213       !albedo for diffuse radiation is taken the same as for direct radiation
214       alb_dif_new(1:knon,:)=alb_dir_new(1:knon,:)
215       !IM 09122015 end
216       !
217    ELSE IF (iflag_albedo==1) THEN
218       !--new parametrization of ocean surface albedo by Sunghye Baek
219       !--albedo for direct and diffuse radiation are different
220       !
221       CALL ocean_albedo(knon,rmu0,knindex,windsp,SFRWL,alb_dir_new,alb_dif_new)
222       !
223       !--ad-hoc correction for model radiative balance tuning
224       alb_dir_new(1:knon,:) = fmagic*alb_dir_new(1:knon,:) + pmagic
225       alb_dif_new(1:knon,:) = fmagic*alb_dif_new(1:knon,:) + pmagic
226       alb_dir_new(1:knon,:)=MIN(MAX(alb_dir_new(1:knon,:),0.0),1.0)
227       alb_dif_new(1:knon,:)=MIN(MAX(alb_dif_new(1:knon,:),0.0),1.0)
228       !
229    ELSE IF (iflag_albedo==2) THEN
230       ! F. Codron albedo read from limit.nc
231       CALL limit_read_rug_alb(itime, dtime, jour,&
232            knon, knindex, z0_lim, alb_eau)
233       DO i =1, knon
234          DO  k=1,nsw
235             alb_dir_new(i,k) = alb_eau(i)
236          ENDDO
237       ENDDO
238       alb_dif_new=alb_dir_new
239    ENDIF
240    !albedo SB <<<
241
242    !******************************************************************************
243    ! Calculate the rugosity
244    !******************************************************************************
245    IF (iflag_z0_oce==0) THEN
246       DO i = 1, knon
247          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
248          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
249               +  0.11*14e-6 / SQRT(cdragm(i) * tmp)
250          z0m(i) = MAX(1.5e-05,z0m(i))
251       ENDDO
252       z0h(1:knon)=z0m(1:knon) ! En attendant mieux
253
254    ELSE IF (iflag_z0_oce==1) THEN
255       DO i = 1, knon
256          tmp = MAX(cepdu2,gustiness(i)+u1(i)**2+v1(i)**2)
257          z0m(i) = 0.018*cdragm(i) * (gustiness(i)+u1(i)**2+v1(i)**2)/RG  &
258               + 0.11*14e-6 / SQRT(cdragm(i) * tmp)
259          z0m(i) = MAX(1.5e-05,z0m(i))
260          z0h(i)=0.4*14e-6 / SQRT(cdragm(i) * tmp)
261       ENDDO
262    ELSE IF (iflag_z0_oce==-1) THEN
263       DO i = 1, knon
264          z0m(i) = z0min
265          z0h(i) = z0min
266       ENDDO
267    ELSE
268       CALL abort_physic(modname,'version non prevue',1)
269    ENDIF
270    !
271    !******************************************************************************
272  END SUBROUTINE surf_ocean
273  !******************************************************************************
274  !
275END MODULE surf_ocean_mod
Note: See TracBrowser for help on using the repository browser.