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

Last change on this file since 3233 was 3194, checked in by emillour, 2 years ago

Generic PCM-WRF:
Small corrections to enable compilation with WRF4.
NC

  • Property svn:executable set to *
File size: 14.7 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
18  USE ocean_slab_mod, ONLY: alb_ice_min
19  USE ocean_slab_mod, ONLY: h_alb_ice
20  use callkeys_mod, only: albedosnow,alb_ocean,albedoco2ice,ok_slab_ocean,Tsaldiff,maxicethick,co2cond
21  use radinc_h, only : L_NSPECTV
22
23  implicit none
24
25!==================================================================
26!     
27!     Purpose
28!     -------
29!     Calculate the surface hydrology and albedo changes.
30!     
31!     Authors
32!     -------
33!     Adapted from LMDTERRE by B. Charnay (2010). Further
34!     Modifications by R. Wordsworth (2010).
35!     Spectral albedo by M. Turbet (2015).
36!     
37!     Called by
38!     ---------
39!     physiq.F
40!     
41!     Calls
42!     -----
43!     none
44!
45!     Notes
46!     -----
47!     rnat is terrain type: 0-ocean; 1-continent
48!     
49!==================================================================
50
51      integer ngrid,nq
52
53!     Inputs
54!     ------
55      real snowlayer
56      parameter (snowlayer=33.0)        ! 33 kg/m^2 of snow, equal to a layer of 3.3 cm
57      real oceantime
58      parameter (oceantime=10*24*3600)
59
60      logical,save :: oceanbulkavg ! relax ocean temperatures to a GLOBAL mean value?
61      logical,save :: activerunoff ! enable simple runoff scheme?
62      logical,save :: oceanalbvary ! ocean albedo varies with the diurnal cycle?
63!$OMP THREADPRIVATE(oceanbulkavg,activerunoff,oceanalbvary)
64
65!     Arguments
66!     ---------
67      real rnat(ngrid) ! I changed this to integer (RW)
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 zfra
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         
133         
134         if (activerunoff) then
135           ALLOCATE(runoff(ngrid))
136           runoff(1:ngrid)=0
137         endif
138
139         ivap=igcm_h2o_vap
140         iliq=igcm_h2o_vap
141         iice=igcm_h2o_ice
142       
143         write(*,*) "hydrol: ivap=",ivap
144         write(*,*) "        iliq=",iliq
145         write(*,*) "        iice=",iice
146
147!     Here's the deal: iice is used in place of igcm_h2o_ice both on the
148!                      surface and in the atmosphere. ivap is used in
149!                      place of igcm_h2o_vap ONLY in the atmosphere, while
150!                      iliq is used in place of igcm_h2o_vap ONLY on the
151!                      surface.
152!                      Soon to be extended to the entire water cycle...
153
154!     LOCAL ocean surface area
155         oceanarea=0.
156         do ig=1,ngrid
157            if(nint(rnat(ig)).eq.0)then
158               oceanarea=oceanarea+cell_area(ig)
159            endif
160         enddo
161
162         if(oceanbulkavg.and.(oceanarea.le.0.))then
163            print*,'How are we supposed to average the ocean'
164            print*,'temperature, when there are no oceans?'
165            call abort
166         endif
167
168         if(activerunoff.and.(oceanarea.le.0.))then
169            print*,'You have enabled runoff, but you have no oceans.'
170            print*,'Where did you think the water was going to go?'
171            call abort
172         endif
173         
174         firstcall = .false.
175      endif
176
177!       call writediagfi(ngrid,"hydrol_rnat"," "," ",2,rnat)
178!!       write (*,*) "oceanarea", oceanarea
179
180!     add physical tendencies already calculated
181!     ------------------------------------------
182
183      do ig=1,ngrid
184         ztsurf(ig) = tsurf(ig) + ptimestep*pdtsurf(ig)
185         pdtsurf_hyd(ig)=0.0
186         do iq=1,nq
187            zqsurf(ig,iq) = qsurf(ig,iq) + ptimestep*dqsurf(ig,iq)
188         enddo   
189      enddo
190
191!      call writediagfi(ngrid,"hydrol_zqsurf_iliq_1"," "," ",2,zqsurf(1:ngrid,iliq))
192 
193      do ig=1,ngrid
194         do iq=1,nq
195            dqs_hyd(ig,iq) = 0.0
196         enddo
197      enddo
198
199      do ig = 1, ngrid
200
201!     Ocean
202!     -----
203         if(nint(rnat(ig)).eq.0)then
204
205!     re-calculate oceanic albedo
206!            if(diurnal.and.oceanalbvary)then
207!               fauxo      = ( 1.47 - ACOS( mu0(ig) ) )/0.15 ! where does this come from (Benjamin)?
208!               albedo(ig) = 1.1*( .03 + .630/( 1. + fauxo*fauxo))
209!               albedo(ig) = MAX(MIN(albedo(ig),0.60),0.04)
210!            else
211               do nw=1,L_NSPECTV
212                  albedo(ig,nw) = alb_ocean ! For now, alb_ocean is defined in inifis_mod.F90. Later we could introduce spectral dependency for alb_ocean.
213               enddo
214!            end if
215
216
217            if(ok_slab_ocean) then
218         
219               zfra = MAX(0.0,MIN(1.0,zqsurf(ig,iice)/45.0))     ! Snow Fraction (Critical height 45kg/m2~15cm)
220               alb_ice=alb_ice_max-(alb_ice_max-alb_ice_min) &   ! Ice Albedo
221               *exp(-sea_ice(ig)/h_alb_ice)
222               ! Albedo final calculation :
223               do nw=1,L_NSPECTV
224                  albedo(ig,nw) = pctsrf_sic(ig)*                                        &
225                                 (albedo_snow_SPECTV(nw)*zfra + alb_ice*(1.0-zfra))      &
226                               + (1.-pctsrf_sic(ig))*alb_ocean
227               enddo
228
229               ! Oceanic ice height, just for diagnostics
230               hice(ig)    = MIN(10.,sea_ice(ig)/rhowater)
231            else !ok_slab_ocean
232
233
234!     calculate oceanic ice height including the latent heat of ice formation
235!     hice is the height of oceanic ice with a maximum of maxicethick.
236               hice(ig)    = zqsurf(ig,iice)/rhowater ! update hice to include recent snowfall
237
238!              twater(ig)  = tsurf(ig) + ptimestep*zdtsurf(ig) &
239               twater(ig)  = ztsurf(ig) - hice(ig)*RLFTT*rhowater/pcapcal(ig)
240               ! this is temperature water would have if we melted entire ocean ice layer
241               hicebis(ig) = hice(ig)
242               hice(ig)    = 0.
243
244               if(twater(ig) .lt. T_h2O_ice_liq)then
245                  E=min((T_h2O_ice_liq+Tsaldiff-twater(ig))*pcapcal(ig),RLFTT*rhowater*maxicethick)
246                  hice(ig)        = E/(RLFTT*rhowater)
247                  hice(ig)        = max(hice(ig),0.0)
248                  hice(ig)        = min(hice(ig),maxicethick)
249                  pdtsurf_hyd(ig) = (hice(ig) - hicebis(ig))*RLFTT*rhowater/pcapcal(ig)/ptimestep             
250                  do nw=1,L_NSPECTV
251                     albedo(ig,nw) = albedo_snow_SPECTV(nw) ! Albedo of ice has been replaced by albedo of snow here. MT2015.
252                  enddo 
253
254!                 if (zqsurf(ig,iice).ge.snowlayer) then
255!                    albedo(ig) = albedoice
256!                 else
257!                    albedo(ig) = albedoocean &
258!                    + (albedosnow - albedoocean)*zqsurf(ig,iice)/snowlayer
259!                 endif
260
261               else
262
263                  pdtsurf_hyd(ig) = -hicebis(ig)*RLFTT*rhowater/pcapcal(ig)/ptimestep
264                  DO nw=1,L_NSPECTV
265                     albedo(ig,nw) = alb_ocean
266                  ENDDO               
267
268               endif
269
270               zqsurf(ig,iliq) = zqsurf(ig,iliq)-(hice(ig)*rhowater-zqsurf(ig,iice))
271               zqsurf(ig,iice) = hice(ig)*rhowater
272
273            endif!(ok_slab_ocean)
274
275
276!     Continent
277!     ---------
278         elseif (nint(rnat(ig)).eq.1) then
279
280!     melt the snow
281            if(ztsurf(ig).gt.T_h2O_ice_liq)then
282               if(zqsurf(ig,iice).gt.1.0e-8)then
283
284                  a     = (ztsurf(ig)-T_h2O_ice_liq)*pcapcal(ig)/RLFTT
285                  b     = zqsurf(ig,iice)
286                  fsnoi = min(a,b)
287
288                  zqsurf(ig,iice) = zqsurf(ig,iice) - fsnoi
289                  zqsurf(ig,iliq) = zqsurf(ig,iliq) + fsnoi
290
291!                 thermal effects
292                  pdtsurf_hyd(ig) = -fsnoi*RLFTT/pcapcal(ig)/ptimestep 
293
294               endif
295            else
296
297!     freeze the water
298               if(zqsurf(ig,iliq).gt.1.0e-8)then
299
300                  a     = -(ztsurf(ig)-T_h2O_ice_liq)*pcapcal(ig)/RLFTT
301                  b     = zqsurf(ig,iliq)
302                 
303                  fsnoi = min(a,b)
304
305                  zqsurf(ig,iice) = zqsurf(ig,iice) + fsnoi
306                  zqsurf(ig,iliq) = zqsurf(ig,iliq) - fsnoi
307
308!                 thermal effects
309                  pdtsurf_hyd(ig) = +fsnoi*RLFTT/pcapcal(ig)/ptimestep 
310
311               endif
312            endif
313           
314!     deal with runoff
315            if(activerunoff)then
316               runoff(ig) = max(zqsurf(ig,iliq) - mx_eau_sol, 0.0)
317               if(ngrid.gt.1)then ! runoff only exists in 3D
318                  if(runoff(ig).ne.0.0)then
319                     zqsurf(ig,iliq) = mx_eau_sol
320!                    runoff is added to ocean at end
321                  endif
322               end if
323
324            endif
325
326!     re-calculate continental albedo
327            DO nw=1,L_NSPECTV
328               albedo(ig,nw) = albedo_bareground(ig)
329            ENDDO
330            if (zqsurf(ig,iice).ge.snowlayer) then
331               DO nw=1,L_NSPECTV
332                  albedo(ig,nw) = albedo_snow_SPECTV(nw)
333               ENDDO
334            else
335               DO nw=1,L_NSPECTV
336                  albedo(ig,nw) = albedo_bareground(ig)                            &
337                               + (albedo_snow_SPECTV(nw) - albedo_bareground(ig))  &
338                                 *zqsurf(ig,iice)/snowlayer
339               ENDDO
340            endif
341
342         else
343
344            print*,'Surface type not recognised in hydrol.F!'
345            print*,'Exiting...'
346            call abort
347
348         endif
349
350      end do ! ig=1,ngrid
351
352!      call writediagfi(ngrid,"hydrol_zqsurf_iliq_2"," "," ",2,zqsurf(1:ngrid,iliq))
353
354!     perform crude bulk averaging of temperature in ocean
355!     ----------------------------------------------------
356      if(oceanbulkavg)then
357
358         oceanarea2=0.       
359         DO ig=1,ngrid
360            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0.))then
361               oceanarea2=oceanarea2+cell_area(ig)*pcapcal(ig)
362            end if
363         END DO
364       
365         tsea=0.
366         DO ig=1,ngrid
367            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0.))then       
368               tsea=tsea+ztsurf(ig)*cell_area(ig)*pcapcal(ig)/oceanarea2
369            end if
370         END DO
371
372         DO ig=1,ngrid
373            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0))then
374               pdtsurf_hyd(ig) = pdtsurf_hyd(ig) + (tsea-ztsurf(ig))/oceantime
375            end if
376         END DO
377
378         print*,'Mean ocean temperature               = ',tsea,' K'
379
380      endif
381
382!     shove all the runoff water into the ocean
383!     -----------------------------------------
384      if(activerunoff)then
385
386!         totalrunoff=0.
387         do ig=1,ngrid
388            runoffamount(ig) = cell_area(ig)*runoff(ig)
389!            if (nint(rnat(ig)).eq.1) then
390!               totalrunoff = totalrunoff + cell_area(ig)*runoff(ig)
391!            endif
392         enddo
393
394!    collect on the full grid
395        call gather(runoffamount,runoffamount_glo)
396        call gather(zqsurf(1:ngrid,iliq),zqsurf_iliq_glo)
397        call gather(rnat,rnat_glo)
398        call gather(cell_area,cell_area_glo)
399
400        if (is_master) then
401                totalrunoff=0.
402                oceanarea_glo=0.
403                do ig=1,klon_glo
404                   if (nint(rnat_glo(ig)).eq.1) then
405                        totalrunoff = totalrunoff + runoffamount_glo(ig)
406                   endif
407                   if (nint(rnat_glo(ig)).eq.0) then
408                        oceanarea_glo = oceanarea_glo + cell_area_glo(ig)
409                   endif
410                enddo
411
412                do ig=1,klon_glo
413                   if (nint(rnat_glo(ig)).eq.0) then
414                        zqsurf_iliq_glo(ig) = zqsurf_iliq_glo(ig) + &
415                    totalrunoff/oceanarea_glo
416                   endif
417                enddo
418
419        endif! is_master
420
421!    scatter the field back on all processes
422        call scatter(zqsurf_iliq_glo,zqsurf(1:ngrid,iliq))
423
424       
425         
426!         do ig=1,ngrid
427!            if (nint(rnat(ig)).eq.0) then
428!               zqsurf(ig,iliq) = zqsurf(ig,iliq) + &
429!                    totalrunoff/oceanarea
430!            endif
431!         enddo
432
433      endif !activerunoff         
434
435!      call writediagfi(ngrid,"hydrol_zqsurf_iliq_3"," "," ",2,zqsurf(1:ngrid,iliq))
436
437!     Re-add the albedo effects of CO2 ice if necessary
438!     -------------------------------------------------
439      if(co2cond)then
440
441         do ig=1,ngrid
442            if (qsurf(ig,igcm_co2_ice).gt.1.) then ! Condition changed - Need now ~1 mm CO2 ice coverage. MT2015
443               DO nw=1,L_NSPECTV
444                  albedo(ig,nw) = albedo_co2_ice_SPECTV(nw)
445               ENDDO
446            endif
447         enddo ! ngrid
448         
449      endif ! co2cond
450
451
452      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 !
453         dqs_hyd(ig,iliq)=(zqsurf(ig,iliq) - qsurf(ig,iliq))/ptimestep - dqsurf(ig,iliq)
454         dqs_hyd(ig,iice)=(zqsurf(ig,iice) - qsurf(ig,iice))/ptimestep - dqsurf(ig,iice)
455      enddo
456
457!      call writediagfi(ngrid,"hydrol_dqs_hyd_iliq"," "," ",2,dqs_hyd(1:ngrid,iliq)) 
458
459      if (activerunoff) then
460        call writediagfi(ngrid,'runoff','Runoff amount',' ',2,runoff)
461      endif
462
463      return
464    end subroutine hydrol
Note: See TracBrowser for help on using the repository browser.