source: trunk/LMDZ.GENERIC/libf/aeronostd/photolysis_online.F @ 3529

Last change on this file since 3529 was 3309, checked in by yjaziri, 7 months ago

GENERIC PCM:

  • Cosmetic + clarifying some variables and comments in chemistry
  • add controle variable for actinique UV flux in photochemistry and output surface UV flux in diagspecUV.nc

YJ

File size: 44.4 KB
Line 
1!==============================================================================
2
3      subroutine photolysis_online(nlayer, alt, press, temp,
4     $           mmean, rm,
5     $           tau, sza, dist_sol, v_phot, e_phot, ig, ngrid,
6     $           nreact)
7
8!***********************************************************************
9!
10!   subject:
11!   --------
12!
13!   photolysis online
14!
15!   VERSION: Extracted from LMDZ.MARS work of Franck Lefevre (Yassin Jaziri)
16!   April 2019 - Yassin Jaziri add updates generic input (Yassin Jaziri)
17!
18!***********************************************************************
19
20      use photolysis_mod
21      use tracer_h
22      use chimiedata_h, only: indexchim, fluxUV
23      use types_asis,   only: nb_phot_hv_max, nb_phot_max, jlabel,
24     $                        reactions
25
26      implicit none
27
28!     input
29
30      integer, intent(in) :: nlayer ! number of atmospheric layers
31      integer, intent(in) :: ngrid  ! number of atmospheric columns
32
33      real,    intent(in), dimension(nlayer) :: press, temp, mmean  ! pressure (hpa)/temperature (k)/mean molecular mass (g.mol-1)
34      real,    intent(in), dimension(nlayer) :: alt                 ! altitude (km)
35      real,    intent(in), dimension(nlayer,nesp) :: rm             ! mixing ratios
36      real,    intent(in) :: tau                                    ! integrated aerosol optical depth at the surface
37      real,    intent(in) :: sza                                    ! solar zenith angle (degrees)
38      real,    intent(in) :: dist_sol                               ! solar distance (au)
39      integer, intent(in) :: ig                                     ! grid point index
40      integer, intent(in) :: nreact                                 ! number of reactions in reactions files
41
42!     output
43
44      real (kind = 8), dimension(nlayer,nb_phot_max) :: v_phot      ! photolysis rates (s-1)
45      real (kind = 8), dimension(nlayer,nb_phot_max) :: e_phot      ! photolysis rates by energie (J.mol-1.s-1)
46 
47!     stellar flux at planet
48
49      real, dimension(nw) :: fplanet                                ! stellar flux (photon.s-1.nm-1.cm-2) at planet
50      real :: factor
51
52!     atmosphere
53
54      real, dimension(nw)               :: albedo_chim              ! Surface albedo calculated on chemistry wavelenght grid
55      real, dimension(nlayer+1)         :: colinc                   ! air column increment (molecule.cm-2)
56      real, dimension(nlayer+1)         :: airlev                   ! air density at each specified altitude (molec/cc)
57      real, dimension(nlayer+1)         :: edir, edn, eup           ! normalised irradiances
58      real, dimension(nlayer+1)         :: fdir, fdn, fup           ! normalised actinic fluxes
59      real, dimension(nlayer+1)         :: saflux                   ! total actinic flux
60      real, dimension(nlayer+1,nw)      :: dtrl                     ! rayleigh optical depth
61      real, dimension(nlayer+1,nw)      :: dtaer                    ! aerosol optical depth
62      real, dimension(nlayer+1,nw)      :: omaer                    ! aerosol single scattering albedo
63      real, dimension(nlayer+1,nw)      :: gaer                     ! aerosol asymmetry parameter
64      real, dimension(nlayer+1,nw)      :: dtcld                    ! cloud optical depth
65      real, dimension(nlayer+1,nw)      :: omcld                    ! cloud single scattering albedo
66      real, dimension(nlayer+1,nw)      :: gcld                     ! cloud asymmetry parameter
67      real, dimension(nlayer+1,nw)      :: dagas                    ! total gas optical depth
68      real, dimension(nlayer+1,nw,nabs) :: dtgas                    ! optical depth for each gas
69      real, dimension(nlayer+1)         :: zpress                   ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1)
70      real, dimension(nlayer+1)         :: zalt                     ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1)
71      real, dimension(nlayer+1)         :: ztemp                    ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1)
72      real, dimension(nlayer+1)         :: zmmean                   ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1)
73
74      integer, dimension(0:nlayer+1)          :: nid
75      real,    dimension(0:nlayer+1,nlayer+1) :: dsdh
76
77      integer :: i, ilay, iw, ialt, iphot, ispe, ij, igas, ireact
78      integer :: ilev, nlev
79      real    :: deltaj
80
81!==== define working vertical grid for the uv radiative code
82
83      nlev = nlayer + 1
84
85      do ilev = 1,nlev-1
86         zpress(ilev) = press(ilev)
87         zalt(ilev)   = alt(ilev)
88         ztemp(ilev)  = temp(ilev)
89         zmmean(ilev) = mmean(ilev)
90      end do
91
92      zpress(nlev) = 0.                  ! top of the atmosphere
93      zalt(nlev)   = zalt(nlev-1) + (zalt(nlev-1) - zalt(nlev-2))
94      ztemp(nlev)  = ztemp(nlev-1)
95      zmmean(nlev) = zmmean(nlev-1)
96
97!==== air column increments and rayleigh optical depth
98
99      call setair(nlev, nw, wl, wc, zpress, ztemp, zmmean, colinc, dtrl,
100     $            airlev)
101
102!====  set surface albedo
103
104      call setalb(nw,wl,ig,ngrid,albedo_chim)
105
106!==== set temperature-dependent cross-sections and optical depths
107
108      iphot  = 0
109      ij     = 0
110      igas   = 0
111      ireact = 0
112      dtgas(:,:,:) = 0.
113
114      do while(iphot<nb_phot_hv_max)
115        ij     = ij     + 1
116        iphot  = iphot  + 1
117        ireact = ireact + 1
118
119        if (tdim(ij).eq.1) then
120          ! Avoid to calculate several times dtgas for a same specie
121          if (jlabelbis(iphot)) then
122            ispe = indexchim(trim(jlabel(iphot,2)))
123            do iw = 1,nw-1
124              do ilay = 1,nlayer
125                dtgas(ilay,iw,ij-igas) = colinc(ilay)*rm(ilay,ispe)
126     $                                *xs(1,iw,ij)
127              end do
128            end do
129          else
130            igas = igas + 1
131          end if
132          do while(reactions(ireact)%rtype/=0)
133            ireact = ireact + 1
134          end do
135          if (reactions(ireact)%three_prod) iphot = iphot + 1
136        else
137          call setsj(nb_phot_hv_max,nlayer,nw,temp,tdim(ij),
138     $               xs(:tdim(ij),:,ij),
139     $               xs_temp(:,ij),sj(:,:,iphot))
140 
141          ! Avoid to calculate several times dtgas for a same specie
142          if (jlabelbis(iphot)) then
143            ispe = indexchim(trim(jlabel(iphot,2)))
144            do iw = 1,nw-1
145              do ilay = 1,nlayer
146                dtgas(ilay,iw,ij-igas) = colinc(ilay)*rm(ilay,ispe)
147     $                                *sj(ilay,iw,iphot)
148              end do
149            end do
150          else
151            igas = igas + 1
152          end if
153 
154          do iw = 1,nw-1
155            do ilay = 1,nlayer
156              sj(ilay,iw,iphot) = sj(ilay,iw,iphot)*qy(iw,ij)
157            end do
158          end do
159
160          do while(reactions(ireact)%rtype/=0)
161            ireact = ireact + 1
162          end do
163          if (reactions(ireact)%three_prod) then
164            iphot = iphot + 1
165            do iw = 1,nw-1
166              do ilay = 1,nlayer
167                sj(ilay,iw,iphot) = sj(ilay,iw,iphot-1)
168              end do
169            end do
170          end if
171
172        endif ! end if tdim .eq. 1
173
174      end do ! end while
175
176!     total gas optical depth
177
178      dagas(:,:) = 0.
179
180      do i = 1,nabs
181         do iw = 1,nw-1
182            do ilay = 1,nlayer
183               dagas(ilay,iw) = dagas(ilay,iw) + dtgas(ilay,iw,i)
184            end do
185         end do
186      end do
187
188!==== set aerosol properties and optical depth
189
190      call setaer(nlev,zalt,tau,nw,dtaer,omaer,gaer)
191
192!==== set cloud properties and optical depth
193
194      call setcld(nlev,nw,dtcld,omcld,gcld)
195
196!==== slant path lengths in spherical geometry
197
198      call sphers(nlev,zalt,sza,dsdh,nid)
199
200!==== stellar flux at planet
201
202      factor = (1./dist_sol)**2.
203      do iw = 1,nw-1
204         fplanet(iw) = fstar1AU(iw)*factor
205      end do
206
207!==== initialise photolysis rates
208
209      v_phot(:,1:nb_phot_hv_max) = 0.
210      e_phot(:,1:nb_phot_hv_max) = 0.
211
212!==== start of wavelength lopp
213
214      do iw = 1,nw-1
215
216!     monochromatic radiative transfer. outputs are:
217!     normalized irradiances     edir(nlayer), edn(nlayer), eup(nlayer)
218!     normalized actinic fluxes  fdir(nlayer), fdn(nlayer), fup(nlayer)
219!     where
220!     dir = direct beam, dn = down-welling diffuse, up = up-welling diffuse
221
222         call rtlink(nlev, nw, iw, albedo_chim(iw),
223     $               sza, dsdh, nid, dtrl,
224     $               dagas, dtcld, omcld, gcld, dtaer, omaer, gaer,
225     $               edir, edn, eup, fdir, fdn, fup)
226
227
228!     spherical actinic flux
229
230         do ilay = 1,nlayer
231           saflux(ilay) = fplanet(iw)*
232     $                     (fdir(ilay) + fdn(ilay) + fup(ilay))
233           fluxUV(ig,iw,ilay) = saflux(ilay)
234         end do
235
236!     photolysis rate integration
237!     (0.12/(wc(iw)*1e-9)) E(wc) en J.mol-1
238
239         do i = 1,nb_phot_hv_max
240            do ilay = 1,nlayer
241               deltaj = saflux(ilay)*sj(ilay,iw,i)
242               v_phot(ilay,i) = v_phot(ilay,i) + deltaj*(wu(iw)-wl(iw))
243               if (wc(iw).le.photoheat_lmax) then
244                  e_phot(ilay,i) = e_phot(ilay,i)
245     $                    + deltaj*(wu(iw)-wl(iw))*(0.12/(wc(iw)*1e-9))
246               end if
247            end do
248         end do
249
250!     eliminate small values
251
252         where (v_phot(:,1:nb_phot_hv_max) < 1.e-30)
253            v_phot(:,1:nb_phot_hv_max) = 0.
254         end where
255         where (e_phot(:,1:nb_phot_hv_max) < 1.e-30)
256            e_phot(:,1:nb_phot_hv_max) = 0.
257         end where
258
259      end do ! iw
260
261      contains
262
263!==============================================================================
264
265      subroutine setair(nlev, nw, wl, wc, press, temp, zmmean,
266     $                  colinc, dtrl, airlev)
267
268*-----------------------------------------------------------------------------*
269*=  PURPOSE:                                                                 =*
270*=  computes air column increments and rayleigh optical depth                =*
271*-----------------------------------------------------------------------------*
272
273      use comcstfi_mod, only:  g, avocado
274
275      implicit none
276
277!     input:
278
279      integer :: nlev, nw
280
281      real, dimension(nw)   :: wl, wc              ! lower and central wavelength grid (nm)
282      real, dimension(nlev) :: press, temp, zmmean ! pressure (hpa), temperature (k), molecular mass (g.mol-1)
283
284!     output:
285
286      real, dimension(nlev)    :: colinc ! air column increments (molecule.cm-2)
287      real, dimension(nlev)    :: airlev ! air density at each specified altitude (molec/cc)
288      real, dimension(nlev,nw) :: dtrl   ! rayleigh optical depth
289
290!     local:
291
292      real :: dp, nu
293      real, dimension(nw) :: srayl
294      integer :: ilev, iw
295
296!     compute column increments
297
298      do ilev = 1, nlev-1
299         dp = (press(ilev) - press(ilev+1))*100.
300         colinc(ilev) = avocado*0.1*dp/(zmmean(ilev)*g)
301      end do
302      colinc(nlev) = 0.0
303
304      do iw = 1, nw - 1
305
306!        co2 rayleigh cross-section
307!        ityaksov et al., chem. phys. lett., 462, 31-34, 2008
308
309!         nu = 1./(wc(iw)*1.e-7)
310!         srayl(iw) = 1.78e-26*nu**(4. + 0.625)
311!         srayl(iw) = srayl(iw)*1.e-20 ! cm2
312
313!        calcul Ityaksov et al., 2008 pour N2
314
315         nu = 1./(wc(iw)*1.e-7)                     ! cm-1
316         srayl(iw) = 1.8e-26*nu**(4. + 0.534)
317         srayl(iw) = srayl(iw)*1.e-20
318
319         do ilev = 1, nlev
320            dtrl(ilev,iw) = colinc(ilev)*srayl(iw)   ! cm2
321         end do
322      end do
323
324!     compute density
325
326      do ilev = 1, nlev
327         airlev(ilev) = press(ilev)/(1.38e-19*temp(ilev))
328      end do
329
330      end subroutine setair
331
332!==============================================================================
333
334      subroutine setaer(nlev,zalt,tau,nw,dtaer,omaer,gaer)
335
336!-----------------------------------------------------------------------------*
337!=  PURPOSE:                                                                 =*
338!=  Set aerosol properties for each specified altitude layer.  Properties    =*
339!=  may be wavelength dependent.                                             =*
340!-----------------------------------------------------------------------------*
341
342      implicit none
343
344!     input
345
346      integer :: nlev                              ! number of vertical layers
347      integer :: nw                                ! number of wavelength grid points
348      real, dimension(nlev) :: zalt                ! altitude (km)
349      real :: tau                                  ! integrated aerosol optical depth at the surface
350
351!     output
352
353      real, dimension(nlev,nw) :: dtaer            ! aerosol optical depth
354      real, dimension(nlev,nw) :: omaer            ! aerosol single scattering albedo
355      real, dimension(nlev,nw) :: gaer             ! aerosol asymmetry parameter
356
357!     local
358
359      integer :: ilev, iw
360      real, dimension(nlev) :: aer                 ! dust extinction
361      real :: omega, g, scaleh, gamma
362      real :: dz, tautot, q0
363
364      omega  = 0.622 ! single scattering albedo : wolff et al.(2010) at 258 nm
365      g      = 0.88  ! asymmetry factor : mateshvili et al. (2007) at 210 nm
366      scaleh = 10.   ! scale height (km)
367      gamma  = 0.03  ! conrath parameter
368
369      dtaer(:,:) = 0.
370      omaer(:,:) = 0.
371      gaer(:,:)  = 0.
372
373!     optical depth profile:
374
375      tautot = 0.
376      do ilev = 1, nlev-1
377         dz = zalt(ilev+1) - zalt(ilev)
378         tautot = tautot + exp(gamma*(1. - exp(zalt(ilev)/scaleh)))*dz
379      end do
380     
381      q0 = tau/tautot
382      do ilev = 1, nlev-1
383         dz = zalt(ilev+1) - zalt(ilev)
384         dtaer(ilev,:) = q0*exp(gamma*(1. - exp(zalt(ilev)/scaleh)))*dz
385         omaer(ilev,:) = omega
386         gaer(ilev,:)  = g
387      end do
388
389      end subroutine setaer
390
391!==============================================================================
392
393      subroutine setcld(nlev,nw,dtcld,omcld,gcld)
394
395!-----------------------------------------------------------------------------*
396!=  PURPOSE:                                                                 =*
397!=  Set cloud properties for each specified altitude layer.  Properties      =*
398!=  may be wavelength dependent.                                             =*
399!-----------------------------------------------------------------------------*
400
401      implicit none
402
403!     input
404
405      integer :: nlev                              ! number of vertical layers
406      integer :: nw                                ! number of wavelength grid points
407
408!     output
409
410      real, dimension(nlev,nw) :: dtcld            ! cloud optical depth
411      real, dimension(nlev,nw) :: omcld            ! cloud single scattering albedo
412      real, dimension(nlev,nw) :: gcld             ! cloud asymmetry parameter
413
414!     local
415
416      integer :: ilev, iw
417
418!     dtcld : optical depth
419!     omcld : single scattering albedo
420!     gcld  : asymmetry factor
421
422      do ilev = 1, nlev - 1
423         do iw = 1, nw - 1
424            dtcld(ilev,iw) = 0.       ! no clouds for the moment
425            omcld(ilev,iw) = 0.99
426            gcld(ilev,iw)  = 0.85
427         end do
428      end do
429
430      end subroutine setcld
431
432!==============================================================================
433
434      subroutine sphers(nlev, z, zen, dsdh, nid)
435
436!-----------------------------------------------------------------------------*
437!=  PURPOSE:                                                                 =*
438!=  Calculate slant path over vertical depth ds/dh in spherical geometry.    =*
439!=  Calculation is based on:  A.Dahlback, and K.Stamnes, A new spheric model =*
440!=  for computing the radiation field available for photolysis and heating   =*
441!=  at twilight, Planet.Space Sci., v39, n5, pp. 671-683, 1991 (Appendix B)  =*
442!-----------------------------------------------------------------------------*
443!=  PARAMETERS:                                                              =*
444!=  NZ      - INTEGER, number of specified altitude levels in the working (I)=*
445!=            grid                                                           =*
446!=  Z       - REAL, specified altitude working grid (km)                  (I)=*
447!=  ZEN     - REAL, solar zenith angle (degrees)                          (I)=*
448!=  DSDH    - REAL, slant path of direct beam through each layer crossed  (O)=*
449!=            when travelling from the top of the atmosphere to layer i;     =*
450!=            DSDH(i,j), i = 0..NZ-1, j = 1..NZ-1                            =*
451!=  NID     - INTEGER, number of layers crossed by the direct beam when   (O)=*
452!=            travelling from the top of the atmosphere to layer i;          =*
453!=            NID(i), i = 0..NZ-1                                            =*
454!-----------------------------------------------------------------------------*
455
456      use comcstfi_mod, only: rad, pi
457
458      implicit none
459
460! input
461
462      integer, intent(in) :: nlev
463      real, dimension(nlev), intent(in) :: z
464      real, intent(in) :: zen
465
466! output
467
468      INTEGER nid(0:nlev)
469      REAL dsdh(0:nlev,nlev)
470
471! more program constants
472
473      REAL re, ze(nlev)
474      REAL  dr
475      real radius ! km
476
477! local
478
479      real :: zenrad, rpsinz, rj, rjp1, dsj, dhj, ga, gb, sm
480      integer :: i, j, k, id, nlay
481
482      REAL zd(0:nlev-1)
483
484!-----------------------------------------------------------------------------
485
486      radius = rad*1.e-3 ! rad [m] -> radius [km]
487      dr = pi/180.
488      zenrad = zen*dr
489
490! number of layers:
491
492      nlay = nlev - 1
493
494! include the elevation above sea level to the radius of Mars:
495
496      re = radius + z(1)
497
498! correspondingly z changed to the elevation above Mars surface:
499
500      DO k = 1, nlev
501         ze(k) = z(k) - z(1)
502      END DO
503
504! inverse coordinate of z
505
506      zd(0) = ze(nlev)
507      DO k = 1, nlay
508        zd(k) = ze(nlev - k)
509      END DO
510
511! initialise dsdh(i,j), nid(i)
512
513      nid(:) = 0.
514      dsdh(:,:) = 0.
515
516! calculate ds/dh of every layer
517
518      do i = 0,nlay
519         rpsinz = (re + zd(i))*sin(zenrad)
520 
521        IF ( (zen .GT. 90.0) .AND. (rpsinz .LT. re) ) THEN
522           nid(i) = -1
523        ELSE
524
525! Find index of layer in which the screening height lies
526
527           id = i
528           if (zen > 90.) then
529              do j = 1,nlay
530                 IF( (rpsinz .LT. ( zd(j-1) + re ) ) .AND.
531     $               (rpsinz .GE. ( zd(j) + re )) ) id = j
532              end do
533           end if
534 
535           do j = 1,id
536              sm = 1.0
537              IF (j .EQ. id .AND. id .EQ. i .AND. zen .GT. 90.0)
538     $            sm = -1.0
539 
540              rj = re + zd(j-1)
541              rjp1 = re + zd(j)
542 
543              dhj = zd(j-1) - zd(j)
544 
545              ga = rj*rj - rpsinz*rpsinz
546              gb = rjp1*rjp1 - rpsinz*rpsinz
547
548              ga = max(ga, 0.)
549              gb = max(gb, 0.)
550 
551              IF (id.GT.i .AND. j.EQ.id) THEN
552                 dsj = sqrt(ga)
553              ELSE
554                 dsj = sqrt(ga) - sm*sqrt(gb)
555              END IF
556              dsdh(i,j) = dsj/dhj
557            end do
558            nid(i) = id
559         end if
560      end do ! i = 0,nlay
561
562      end subroutine sphers
563
564!==============================================================================
565
566      SUBROUTINE rtlink(nlev, nw, iw, ag, zen, dsdh, nid, dtrl,
567     $                  dagas, dtcld, omcld, gcld, dtaer, omaer, gaer,
568     $                  edir, edn, eup, fdir, fdn, fup)
569
570      implicit none
571
572! input
573
574      integer, intent(in) :: nlev, nw, iw                       ! number of wavelength grid points
575      REAL ag
576      REAL zen
577      REAL dsdh(0:nlev,nlev)
578      INTEGER nid(0:nlev)
579
580      REAL dtrl(nlev,nw)
581      REAL dagas(nlev,nw)
582      REAL dtcld(nlev,nw), omcld(nlev,nw), gcld(nlev,nw)
583      REAL dtaer(nlev,nw), omaer(nlev,nw), gaer(nlev,nw)
584
585! output
586
587      REAL edir(nlev), edn(nlev), eup(nlev)
588      REAL fdir(nlev), fdn(nlev), fup(nlev)
589
590! local:
591
592      REAL dt(nlev), om(nlev), g(nlev)
593      REAL dtabs,dtsct,dscld,dsaer,dacld,daaer
594      INTEGER i, ii
595      real, parameter :: largest = 1.e+36
596
597! specific two ps2str
598
599      REAL ediri(nlev), edni(nlev), eupi(nlev)
600      REAL fdiri(nlev), fdni(nlev), fupi(nlev)
601
602      logical, save :: delta = .true.
603
604!_______________________________________________________________________
605
606! initialize:
607
608      do i = 1, nlev
609         fdir(i) = 0.
610         fup(i)  = 0.
611         fdn(i)  = 0.
612         edir(i) = 0.
613         eup(i)  = 0.
614         edn(i)  = 0.
615      end do
616
617      do i = 1, nlev - 1
618         dscld = dtcld(i,iw)*omcld(i,iw)
619         dacld = dtcld(i,iw)*(1.-omcld(i,iw))
620
621         dsaer = dtaer(i,iw)*omaer(i,iw)
622         daaer = dtaer(i,iw)*(1.-omaer(i,iw))
623
624         dtsct = dtrl(i,iw) + dscld + dsaer
625         dtabs = dagas(i,iw) + dacld + daaer
626
627         dtabs = amax1(dtabs,1./largest)
628         dtsct = amax1(dtsct,1./largest)
629
630!     invert z-coordinate:
631
632         ii = nlev - i
633         dt(ii) = dtsct + dtabs
634         om(ii) = dtsct/(dtsct + dtabs)
635         IF(dtsct .EQ. 1./largest) om(ii) = 1./largest
636         g(ii) = (gcld(i,iw)*dscld +
637     $            gaer(i,iw)*dsaer)/dtsct
638      end do
639
640!     call rt routine:
641
642      call ps2str(nlev, zen, ag, dt, om, g,
643     $            dsdh, nid, delta,
644     $            fdiri, fupi, fdni, ediri, eupi, edni)
645
646!     output (invert z-coordinate)
647
648      do i = 1, nlev
649         ii = nlev - i + 1
650         fdir(i) = fdiri(ii)
651         fup(i) = fupi(ii)
652         fdn(i) = fdni(ii)
653         edir(i) = ediri(ii)
654         eup(i) = eupi(ii)
655         edn(i) = edni(ii)
656      end do
657
658      end subroutine rtlink
659
660*=============================================================================*
661
662      subroutine ps2str(nlev,zen,rsfc,tauu,omu,gu,
663     $                  dsdh, nid, delta,
664     $                  fdr, fup, fdn, edr, eup, edn)
665
666!-----------------------------------------------------------------------------*
667!=  PURPOSE:                                                                 =*
668!=  Solve two-stream equations for multiple layers.  The subroutine is based =*
669!=  on equations from:  Toon et al., J.Geophys.Res., v94 (D13), Nov 20, 1989.=*
670!=  It contains 9 two-stream methods to choose from.  A pseudo-spherical     =*
671!=  correction has also been added.                                          =*
672!-----------------------------------------------------------------------------*
673!=  PARAMETERS:                                                              =*
674!=  NLEVEL  - INTEGER, number of specified altitude levels in the working (I)=*
675!=            grid                                                           =*
676!=  ZEN     - REAL, solar zenith angle (degrees)                          (I)=*
677!=  RSFC    - REAL, surface albedo at current wavelength                  (I)=*
678!=  TAUU    - REAL, unscaled optical depth of each layer                  (I)=*
679!=  OMU     - REAL, unscaled single scattering albedo of each layer       (I)=*
680!=  GU      - REAL, unscaled asymmetry parameter of each layer            (I)=*
681!=  DSDH    - REAL, slant path of direct beam through each layer crossed  (I)=*
682!=            when travelling from the top of the atmosphere to layer i;     =*
683!=            DSDH(i,j), i = 0..NZ-1, j = 1..NZ-1                            =*
684!=  NID     - INTEGER, number of layers crossed by the direct beam when   (I)=*
685!=            travelling from the top of the atmosphere to layer i;          =*
686!=            NID(i), i = 0..NZ-1                                            =*
687!=  DELTA   - LOGICAL, switch to use delta-scaling                        (I)=*
688!=            .TRUE. -> apply delta-scaling                                  =*
689!=            .FALSE.-> do not apply delta-scaling                           =*
690!=  FDR     - REAL, contribution of the direct component to the total     (O)=*
691!=            actinic flux at each altitude level                            =*
692!=  FUP     - REAL, contribution of the diffuse upwelling component to    (O)=*
693!=            the total actinic flux at each altitude level                  =*
694!=  FDN     - REAL, contribution of the diffuse downwelling component to  (O)=*
695!=            the total actinic flux at each altitude level                  =*
696!=  EDR     - REAL, contribution of the direct component to the total     (O)=*
697!=            spectral irradiance at each altitude level                     =*
698!=  EUP     - REAL, contribution of the diffuse upwelling component to    (O)=*
699!=            the total spectral irradiance at each altitude level           =*
700!=  EDN     - REAL, contribution of the diffuse downwelling component to  (O)=*
701*=            the total spectral irradiance at each altitude level           =*
702!-----------------------------------------------------------------------------*
703
704      implicit none
705
706! input:
707
708      INTEGER nlev
709      REAL zen, rsfc
710      REAL tauu(nlev), omu(nlev), gu(nlev)
711      REAL dsdh(0:nlev,nlev)
712      INTEGER nid(0:nlev)
713      LOGICAL delta
714
715! output:
716
717      REAL fup(nlev),fdn(nlev),fdr(nlev)
718      REAL eup(nlev),edn(nlev),edr(nlev)
719
720! local:
721
722      REAL tausla(0:nlev), tauc(0:nlev)
723      REAL mu2(0:nlev), mu, sum
724
725! internal coefficients and matrix
726
727      REAL lam(nlev),taun(nlev),bgam(nlev)
728      REAL e1(nlev),e2(nlev),e3(nlev),e4(nlev)
729      REAL cup(nlev),cdn(nlev),cuptn(nlev),cdntn(nlev)
730      REAL mu1(nlev)
731      INTEGER row
732      REAL a(2*nlev),b(2*nlev),d(2*nlev),e(2*nlev),y(2*nlev)
733
734! other:
735
736      REAL pifs, fdn0
737      REAL gi(nlev), omi(nlev), tempg
738      REAL f, g, om
739      REAL gam1, gam2, gam3, gam4
740      real, parameter :: largest = 1.e+36
741      real, parameter :: precis = 1.e-7
742
743! For calculations of Associated Legendre Polynomials for GAMA1,2,3,4
744! in delta-function, modified quadrature, hemispheric constant,
745! Hybrid modified Eddington-delta function metods, p633,Table1.
746! W.E.Meador and W.R.Weaver, GAS,1980,v37,p.630
747! W.J.Wiscombe and G.W. Grams, GAS,1976,v33,p2440,
748! uncomment the following two lines and the appropriate statements further
749! down.
750!     REAL YLM0, YLM2, YLM4, YLM6, YLM8, YLM10, YLM12, YLMS, BETA0,
751!    >     BETA1, BETAn, amu1, subd
752
753      REAL expon, expon0, expon1, divisr, temp, up, dn
754      REAL ssfc
755      INTEGER nlayer, mrows, lev
756
757      INTEGER i, j
758
759! Some additional program constants:
760
761      real pi, dr
762      REAL eps
763      PARAMETER (eps = 1.E-3)
764!_______________________________________________________________________
765
766! MU = cosine of solar zenith angle
767! RSFC = surface albedo
768! TAUU =  unscaled optical depth of each layer
769! OMU  =  unscaled single scattering albedo
770! GU   =  unscaled asymmetry factor
771! KLEV = max dimension of number of layers in atmosphere
772! NLAYER = number of layers in the atmosphere
773! NLEVEL = nlayer + 1 = number of levels
774
775! initial conditions:  pi*solar flux = 1;  diffuse incidence = 0
776
777      pifs = 1.     
778      fdn0 = 0.
779
780      nlayer = nlev - 1
781
782      pi = acos(-1.)
783      dr = pi/180.
784      mu = COS(zen*dr)
785
786!************* compute coefficients for each layer:
787! GAM1 - GAM4 = 2-stream coefficients, different for different approximations
788! EXPON0 = calculation of e when TAU is zero
789! EXPON1 = calculation of e when TAU is TAUN
790! CUP and CDN = calculation when TAU is zero
791! CUPTN and CDNTN = calc. when TAU is TAUN
792! DIVISR = prevents division by zero
793
794      do j = 0, nlev
795         tauc(j) = 0.
796         tausla(j) = 0.
797         mu2(j) = 1./SQRT(largest)
798      end do
799
800      IF (.NOT. delta) THEN
801         DO i = 1, nlayer
802            gi(i) = gu(i)
803            omi(i) = omu(i)
804            taun(i) = tauu(i)
805         END DO
806      ELSE
807
808! delta-scaling. Have to be done for delta-Eddington approximation,
809! delta discrete ordinate, Practical Improved Flux Method, delta function,
810! and Hybrid modified Eddington-delta function methods approximations
811
812         DO i = 1, nlayer
813            f = gu(i)*gu(i)
814            gi(i) = (gu(i) - f)/(1 - f)
815            omi(i) = (1 - f)*omu(i)/(1 - omu(i)*f)       
816            taun(i) = (1 - omu(i)*f)*tauu(i)
817         END DO
818      END IF
819
820! calculate slant optical depth at the top of the atmosphere when zen>90.
821! in this case, higher altitude of the top layer is recommended.
822
823      IF (zen .GT. 90.0) THEN
824         IF (nid(0) .LT. 0) THEN
825            tausla(0) = largest
826         ELSE
827            sum = 0.0
828            DO j = 1, nid(0)
829               sum = sum + 2.*taun(j)*dsdh(0,j)
830            END DO
831            tausla(0) = sum
832         END IF
833      END IF
834
835      DO 11, i = 1, nlayer
836          g = gi(i)
837          om = omi(i)
838          tauc(i) = tauc(i-1) + taun(i)
839
840! stay away from 1 by precision.  For g, also stay away from -1
841
842          tempg = AMIN1(abs(g),1. - precis)
843          g = SIGN(tempg,g)
844          om = AMIN1(om,1.-precis)
845
846! calculate slant optical depth
847             
848          IF (nid(i) .LT. 0) THEN
849             tausla(i) = largest
850          ELSE
851             sum = 0.0
852             DO j = 1, MIN(nid(i),i)
853                sum = sum + taun(j)*dsdh(i,j)
854             END DO
855             DO j = MIN(nid(i),i)+1,nid(i)
856                sum = sum + 2.*taun(j)*dsdh(i,j)
857             END DO
858             tausla(i) = sum
859             IF (tausla(i) .EQ. tausla(i-1)) THEN
860                mu2(i) = SQRT(largest)
861             ELSE
862                mu2(i) = (tauc(i)-tauc(i-1))/(tausla(i)-tausla(i-1))
863                mu2(i) = SIGN( AMAX1(ABS(mu2(i)),1./SQRT(largest)),
864     $                     mu2(i) )
865             END IF
866          END IF
867
868!** the following gamma equations are from pg 16,289, Table 1
869!** save mu1 for each approx. for use in converting irradiance to actinic flux
870
871! Eddington approximation(Joseph et al., 1976, JAS, 33, 2452):
872
873c       gam1 =  (7. - om*(4. + 3.*g))/4.
874c       gam2 = -(1. - om*(4. - 3.*g))/4.
875c       gam3 = (2. - 3.*g*mu)/4.
876c       gam4 = 1. - gam3
877c       mu1(i) = 0.5
878
879* quadrature (Liou, 1973, JAS, 30, 1303-1326; 1974, JAS, 31, 1473-1475):
880
881c          gam1 = 1.7320508*(2. - om*(1. + g))/2.
882c          gam2 = 1.7320508*om*(1. - g)/2.
883c          gam3 = (1. - 1.7320508*g*mu)/2.
884c          gam4 = 1. - gam3
885c          mu1(i) = 1./sqrt(3.)
886         
887* hemispheric mean (Toon et al., 1089, JGR, 94, 16287):
888
889           gam1 = 2. - om*(1. + g)
890           gam2 = om*(1. - g)
891           gam3 = (2. - g*mu)/4.
892           gam4 = 1. - gam3
893           mu1(i) = 0.5
894
895* PIFM  (Zdunkovski et al.,1980, Conrib.Atmos.Phys., 53, 147-166):
896c         GAM1 = 0.25*(8. - OM*(5. + 3.*G))
897c         GAM2 = 0.75*OM*(1.-G)
898c         GAM3 = 0.25*(2.-3.*G*MU)
899c         GAM4 = 1. - GAM3
900c         mu1(i) = 0.5
901
902* delta discrete ordinates  (Schaller, 1979, Contrib.Atmos.Phys, 52, 17-26):
903c         GAM1 = 0.5*1.7320508*(2. - OM*(1. + G))
904c         GAM2 = 0.5*1.7320508*OM*(1.-G)
905c         GAM3 = 0.5*(1.-1.7320508*G*MU)
906c         GAM4 = 1. - GAM3
907c         mu1(i) = 1./sqrt(3.)
908
909* Calculations of Associated Legendre Polynomials for GAMA1,2,3,4
910* in delta-function, modified quadrature, hemispheric constant,
911* Hybrid modified Eddington-delta function metods, p633,Table1.
912* W.E.Meador and W.R.Weaver, GAS,1980,v37,p.630
913* W.J.Wiscombe and G.W. Grams, GAS,1976,v33,p2440
914c      YLM0 = 2.
915c      YLM2 = -3.*G*MU
916c      YLM4 = 0.875*G**3*MU*(5.*MU**2-3.)
917c      YLM6=-0.171875*G**5*MU*(15.-70.*MU**2+63.*MU**4)
918c     YLM8=+0.073242*G**7*MU*(-35.+315.*MU**2-693.*MU**4
919c    *+429.*MU**6)
920c     YLM10=-0.008118*G**9*MU*(315.-4620.*MU**2+18018.*MU**4
921c    *-25740.*MU**6+12155.*MU**8)
922c     YLM12=0.003685*G**11*MU*(-693.+15015.*MU**2-90090.*MU**4
923c    *+218790.*MU**6-230945.*MU**8+88179.*MU**10)
924c      YLMS=YLM0+YLM2+YLM4+YLM6+YLM8+YLM10+YLM12
925c      YLMS=0.25*YLMS
926c      BETA0 = YLMS
927c
928c         amu1=1./1.7320508
929c      YLM0 = 2.
930c      YLM2 = -3.*G*amu1
931c      YLM4 = 0.875*G**3*amu1*(5.*amu1**2-3.)
932c      YLM6=-0.171875*G**5*amu1*(15.-70.*amu1**2+63.*amu1**4)
933c     YLM8=+0.073242*G**7*amu1*(-35.+315.*amu1**2-693.*amu1**4
934c    *+429.*amu1**6)
935c     YLM10=-0.008118*G**9*amu1*(315.-4620.*amu1**2+18018.*amu1**4
936c    *-25740.*amu1**6+12155.*amu1**8)
937c     YLM12=0.003685*G**11*amu1*(-693.+15015.*amu1**2-90090.*amu1**4
938c    *+218790.*amu1**6-230945.*amu1**8+88179.*amu1**10)
939c      YLMS=YLM0+YLM2+YLM4+YLM6+YLM8+YLM10+YLM12
940c      YLMS=0.25*YLMS
941c      BETA1 = YLMS
942c
943c         BETAn = 0.25*(2. - 1.5*G-0.21875*G**3-0.085938*G**5
944c    *-0.045776*G**7)
945
946
947* Hybrid modified Eddington-delta function(Meador and Weaver,1980,JAS,37,630):
948c         subd=4.*(1.-G*G*(1.-MU))
949c         GAM1 = (7.-3.*G*G-OM*(4.+3.*G)+OM*G*G*(4.*BETA0+3.*G))/subd
950c         GAM2 =-(1.-G*G-OM*(4.-3.*G)-OM*G*G*(4.*BETA0+3.*G-4.))/subd
951c         GAM3 = BETA0
952c         GAM4 = 1. - GAM3
953c         mu1(i) = (1. - g*g*(1.- mu) )/(2. - g*g)
954
955*****
956* delta function  (Meador, and Weaver, 1980, JAS, 37, 630):
957c         GAM1 = (1. - OM*(1. - beta0))/MU
958c         GAM2 = OM*BETA0/MU
959c         GAM3 = BETA0
960c         GAM4 = 1. - GAM3
961c         mu1(i) = mu
962*****
963* modified quadrature (Meador, and Weaver, 1980, JAS, 37, 630):
964c         GAM1 = 1.7320508*(1. - OM*(1. - beta1))
965c         GAM2 = 1.7320508*OM*beta1
966c         GAM3 = BETA0
967c         GAM4 = 1. - GAM3
968c         mu1(i) = 1./sqrt(3.)
969
970* hemispheric constant (Toon et al., 1989, JGR, 94, 16287):
971c         GAM1 = 2.*(1. - OM*(1. - betan))
972c         GAM2 = 2.*OM*BETAn
973c         GAM3 = BETA0
974c         GAM4 = 1. - GAM3
975c         mu1(i) = 0.5
976
977*****
978
979* lambda = pg 16,290 equation 21
980* big gamma = pg 16,290 equation 22
981* if gam2 = 0., then bgam = 0.
982
983         lam(i) = sqrt(gam1*gam1 - gam2*gam2)
984
985         IF (gam2 .NE. 0.) THEN
986            bgam(i) = (gam1 - lam(i))/gam2
987         ELSE
988            bgam(i) = 0.
989         END IF
990
991         expon = EXP(-lam(i)*taun(i))
992
993* e1 - e4 = pg 16,292 equation 44
994         
995         e1(i) = 1. + bgam(i)*expon
996         e2(i) = 1. - bgam(i)*expon
997         e3(i) = bgam(i) + expon
998         e4(i) = bgam(i) - expon
999
1000* the following sets up for the C equations 23, and 24
1001* found on page 16,290
1002* prevent division by zero (if LAMBDA=1/MU, shift 1/MU^2 by EPS = 1.E-3
1003* which is approx equiv to shifting MU by 0.5*EPS* (MU)**3
1004
1005         expon0 = EXP(-tausla(i-1))
1006         expon1 = EXP(-tausla(i))
1007         
1008         divisr = lam(i)*lam(i) - 1./(mu2(i)*mu2(i))
1009         temp = AMAX1(eps,abs(divisr))
1010         divisr = SIGN(temp,divisr)
1011
1012         up = om*pifs*((gam1 - 1./mu2(i))*gam3 + gam4*gam2)/divisr
1013         dn = om*pifs*((gam1 + 1./mu2(i))*gam4 + gam2*gam3)/divisr
1014         
1015* cup and cdn are when tau is equal to zero
1016* cuptn and cdntn are when tau is equal to taun
1017
1018         cup(i) = up*expon0
1019         cdn(i) = dn*expon0
1020         cuptn(i) = up*expon1
1021         cdntn(i) = dn*expon1
1022 
1023   11 CONTINUE
1024
1025***************** set up matrix ******
1026* ssfc = pg 16,292 equation 37  where pi Fs is one (unity).
1027
1028      ssfc = rsfc*mu*EXP(-tausla(nlayer))*pifs
1029
1030* MROWS = the number of rows in the matrix
1031
1032      mrows = 2*nlayer     
1033     
1034* the following are from pg 16,292  equations 39 - 43.
1035* set up first row of matrix:
1036
1037      i = 1
1038      a(1) = 0.
1039      b(1) = e1(i)
1040      d(1) = -e2(i)
1041      e(1) = fdn0 - cdn(i)
1042
1043      row=1
1044
1045* set up odd rows 3 thru (MROWS - 1):
1046
1047      i = 0
1048      DO 20, row = 3, mrows - 1, 2
1049         i = i + 1
1050         a(row) = e2(i)*e3(i) - e4(i)*e1(i)
1051         b(row) = e1(i)*e1(i + 1) - e3(i)*e3(i + 1)
1052         d(row) = e3(i)*e4(i + 1) - e1(i)*e2(i + 1)
1053         e(row) = e3(i)*(cup(i + 1) - cuptn(i)) +
1054     $        e1(i)*(cdntn(i) - cdn(i + 1))
1055   20 CONTINUE
1056
1057* set up even rows 2 thru (MROWS - 2):
1058
1059      i = 0
1060      DO 30, row = 2, mrows - 2, 2
1061         i = i + 1
1062         a(row) = e2(i + 1)*e1(i) - e3(i)*e4(i + 1)
1063         b(row) = e2(i)*e2(i + 1) - e4(i)*e4(i + 1)
1064         d(row) = e1(i + 1)*e4(i + 1) - e2(i + 1)*e3(i + 1)
1065         e(row) = (cup(i + 1) - cuptn(i))*e2(i + 1) -
1066     $        (cdn(i + 1) - cdntn(i))*e4(i + 1)
1067   30 CONTINUE
1068
1069* set up last row of matrix at MROWS:
1070
1071      row = mrows
1072      i = nlayer
1073     
1074      a(row) = e1(i) - rsfc*e3(i)
1075      b(row) = e2(i) - rsfc*e4(i)
1076      d(row) = 0.
1077      e(row) = ssfc - cuptn(i) + rsfc*cdntn(i)
1078
1079* solve tri-diagonal matrix:
1080
1081      CALL tridiag(a, b, d, e, y, mrows)
1082
1083**** unfold solution of matrix, compute output fluxes:
1084
1085      row = 1
1086      lev = 1
1087      j = 1
1088     
1089* the following equations are from pg 16,291  equations 31 & 32
1090
1091      fdr(lev) = EXP( -tausla(0) )
1092      edr(lev) = mu * fdr(lev)
1093      edn(lev) = fdn0
1094      eup(lev) =  y(row)*e3(j) - y(row + 1)*e4(j) + cup(j)
1095      fdn(lev) = edn(lev)/mu1(lev)
1096      fup(lev) = eup(lev)/mu1(lev)
1097
1098      DO 60, lev = 2, nlayer + 1
1099         fdr(lev) = EXP(-tausla(lev-1))
1100         edr(lev) =  mu *fdr(lev)
1101         edn(lev) =  y(row)*e3(j) + y(row + 1)*e4(j) + cdntn(j)
1102         eup(lev) =  y(row)*e1(j) + y(row + 1)*e2(j) + cuptn(j)
1103         fdn(lev) = edn(lev)/mu1(j)
1104         fup(lev) = eup(lev)/mu1(j)
1105
1106         row = row + 2
1107         j = j + 1
1108   60 CONTINUE
1109
1110      end subroutine ps2str
1111
1112*=============================================================================*
1113
1114      subroutine tridiag(a,b,c,r,u,n)
1115
1116!_______________________________________________________________________
1117! solves tridiagonal system.  From Numerical Recipies, p. 40
1118!_______________________________________________________________________
1119
1120      IMPLICIT NONE
1121
1122! input:
1123
1124      INTEGER n
1125      REAL a, b, c, r
1126      DIMENSION a(n),b(n),c(n),r(n)
1127
1128! output:
1129
1130      REAL u
1131      DIMENSION u(n)
1132
1133! local:
1134
1135      INTEGER j
1136
1137      REAL bet, gam
1138      DIMENSION gam(n)
1139!_______________________________________________________________________
1140
1141      IF (b(1) .EQ. 0.) STOP 1001
1142      bet   = b(1)
1143      u(1) = r(1)/bet
1144      DO 11, j = 2, n   
1145         gam(j) = c(j - 1)/bet
1146         bet = b(j) - a(j)*gam(j)
1147         IF (bet .EQ. 0.) STOP 2002
1148         u(j) = (r(j) - a(j)*u(j - 1))/bet
1149   11 CONTINUE
1150      DO 12, j = n - 1, 1, -1 
1151         u(j) = u(j) - gam(j + 1)*u(j + 1)
1152   12 CONTINUE
1153!_______________________________________________________________________
1154
1155      end subroutine tridiag
1156
1157!==============================================================================
1158 
1159      subroutine setalb(nw,wl,ig,ngrid,albedo_chim)
1160 
1161!-----------------------------------------------------------------------------*
1162!=  PURPOSE:                                                                 =*
1163!=  Set the albedo of the surface.  The albedo is assumed to be Lambertian,  =*
1164!=  i.e., the reflected light is isotropic, and idependt of the direction    =*
1165!=  of incidence of light.  Albedo can be chosen to be wavelength dependent. =*
1166!-----------------------------------------------------------------------------*
1167!=  PARAMETERS:                                                              =*
1168!=  NW      - INTEGER, number of specified intervals + 1 in working       (I)=*
1169!=            wavelength grid                                                =*
1170!=  WL     - REAL, vector of lower limits of wavelength intervals in      (I)=*
1171!=           working wavelength grid                                         =*
1172!=  ALBEDO  - REAL, surface albedo at each specified wavelength           (O)=*
1173!-----------------------------------------------------------------------------*
1174
1175      use chimiedata_h,  only: albedo_snow_chim, albedo_co2_ice_chim
1176!      use slab_ice_h,    only: h_alb_ice, alb_ice_min, alb_ice_max
1177      use ocean_slab_mod, only: h_alb_ice,alb_ice_min, snow_min
1178      use tracer_h,      only: igcm_h2o_ice, igcm_co2_ice
1179      use callkeys_mod,  only: ok_slab_ocean, co2cond, alb_ocean,
1180     &                         hydrology
1181      use phys_state_var_mod, only: albedo_bareground,
1182     &                              rnat, qsurf, sea_ice,
1183     &                              pctsrf_sic, tsurf, capcal
1184      use watercommon_h, only: T_h2O_ice_liq, RLFTT, rhowater
1185
1186      implicit none
1187
1188! input: (wavelength working grid data)
1189
1190      integer, intent(in) :: ngrid  ! number of atmospheric columns
1191      INTEGER nw
1192      INTEGER ig             ! grid point index
1193      REAL wl(nw)
1194
1195! output:
1196
1197      REAL albedo_chim(nw)
1198
1199! local:
1200
1201      INTEGER iw
1202!      REAL alb
1203      real zfra, alb_ice, twater, hice
1204      real snowlayer
1205      parameter (snowlayer=33.0)        ! 33 kg/m^2 of snow, equal to a layer of 3.3 cm
1206
1207!     0.015: mean value from clancy et al., icarus, 49-63, 1999.
1208
1209!      alb = 0.015
1210!      alb = albedo_phys
1211
1212!      do iw = 1, nw - 1
1213!         albedo_chim(iw) = alb
1214!      end do
1215
1216!     See hydrol.F90 for taking into account new tendencies
1217
1218      if(hydrology)then
1219
1220        if(nint(rnat(ig)).eq.0)then
1221
1222          if(ok_slab_ocean) then
1223       
1224            zfra = MAX(0.0,MIN(1.0,qsurf(ig,igcm_h2o_ice)/snow_min)) ! Critical snow height (in kg/m2) from ocean_slab_ice routine.
1225            ! 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).
1226
1227            ! Albedo final calculation :
1228            do iw=1,nw - 1
1229               alb_ice=albedo_snow_chim(iw)
1230     &         -(albedo_snow_chim(iw)-alb_ice_min)
1231     &         *exp(-sea_ice(ig)/h_alb_ice) ! this replaces the formulation from BC2014
1232               ! More details on the parameterization of sea ice albedo vs thickness is provided in the wiki :
1233               ! https://lmdz-forge.lmd.jussieu.fr/mediawiki/Planets/index.php/Slab_ocean_model
1234               ! sea_ice is the ice thickness (calculated in ocean_slab routine) in kg/m2 ; h_alb_ice is fixed to 275.1kg/m2 i.e. 30cm based on comparisons with Brandt et al. 2005
1235               albedo_chim(iw) = pctsrf_sic(ig)*
1236     &                          (albedo_snow_chim(iw)*zfra
1237     &                          + alb_ice*(1.0-zfra))
1238     &                          + (1.-pctsrf_sic(ig))*alb_ocean
1239            enddo 
1240
1241          else !ok_slab_ocean
1242           
1243           
1244!       calculate oceanic ice height including the latent heat of ice formation
1245!       hice is the height of oceanic ice with a maximum of maxicethick.
1246            hice    = qsurf(ig,igcm_h2o_ice)/rhowater ! update hice to include recent snowfall
1247            twater  = tsurf(ig) - hice*RLFTT*rhowater/capcal(ig)
1248            ! this is temperature water would have if we melted entire ocean ice layer
1249
1250            if(twater .lt. T_h2O_ice_liq)then
1251
1252              do iw=1,nw - 1
1253                albedo_chim(iw) = albedo_snow_chim(iw) ! Albedo of ice has been replaced by albedo of snow here. MT2015.
1254              enddo
1255
1256            else
1257
1258              DO iw=1,nw - 1
1259                albedo_chim(iw) = alb_ocean
1260                if(ngrid.eq.1) then
1261                  albedo_chim(iw) = albedo_bareground(ig)
1262                endif
1263              ENDDO
1264
1265            endif
1266
1267          endif!(ok_slab_ocean)
1268
1269
1270!       Continent
1271!       ---------
1272        elseif (nint(rnat(ig)).eq.1) then
1273
1274!       re-calculate continental albedo
1275          if (qsurf(ig,igcm_h2o_ice).ge.snowlayer) then
1276            DO iw=1,nw - 1
1277              albedo_chim(iw) = albedo_snow_chim(iw)
1278            ENDDO
1279          else
1280            DO iw=1,nw - 1
1281              albedo_chim(iw) = albedo_bareground(ig)
1282     &                        + (albedo_snow_chim(iw)
1283     &                           - albedo_bareground(ig))
1284     &                        *qsurf(ig,igcm_h2o_ice)/snowlayer
1285            ENDDO
1286          endif
1287
1288        else
1289
1290          print*,'Surface type not recognised in photolysis_online.F!'
1291          print*,'Exiting...'
1292          call abort
1293
1294        endif
1295      else
1296        albedo_chim(:) = albedo_bareground(ig)
1297      endif ! end if hydrology
1298
1299!     Re-add the albedo effects of CO2 ice if necessary
1300!     -------------------------------------------------
1301      if(co2cond)then
1302
1303          if (qsurf(ig,igcm_co2_ice).gt.1.) then ! Condition changed - Need now ~1 mm CO2 ice coverage. MT2015
1304            DO iw=1,nw - 1
1305              albedo_chim(iw) = albedo_co2_ice_chim(iw)
1306            ENDDO
1307          endif
1308           
1309      endif ! co2cond
1310
1311      end subroutine setalb
1312
1313!==============================================================================
1314 
1315      subroutine setsj(nd,nlayer,nw,tlay,tdiml,xsl,xs_templ,sjl)
1316
1317
1318      implicit none
1319
1320! input: (wavelength working grid data)
1321
1322      integer :: nd                                      ! number of photolysis rates
1323      integer :: nlayer                                  ! number of vertical layers
1324      integer :: nw                                      ! number of wavelength grid points
1325      integer :: tdiml                                   ! number of different temperature cross section file
1326      real    :: xsl(tdiml,nw)                           ! cross section (cm2) from reading files
1327      real    :: xs_templ(tdiml)                         ! temperature of the cross section (cm2)
1328      real    :: tlay(nlayer)                            ! temperature (K)
1329
1330! output:
1331
1332      real    :: sjl(nlayer,nw)                          ! output cross section (cm2)
1333
1334
1335! local:
1336
1337      INTEGER ilay,iw,tpos
1338
1339        do iw = 1,nw-1
1340          do ilay = 1,nlayer
1341            tpos = 1
1342            do while(tlay(ilay)>xs_templ(tpos) .and. tpos<tdiml)
1343              tpos = tpos + 1
1344            end do
1345            if (tpos.eq.1 .or. tpos.eq.tdiml) then
1346              sjl(ilay,iw) = xsl(tpos,iw)
1347            else
1348              sjl(ilay,iw) = ( (xsl(tpos,iw)-xsl(tpos-1,iw))*
1349     $                       (tlay(ilay)-xs_templ(tpos-1))
1350     $                       /(xs_templ(tpos)-xs_templ(tpos-1))
1351     $                       + xsl(tpos-1,iw) )
1352            end if
1353          end do
1354        end do
1355
1356      end subroutine setsj
1357
1358      end subroutine photolysis_online
Note: See TracBrowser for help on using the repository browser.