source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/hydrol.F90 @ 222

Last change on this file since 222 was 222, checked in by ymipsl, 10 years ago

Creating temporary dynamico/lmdz/saturn branche

YM

  • Property svn:executable set to *
File size: 11.1 KB
Line 
1subroutine hydrol(ngrid,nq,ptimestep,rnat,tsurf,          &
2     qsurf,dqsurf,dqs_hyd,pcapcal,                &
3     albedo0,albedo,mu0,pdtsurf,pdtsurf_hyd,hice, &
4     pctsrf_sic,sea_ice)
5
6  use ioipsl_getincom 
7  use watercommon_h, only: T_h2O_ice_liq, RLFTT, rhowater, mx_eau_sol
8  USE surfdat_h
9  use comdiurn_h
10  USE comgeomfi_h
11  USE tracer_h
12  use slab_ice_h
13
14  implicit none
15
16!==================================================================
17!     
18!     Purpose
19!     -------
20!     Calculate the surface hydrology and albedo changes.
21!     
22!     Authors
23!     -------
24!     Adapted from LMDTERRE by B. Charnay (2010). Further
25!     modifications by R. Wordsworth (2010).
26!     
27!     Called by
28!     ---------
29!     physiq.F
30!     
31!     Calls
32!     -----
33!     none
34!
35!     Notes
36!     -----
37!     rnat is terrain type: 0-ocean; 1-continent
38!     
39!==================================================================
40
41#include "dimensions.h"
42#include "dimphys.h"
43#include "comcstfi.h"
44#include "callkeys.h"
45
46      integer ngrid,nq
47
48!     Inputs
49!     ------
50      real albedoice
51      save albedoice
52
53      real snowlayer
54      parameter (snowlayer=33.0)        ! 33 kg/m^2 of snow, equal to a layer of 3.3 cm
55      real oceantime
56      parameter (oceantime=10*24*3600)
57
58      logical,save :: oceanbulkavg ! relax ocean temperatures to a GLOBAL mean value?
59      logical,save :: activerunoff ! enable simple runoff scheme?
60      logical,save :: oceanalbvary ! ocean albedo varies with the diurnal cycle?
61
62!     Arguments
63!     ---------
64      real rnat(ngrid) ! I changed this to integer (RW)
65      real,dimension(:),allocatable,save :: runoff
66      real totalrunoff, tsea, oceanarea
67      save oceanarea
68
69      real ptimestep
70      real mu0(ngrid)
71      real qsurf(ngrid,nq), tsurf(ngrid)
72      real dqsurf(ngrid,nq), pdtsurf(ngrid)
73      real hice(ngrid)
74      real albedo0(ngrid), albedo(ngrid)
75      real pctsrf_sic(ngrid), sea_ice(ngrid)
76
77      real oceanarea2
78
79!     Output
80!     ------
81      real dqs_hyd(ngrid,nq)
82      real pdtsurf_hyd(ngrid)
83
84!     Local
85!     -----
86      real a,b,E
87      integer ig,iq, icap ! wld like to remove icap
88      real fsnoi, subli, fauxo
89      real twater(ngrid)
90      real pcapcal(ngrid)
91      real hicebis(ngrid)
92      real zqsurf(ngrid,nq)
93      real ztsurf(ngrid)
94      real albedo_sic, alb_ice
95      real zfra
96
97      integer, save :: ivap, iliq, iice
98
99      logical, save :: firstcall
100
101      data firstcall /.true./
102
103
104      if(firstcall)then
105
106         oceanbulkavg=.false.
107         oceanalbvary=.false.
108         write(*,*)"Activate runnoff into oceans?"
109         activerunoff=.false.
110         call getin("activerunoff",activerunoff)
111         write(*,*)" activerunoff = ",activerunoff
112         
113         
114         
115         if (activerunoff) ALLOCATE(runoff(ngrid))
116
117         ivap=igcm_h2o_vap
118         iliq=igcm_h2o_vap
119         iice=igcm_h2o_ice
120       
121         write(*,*) "hydrol: ivap=",ivap
122         write(*,*) "        iliq=",iliq
123         write(*,*) "        iice=",iice
124
125!     Here's the deal: iice is used in place of igcm_h2o_ice both on the
126!                      surface and in the atmosphere. ivap is used in
127!                      place of igcm_h2o_vap ONLY in the atmosphere, while
128!                      iliq is used in place of igcm_h2o_vap ONLY on the
129!                      surface.
130
131!                      Soon to be extended to the entire water cycle...
132
133!     Ice albedo = snow albedo for now
134         albedoice=albedosnow
135
136!     Total ocean surface area
137         oceanarea=0.
138         do ig=1,ngrid
139            if(nint(rnat(ig)).eq.0)then
140               oceanarea=oceanarea+area(ig)
141            endif
142         enddo
143
144         if(oceanbulkavg.and.(oceanarea.le.0.))then
145            print*,'How are we supposed to average the ocean'
146            print*,'temperature, when there are no oceans?'
147            call abort
148         endif
149
150         if(activerunoff.and.(oceanarea.le.0.))then
151            print*,'You have enabled runoff, but you have no oceans.'
152            print*,'Where did you think the water was going to go?'
153            call abort
154         endif
155         
156         firstcall = .false.
157      endif
158
159!     add physical tendencies already calculated
160!     ------------------------------------------
161
162      do ig=1,ngrid
163         ztsurf(ig) = tsurf(ig) + ptimestep*pdtsurf(ig)
164         pdtsurf_hyd(ig)=0.0
165         do iq=1,nq
166            zqsurf(ig,iq) = qsurf(ig,iq) + ptimestep*dqsurf(ig,iq)
167         enddo   
168      enddo
169 
170      do ig=1,ngrid
171         do iq=1,nq
172            dqs_hyd(ig,iq) = 0.0
173         enddo
174      enddo
175
176      do ig = 1, ngrid
177
178!     Ocean
179!     -----
180         if(nint(rnat(ig)).eq.0)then
181
182!     re-calculate oceanic albedo
183!            if(diurnal.and.oceanalbvary)then
184!               fauxo      = ( 1.47 - ACOS( mu0(ig) ) )/0.15 ! where does this come from (Benjamin)?
185!               albedo(ig) = 1.1*( .03 + .630/( 1. + fauxo*fauxo))
186!               albedo(ig) = MAX(MIN(albedo(ig),0.60),0.04)
187!            else
188               albedo(ig) = alb_ocean ! modif Benjamin
189!            end if
190
191
192          if(ok_slab_ocean) then
193
194! Fraction neige (hauteur critique 45kg/m2~15cm)
195            zfra = MAX(0.0,MIN(1.0,zqsurf(ig,iice)/45.0))
196! Albedo glace
197            alb_ice=alb_ice_max-(alb_ice_max-alb_ice_min) &
198            *exp(-sea_ice(ig)/h_alb_ice)
199! Albedo glace+neige
200            albedo_sic= albedosnow*zfra + alb_ice*(1.0-zfra)
201
202! Albedo final
203            albedo(ig) = pctsrf_sic(ig)*albedo_sic + (1.-pctsrf_sic(ig))*alb_ocean
204!     oceanic ice height, just for diagnostics
205            hice(ig)    = MIN(10.,sea_ice(ig)/rhowater)
206          else !ok_slab_ocean
207
208
209!     calculate oceanic ice height including the latent heat of ice formation
210!     hice is the height of oceanic ice with a maximum of maxicethick.
211            hice(ig)    = zqsurf(ig,iice)/rhowater ! update hice to include recent snowfall
212
213!            twater(ig)  = tsurf(ig) + ptimestep*zdtsurf(ig) &
214            twater(ig)  = ztsurf(ig) - hice(ig)*RLFTT*rhowater/pcapcal(ig)
215            ! this is temperature water would have if we melted entire ocean ice layer
216            hicebis(ig) = hice(ig)
217            hice(ig)    = 0.
218
219            if(twater(ig) .lt. T_h2O_ice_liq)then
220               E=min((T_h2O_ice_liq+Tsaldiff-twater(ig))*pcapcal(ig),RLFTT*rhowater*maxicethick)
221               hice(ig)        = E/(RLFTT*rhowater)
222               hice(ig)        = max(hice(ig),0.0)
223               hice(ig)        = min(hice(ig),maxicethick)
224               pdtsurf_hyd(ig) = (hice(ig) - hicebis(ig))*RLFTT*rhowater/pcapcal(ig)/ptimestep             
225               albedo(ig)      = albedoice
226
227!               if (zqsurf(ig,iice).ge.snowlayer) then
228!                  albedo(ig) = albedoice
229!               else
230!                  albedo(ig) = albedoocean &
231!                       + (albedosnow - albedoocean)*zqsurf(ig,iice)/snowlayer
232!               endif
233
234            else
235
236               pdtsurf_hyd(ig) = -hicebis(ig)*RLFTT*rhowater/pcapcal(ig)/ptimestep               
237               albedo(ig)      = alb_ocean
238
239            endif
240
241            zqsurf(ig,iliq) = zqsurf(ig,iliq)-(hice(ig)*rhowater-zqsurf(ig,iice))
242            zqsurf(ig,iice) = hice(ig)*rhowater
243
244          endif!(ok_slab_ocean)
245
246
247!     Continent
248!     ---------
249         elseif (nint(rnat(ig)).eq.1) then
250
251!     melt the snow
252            if(ztsurf(ig).gt.T_h2O_ice_liq)then
253               if(zqsurf(ig,iice).gt.1.0e-8)then
254
255                  a     = (ztsurf(ig)-T_h2O_ice_liq)*pcapcal(ig)/RLFTT
256                  b     = zqsurf(ig,iice)
257                  fsnoi = min(a,b)
258
259                  zqsurf(ig,iice) = zqsurf(ig,iice) - fsnoi
260                  zqsurf(ig,iliq) = zqsurf(ig,iliq) + fsnoi
261
262!                 thermal effects
263                  pdtsurf_hyd(ig) = -fsnoi*RLFTT/pcapcal(ig)/ptimestep 
264
265               endif
266            else
267
268!     freeze the water
269               if(zqsurf(ig,iliq).gt.1.0e-8)then
270
271                  a     = -(ztsurf(ig)-T_h2O_ice_liq)*pcapcal(ig)/RLFTT
272                  b     = zqsurf(ig,iliq)
273                 
274                  fsnoi = min(a,b)
275
276                  zqsurf(ig,iice) = zqsurf(ig,iice) + fsnoi
277                  zqsurf(ig,iliq) = zqsurf(ig,iliq) - fsnoi
278
279!                 thermal effects
280                  pdtsurf_hyd(ig) = +fsnoi*RLFTT/pcapcal(ig)/ptimestep 
281
282               endif
283            endif
284           
285!     deal with runoff
286            if(activerunoff)then
287
288               runoff(ig) = max(zqsurf(ig,iliq) - mx_eau_sol, 0.0)
289               if(ngrid.gt.1)then ! runoff only exists in 3D
290                  if(runoff(ig).ne.0.0)then
291                     zqsurf(ig,iliq) = mx_eau_sol
292!                    runoff is added to ocean at end
293                  endif
294               end if
295
296            endif
297
298!     re-calculate continental albedo
299            albedo(ig) = albedo0(ig) ! albedo0 = base values
300            if (zqsurf(ig,iice).ge.snowlayer) then
301               albedo(ig) = albedosnow
302            else
303               albedo(ig) = albedo0(ig) &
304                    + (albedosnow - albedo0(ig))*zqsurf(ig,iice)/snowlayer
305            endif
306
307         else
308
309            print*,'Surface type not recognised in hydrol.F!'
310            print*,'Exiting...'
311            call abort
312
313         endif
314
315      end do ! ig=1,ngrid
316
317!     perform crude bulk averaging of temperature in ocean
318!     ----------------------------------------------------
319      if(oceanbulkavg)then
320
321         oceanarea2=0.       
322         DO ig=1,ngrid
323            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0.))then
324               oceanarea2=oceanarea2+area(ig)*pcapcal(ig)
325            end if
326         END DO
327       
328         tsea=0.
329         DO ig=1,ngrid
330            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0.))then       
331               tsea=tsea+ztsurf(ig)*area(ig)*pcapcal(ig)/oceanarea2
332            end if
333         END DO
334
335         DO ig=1,ngrid
336            if((nint(rnat(ig)).eq.0).and.(hice(ig).eq.0))then
337               pdtsurf_hyd(ig) = pdtsurf_hyd(ig) + (tsea-ztsurf(ig))/oceantime
338            end if
339         END DO
340
341         print*,'Mean ocean temperature               = ',tsea,' K'
342
343      endif
344
345!     shove all the runoff water into the ocean
346!     -----------------------------------------
347      if(activerunoff)then
348
349         totalrunoff=0.
350         do ig=1,ngrid
351            if (nint(rnat(ig)).eq.1) then
352               totalrunoff = totalrunoff + area(ig)*runoff(ig)
353            endif
354         enddo
355         
356         do ig=1,ngrid
357            if (nint(rnat(ig)).eq.0) then
358               zqsurf(ig,iliq) = zqsurf(ig,iliq) + &
359                    totalrunoff/oceanarea
360            endif
361         enddo
362
363      endif         
364
365
366!     Re-add the albedo effects of CO2 ice if necessary
367!     -------------------------------------------------
368      if(co2cond)then
369
370         icap=1
371         do ig=1,ngrid
372            if (qsurf(ig,igcm_co2_ice).gt.0) then
373               albedo(ig) = albedice(icap)
374            endif
375         enddo
376         
377      endif
378
379
380      do ig=1,ngrid
381         dqs_hyd(ig,iliq)=(zqsurf(ig,iliq) - qsurf(ig,iliq))/ptimestep
382         dqs_hyd(ig,iice)=(zqsurf(ig,iice) - qsurf(ig,iice))/ptimestep
383      enddo
384
385      if (activerunoff) then
386        call writediagfi(ngrid,'runoff','Runoff amount',' ',2,runoff)
387      endif
388
389      return
390    end subroutine hydrol
Note: See TracBrowser for help on using the repository browser.