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

Last change on this file since 255 was 253, checked in by emillour, 14 years ago

Generic GCM

  • Massive update to version 0.7

EM+RW

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