source: LMDZ6/trunk/libf/phylmdiso/ocean_forced_mod.F90 @ 4143

Last change on this file since 4143 was 4143, checked in by dcugnet, 2 years ago
  • Some variables are renamed or replaced by direct equivalents:
    • iso_indnum -> tracers(:)%iso_iName
    • niso_possibles -> niso
    • iqiso -> iqIsoPha ; index_trac -> itZonIso
    • ok_iso_verif -> isoCheck
    • ntraceurs_zone -> nzone ; ntraciso -> ntiso
    • qperemin -> min_qparent ; masseqmin -> min_qmass ; ratiomin -> min_ratio
  • Some renamed variables are only aliased with the older name (using USE <module>, ONLY: <oldName> => <newName>) in routines where they are repeated many times.
  • Few hard-coded indexes are now computed (examples: ilic, iso, ivap, irneb, iq_vap, iq_liq, iso_H2O, iso_HDO, iso_HTO, iso_O17, iso_O18).
  • The IF(isoCheck) test is now embedded in the check_isotopes_seq and check_isotopes_loc routines (lighter calling).
File size: 19.2 KB
Line 
1!
2! $Id: ocean_forced_mod.F90 3327 2018-05-16 13:58:59Z musat $
3!
4MODULE ocean_forced_mod
5!
6! This module is used for both the sub-surfaces ocean and sea-ice for the case of a
7! forced ocean,  "ocean=force".
8!
9  IMPLICIT NONE
10
11CONTAINS
12!
13!****************************************************************************************
14!
15  SUBROUTINE ocean_forced_noice( &
16       itime, dtime, jour, knon, knindex, &
17       p1lay, cdragh, cdragq, cdragm, precip_rain, precip_snow, &
18       temp_air, spechum, &
19       AcoefH, AcoefQ, BcoefH, BcoefQ, &
20       AcoefU, AcoefV, BcoefU, BcoefV, &
21       ps, u1, v1, gustiness, tsurf_in, &
22       radsol, snow, agesno, &
23       qsurf, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
24       tsurf_new, dflux_s, dflux_l, sens_prec_liq, rhoa &
25#ifdef ISO
26            ,xtprecip_rain, xtprecip_snow, xtspechum,Roce,rlat, &
27            xtsnow,xtevap,h1 & 
28#endif           
29            )
30!
31! This subroutine treats the "open ocean", all grid points that are not entierly covered
32! by ice.
33! The routine receives data from climatologie file limit.nc and does some calculations at the
34! surface.
35!
36    USE dimphy
37    USE calcul_fluxs_mod
38    USE limit_read_mod
39    USE mod_grid_phy_lmdz
40    USE indice_sol_mod
41    USE phys_output_var_mod, ONLY : sens_prec_liq_o, sens_prec_sol_o, lat_prec_liq_o, lat_prec_sol_o
42    use config_ocean_skin_m, only: activate_ocean_skin
43#ifdef ISO
44  USE infotrac_phy, ONLY: ntiso,niso
45    USE isotopes_routines_mod, ONLY: calcul_iso_surf_oce_vectall, &
46&       calcul_iso_surf_sic_vectall   
47#ifdef ISOVERIF
48    USE isotopes_mod, ONLY: iso_eau,ridicule
49    !USE isotopes_verif_mod, ONLY: errmax,errmaxrel,iso_verif_egalite_choix
50    USE isotopes_verif_mod
51#endif
52#endif
53
54    INCLUDE "YOMCST.h"
55    INCLUDE "clesphys.h"
56    INCLUDE "flux_arp.h"
57
58! Input arguments
59!****************************************************************************************
60    INTEGER, INTENT(IN)                      :: itime, jour, knon
61    INTEGER, DIMENSION(klon), INTENT(IN)     :: knindex
62    REAL, INTENT(IN)                         :: dtime
63    REAL, DIMENSION(klon), INTENT(IN)        :: p1lay
64    REAL, DIMENSION(klon), INTENT(IN)        :: cdragh, cdragq, cdragm
65    REAL, DIMENSION(klon), INTENT(IN)        :: precip_rain, precip_snow
66    REAL, DIMENSION(klon), INTENT(IN)        :: temp_air, spechum
67    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefH, AcoefQ, BcoefH, BcoefQ
68    REAL, DIMENSION(klon), INTENT(IN)        :: AcoefU, AcoefV, BcoefU, BcoefV
69    REAL, DIMENSION(klon), INTENT(IN)        :: ps
70    REAL, DIMENSION(klon), INTENT(IN)        :: u1, v1, gustiness
71    REAL, DIMENSION(klon), INTENT(IN)        :: tsurf_in
72    real, intent(in):: rhoa(:) ! (knon) density of moist air  (kg / m3)
73
74#ifdef ISO
75    REAL, DIMENSION(ntiso,klon), INTENT(IN)  :: xtprecip_rain, xtprecip_snow
76    REAL, DIMENSION(ntiso,klon), INTENT(IN)  :: xtspechum
77    real, dimension(klon), intent(IN) :: rlat
78#endif
79
80! In/Output arguments
81!****************************************************************************************
82    REAL, DIMENSION(klon), INTENT(INOUT)     :: radsol
83    REAL, DIMENSION(klon), INTENT(INOUT)     :: snow
84    REAL, DIMENSION(klon), INTENT(INOUT)     :: agesno !? put to 0 in ocean
85#ifdef ISO     
86    REAL, DIMENSION(niso,klon), INTENT(IN)    :: xtsnow
87    REAL, DIMENSION(niso,klon), INTENT(INOUT)    :: Roce
88#endif
89
90! Output arguments
91!****************************************************************************************
92    REAL, DIMENSION(klon), INTENT(OUT)       :: qsurf
93    REAL, DIMENSION(klon), INTENT(OUT)       :: evap, fluxsens, fluxlat
94    REAL, DIMENSION(klon), INTENT(OUT)       :: flux_u1, flux_v1
95    REAL, DIMENSION(klon), INTENT(OUT)       :: tsurf_new
96    REAL, DIMENSION(klon), INTENT(OUT)       :: dflux_s, dflux_l
97    REAL, intent(out):: sens_prec_liq(:) ! (knon)
98
99#ifdef ISO     
100    REAL, DIMENSION(ntiso,klon), INTENT(OUT) :: xtevap ! isotopes in evaporation flux
101    REAL, DIMENSION(klon), INTENT(out)    :: h1 ! just a diagnostic, not useful for the simulation
102#endif
103
104! Local variables
105!****************************************************************************************
106    INTEGER                     :: i, j
107    REAL, DIMENSION(klon)       :: cal, beta, dif_grnd
108    REAL, DIMENSION(klon)       :: alb_neig, tsurf_lim, zx_sl
109    REAL, DIMENSION(klon)       :: u0, v0
110    REAL, DIMENSION(klon)       :: u1_lay, v1_lay
111    LOGICAL                     :: check=.FALSE.
112    REAL sens_prec_sol(knon)
113    REAL, DIMENSION(klon) :: lat_prec_liq, lat_prec_sol   
114
115#ifdef ISO   
116      real, parameter :: t_coup = 273.15     
117#endif
118
119!****************************************************************************************
120! Start calculation
121!****************************************************************************************
122    IF (check) WRITE(*,*)' Entering ocean_forced_noice'
123   
124#ifdef ISO
125#ifdef ISOVERIF
126        do i=1,knon
127          if (iso_eau.gt.0) then         
128              call iso_verif_egalite_choix(xtspechum(iso_eau,i), &
129     &                  spechum(i),'ocean_forced_mod 111', &
130     &                  errmax,errmaxrel)     
131              call iso_verif_egalite_choix(snow(i), &
132     &          xtsnow(iso_eau,i),'ocean_forced_mod 117', &
133     &          errmax,errmaxrel)
134           endif !if (iso_eau.gt.0) then
135        enddo !do i=1,knon
136#endif     
137#endif
138
139!****************************************************************************************
140! 1)   
141! Read sea-surface temperature from file limit.nc
142!
143!****************************************************************************************
144!--sb:
145!!jyg    if (knon.eq.1) then ! single-column model
146    if (klon_glo.eq.1) then ! single-column model
147      ! EV: now surface Tin flux_arp.h
148      !CALL read_tsurf1d(knon,tsurf_lim) ! new
149       DO i = 1, knon
150        tsurf_lim(i) = tg
151       ENDDO
152
153    else ! GCM
154      CALL limit_read_sst(knon,knindex,tsurf_lim &
155#ifdef ISO
156     &  ,Roce,rlat   &
157#endif     
158     &   )
159    endif ! knon
160!sb--
161
162!****************************************************************************************
163! 2)
164! Flux calculation
165!
166!****************************************************************************************
167! Set some variables for calcul_fluxs
168    !cal = 0.
169    !beta = 1.
170    !dif_grnd = 0.
171   
172   
173    ! EV: use calbeta to calculate beta
174    ! Need to initialize qsurf for calbeta but it is not modified by this routine
175    qsurf(:)=0.
176    CALL calbeta(dtime, is_oce, knon, snow, qsurf, beta, cal, dif_grnd)
177
178
179    alb_neig(:) = 0.
180    agesno(:) = 0.
181    lat_prec_liq = 0.; lat_prec_sol = 0.
182
183! Suppose zero surface speed
184    u0(:)=0.0
185    v0(:)=0.0
186    u1_lay(:) = u1(:) - u0(:)
187    v1_lay(:) = v1(:) - v0(:)
188
189! Calcul de tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l and qsurf
190    CALL calcul_fluxs(knon, is_oce, dtime, &
191         merge(tsurf_in, tsurf_lim, activate_ocean_skin == 2), p1lay, cal, &
192         beta, cdragh, cdragq, ps, &
193         precip_rain, precip_snow, snow, qsurf,  &
194         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
195         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
196         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
197         sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, rhoa)
198    if (activate_ocean_skin == 2) tsurf_new = tsurf_lim
199
200    do j = 1, knon
201      i = knindex(j)
202      sens_prec_liq_o(i,1) = sens_prec_liq(j)
203      sens_prec_sol_o(i,1) = sens_prec_sol(j)
204      lat_prec_liq_o(i,1) = lat_prec_liq(j)
205      lat_prec_sol_o(i,1) = lat_prec_sol(j)
206    enddo
207
208
209! - Flux calculation at first modele level for U and V
210    CALL calcul_flux_wind(knon, dtime, &
211         u0, v0, u1, v1, gustiness, cdragm, &
212         AcoefU, AcoefV, BcoefU, BcoefV, &
213         p1lay, temp_air, &
214         flux_u1, flux_v1) 
215 
216
217#ifdef ISO     
218        call calcul_iso_surf_oce_vectall(klon, knon,t_coup, &
219     &    ps,tsurf_new,spechum,u1_lay, v1_lay, xtspechum, &
220     &    evap, Roce,xtevap,h1 &
221#ifdef ISOTRAC
222     &    ,knindex &
223#endif
224     &   )
225#endif         
226
227#ifdef ISO
228#ifdef ISOVERIF
229!          write(*,*) 'ocean_forced_mod 176: sortie de ocean_forced_noice'
230          if (iso_eau.gt.0) then
231           do i=1,knon               
232              call iso_verif_egalite_choix(snow(i), &
233     &          xtsnow(iso_eau,i),'ocean_forced_mod 180', &
234     &          errmax,errmaxrel)
235           enddo ! do j=1,knon
236          endif !if (iso_eau.gt.0) then
237#endif
238#endif   
239
240  END SUBROUTINE ocean_forced_noice
241!
242!***************************************************************************************
243!
244  SUBROUTINE ocean_forced_ice( &
245       itime, dtime, jour, knon, knindex, &
246       tsurf_in, p1lay, cdragh, cdragm, precip_rain, precip_snow, temp_air, spechum, &
247       AcoefH, AcoefQ, BcoefH, BcoefQ, &
248       AcoefU, AcoefV, BcoefU, BcoefV, &
249       ps, u1, v1, gustiness, &
250       radsol, snow, qsol, agesno, tsoil, &
251       qsurf, alb1_new, alb2_new, evap, fluxsens, fluxlat, flux_u1, flux_v1, &
252       tsurf_new, dflux_s, dflux_l, rhoa &
253#ifdef ISO
254            ,xtprecip_rain, xtprecip_snow, xtspechum,Roce, &
255            xtsnow, xtsol,xtevap,Rland_ice & 
256#endif           
257            )
258!
259! This subroutine treats the ocean where there is ice.
260! The routine reads data from climatologie file and does flux calculations at the
261! surface.
262!
263    USE dimphy
264    USE geometry_mod, ONLY: longitude,latitude
265
266    USE calcul_fluxs_mod
267    USE surface_data,     ONLY : calice, calsno
268    USE limit_read_mod
269    USE fonte_neige_mod,  ONLY : fonte_neige
270    USE indice_sol_mod
271    USE phys_output_var_mod, ONLY : sens_prec_liq_o, sens_prec_sol_o, lat_prec_liq_o, lat_prec_sol_o
272#ifdef ISO
273  USE infotrac_phy, ONLY: niso, ntiso
274    USE isotopes_routines_mod, ONLY: calcul_iso_surf_oce_vectall, &
275&       calcul_iso_surf_sic_vectall
276#ifdef ISOVERIF
277    USE isotopes_mod, ONLY: iso_eau,ridicule
278    !USE isotopes_verif_mod, ONLY: errmax,errmaxrel,iso_verif_egalite_choix
279    USE isotopes_verif_mod
280#endif
281#endif
282
283!   INCLUDE "indicesol.h"
284    INCLUDE "dimsoil.h"
285    INCLUDE "YOMCST.h"
286    INCLUDE "clesphys.h"
287    INCLUDE "flux_arp.h"
288
289! Input arguments
290!****************************************************************************************
291    INTEGER, INTENT(IN)                  :: itime, jour, knon
292    INTEGER, DIMENSION(klon), INTENT(IN) :: knindex
293    REAL, INTENT(IN)                     :: dtime
294    REAL, DIMENSION(klon), INTENT(IN)    :: tsurf_in
295    REAL, DIMENSION(klon), INTENT(IN)    :: p1lay
296    REAL, DIMENSION(klon), INTENT(IN)    :: cdragh, cdragm
297    REAL, DIMENSION(klon), INTENT(IN)    :: precip_rain, precip_snow
298    REAL, DIMENSION(klon), INTENT(IN)    :: temp_air, spechum
299    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefH, AcoefQ, BcoefH, BcoefQ
300    REAL, DIMENSION(klon), INTENT(IN)    :: AcoefU, AcoefV, BcoefU, BcoefV
301    REAL, DIMENSION(klon), INTENT(IN)    :: ps
302    REAL, DIMENSION(klon), INTENT(IN)    :: u1, v1, gustiness
303    real, intent(in):: rhoa(:) ! (knon) density of moist air  (kg / m3)
304#ifdef ISO
305    REAL, DIMENSION(ntiso,klon), INTENT(IN)    :: xtprecip_rain, xtprecip_snow
306    REAL, DIMENSION(ntiso,klon), INTENT(IN)    :: xtspechum
307    REAL, DIMENSION(niso,klon),  INTENT(IN)    :: Roce
308    REAL, DIMENSION(niso,klon),  INTENT(IN)    :: Rland_ice
309#endif
310
311! In/Output arguments
312!****************************************************************************************
313    REAL, DIMENSION(klon), INTENT(INOUT)          :: radsol
314    REAL, DIMENSION(klon), INTENT(INOUT)          :: snow, qsol
315    REAL, DIMENSION(klon), INTENT(INOUT)          :: agesno
316    REAL, DIMENSION(klon, nsoilmx), INTENT(INOUT) :: tsoil
317#ifdef ISO     
318    REAL, DIMENSION(niso,klon), INTENT(INOUT)    :: xtsnow
319    REAL, DIMENSION(niso,klon), INTENT(IN)    :: xtsol
320#endif
321
322! Output arguments
323!****************************************************************************************
324    REAL, DIMENSION(klon), INTENT(OUT)            :: qsurf
325    REAL, DIMENSION(klon), INTENT(OUT)            :: alb1_new  ! new albedo in visible SW interval
326    REAL, DIMENSION(klon), INTENT(OUT)            :: alb2_new  ! new albedo in near IR interval
327    REAL, DIMENSION(klon), INTENT(OUT)            :: evap, fluxsens, fluxlat
328    REAL, DIMENSION(klon), INTENT(OUT)            :: flux_u1, flux_v1
329    REAL, DIMENSION(klon), INTENT(OUT)            :: tsurf_new
330    REAL, DIMENSION(klon), INTENT(OUT)            :: dflux_s, dflux_l     
331#ifdef ISO     
332    REAL, DIMENSION(ntiso,klon), INTENT(OUT)    :: xtevap
333#endif     
334
335! Local variables
336!****************************************************************************************
337    LOGICAL                     :: check=.FALSE.
338    INTEGER                     :: i, j
339    REAL                        :: zfra
340    REAL, PARAMETER             :: t_grnd=271.35
341    REAL, DIMENSION(klon)       :: cal, beta, dif_grnd, capsol
342    REAL, DIMENSION(klon)       :: alb_neig, tsurf_tmp
343    REAL, DIMENSION(klon)       :: soilcap, soilflux
344    REAL, DIMENSION(klon)       :: u0, v0
345    REAL, DIMENSION(klon)       :: u1_lay, v1_lay
346    REAL sens_prec_liq(knon), sens_prec_sol (knon)
347    REAL, DIMENSION(klon)       :: lat_prec_liq, lat_prec_sol   
348
349#ifdef ISO       
350      real, parameter :: t_coup = 273.15
351      real, dimension(klon) :: fq_fonte_diag
352      real, dimension(klon) :: fqfonte_diag
353      real, dimension(klon) ::  snow_evap_diag
354      real, dimension(klon) ::  fqcalving_diag
355      real, DIMENSION(klon) :: run_off_lic_diag
356      real :: coeff_rel_diag
357      real  :: max_eau_sol_diag 
358      real, dimension(klon)  ::  runoff_diag   
359      integer ixt
360      REAL, DIMENSION(niso,klon) :: xtsnow_prec,xtsol_prec
361      REAL, DIMENSION(klon) :: snow_prec,qsol_prec 
362#endif     
363
364!****************************************************************************************
365! Start calculation
366!****************************************************************************************
367    IF (check) WRITE(*,*)'Entering surface_seaice, knon=',knon
368
369!****************************************************************************************
370! 1)
371! Flux calculation : tsurf_new, evap, fluxlat, fluxsens, flux_u1, flux_v1
372!                    dflux_s, dflux_l and qsurf
373!****************************************************************************************
374
375    tsurf_tmp(:) = tsurf_in(:)
376
377! calculate the parameters cal, beta, capsol and dif_grnd and then recalculate cal
378    CALL calbeta(dtime, is_sic, knon, snow, qsol, beta, capsol, dif_grnd)
379
380   
381    IF (soil_model) THEN
382! update tsoil and calculate soilcap and soilflux
383       CALL soil(dtime, is_sic, knon, snow, tsurf_tmp, qsol, &
384        & longitude(knindex(1:knon)), latitude(knindex(1:knon)), tsoil,soilcap, soilflux)
385
386       cal(1:knon) = RCPD / soilcap(1:knon)
387       radsol(1:knon) = radsol(1:knon)  + soilflux(1:knon)
388       dif_grnd = 1.0 / tau_gl
389    ELSE
390       dif_grnd = 1.0 / tau_gl
391       cal = RCPD * calice
392       WHERE (snow > 0.0) cal = RCPD * calsno
393    ENDIF
394
395!    beta = 1.0
396    lat_prec_liq = 0.; lat_prec_sol = 0.
397
398! Suppose zero surface speed
399    u0(:)=0.0
400    v0(:)=0.0
401    u1_lay(:) = u1(:) - u0(:)
402    v1_lay(:) = v1(:) - v0(:)
403    CALL calcul_fluxs(knon, is_sic, dtime, &
404         tsurf_tmp, p1lay, cal, beta, cdragh, cdragh, ps, &
405         precip_rain, precip_snow, snow, qsurf,  &
406         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
407         f_qsat_oce,AcoefH, AcoefQ, BcoefH, BcoefQ, &
408         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
409         sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, rhoa)
410    do j = 1, knon
411      i = knindex(j)
412      sens_prec_liq_o(i,2) = sens_prec_liq(j)
413      sens_prec_sol_o(i,2) = sens_prec_sol(j)
414      lat_prec_liq_o(i,2) = lat_prec_liq(j)
415      lat_prec_sol_o(i,2) = lat_prec_sol(j)
416    enddo
417
418! - Flux calculation at first modele level for U and V
419    CALL calcul_flux_wind(knon, dtime, &
420         u0, v0, u1, v1, gustiness, cdragm, &
421         AcoefU, AcoefV, BcoefU, BcoefV, &
422         p1lay, temp_air, &
423         flux_u1, flux_v1) 
424
425!****************************************************************************************
426! 2)
427! Calculations due to snow and runoff
428!
429!****************************************************************************************
430
431#ifdef ISO
432   ! verif
433#ifdef ISOVERIF
434     do i=1,knon
435       if (iso_eau.gt.0) then
436           if (snow(i).gt.ridicule) then
437             call iso_verif_egalite_choix(xtsnow(iso_eau,i),snow(i), &
438     &           'interfsurf 964',errmax,errmaxrel)
439            endif !if ((snow(i).gt.ridicule)) then
440        endif !if (iso_eau.gt.0) then     
441      enddo !do i=1,knon 
442#endif
443   ! end verif
444
445    do i=1,knon
446      snow_prec(i)=snow(i)
447      do ixt=1,niso
448      xtsnow_prec(ixt,i)=xtsnow(ixt,i)
449      enddo !do ixt=1,niso
450      ! initialisation:
451      fq_fonte_diag(i)=0.0
452      fqfonte_diag(i)=0.0
453      snow_evap_diag(i)=0.0
454    enddo !do i=1,knon
455#endif
456
457    CALL fonte_neige( knon, is_sic, knindex, dtime, &
458         tsurf_tmp, precip_rain, precip_snow, &
459         snow, qsol, tsurf_new, evap &
460#ifdef ISO   
461     & ,fq_fonte_diag,fqfonte_diag,snow_evap_diag,fqcalving_diag   &
462     & ,max_eau_sol_diag,runoff_diag,run_off_lic_diag,coeff_rel_diag   &
463#endif
464     &   )
465
466
467#ifdef ISO
468! isotopes: tout est externalisé
469!#ifdef ISOVERIF
470!        write(*,*) 'ocean_forced_mod 377: call calcul_iso_surf_sic_vectall'
471!        write(*,*) 'klon,knon=',klon,knon
472!#endif
473         call calcul_iso_surf_sic_vectall(klon,knon, &
474             &   evap,snow_evap_diag,Tsurf_new,Roce,snow, &
475             &   fq_fonte_diag,fqfonte_diag,dtime,t_coup, &
476             &   precip_snow,xtprecip_snow,xtprecip_rain, snow_prec,xtsnow_prec, &
477             &   xtspechum,spechum,ps, &
478             &   xtevap,xtsnow,fqcalving_diag, &
479             &   knindex,is_sic,run_off_lic_diag,coeff_rel_diag,Rland_ice &
480     &  )   
481#ifdef ISOVERIF
482        !write(*,*) 'ocean_forced_mod 391: sortie calcul_iso_surf_sic_vectall'
483          if (iso_eau.gt.0) then
484           do i=1,knon 
485              call iso_verif_egalite_choix(snow(i), &
486     &          xtsnow(iso_eau,i),'ocean_forced_mod 396', &
487     &          errmax,errmaxrel)
488           enddo ! do j=1,knon
489          endif !if (iso_eau.gt.0) then
490#endif
491!#ifdef ISOVERIF
492#endif   
493!#ifdef ISO
494   
495! Calculation of albedo at snow (alb_neig) and update the age of snow (agesno)
496!
497    CALL albsno(klon, knon, dtime, agesno(:), alb_neig(:), precip_snow(:)) 
498
499    WHERE (snow(1:knon) .LT. 0.0001) agesno(1:knon) = 0.
500
501    alb1_new(:) = 0.0
502    DO i=1, knon
503       zfra = MAX(0.0,MIN(1.0,snow(i)/(snow(i)+10.0)))
504       alb1_new(i) = alb_neig(i) * zfra +  0.6 * (1.0-zfra)
505    ENDDO
506
507    alb2_new(:) = alb1_new(:)
508
509  END SUBROUTINE ocean_forced_ice
510
511!************************************************************************
512! 1D case
513!************************************************************************
514!  SUBROUTINE read_tsurf1d(knon,sst_out)
515!
516! This subroutine specifies the surface temperature to be used in 1D simulations
517!
518!      USE dimphy, ONLY : klon
519!
520!      INTEGER, INTENT(IN)                  :: knon     ! nomber of points on compressed grid
521!      REAL, DIMENSION(klon), INTENT(OUT)   :: sst_out  ! tsurf used to force the single-column model
522!
523!       INTEGER :: i
524! COMMON defined in lmdz1d.F:
525!       real ts_cur
526!       common /sst_forcing/ts_cur
527!
528!       DO i = 1, knon
529!        sst_out(i) = ts_cur
530!       ENDDO
531!
532!      END SUBROUTINE read_tsurf1d
533!
534!
535!************************************************************************
536END MODULE ocean_forced_mod
537
538
539
540
541
542
Note: See TracBrowser for help on using the repository browser.