source: trunk/LMDZ.GENERIC/libf/phystd/hydrol.F90 @ 3279

Last change on this file since 3279 was 3270, checked in by mturbet, 23 months ago

physiq_mod.F90 cleaning + add moist and dry adjustment heating tendencies

  • Property svn:executable set to *
File size: 15.6 KB
Line 
1subroutine hydrol(ngrid,nq,ptimestep,rnat,tsurf,  &
2     qsurf,dqsurf,dqs_hyd,pcapcal,                &
3     albedo,albedo_bareground,                    &
4     albedo_snow_SPECTV,albedo_co2_ice_SPECTV,    &
5     mu0,pdtsurf,pdtsurf_hyd,hice,                &
6     pctsrf_sic,sea_ice)
7
8  use ioipsl_getin_p_mod, only: getin_p
9  use mod_grid_phy_lmdz, only : klon_glo ! # of physics point on full grid
10  use mod_phys_lmdz_para, only : is_master, gather, scatter
11  use watercommon_h, only: T_h2O_ice_liq, RLFTT, rhowater, mx_eau_sol
12  USE surfdat_h
13  use comdiurn_h
14  USE geometry_mod, only: cell_area
15  USE tracer_h
16!  use slab_ice_h
17  USE ocean_slab_mod, ONLY: alb_ice_max,alb_ice_min,h_alb_ice,snow_min
18  use callkeys_mod, only: albedosnow,alb_ocean,albedoco2ice,ok_slab_ocean,Tsaldiff,maxicethick,co2cond
19  use radinc_h, only : L_NSPECTV
20
21  implicit none
22
23!==================================================================
24!     
25!     Purpose
26!     -------
27!     Calculate the surface hydrology and albedo changes.
28!     Both for oceanic and continental regions
29!     
30!     Authors
31!     -------
32!     Adapted from LMDTERRE by B. Charnay (2010). Further
33!     Modifications by R. Wordsworth (2010).
34!     Spectral albedo by M. Turbet (2015).
35!     
36!     Called by
37!     ---------
38!     physiq.F
39!     
40!     Calls
41!     -----
42!     none
43!
44!     Notes
45!     -----
46!     rnat is terrain type: 0-ocean; 1-continent
47!     
48!==================================================================
49
50      integer ngrid,nq
51
52!     Inputs
53!     ------
54      real snowlayer
55      parameter (snowlayer=33.0)        ! 33 kg/m^2 of snow, equal to a layer of 3.3 cm
56     
57      real oceantime ! this is a relaxation timescale for the oceanbulkavg parameterization
58      parameter (oceantime=10*24*3600)
59      logical,save :: oceanbulkavg ! simple parameterization to relax ocean temperatures to the global mean value (crude, 0th order parameterization to mimick ocean heat transport)
60     
61      logical,save :: activerunoff ! enable simple runoff scheme?
62      logical,save :: oceanalbvary ! simple parameterization to account for the effect of solar zenith angle on the ocean albedo (for the moment it is not used, but to be included in the future)
63!$OMP THREADPRIVATE(oceanbulkavg,activerunoff,oceanalbvary)
64
65!     Arguments
66!     ---------
67      real rnat(ngrid) ! rnat is terrain type: 0-ocean; 1-continent
68      real,dimension(:),allocatable,save :: runoff
69      real totalrunoff, tsea, oceanarea
70      save oceanarea
71!$OMP THREADPRIVATE(runoff,oceanarea)
72
73      real ptimestep
74      real mu0(ngrid)
75      real qsurf(ngrid,nq), tsurf(ngrid)
76      real dqsurf(ngrid,nq), pdtsurf(ngrid)
77      real hice(ngrid)
78      real albedo(ngrid,L_NSPECTV)
79      real albedo_bareground(ngrid)
80      real albedo_snow_SPECTV(L_NSPECTV)
81      real albedo_co2_ice_SPECTV(L_NSPECTV)
82      real pctsrf_sic(ngrid), sea_ice(ngrid)
83
84      real oceanarea2
85
86!     Output
87!     ------
88      real dqs_hyd(ngrid,nq)
89      real pdtsurf_hyd(ngrid)
90
91!     Local
92!     -----
93      real a,b,E
94      integer ig,iq, nw
95      real fsnoi, subli, fauxo
96      real twater(ngrid)
97      real pcapcal(ngrid)
98      real hicebis(ngrid)
99      real zqsurf(ngrid,nq)
100      real ztsurf(ngrid)
101      real albedo_sic, alb_ice
102      real frac_snow
103
104      integer, save :: ivap, iliq, iice
105!$OMP THREADPRIVATE(ivap,iliq,iice)
106
107      logical, save :: firstcall=.true.
108!$OMP THREADPRIVATE(firstcall)
109
110      real :: runoffamount(ngrid)
111!#ifdef CPP_PARA
112      real :: runoffamount_glo(klon_glo)
113      real :: zqsurf_iliq_glo(klon_glo)
114      real :: rnat_glo(klon_glo)
115      real :: oceanarea_glo
116      real :: cell_area_glo(klon_glo)
117!#else
118!      real :: runoffamount_glo(ngrid)
119!      real :: zqsurf_iliq_glo(ngrid)
120!#endif
121
122
123      if(firstcall)then
124
125         oceanbulkavg=.false.
126         oceanalbvary=.false.
127         write(*,*)"Activate runnoff into oceans?"
128         activerunoff=.false.
129         call getin_p("activerunoff",activerunoff)
130         write(*,*)" activerunoff = ",activerunoff
131         
132         if (activerunoff) then
133           ALLOCATE(runoff(ngrid))
134           runoff(1:ngrid)=0
135         endif
136
137         ivap=igcm_h2o_vap
138         iliq=igcm_h2o_vap
139         iice=igcm_h2o_ice
140       
141         write(*,*) "hydrol: ivap=",ivap
142         write(*,*) "        iliq=",iliq
143         write(*,*) "        iice=",iice
144
145!     Here's the deal: iice is used in place of igcm_h2o_ice both on the
146!                      surface and in the atmosphere. ivap is used in
147!                      place of igcm_h2o_vap ONLY in the atmosphere, while
148!                      iliq is used in place of igcm_h2o_vap ONLY on the
149!                      surface.
150!                      Soon to be extended to the entire water cycle...
151
152!     LOCAL ocean surface area
153         oceanarea=0.
154         do ig=1,ngrid
155            if(nint(rnat(ig)).eq.0)then
156               oceanarea=oceanarea+cell_area(ig)
157            endif
158         enddo
159
160         if(oceanbulkavg.and.(oceanarea.le.0.))then
161            print*,'How are we supposed to average the ocean'
162            print*,'temperature, when there are no oceans?'
163            call abort
164         endif
165
166         if(activerunoff.and.(oceanarea.le.0.))then
167            print*,'You have enabled runoff, but you have no oceans.'
168            print*,'Where did you think the water was going to go?'
169            call abort
170         endif
171         
172         firstcall = .false.
173      endif
174
175!      write (*,*) "oceanarea", oceanarea
176
177!     add physical tendencies already calculated
178!     ------------------------------------------
179
180      do ig=1,ngrid
181         ztsurf(ig) = tsurf(ig) + ptimestep*pdtsurf(ig)
182         pdtsurf_hyd(ig)=0.0
183         do iq=1,nq
184            zqsurf(ig,iq) = qsurf(ig,iq) + ptimestep*dqsurf(ig,iq)
185         enddo   
186      enddo
187 
188      do ig=1,ngrid
189         do iq=1,nq
190            dqs_hyd(ig,iq) = 0.0
191         enddo
192      enddo
193
194      do ig = 1, ngrid
195
196!     Ocean regions (rnat = 0)
197!     -----------------------
198         if(nint(rnat(ig)).eq.0)then
199
200!     Parameterization (not used for the moment) to compute the effect of solar zenith angle on the albedo
201!     --------------------------
202!
203!            if(diurnal.and.oceanalbvary)then
204!               fauxo      = ( 1.47 - ACOS( mu0(ig) ) )/0.15 ! where does this come from (Benjamin)?
205!               albedo(ig) = 1.1*( .03 + .630/( 1. + fauxo*fauxo))
206!               albedo(ig) = MAX(MIN(albedo(ig),0.60),0.04)
207!            else
208!
209!               do nw=1,L_NSPECTV
210!                  albedo(ig,nw) = alb_ocean ! For now, alb_ocean is defined in inifis_mod.F90. Later we could introduce spectral dependency for alb_ocean.
211!              enddo
212!            end if
213
214            ! we first start by fixing the albedo of oceanic grid to that of the ocean
215            do nw=1,L_NSPECTV
216               albedo(ig,nw) = alb_ocean ! For now, alb_ocean is defined in inifis_mod.F90. Later we could introduce spectral dependency for alb_ocean.
217            enddo
218
219            if(ok_slab_ocean) then ! if ocean heat transport param activated
220         
221               frac_snow = MAX(0.0,MIN(1.0,zqsurf(ig,iice)/snow_min)) ! Critical snow height (in kg/m2) from ocean_slab_ice routine.
222               ! Standard value should be 15kg/m2 (i.e. about 5 cm). Note that in the previous ocean param. (from BC2014), this value was 45kg/m2 (i.e. about 15cm).
223               alb_ice=alb_ice_max-(alb_ice_max-alb_ice_min)*exp(-sea_ice(ig)/h_alb_ice) ! this is the old formulation from BC2014 (earth-centric so will be replaced)
224               ! Albedo final calculation :
225               do nw=1,L_NSPECTV
226                  albedo(ig,nw) = pctsrf_sic(ig)*                                        &
227                                 (albedo_snow_SPECTV(nw)*frac_snow + alb_ice*(1.0-frac_snow))      &
228                               + (1.-pctsrf_sic(ig))*alb_ocean
229               enddo
230
231               ! Oceanic ice height, just for diagnostics
232               hice(ig)    = MIN(10.,sea_ice(ig)/rhowater)
233            else !ok_slab_ocean ; here this is the case where we are dealing with a static ocean
234
235
236               ! calculate oceanic ice height including the latent heat of ice formation
237               ! hice is the height of oceanic ice with a maximum of maxicethick.
238               hice(ig)    = zqsurf(ig,iice)/rhowater ! update hice to include recent snowfall
239               twater(ig)  = ztsurf(ig) - hice(ig)*RLFTT*rhowater/pcapcal(ig) ! this is the temperature water would have if we melted the entire ocean ice layer
240               hicebis(ig) = hice(ig)
241               hice(ig)    = 0.
242
243               if(twater(ig) .lt. T_h2O_ice_liq)then
244                  E=min((T_h2O_ice_liq+Tsaldiff-twater(ig))*pcapcal(ig),RLFTT*rhowater*maxicethick)
245                  hice(ig)        = E/(RLFTT*rhowater)
246                  hice(ig)        = max(hice(ig),0.0)
247                  hice(ig)        = min(hice(ig),maxicethick)
248                  pdtsurf_hyd(ig) = (hice(ig) - hicebis(ig))*RLFTT*rhowater/pcapcal(ig)/ptimestep             
249                  do nw=1,L_NSPECTV
250                     albedo(ig,nw) = albedo_snow_SPECTV(nw) ! Albedo of ice has been replaced by albedo of snow here. MT2015.
251                  enddo 
252
253!                 if (zqsurf(ig,iice).ge.snowlayer) then
254!                    albedo(ig) = albedoice
255!                 else
256!                    albedo(ig) = albedoocean &
257!                    + (albedosnow - albedoocean)*zqsurf(ig,iice)/snowlayer
258!                 endif
259
260               else
261
262                  pdtsurf_hyd(ig) = -hicebis(ig)*RLFTT*rhowater/pcapcal(ig)/ptimestep
263                  DO nw=1,L_NSPECTV
264                     albedo(ig,nw) = alb_ocean
265                  ENDDO               
266
267               endif
268
269               zqsurf(ig,iliq) = zqsurf(ig,iliq)-(hice(ig)*rhowater-zqsurf(ig,iice))
270               zqsurf(ig,iice) = hice(ig)*rhowater
271
272            endif!(ok_slab_ocean)
273
274
275!     Continental regions (rnat = 1)
276!     -----------------------
277         elseif (nint(rnat(ig)).eq.1) then
278
279            ! melt the snow
280            if(ztsurf(ig).gt.T_h2O_ice_liq)then
281               if(zqsurf(ig,iice).gt.1.0e-8)then
282
283                  a     = (ztsurf(ig)-T_h2O_ice_liq)*pcapcal(ig)/RLFTT
284                  b     = zqsurf(ig,iice)
285                  fsnoi = min(a,b)
286
287                  zqsurf(ig,iice) = zqsurf(ig,iice) - fsnoi
288                  zqsurf(ig,iliq) = zqsurf(ig,iliq) + fsnoi
289
290                  ! thermal effects
291                  pdtsurf_hyd(ig) = -fsnoi*RLFTT/pcapcal(ig)/ptimestep 
292
293               endif
294            else
295
296               ! freeze the water
297               if(zqsurf(ig,iliq).gt.1.0e-8)then
298
299                  a     = -(ztsurf(ig)-T_h2O_ice_liq)*pcapcal(ig)/RLFTT
300                  b     = zqsurf(ig,iliq)
301                 
302                  fsnoi = min(a,b)
303
304                  zqsurf(ig,iice) = zqsurf(ig,iice) + fsnoi
305                  zqsurf(ig,iliq) = zqsurf(ig,iliq) - fsnoi
306
307                  ! thermal effects
308                  pdtsurf_hyd(ig) = +fsnoi*RLFTT/pcapcal(ig)/ptimestep 
309
310               endif
311            endif
312           
313            ! add runoff (to simulate transport of water from continental regions to oceanic regions ; in practice, this prevents liquid water to build up too much on continental regions)
314            if(activerunoff)then
315               runoff(ig) = max(zqsurf(ig,iliq) - mx_eau_sol, 0.0)
316               if(ngrid.gt.1)then ! runoff only exists in 3D
317                  if(runoff(ig).ne.0.0)then
318                     zqsurf(ig,iliq) = mx_eau_sol
319                    ! note: runoff is added to ocean at end
320                  endif
321               end if
322
323            endif
324
325            ! re-calculate continental albedo
326            DO nw=1,L_NSPECTV
327               albedo(ig,nw) = albedo_bareground(ig)
328            ENDDO
329            if (zqsurf(ig,iice).ge.snowlayer) then
330               DO nw=1,L_NSPECTV
331                  albedo(ig,nw) = albedo_snow_SPECTV(nw)
332               ENDDO
333            else
334               DO nw=1,L_NSPECTV
335                  albedo(ig,nw) = albedo_bareground(ig)                            &
336                               + (albedo_snow_SPECTV(nw) - albedo_bareground(ig))  &
337                                 *zqsurf(ig,iice)/snowlayer
338               ENDDO
339            endif
340
341         else
342
343            print*,'Surface type not recognised in hydrol.F!'
344            print*,'Exiting...'
345            call abort
346
347         endif
348
349      end do ! ig=1,ngrid
350
351
352!     simple parameterization to perform crude bulk averaging of temperature in ocean
353!     ----------------------------------------------------
354      if(oceanbulkavg)then
355
356         oceanarea2=0.       
357         DO ig=1,ngrid
358            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0.))then
359               oceanarea2=oceanarea2+cell_area(ig)*pcapcal(ig)
360            end if
361         END DO
362       
363         tsea=0.
364         DO ig=1,ngrid
365            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0.))then       
366               tsea=tsea+ztsurf(ig)*cell_area(ig)*pcapcal(ig)/oceanarea2
367            end if
368         END DO
369
370         DO ig=1,ngrid
371            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0))then
372               pdtsurf_hyd(ig) = pdtsurf_hyd(ig) + (tsea-ztsurf(ig))/oceantime
373            end if
374         END DO
375
376         print*,'Mean ocean temperature               = ',tsea,' K'
377
378      endif
379
380!     shove all the runoff water into the ocean
381!     -----------------------------------------
382      if(activerunoff)then
383
384!         totalrunoff=0.
385         do ig=1,ngrid
386            runoffamount(ig) = cell_area(ig)*runoff(ig)
387!            if (nint(rnat(ig)).eq.1) then
388!               totalrunoff = totalrunoff + cell_area(ig)*runoff(ig)
389!            endif
390         enddo
391
392        ! collect on the full grid
393        call gather(runoffamount,runoffamount_glo)
394        call gather(zqsurf(1:ngrid,iliq),zqsurf_iliq_glo)
395        call gather(rnat,rnat_glo)
396        call gather(cell_area,cell_area_glo)
397
398        if (is_master) then
399                totalrunoff=0.
400                oceanarea_glo=0.
401                do ig=1,klon_glo
402                   if (nint(rnat_glo(ig)).eq.1) then
403                        totalrunoff = totalrunoff + runoffamount_glo(ig)
404                   endif
405                   if (nint(rnat_glo(ig)).eq.0) then
406                        oceanarea_glo = oceanarea_glo + cell_area_glo(ig)
407                   endif
408                enddo
409
410                do ig=1,klon_glo
411                   if (nint(rnat_glo(ig)).eq.0) then
412                        zqsurf_iliq_glo(ig) = zqsurf_iliq_glo(ig) + &
413                    totalrunoff/oceanarea_glo
414                   endif
415                enddo
416
417        endif! is_master
418
419        ! scatter the field back on all processes
420        call scatter(zqsurf_iliq_glo,zqsurf(1:ngrid,iliq))
421
422       
423         
424!         do ig=1,ngrid
425!            if (nint(rnat(ig)).eq.0) then
426!               zqsurf(ig,iliq) = zqsurf(ig,iliq) + &
427!                    totalrunoff/oceanarea
428!            endif
429!         enddo
430
431      endif !activerunoff         
432
433!     Re-add the albedo effects of CO2 ice if necessary
434!     -------------------------------------------------
435      if(co2cond)then
436
437         do ig=1,ngrid
438            if (qsurf(ig,igcm_co2_ice).gt.1.) then ! Condition changed - Need now ~1 mm CO2 ice coverage. MT2015
439               DO nw=1,L_NSPECTV
440                  albedo(ig,nw) = albedo_co2_ice_SPECTV(nw)
441               ENDDO
442            endif
443         enddo ! ngrid
444         
445      endif ! co2cond
446
447
448      do ig=1,ngrid ! We calculate here the tracer tendencies. Don't forget that we have to retrieve the dqsurf tendencies we added at the beginning of the routine !
449         dqs_hyd(ig,iliq)=(zqsurf(ig,iliq) - qsurf(ig,iliq))/ptimestep - dqsurf(ig,iliq)
450         dqs_hyd(ig,iice)=(zqsurf(ig,iice) - qsurf(ig,iice))/ptimestep - dqsurf(ig,iice)
451      enddo
452
453      if (activerunoff) then
454        call writediagfi(ngrid,'runoff','Runoff amount',' ',2,runoff)
455      endif
456
457      return
458    end subroutine hydrol
Note: See TracBrowser for help on using the repository browser.