Changeset 3309 for trunk/LMDZ.GENERIC
- Timestamp:
- Apr 22, 2024, 3:03:22 AM (7 months ago)
- Location:
- trunk/LMDZ.GENERIC
- Files:
-
- 2 added
- 8 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/changelog.txt
r3299 r3309 1906 1906 need be provided in a "volcano.def" ASCII file read at runtime. 1907 1907 (see routine "read_volcano" in volcano.F90) 1908 1909 == 20/11/2023 == YJ 1910 - Cosmetic + clarifying some variables and comments in chemistry 1911 - add controle variable for actinique UV flux in photochemistry 1912 and output surface UV flux in diagspecUV.nc -
trunk/LMDZ.GENERIC/libf/aeronostd/calchim_asis.F90
r3155 r3309 5 5 pu,pdu,pv,pdv,surfdust,surfice,icount,zdtchim) 6 6 7 use tracer_h, only: mmol, noms, nesp, is_chim8 9 use conc_mod, only: mmean ! mean molecular mass of the atmosphere7 use tracer_h, only: mmol, noms, nesp, is_chim 8 9 use conc_mod, only: mmean ! mean molecular mass of the atmosphere 10 10 USE comcstfi_mod 11 11 use callkeys_mod … … 14 14 use chimiedata_h 15 15 use radcommon_h, only: glat 16 use wstats_mod, only: wstats16 use wstats_mod, only: wstats 17 17 18 18 implicit none … … 256 256 ! chemistry in lower atmosphere 257 257 258 ! if (photochem) then259 260 258 call photochemistry_asis(nlayer,ngrid, & 261 259 ig,lswitch,zycol,szacol,fractcol,ptimestep,& … … 338 336 ! value of parameter 'output' to trigger writting of outputs 339 337 ! is set above at the declaration of the variable. 340 if ( photochem .and.output) then338 if (output) then 341 339 342 340 ! photoloysis … … 417 415 endif 418 416 endif ! end depos 417 418 if (jonline .and. output_writediagspecUV) then 419 420 !flux at the surface 421 call writediagspecUV(ngrid,"flux","flux(lon,lat,band)","photon.s-1.nm-1.cm-2",3,fluxUV(:,:,1)) 422 423 endif ! end jonline 419 424 420 425 end if ! of if (output) -
trunk/LMDZ.GENERIC/libf/aeronostd/chimiedata_h.F90
r2542 r3309 56 56 real, save, allocatable :: albedo_snow_chim(:) ! albedo snow on chemistry wavelenght grid 57 57 real, save, allocatable :: albedo_co2_ice_chim(:) ! albedo co2 ice on chemistry wavelenght grid 58 real, save, allocatable :: fluxUV(:,:,:) ! total actinic flux at each ngrid, wavelength and altitude (photon.s-1.nm-1.cm-2) 58 59 59 60 !-------------------------------------------- -
trunk/LMDZ.GENERIC/libf/aeronostd/photochemistry_asis.F90
r2553 r3309 137 137 allocate(albedo_snow_chim(nw)) 138 138 allocate(albedo_co2_ice_chim(nw)) 139 allocate(fluxUV(ngrid,nw,nlayer)) 140 fluxUV(:,:,:) = 0. 139 141 ! Step 1 : Initialisation of the Spectral Albedos. 140 142 DO iw=1,nw -
trunk/LMDZ.GENERIC/libf/aeronostd/photolysis_mod.F90
r2553 r3309 32 32 ! solar flux 33 33 34 real, allocatable, save :: f (:) ! solar flux (photon.s-1.nm-1.cm-2) at 1 au34 real, allocatable, save :: fstar1AU(:) ! stellar flux (photon.s-1.nm-1.cm-2) at 1 au 35 35 36 36 ! cross-sections and quantum yields … … 69 69 call gridw(nw,wl,wc,wu,mopt) 70 70 71 allocate(f (nw))71 allocate(fstar1AU(nw)) 72 72 73 73 ! read and grid solar flux data 74 74 75 call rdsolarflux(nw,wl,wc,f )75 call rdsolarflux(nw,wl,wc,fstar1AU) 76 76 77 77 ! read maximum wavelength until photochemical heating is calculated … … 474 474 !============================================================================== 475 475 476 subroutine rdsolarflux(nw,wl,wc,f )476 subroutine rdsolarflux(nw,wl,wc,fstar1AU) 477 477 478 478 ! Read and re-grid solar flux data. … … 490 490 ! output 491 491 492 real, dimension(nw) :: f ! solar flux (photon.s-1.nm-1.cm-2)492 real, dimension(nw) :: fstar1AU ! stellar flux (photon.s-1.nm-1.cm-2) 493 493 494 494 ! local … … 568 568 ! 5.039e11 = 1.e-4*1e-9/(hc = 6.62e-34*2.998e8) 569 569 570 ! f need to be in photon.s-1.nm-1.cm-2570 ! fstar1AU need to be in photon.s-1.nm-1.cm-2 571 571 DO iw = 1, nw-1 572 572 !f(iw) = yg1(iw)*wc(iw)*5.039e11 ! If yg1 in w.m-2.nm-1 573 f (iw) = yg1(iw)573 fstar1AU(iw) = yg1(iw) 574 574 ENDDO 575 575 -
trunk/LMDZ.GENERIC/libf/aeronostd/photolysis_online.F
r3291 r3309 20 20 use photolysis_mod 21 21 use tracer_h 22 use chimiedata_h, only: indexchim 22 use chimiedata_h, only: indexchim, fluxUV 23 23 use types_asis, only: nb_phot_hv_max, nb_phot_max, jlabel, 24 24 $ reactions … … 45 45 real (kind = 8), dimension(nlayer,nb_phot_max) :: e_phot ! photolysis rates by energie (J.mol-1.s-1) 46 46 47 ! s olar flux at mars48 49 real, dimension(nw) :: f mars ! solar flux (w.m-2.nm-1)47 ! stellar flux at planet 48 49 real, dimension(nw) :: fplanet ! stellar flux (photon.s-1.nm-1.cm-2) at planet 50 50 real :: factor 51 51 … … 198 198 call sphers(nlev,zalt,sza,dsdh,nid) 199 199 200 !==== s olar flux at mars200 !==== stellar flux at planet 201 201 202 202 factor = (1./dist_sol)**2. 203 203 do iw = 1,nw-1 204 f mars(iw) = f(iw)*factor204 fplanet(iw) = fstar1AU(iw)*factor 205 205 end do 206 206 … … 229 229 230 230 do ilay = 1,nlayer 231 saflux(ilay) = fmars(iw)*(fdir(ilay) + fdn(ilay) + fup(ilay)) 231 saflux(ilay) = fplanet(iw)* 232 $ (fdir(ilay) + fdn(ilay) + fup(ilay)) 233 fluxUV(ig,iw,ilay) = saflux(ilay) 232 234 end do 233 235 -
trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90
r3299 r3309 87 87 logical,save :: depos 88 88 logical,save :: haze 89 !$OMP THREADPRIVATE(photochem,photoheat,jonline,depos) 89 logical,save :: output_writediagspecUV 90 !$OMP THREADPRIVATE(photochem,photoheat,jonline,depos,output_writediagspecUV) 90 91 logical,save :: calllott_nonoro 91 92 !$OMP THREADPRIVATE(calllott_nonoro) -
trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90
r3299 r3309 543 543 if (is_master) write(*,*)trim(rname)//": haze = ",haze 544 544 545 if (is_master) write(*,*)trim(rname)//": Output photochemical surface UV flux ?" 546 output_writediagspecUV=.false. ! default value 547 call getin_p("output_writediagspecUV",output_writediagspecUV) 548 if (is_master) write(*,*)trim(rname)//": output_writediagspecUV = ",output_writediagspecUV 549 545 550 ! Global1D mean and solar zenith angle 546 551
Note: See TracChangeset
for help on using the changeset viewer.