Changeset 3755
- Timestamp:
- May 6, 2025, 9:45:32 AM (5 weeks ago)
- Location:
- trunk/LMDZ.VENUS/libf/phyvenus
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.VENUS/libf/phyvenus/chemparam_mod.F90
r3689 r3755 21 21 i_osso_cis, i_osso_trans, i_s2o2_cyc, & 22 22 i_ocs, i_hso3, i_h2so4, i_s2, & 23 i_clso2, i_ oscl, i_n2, i_he, i_n, i_no,&24 i_n o2, i_n2d,&23 i_clso2, i_cl2so2, i_oscl, i_n2, i_he, & 24 i_n, i_no, i_no2, i_n2d, & 25 25 i_co2plus, i_coplus, i_oplus, i_o2plus, & 26 26 i_n2plus, i_hplus, i_h2oplus, i_nplus, & … … 1542 1542 i_s2 = 0 1543 1543 i_clso2 = 0 1544 i_cl2so2 = 0 1544 1545 i_oscl = 0 1545 1546 i_n2 = 0 … … 1744 1745 m_tr(i_clso2) = 99.517 1745 1746 type_tr(i_clso2) = 1 1747 case('cl2so2') 1748 i_cl2so2 = i 1749 print*,'cl2so2',i_cl2so2 1750 m_tr(i_cl2so2) = 134.97 1751 type_tr(i_cl2so2) = 1 1746 1752 case('oscl') 1747 1753 i_oscl = i -
trunk/LMDZ.VENUS/libf/phyvenus/photochemistry_venus.F90
r3722 r3755 212 212 i_cl, i_clo, i_cl2, i_hocl, i_so2, i_so, i_so3, i_s2, & 213 213 i_osso_cis, i_osso_trans, i_s2o2_cyc, i_clso2, & 214 i_ ocs, i_cocl2, i_h2so4,&214 i_cl2so2, i_ocs, i_cocl2, i_h2so4, & 215 215 i_no2, i_no, i_n2, i_n2d, & 216 216 nesp, tr, sza_input, dist_sol, v_phot) … … 652 652 653 653 !=========================================================== 654 ! Cl2SO2 + hv -> Cl + ClSO2 655 !=========================================================== 656 657 nb_phot = nb_phot + 1 658 659 indice_phot(nb_phot) = z3spec(1.0, i_cl2so2, 1.0, i_cl, 1.0, i_clso2) 660 661 !=========================================================== 654 662 ! OCS + hv -> CO + S 655 663 !=========================================================== … … 1429 1437 1430 1438 !=========================================================== 1431 ! f030 : ClSO2 + ClSO2 -> Cl2 +SO2 + SO21439 ! f030 : ClSO2 + ClSO2 -> Cl2SO2 + SO2 1432 1440 !=========================================================== 1433 1441 1434 1442 nb_reaction_3 = nb_reaction_3 + 1 1435 1443 1436 indice_3(nb_reaction_3) = z3spec(2.0, i_clso2, 1.0, i_cl2 , 2.0, i_so2)1444 indice_3(nb_reaction_3) = z3spec(2.0, i_clso2, 1.0, i_cl2so2, 1.0, i_so2) 1437 1445 1438 1446 !=========================================================== … … 3540 3548 v_4(:,nb_reaction_4) = f029(:) 3541 3549 3542 !--- f030: clso2 + clso2 -> cl2 +so2 + so23550 !--- f030: clso2 + clso2 -> cl2so2 + so2 3543 3551 3544 3552 ! moses et al. 2002 3545 3553 3546 3554 ! f030(:) = 5.0E-13 3555 3556 ! croce and cobos, 2018 3547 3557 3548 3558 do iz = 1,nz -
trunk/LMDZ.VENUS/libf/phyvenus/photolysis_mod.F90
r3722 r3755 5 5 ! photolysis 6 6 7 integer, save :: nphot = 2 8! number of photolysis7 integer, save :: nphot = 29 ! number of photolysis 8 8 9 9 !$OMP THREADPRIVATE(nphot) 10 10 11 integer, parameter :: nabs = 2 5! number of absorbing gases11 integer, parameter :: nabs = 26 ! number of absorbing gases 12 12 13 13 ! spectral grid … … 48 48 real, dimension(nw), save :: xss2o2_cyc ! s2o2_cyc absorption cross-section (cm2) 49 49 real, dimension(nw), save :: xsclso2 ! clso2 absorption cross-section (cm2) 50 real, dimension(nw), save :: xscl2so2 ! cl2so2 absorption cross-section (cm2) 50 51 real, dimension(nw), save :: xsocs ! cos absorption cross-section (cm2) 51 52 real, dimension(nw), save :: xscocl2 ! cocl2 absorption cross-section (cm2) … … 62 63 real, dimension(nw), save :: albedo ! surface albedo 63 64 64 !$OMP THREADPRIVATE(f, xsco2_195, xsco2_295, xsco2_370, yieldco2, xso2_150, xso2_200, xso2_250, xso2_300, yieldo2, xso3_218, xso3_298, xs h2o, xsh2o2, xsho2, xsh2, yieldh2, xsno2, xsno2_220, xsno2_294, yldno2_248, yldno2_298)65 !$OMP THREADPRIVATE(f, xsco2_195, xsco2_295, xsco2_370, yieldco2, xso2_150, xso2_200, xso2_250, xso2_300, yieldo2, xso3_218, xso3_298, xsclso2, xscl2so2, xsh2o, xsh2o2, xsho2, xsh2, yieldh2, xsno2, xsno2_220, xsno2_294, yldno2_248, yldno2_298) 65 66 !$OMP THREADPRIVATE(xsno, yieldno, xsn2, yieldn2, xshdo, albedo) 66 67 … … 148 149 149 150 call rdxsclso2(nw,wl,wc,xsclso2) 151 152 ! read and grid cl2so2 cross-sections 153 154 call rdxscl2so2(nw,wl,wc,xscl2so2) 150 155 151 156 ! read and grid ocs cross-sections … … 2795 2800 !============================================================================== 2796 2801 2802 subroutine rdxscl2so2(nw, wl, wc, yg) 2803 2804 !-----------------------------------------------------------------------------* 2805 != PURPOSE: =* 2806 != Read cl2so2 cross-sections =* 2807 != From Uthman et al., J. Phys. Chem., 1978 =* 2808 !-----------------------------------------------------------------------------* 2809 != PARAMETERS: =* 2810 != NW - INTEGER, number of specified intervals + 1 in working (I)=* 2811 != wavelength grid =* 2812 !-----------------------------------------------------------------------------* 2813 2814 USE mod_phys_lmdz_para, ONLY: is_master 2815 USE mod_phys_lmdz_transfert_para, ONLY: bcast 2816 2817 IMPLICIT NONE 2818 2819 ! input 2820 2821 integer :: nw ! number of wavelength grid points 2822 real, dimension(nw) :: wl, wc ! lower and central wavelength for each interval 2823 2824 ! output 2825 2826 real, dimension(nw) :: yg ! clso2 cross-sections (cm2) 2827 2828 ! local 2829 2830 real, parameter :: deltax = 1.e-4 2831 integer, parameter :: kdata = 1500 2832 2833 real, dimension(kdata) :: x1, y1 2834 real :: qy, lambda 2835 integer :: i, iw, n, ierr 2836 integer :: kin, kout ! input/output logical units 2837 character*100 fil 2838 2839 kin = 10 2840 2841 !*** cross sections from Uthman et al., J. Phys. Chem., 1978 2842 2843 fil = 'cross_sections/cl2so2_Uthman1978_295K_190_300nm.txt' 2844 print*, 'section efficace Cl2SO2: ', fil 2845 2846 if (is_master) then 2847 2848 n = 56 2849 2850 OPEN(kin,FILE=fil,STATUS='OLD') 2851 2852 do i = 1,7 2853 read(kin,*) 2854 end do 2855 do i = 1,n 2856 read(kin,*) x1(i), y1(i) 2857 end do 2858 close(kin) 2859 2860 call addpnt(x1,y1,kdata,n,x1(1)*(1.-deltax),0.) 2861 call addpnt(x1,y1,kdata,n, 0.,0.) 2862 call addpnt(x1,y1,kdata,n,x1(n)*(1.+deltax),0.) 2863 call addpnt(x1,y1,kdata,n, 1e38,0.) 2864 2865 call inter2(nw,wl,yg,n,x1,y1,ierr) 2866 2867 if (ierr /= 0) then 2868 write(*,*) ierr, fil 2869 stop 2870 end if 2871 2872 do iw = 1,nw - 1 2873 lambda = wc(iw) 2874 ! write(55,'(f8.3,3e12.4)') lambda, yg(iw) 2875 end do 2876 2877 end if !is_master 2878 2879 call bcast(yg) 2880 2881 end subroutine rdxscl2so2 2882 2883 !============================================================================== 2884 2797 2885 subroutine rdxsclo(nw, wl, yg) 2798 2886 -
trunk/LMDZ.VENUS/libf/phyvenus/photolysis_online.F
r3722 r3755 6 6 $ i_oh, i_ho2, i_h2o2, i_h2o, i_h, i_hcl, 7 7 $ i_cl, i_clo, i_cl2, i_hocl, i_so2, i_so, i_so3, i_s2, 8 $ i_osso_cis, i_osso_trans, i_s2o2_cyc, i_clso2, 9 $ i_ ocs, i_cocl2, i_h2so4,8 $ i_osso_cis, i_osso_trans, i_s2o2_cyc, i_clso2, 9 $ i_cl2so2, i_ocs, i_cocl2, i_h2so4, 10 10 $ i_no2, i_no, i_n2, i_n2d, 11 11 & nesp, rm, sza, dist_sol, v_phot) … … 23 23 $ i_oh, i_ho2, i_h2o2, i_h2o, i_h, i_hcl, 24 24 $ i_cl, i_clo, i_cl2, i_hocl, i_so2, i_so, 25 $ i_so3, i_s2, 26 $ i_ osso_cis, i_osso_trans, i_s2o2_cyc,27 $ i_clso2, i_ ocs, i_cocl2, i_h2so4,25 $ i_so3, i_s2, i_osso_cis, i_osso_trans, 26 $ i_s2o2_cyc, 27 $ i_clso2, i_cl2so2, i_ocs, i_cocl2, i_h2so4, 28 28 $ i_no2, i_no, i_n2, i_n2d 29 29 … … 72 72 $ j_h2o, j_h2o2, j_ho2, j_h, j_hcl, j_cl2, j_hocl, j_clo, 73 73 $ j_so2, j_so, j_so3, j_s2, j_osso_cis, j_osso_trans, 74 $ j_s2o2_cyc, j_clso2, j_ ocs, j_cocl2, j_h2so4,74 $ j_s2o2_cyc, j_clso2, j_cl2so2, j_ocs, j_cocl2, j_h2so4, 75 75 $ j_no2, j_no, j_n2, j_h2 76 76 77 77 integer :: a_o2, a_co2, a_o3, a_h2o, a_h2o2, a_ho2, a_hcl, a_cl2, 78 78 $ a_hocl, a_clo, a_so2, a_so, a_so3, a_s2, a_osso_cis, 79 $ a_osso_trans, a_s2o2_cyc, a_clso2, a_ ocs,79 $ a_osso_trans, a_s2o2_cyc, a_clso2, a_cl2so2, a_ocs, 80 80 $ a_cocl2, a_h2so4, a_no2, a_no, a_n2, a_h2 81 81 … … 105 105 a_s2o2_cyc = 18 ! s2o2_cyc 106 106 a_clso2 = 19 ! clso2 107 a_ocs = 20 ! ocs 108 a_cocl2 = 21 ! cocl2 109 a_h2so4 = 22 ! h2so4 110 a_no2 = 23 ! no2 111 a_no = 24 ! no 112 a_n2 = 25 ! n2 107 a_cl2so2 = 20 ! cl2so2 108 a_ocs = 21 ! ocs 109 a_cocl2 = 22 ! cocl2 110 a_h2so4 = 23 ! h2so4 111 a_no2 = 24 ! no2 112 a_no = 25 ! no 113 a_n2 = 26 ! n2 113 114 114 115 ! photodissociation rates numbering. … … 137 138 j_s2o2_cyc = 21 ! s2o2_cyc + hv -> so + so 138 139 j_clso2 = 22 ! clso2 + hv -> cl + so2 139 j_ocs = 23 ! ocs + hv -> co + s 140 j_cocl2 = 24 ! cocl2 + hv -> 2cl + co 141 j_h2so4 = 25 ! h2so4 + hv -> so3 + h2o 142 j_no2 = 26 ! no2 + hv -> no + o 143 j_no = 27 ! no + hv -> n + o 144 j_n2 = 28 ! n2 + hv -> n(2d) + n 140 j_cl2so2 = 23 ! cl2so2 + hv -> cl + clso2 141 j_ocs = 24 ! ocs + hv -> co + s 142 j_cocl2 = 25 ! cocl2 + hv -> 2cl + co 143 j_h2so4 = 26 ! h2so4 + hv -> so3 + h2o 144 j_no2 = 27 ! no2 + hv -> no + o 145 j_no = 28 ! no + hv -> n + o 146 j_n2 = 29 ! n2 + hv -> n(2d) + n 145 147 146 148 ! j_hdo_od = ! hdo + hv -> od + h … … 240 242 dtgas(ilay,iw,a_clso2) = 241 243 $ colinc(ilay)*rm(ilay,i_clso2)*xsclso2(iw) 244 dtgas(ilay,iw,a_cl2so2) = 245 $ colinc(ilay)*rm(ilay,i_cl2so2)*xscl2so2(iw) 242 246 dtgas(ilay,iw,a_ocs) = colinc(ilay)*rm(ilay,i_ocs)*xsocs(iw) 243 247 dtgas(ilay,iw,a_cocl2) = colinc(ilay)*rm(ilay,i_cocl2) … … 279 283 sj(ilay,iw,j_s2o2_cyc) = xss2o2_cyc(iw) ! s2o2_cyc 280 284 sj(ilay,iw,j_clso2) = xsclso2(iw) ! clso2 285 sj(ilay,iw,j_cl2so2) = xscl2so2(iw) ! cl2so2 281 286 sj(ilay,iw,j_ocs) = xsocs(iw) ! ocs 282 287 sj(ilay,iw,j_cocl2) = xscocl2(iw) ! cocl2
Note: See TracChangeset
for help on using the changeset viewer.