!============================================================================== subroutine photolysis_online(nlayer, alt, press, temp, $ mmean, rm, $ tau, sza, dist_sol, v_phot, e_phot, ig, ngrid, $ nreact) !*********************************************************************** ! ! subject: ! -------- ! ! photolysis online ! ! VERSION: Extracted from LMDZ.MARS work of Franck Lefevre (Yassin Jaziri) ! April 2019 - Yassin Jaziri add updates generic input (Yassin Jaziri) ! !*********************************************************************** use photolysis_mod use tracer_h use chimiedata_h, only: indexchim, fluxUV use types_asis, only: nb_phot_hv_max, nb_phot_max, jlabel, $ reactions implicit none ! input integer, intent(in) :: nlayer ! number of atmospheric layers integer, intent(in) :: ngrid ! number of atmospheric columns real, intent(in), dimension(nlayer) :: press, temp, mmean ! pressure (hpa)/temperature (k)/mean molecular mass (g.mol-1) real, intent(in), dimension(nlayer) :: alt ! altitude (km) real, intent(in), dimension(nlayer,nesp) :: rm ! mixing ratios real, intent(in) :: tau ! integrated aerosol optical depth at the surface real, intent(in) :: sza ! solar zenith angle (degrees) real, intent(in) :: dist_sol ! solar distance (au) integer, intent(in) :: ig ! grid point index integer, intent(in) :: nreact ! number of reactions in reactions files ! output real (kind = 8), dimension(nlayer,nb_phot_max) :: v_phot ! photolysis rates (s-1) real (kind = 8), dimension(nlayer,nb_phot_max) :: e_phot ! photolysis rates by energie (J.mol-1.s-1) ! stellar flux at planet real, dimension(nw) :: fplanet ! stellar flux (photon.s-1.nm-1.cm-2) at planet real :: factor ! atmosphere real, dimension(nw) :: albedo_chim ! Surface albedo calculated on chemistry wavelenght grid real, dimension(nlayer+1) :: colinc ! air column increment (molecule.cm-2) real, dimension(nlayer+1) :: airlev ! air density at each specified altitude (molec/cc) real, dimension(nlayer+1) :: edir, edn, eup ! normalised irradiances real, dimension(nlayer+1) :: fdir, fdn, fup ! normalised actinic fluxes real, dimension(nlayer+1) :: saflux ! total actinic flux real, dimension(nlayer+1,nw) :: dtrl ! rayleigh optical depth real, dimension(nlayer+1,nw) :: dtaer ! aerosol optical depth real, dimension(nlayer+1,nw) :: omaer ! aerosol single scattering albedo real, dimension(nlayer+1,nw) :: gaer ! aerosol asymmetry parameter real, dimension(nlayer+1,nw) :: dtcld ! cloud optical depth real, dimension(nlayer+1,nw) :: omcld ! cloud single scattering albedo real, dimension(nlayer+1,nw) :: gcld ! cloud asymmetry parameter real, dimension(nlayer+1,nw) :: dagas ! total gas optical depth real, dimension(nlayer+1,nw,nabs) :: dtgas ! optical depth for each gas real, dimension(nlayer+1) :: zpress ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1) real, dimension(nlayer+1) :: zalt ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1) real, dimension(nlayer+1) :: ztemp ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1) real, dimension(nlayer+1) :: zmmean ! pressure (hpa)/altitude (km)/temperature (k)/mean molecular mass (g.mol-1) integer, dimension(0:nlayer+1) :: nid real, dimension(0:nlayer+1,nlayer+1) :: dsdh integer :: i, ilay, iw, ialt, iphot, ispe, ij, igas, ireact integer :: ilev, nlev real :: deltaj !==== define working vertical grid for the uv radiative code nlev = nlayer + 1 do ilev = 1,nlev-1 zpress(ilev) = press(ilev) zalt(ilev) = alt(ilev) ztemp(ilev) = temp(ilev) zmmean(ilev) = mmean(ilev) end do zpress(nlev) = 0. ! top of the atmosphere zalt(nlev) = zalt(nlev-1) + (zalt(nlev-1) - zalt(nlev-2)) ztemp(nlev) = ztemp(nlev-1) zmmean(nlev) = zmmean(nlev-1) !==== air column increments and rayleigh optical depth call setair(nlev, nw, wl, wc, zpress, ztemp, zmmean, colinc, dtrl, $ airlev) !==== set surface albedo call setalb(nw,wl,ig,ngrid,albedo_chim) !==== set temperature-dependent cross-sections and optical depths iphot = 0 ij = 0 igas = 0 ireact = 0 dtgas(:,:,:) = 0. do while(iphot radius [km] dr = pi/180. zenrad = zen*dr ! number of layers: nlay = nlev - 1 ! include the elevation above sea level to the radius of Mars: re = radius + z(1) ! correspondingly z changed to the elevation above Mars surface: DO k = 1, nlev ze(k) = z(k) - z(1) END DO ! inverse coordinate of z zd(0) = ze(nlev) DO k = 1, nlay zd(k) = ze(nlev - k) END DO ! initialise dsdh(i,j), nid(i) nid(:) = 0. dsdh(:,:) = 0. ! calculate ds/dh of every layer do i = 0,nlay rpsinz = (re + zd(i))*sin(zenrad) IF ( (zen .GT. 90.0) .AND. (rpsinz .LT. re) ) THEN nid(i) = -1 ELSE ! Find index of layer in which the screening height lies id = i if (zen > 90.) then do j = 1,nlay IF( (rpsinz .LT. ( zd(j-1) + re ) ) .AND. $ (rpsinz .GE. ( zd(j) + re )) ) id = j end do end if do j = 1,id sm = 1.0 IF (j .EQ. id .AND. id .EQ. i .AND. zen .GT. 90.0) $ sm = -1.0 rj = re + zd(j-1) rjp1 = re + zd(j) dhj = zd(j-1) - zd(j) ga = rj*rj - rpsinz*rpsinz gb = rjp1*rjp1 - rpsinz*rpsinz ga = max(ga, 0.) gb = max(gb, 0.) IF (id.GT.i .AND. j.EQ.id) THEN dsj = sqrt(ga) ELSE dsj = sqrt(ga) - sm*sqrt(gb) END IF dsdh(i,j) = dsj/dhj end do nid(i) = id end if end do ! i = 0,nlay end subroutine sphers !============================================================================== SUBROUTINE rtlink(nlev, nw, iw, ag, zen, dsdh, nid, dtrl, $ dagas, dtcld, omcld, gcld, dtaer, omaer, gaer, $ edir, edn, eup, fdir, fdn, fup) implicit none ! input integer, intent(in) :: nlev, nw, iw ! number of wavelength grid points REAL ag REAL zen REAL dsdh(0:nlev,nlev) INTEGER nid(0:nlev) REAL dtrl(nlev,nw) REAL dagas(nlev,nw) REAL dtcld(nlev,nw), omcld(nlev,nw), gcld(nlev,nw) REAL dtaer(nlev,nw), omaer(nlev,nw), gaer(nlev,nw) ! output REAL edir(nlev), edn(nlev), eup(nlev) REAL fdir(nlev), fdn(nlev), fup(nlev) ! local: REAL dt(nlev), om(nlev), g(nlev) REAL dtabs,dtsct,dscld,dsaer,dacld,daaer INTEGER i, ii real, parameter :: largest = 1.e+36 ! specific two ps2str REAL ediri(nlev), edni(nlev), eupi(nlev) REAL fdiri(nlev), fdni(nlev), fupi(nlev) logical, save :: delta = .true. !_______________________________________________________________________ ! initialize: do i = 1, nlev fdir(i) = 0. fup(i) = 0. fdn(i) = 0. edir(i) = 0. eup(i) = 0. edn(i) = 0. end do do i = 1, nlev - 1 dscld = dtcld(i,iw)*omcld(i,iw) dacld = dtcld(i,iw)*(1.-omcld(i,iw)) dsaer = dtaer(i,iw)*omaer(i,iw) daaer = dtaer(i,iw)*(1.-omaer(i,iw)) dtsct = dtrl(i,iw) + dscld + dsaer dtabs = dagas(i,iw) + dacld + daaer dtabs = amax1(dtabs,1./largest) dtsct = amax1(dtsct,1./largest) ! invert z-coordinate: ii = nlev - i dt(ii) = dtsct + dtabs om(ii) = dtsct/(dtsct + dtabs) IF(dtsct .EQ. 1./largest) om(ii) = 1./largest g(ii) = (gcld(i,iw)*dscld + $ gaer(i,iw)*dsaer)/dtsct end do ! call rt routine: call ps2str(nlev, zen, ag, dt, om, g, $ dsdh, nid, delta, $ fdiri, fupi, fdni, ediri, eupi, edni) ! output (invert z-coordinate) do i = 1, nlev ii = nlev - i + 1 fdir(i) = fdiri(ii) fup(i) = fupi(ii) fdn(i) = fdni(ii) edir(i) = ediri(ii) eup(i) = eupi(ii) edn(i) = edni(ii) end do end subroutine rtlink *=============================================================================* subroutine ps2str(nlev,zen,rsfc,tauu,omu,gu, $ dsdh, nid, delta, $ fdr, fup, fdn, edr, eup, edn) !-----------------------------------------------------------------------------* != PURPOSE: =* != Solve two-stream equations for multiple layers. The subroutine is based =* != on equations from: Toon et al., J.Geophys.Res., v94 (D13), Nov 20, 1989.=* != It contains 9 two-stream methods to choose from. A pseudo-spherical =* != correction has also been added. =* !-----------------------------------------------------------------------------* != PARAMETERS: =* != NLEVEL - INTEGER, number of specified altitude levels in the working (I)=* != grid =* != ZEN - REAL, solar zenith angle (degrees) (I)=* != RSFC - REAL, surface albedo at current wavelength (I)=* != TAUU - REAL, unscaled optical depth of each layer (I)=* != OMU - REAL, unscaled single scattering albedo of each layer (I)=* != GU - REAL, unscaled asymmetry parameter of each layer (I)=* != DSDH - REAL, slant path of direct beam through each layer crossed (I)=* != when travelling from the top of the atmosphere to layer i; =* != DSDH(i,j), i = 0..NZ-1, j = 1..NZ-1 =* != NID - INTEGER, number of layers crossed by the direct beam when (I)=* != travelling from the top of the atmosphere to layer i; =* != NID(i), i = 0..NZ-1 =* != DELTA - LOGICAL, switch to use delta-scaling (I)=* != .TRUE. -> apply delta-scaling =* != .FALSE.-> do not apply delta-scaling =* != FDR - REAL, contribution of the direct component to the total (O)=* != actinic flux at each altitude level =* != FUP - REAL, contribution of the diffuse upwelling component to (O)=* != the total actinic flux at each altitude level =* != FDN - REAL, contribution of the diffuse downwelling component to (O)=* != the total actinic flux at each altitude level =* != EDR - REAL, contribution of the direct component to the total (O)=* != spectral irradiance at each altitude level =* != EUP - REAL, contribution of the diffuse upwelling component to (O)=* != the total spectral irradiance at each altitude level =* != EDN - REAL, contribution of the diffuse downwelling component to (O)=* *= the total spectral irradiance at each altitude level =* !-----------------------------------------------------------------------------* implicit none ! input: INTEGER nlev REAL zen, rsfc REAL tauu(nlev), omu(nlev), gu(nlev) REAL dsdh(0:nlev,nlev) INTEGER nid(0:nlev) LOGICAL delta ! output: REAL fup(nlev),fdn(nlev),fdr(nlev) REAL eup(nlev),edn(nlev),edr(nlev) ! local: REAL tausla(0:nlev), tauc(0:nlev) REAL mu2(0:nlev), mu, sum ! internal coefficients and matrix REAL lam(nlev),taun(nlev),bgam(nlev) REAL e1(nlev),e2(nlev),e3(nlev),e4(nlev) REAL cup(nlev),cdn(nlev),cuptn(nlev),cdntn(nlev) REAL mu1(nlev) INTEGER row REAL a(2*nlev),b(2*nlev),d(2*nlev),e(2*nlev),y(2*nlev) ! other: REAL pifs, fdn0 REAL gi(nlev), omi(nlev), tempg REAL f, g, om REAL gam1, gam2, gam3, gam4 real, parameter :: largest = 1.e+36 real, parameter :: precis = 1.e-7 ! For calculations of Associated Legendre Polynomials for GAMA1,2,3,4 ! in delta-function, modified quadrature, hemispheric constant, ! Hybrid modified Eddington-delta function metods, p633,Table1. ! W.E.Meador and W.R.Weaver, GAS,1980,v37,p.630 ! W.J.Wiscombe and G.W. Grams, GAS,1976,v33,p2440, ! uncomment the following two lines and the appropriate statements further ! down. ! REAL YLM0, YLM2, YLM4, YLM6, YLM8, YLM10, YLM12, YLMS, BETA0, ! > BETA1, BETAn, amu1, subd REAL expon, expon0, expon1, divisr, temp, up, dn REAL ssfc INTEGER nlayer, mrows, lev INTEGER i, j ! Some additional program constants: real pi, dr REAL eps PARAMETER (eps = 1.E-3) !_______________________________________________________________________ ! MU = cosine of solar zenith angle ! RSFC = surface albedo ! TAUU = unscaled optical depth of each layer ! OMU = unscaled single scattering albedo ! GU = unscaled asymmetry factor ! KLEV = max dimension of number of layers in atmosphere ! NLAYER = number of layers in the atmosphere ! NLEVEL = nlayer + 1 = number of levels ! initial conditions: pi*solar flux = 1; diffuse incidence = 0 pifs = 1. fdn0 = 0. nlayer = nlev - 1 pi = acos(-1.) dr = pi/180. mu = COS(zen*dr) !************* compute coefficients for each layer: ! GAM1 - GAM4 = 2-stream coefficients, different for different approximations ! EXPON0 = calculation of e when TAU is zero ! EXPON1 = calculation of e when TAU is TAUN ! CUP and CDN = calculation when TAU is zero ! CUPTN and CDNTN = calc. when TAU is TAUN ! DIVISR = prevents division by zero do j = 0, nlev tauc(j) = 0. tausla(j) = 0. mu2(j) = 1./SQRT(largest) end do IF (.NOT. delta) THEN DO i = 1, nlayer gi(i) = gu(i) omi(i) = omu(i) taun(i) = tauu(i) END DO ELSE ! delta-scaling. Have to be done for delta-Eddington approximation, ! delta discrete ordinate, Practical Improved Flux Method, delta function, ! and Hybrid modified Eddington-delta function methods approximations DO i = 1, nlayer f = gu(i)*gu(i) gi(i) = (gu(i) - f)/(1 - f) omi(i) = (1 - f)*omu(i)/(1 - omu(i)*f) taun(i) = (1 - omu(i)*f)*tauu(i) END DO END IF ! calculate slant optical depth at the top of the atmosphere when zen>90. ! in this case, higher altitude of the top layer is recommended. IF (zen .GT. 90.0) THEN IF (nid(0) .LT. 0) THEN tausla(0) = largest ELSE sum = 0.0 DO j = 1, nid(0) sum = sum + 2.*taun(j)*dsdh(0,j) END DO tausla(0) = sum END IF END IF DO 11, i = 1, nlayer g = gi(i) om = omi(i) tauc(i) = tauc(i-1) + taun(i) ! stay away from 1 by precision. For g, also stay away from -1 tempg = AMIN1(abs(g),1. - precis) g = SIGN(tempg,g) om = AMIN1(om,1.-precis) ! calculate slant optical depth IF (nid(i) .LT. 0) THEN tausla(i) = largest ELSE sum = 0.0 DO j = 1, MIN(nid(i),i) sum = sum + taun(j)*dsdh(i,j) END DO DO j = MIN(nid(i),i)+1,nid(i) sum = sum + 2.*taun(j)*dsdh(i,j) END DO tausla(i) = sum IF (tausla(i) .EQ. tausla(i-1)) THEN mu2(i) = SQRT(largest) ELSE mu2(i) = (tauc(i)-tauc(i-1))/(tausla(i)-tausla(i-1)) mu2(i) = SIGN( AMAX1(ABS(mu2(i)),1./SQRT(largest)), $ mu2(i) ) END IF END IF !** the following gamma equations are from pg 16,289, Table 1 !** save mu1 for each approx. for use in converting irradiance to actinic flux ! Eddington approximation(Joseph et al., 1976, JAS, 33, 2452): c gam1 = (7. - om*(4. + 3.*g))/4. c gam2 = -(1. - om*(4. - 3.*g))/4. c gam3 = (2. - 3.*g*mu)/4. c gam4 = 1. - gam3 c mu1(i) = 0.5 * quadrature (Liou, 1973, JAS, 30, 1303-1326; 1974, JAS, 31, 1473-1475): c gam1 = 1.7320508*(2. - om*(1. + g))/2. c gam2 = 1.7320508*om*(1. - g)/2. c gam3 = (1. - 1.7320508*g*mu)/2. c gam4 = 1. - gam3 c mu1(i) = 1./sqrt(3.) * hemispheric mean (Toon et al., 1089, JGR, 94, 16287): gam1 = 2. - om*(1. + g) gam2 = om*(1. - g) gam3 = (2. - g*mu)/4. gam4 = 1. - gam3 mu1(i) = 0.5 * PIFM (Zdunkovski et al.,1980, Conrib.Atmos.Phys., 53, 147-166): c GAM1 = 0.25*(8. - OM*(5. + 3.*G)) c GAM2 = 0.75*OM*(1.-G) c GAM3 = 0.25*(2.-3.*G*MU) c GAM4 = 1. - GAM3 c mu1(i) = 0.5 * delta discrete ordinates (Schaller, 1979, Contrib.Atmos.Phys, 52, 17-26): c GAM1 = 0.5*1.7320508*(2. - OM*(1. + G)) c GAM2 = 0.5*1.7320508*OM*(1.-G) c GAM3 = 0.5*(1.-1.7320508*G*MU) c GAM4 = 1. - GAM3 c mu1(i) = 1./sqrt(3.) * Calculations of Associated Legendre Polynomials for GAMA1,2,3,4 * in delta-function, modified quadrature, hemispheric constant, * Hybrid modified Eddington-delta function metods, p633,Table1. * W.E.Meador and W.R.Weaver, GAS,1980,v37,p.630 * W.J.Wiscombe and G.W. Grams, GAS,1976,v33,p2440 c YLM0 = 2. c YLM2 = -3.*G*MU c YLM4 = 0.875*G**3*MU*(5.*MU**2-3.) c YLM6=-0.171875*G**5*MU*(15.-70.*MU**2+63.*MU**4) c YLM8=+0.073242*G**7*MU*(-35.+315.*MU**2-693.*MU**4 c *+429.*MU**6) c YLM10=-0.008118*G**9*MU*(315.-4620.*MU**2+18018.*MU**4 c *-25740.*MU**6+12155.*MU**8) c YLM12=0.003685*G**11*MU*(-693.+15015.*MU**2-90090.*MU**4 c *+218790.*MU**6-230945.*MU**8+88179.*MU**10) c YLMS=YLM0+YLM2+YLM4+YLM6+YLM8+YLM10+YLM12 c YLMS=0.25*YLMS c BETA0 = YLMS c c amu1=1./1.7320508 c YLM0 = 2. c YLM2 = -3.*G*amu1 c YLM4 = 0.875*G**3*amu1*(5.*amu1**2-3.) c YLM6=-0.171875*G**5*amu1*(15.-70.*amu1**2+63.*amu1**4) c YLM8=+0.073242*G**7*amu1*(-35.+315.*amu1**2-693.*amu1**4 c *+429.*amu1**6) c YLM10=-0.008118*G**9*amu1*(315.-4620.*amu1**2+18018.*amu1**4 c *-25740.*amu1**6+12155.*amu1**8) c YLM12=0.003685*G**11*amu1*(-693.+15015.*amu1**2-90090.*amu1**4 c *+218790.*amu1**6-230945.*amu1**8+88179.*amu1**10) c YLMS=YLM0+YLM2+YLM4+YLM6+YLM8+YLM10+YLM12 c YLMS=0.25*YLMS c BETA1 = YLMS c c BETAn = 0.25*(2. - 1.5*G-0.21875*G**3-0.085938*G**5 c *-0.045776*G**7) * Hybrid modified Eddington-delta function(Meador and Weaver,1980,JAS,37,630): c subd=4.*(1.-G*G*(1.-MU)) c GAM1 = (7.-3.*G*G-OM*(4.+3.*G)+OM*G*G*(4.*BETA0+3.*G))/subd c GAM2 =-(1.-G*G-OM*(4.-3.*G)-OM*G*G*(4.*BETA0+3.*G-4.))/subd c GAM3 = BETA0 c GAM4 = 1. - GAM3 c mu1(i) = (1. - g*g*(1.- mu) )/(2. - g*g) ***** * delta function (Meador, and Weaver, 1980, JAS, 37, 630): c GAM1 = (1. - OM*(1. - beta0))/MU c GAM2 = OM*BETA0/MU c GAM3 = BETA0 c GAM4 = 1. - GAM3 c mu1(i) = mu ***** * modified quadrature (Meador, and Weaver, 1980, JAS, 37, 630): c GAM1 = 1.7320508*(1. - OM*(1. - beta1)) c GAM2 = 1.7320508*OM*beta1 c GAM3 = BETA0 c GAM4 = 1. - GAM3 c mu1(i) = 1./sqrt(3.) * hemispheric constant (Toon et al., 1989, JGR, 94, 16287): c GAM1 = 2.*(1. - OM*(1. - betan)) c GAM2 = 2.*OM*BETAn c GAM3 = BETA0 c GAM4 = 1. - GAM3 c mu1(i) = 0.5 ***** * lambda = pg 16,290 equation 21 * big gamma = pg 16,290 equation 22 * if gam2 = 0., then bgam = 0. lam(i) = sqrt(gam1*gam1 - gam2*gam2) IF (gam2 .NE. 0.) THEN bgam(i) = (gam1 - lam(i))/gam2 ELSE bgam(i) = 0. END IF expon = EXP(-lam(i)*taun(i)) * e1 - e4 = pg 16,292 equation 44 e1(i) = 1. + bgam(i)*expon e2(i) = 1. - bgam(i)*expon e3(i) = bgam(i) + expon e4(i) = bgam(i) - expon * the following sets up for the C equations 23, and 24 * found on page 16,290 * prevent division by zero (if LAMBDA=1/MU, shift 1/MU^2 by EPS = 1.E-3 * which is approx equiv to shifting MU by 0.5*EPS* (MU)**3 expon0 = EXP(-tausla(i-1)) expon1 = EXP(-tausla(i)) divisr = lam(i)*lam(i) - 1./(mu2(i)*mu2(i)) temp = AMAX1(eps,abs(divisr)) divisr = SIGN(temp,divisr) up = om*pifs*((gam1 - 1./mu2(i))*gam3 + gam4*gam2)/divisr dn = om*pifs*((gam1 + 1./mu2(i))*gam4 + gam2*gam3)/divisr * cup and cdn are when tau is equal to zero * cuptn and cdntn are when tau is equal to taun cup(i) = up*expon0 cdn(i) = dn*expon0 cuptn(i) = up*expon1 cdntn(i) = dn*expon1 11 CONTINUE ***************** set up matrix ****** * ssfc = pg 16,292 equation 37 where pi Fs is one (unity). ssfc = rsfc*mu*EXP(-tausla(nlayer))*pifs * MROWS = the number of rows in the matrix mrows = 2*nlayer * the following are from pg 16,292 equations 39 - 43. * set up first row of matrix: i = 1 a(1) = 0. b(1) = e1(i) d(1) = -e2(i) e(1) = fdn0 - cdn(i) row=1 * set up odd rows 3 thru (MROWS - 1): i = 0 DO 20, row = 3, mrows - 1, 2 i = i + 1 a(row) = e2(i)*e3(i) - e4(i)*e1(i) b(row) = e1(i)*e1(i + 1) - e3(i)*e3(i + 1) d(row) = e3(i)*e4(i + 1) - e1(i)*e2(i + 1) e(row) = e3(i)*(cup(i + 1) - cuptn(i)) + $ e1(i)*(cdntn(i) - cdn(i + 1)) 20 CONTINUE * set up even rows 2 thru (MROWS - 2): i = 0 DO 30, row = 2, mrows - 2, 2 i = i + 1 a(row) = e2(i + 1)*e1(i) - e3(i)*e4(i + 1) b(row) = e2(i)*e2(i + 1) - e4(i)*e4(i + 1) d(row) = e1(i + 1)*e4(i + 1) - e2(i + 1)*e3(i + 1) e(row) = (cup(i + 1) - cuptn(i))*e2(i + 1) - $ (cdn(i + 1) - cdntn(i))*e4(i + 1) 30 CONTINUE * set up last row of matrix at MROWS: row = mrows i = nlayer a(row) = e1(i) - rsfc*e3(i) b(row) = e2(i) - rsfc*e4(i) d(row) = 0. e(row) = ssfc - cuptn(i) + rsfc*cdntn(i) * solve tri-diagonal matrix: CALL tridiag(a, b, d, e, y, mrows) **** unfold solution of matrix, compute output fluxes: row = 1 lev = 1 j = 1 * the following equations are from pg 16,291 equations 31 & 32 fdr(lev) = EXP( -tausla(0) ) edr(lev) = mu * fdr(lev) edn(lev) = fdn0 eup(lev) = y(row)*e3(j) - y(row + 1)*e4(j) + cup(j) fdn(lev) = edn(lev)/mu1(lev) fup(lev) = eup(lev)/mu1(lev) DO 60, lev = 2, nlayer + 1 fdr(lev) = EXP(-tausla(lev-1)) edr(lev) = mu *fdr(lev) edn(lev) = y(row)*e3(j) + y(row + 1)*e4(j) + cdntn(j) eup(lev) = y(row)*e1(j) + y(row + 1)*e2(j) + cuptn(j) fdn(lev) = edn(lev)/mu1(j) fup(lev) = eup(lev)/mu1(j) row = row + 2 j = j + 1 60 CONTINUE end subroutine ps2str *=============================================================================* subroutine tridiag(a,b,c,r,u,n) !_______________________________________________________________________ ! solves tridiagonal system. From Numerical Recipies, p. 40 !_______________________________________________________________________ IMPLICIT NONE ! input: INTEGER n REAL a, b, c, r DIMENSION a(n),b(n),c(n),r(n) ! output: REAL u DIMENSION u(n) ! local: INTEGER j REAL bet, gam DIMENSION gam(n) !_______________________________________________________________________ IF (b(1) .EQ. 0.) STOP 1001 bet = b(1) u(1) = r(1)/bet DO 11, j = 2, n gam(j) = c(j - 1)/bet bet = b(j) - a(j)*gam(j) IF (bet .EQ. 0.) STOP 2002 u(j) = (r(j) - a(j)*u(j - 1))/bet 11 CONTINUE DO 12, j = n - 1, 1, -1 u(j) = u(j) - gam(j + 1)*u(j + 1) 12 CONTINUE !_______________________________________________________________________ end subroutine tridiag !============================================================================== subroutine setalb(nw,wl,ig,ngrid,albedo_chim) !-----------------------------------------------------------------------------* != PURPOSE: =* != Set the albedo of the surface. The albedo is assumed to be Lambertian, =* != i.e., the reflected light is isotropic, and idependt of the direction =* != of incidence of light. Albedo can be chosen to be wavelength dependent. =* !-----------------------------------------------------------------------------* != PARAMETERS: =* != NW - INTEGER, number of specified intervals + 1 in working (I)=* != wavelength grid =* != WL - REAL, vector of lower limits of wavelength intervals in (I)=* != working wavelength grid =* != ALBEDO - REAL, surface albedo at each specified wavelength (O)=* !-----------------------------------------------------------------------------* use chimiedata_h, only: albedo_snow_chim, albedo_co2_ice_chim ! use slab_ice_h, only: h_alb_ice, alb_ice_min, alb_ice_max use ocean_slab_mod, only: h_alb_ice,alb_ice_min, snow_min use tracer_h, only: igcm_h2o_ice, igcm_co2_ice use callkeys_mod, only: ok_slab_ocean, co2cond, alb_ocean, & hydrology use phys_state_var_mod, only: albedo_bareground, & rnat, qsurf, sea_ice, & pctsrf_sic, tsurf, capcal use watercommon_h, only: T_h2O_ice_liq, RLFTT, rhowater implicit none ! input: (wavelength working grid data) integer, intent(in) :: ngrid ! number of atmospheric columns INTEGER nw INTEGER ig ! grid point index REAL wl(nw) ! output: REAL albedo_chim(nw) ! local: INTEGER iw ! REAL alb real zfra, alb_ice, twater, hice real snowlayer parameter (snowlayer=33.0) ! 33 kg/m^2 of snow, equal to a layer of 3.3 cm ! 0.015: mean value from clancy et al., icarus, 49-63, 1999. ! alb = 0.015 ! alb = albedo_phys ! do iw = 1, nw - 1 ! albedo_chim(iw) = alb ! end do ! See hydrol.F90 for taking into account new tendencies if(hydrology)then if(nint(rnat(ig)).eq.0)then if(ok_slab_ocean) then 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. ! 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). ! Albedo final calculation : do iw=1,nw - 1 alb_ice=albedo_snow_chim(iw) & -(albedo_snow_chim(iw)-alb_ice_min) & *exp(-sea_ice(ig)/h_alb_ice) ! this replaces the formulation from BC2014 ! More details on the parameterization of sea ice albedo vs thickness is provided in the wiki : ! https://lmdz-forge.lmd.jussieu.fr/mediawiki/Planets/index.php/Slab_ocean_model ! 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 albedo_chim(iw) = pctsrf_sic(ig)* & (albedo_snow_chim(iw)*zfra & + alb_ice*(1.0-zfra)) & + (1.-pctsrf_sic(ig))*alb_ocean enddo else !ok_slab_ocean ! calculate oceanic ice height including the latent heat of ice formation ! hice is the height of oceanic ice with a maximum of maxicethick. hice = qsurf(ig,igcm_h2o_ice)/rhowater ! update hice to include recent snowfall twater = tsurf(ig) - hice*RLFTT*rhowater/capcal(ig) ! this is temperature water would have if we melted entire ocean ice layer if(twater .lt. T_h2O_ice_liq)then do iw=1,nw - 1 albedo_chim(iw) = albedo_snow_chim(iw) ! Albedo of ice has been replaced by albedo of snow here. MT2015. enddo else DO iw=1,nw - 1 albedo_chim(iw) = alb_ocean if(ngrid.eq.1) then albedo_chim(iw) = albedo_bareground(ig) endif ENDDO endif endif!(ok_slab_ocean) ! Continent ! --------- elseif (nint(rnat(ig)).eq.1) then ! re-calculate continental albedo if (qsurf(ig,igcm_h2o_ice).ge.snowlayer) then DO iw=1,nw - 1 albedo_chim(iw) = albedo_snow_chim(iw) ENDDO else DO iw=1,nw - 1 albedo_chim(iw) = albedo_bareground(ig) & + (albedo_snow_chim(iw) & - albedo_bareground(ig)) & *qsurf(ig,igcm_h2o_ice)/snowlayer ENDDO endif else print*,'Surface type not recognised in photolysis_online.F!' print*,'Exiting...' call abort endif else albedo_chim(:) = albedo_bareground(ig) endif ! end if hydrology ! Re-add the albedo effects of CO2 ice if necessary ! ------------------------------------------------- if(co2cond)then if (qsurf(ig,igcm_co2_ice).gt.1.) then ! Condition changed - Need now ~1 mm CO2 ice coverage. MT2015 DO iw=1,nw - 1 albedo_chim(iw) = albedo_co2_ice_chim(iw) ENDDO endif endif ! co2cond end subroutine setalb !============================================================================== subroutine setsj(nd,nlayer,nw,tlay,tdiml,xsl,xs_templ,sjl) implicit none ! input: (wavelength working grid data) integer :: nd ! number of photolysis rates integer :: nlayer ! number of vertical layers integer :: nw ! number of wavelength grid points integer :: tdiml ! number of different temperature cross section file real :: xsl(tdiml,nw) ! cross section (cm2) from reading files real :: xs_templ(tdiml) ! temperature of the cross section (cm2) real :: tlay(nlayer) ! temperature (K) ! output: real :: sjl(nlayer,nw) ! output cross section (cm2) ! local: INTEGER ilay,iw,tpos do iw = 1,nw-1 do ilay = 1,nlayer tpos = 1 do while(tlay(ilay)>xs_templ(tpos) .and. tpos