Changeset 2851 for trunk/LMDZ.VENUS/libf
- Timestamp:
- Dec 21, 2022, 10:55:05 AM (2 years ago)
- Location:
- trunk/LMDZ.VENUS/libf/phyvenus
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.VENUS/libf/phyvenus/dyn1d/rcm1d.F
r2684 r2851 89 89 real :: dummy 90 90 91 character*8 specname( 36)92 real mmol( 36)91 character*8 specname(40) 92 real mmol(40) 93 93 94 94 c======================================================================= … … 103 103 & 99., 32., 48., 64., 80., 104 104 & 96., 60., 81., 98., 64., 105 & 99., 83., 28., 4., 18.,106 & 98./)105 & 99., 83., 28., 4., 46., 106 & 30., 14., 14., 18., 98./) 107 107 108 108 specname = (/"co2", "co", "h2", "h2o", "o1d", … … 112 112 & "cocl2", "s", "so", "so2", "so3", 113 113 & "s2o2", "ocs", "hso3", "h2so4", "s2", 114 & "clso2", "oscl", "n2", "he", " h2oliq",115 & " h2so4liq"/)114 & "clso2", "oscl", "n2", "he", "no2", "no", 115 & "n", "n2d", "h2oliq", "h2so4liq"/) 116 116 117 117 c ------------------------------------------------------ … … 514 514 c ------------------------------------------------------------ 515 515 516 if ((idt==1) .OR. (idt==9600) .OR. (idt==19200) .OR. 517 $ (idt==28800) .OR. (idt==38400) .OR. (idt==48000)) then 516 if (mod(idt,9600)==0) then 518 517 DO ilayer=1,nlayer 519 518 write (5,'(100e12.4)')play(ilayer)/100., q(ilayer,:) … … 523 522 ENDDO ! fin de la boucle temporelle 524 523 525 cclose(5)524 close(5) 526 525 c ======================================================== 527 526 c GESTION DES SORTIE -
trunk/LMDZ.VENUS/libf/phyvenus/photochemistry_venus.F90
r2843 r2851 98 98 ! if prior to jvenus.20211025, set nztable = 201 below 99 99 100 integer, parameter :: nj = 22, nztable = 281, nsza = 27, nso2 = 13 101 !integer, parameter :: nztable = 281, nsza = 27, nso2 = 13 100 integer, parameter :: nj = 23, nztable = 281, nsza = 27, nso2 = 13 102 101 real, dimension(nso2,nsza,nztable,nj), save :: jphot ! nj must be equal to nphot 103 102 real, dimension(nztable), save :: table_colair … … 215 214 call photolysis_online(nz, nb_phot_max, zlocal, p, & 216 215 t, mumean, i_co2,i_co, i_o, i_o1d, & 217 i_o2, i_o3, i_ oh, i_ho2, i_h2o2, i_h2o,&216 i_o2, i_o3, i_h2, i_oh, i_ho2, i_h2o2, i_h2o,& 218 217 i_h, i_hcl, i_cl2, i_hocl, i_so2, i_so, & 219 218 i_so3, i_clo, i_ocs, i_cocl2, i_h2so4, i_cl,& … … 263 262 end if 264 263 else 265 call phot( 264 call phot(nj, nztable, nsza, nso2, sza_input, dist_sol, mumean, tr(:,i_co2), tr(:,i_so2), & 266 265 jphot, table_colair, table_colso2, table_sza, nz, nb_phot_max, t, p, v_phot) 267 266 end if … … 529 528 530 529 !=========================================================== 530 ! H2 + hv -> H + H 531 !=========================================================== 532 533 nb_phot = nb_phot + 1 534 535 indice_phot(nb_phot) = z3spec(1.0, i_h2, 2.0, i_h, 0.0, i_dummy) 536 537 538 !=========================================================== 531 539 ! H2O + hv -> H + OH 532 540 !=========================================================== … … 997 1005 998 1006 nb_reaction_3 = nb_reaction_3 + 1 999 1000 1007 indice_3(nb_reaction_3) = z3spec(2.0, i_ho2, 1.0, i_h2o2, 1.0, i_o2) 1001 1008 … … 2526 2533 !================================================================ 2527 2534 2528 2529 2535 USE chemparam_mod 2530 2536 USE photolysis_mod, only : nphot … … 2564 2570 integer :: iz 2565 2571 real :: ak0, ak1, xpo, rate, rate1, rate2, pi, gam 2566 real :: k1a0, k1b0, k1ainf, k1a, k1b, fc, fx, x, y2572 real :: k1a0, k1b0, k1ainf, k1a, k1b, k0, kinf, kf, kint, kca, fc, fx, x, y 2567 2573 integer :: nb_phot, nb_reaction_3, nb_reaction_4 2568 2574 real, dimension(nz) :: surfice1d, surfdust1d … … 2606 2612 i062, i063, & 2607 2613 j001, j002 2608 2609 2614 !---------------------------------------------------------------------- 2610 2615 ! initialisation … … 2621 2626 !---------------------------------------------------------------------- 2622 2627 2623 !--- a001: o + o2 + (co2 or M) -> o3 + (co2 or M)2624 2625 2628 ! jpl 2003 2626 2629 ! ------ BEFORE VCD 2.1 ------- … … 2633 2636 ( 2.5 * c(:,i_co2) + & 2634 2637 1.0 * (conc(:)-c(:,i_co2)) ) 2635 ! ----------------------------- 2636 2637 nb_reaction_4 = nb_reaction_4 + 1 2638 v_4(:,nb_reaction_4) = a001(:) 2639 2640 !--- a002: o + o + (co2 or M) -> o2(Dg) + (co2 or M) 2638 nb_reaction_4 = nb_reaction_4 + 1 2639 v_4(:,nb_reaction_4) = a001(:) 2640 2641 2641 2642 2642 ! Tsang and Hampson, J. Chem. Phys. Ref. Data, 15, 1087, 1986 2643 2644 ! a002(:) = 2.5*5.2E-35*exp(900./t(:))*conc(:)2645 2646 ! Campbell and Gray, Chem. Phys. Lett., 18, 607, 19732647 2643 ! ------ BEFORE VCD 2.1 ------- 2648 2644 ! a002(:) = 2.5*9.46E-34*exp(485./t(:))*conc(:) ! nist expression … … 2664 2660 nb_reaction_3 = nb_reaction_3 + 1 2665 2661 v_3(:,nb_reaction_3) = a002(:) 2666 2667 2662 ind_orec = nb_reaction_3 2668 2663 2669 ! ---a003: o + o3 -> o2 + o22664 ! a003: o + o3 -> o2 + o2 2670 2665 2671 2666 ! jpl 2003 … … 2692 2687 2693 2688 ! jpl 2006 2694 2689 2695 2690 b002(:) = 1.63E-10*exp(60./t(:)) 2696 2691 2697 2692 nb_reaction_4 = nb_reaction_4 + 1 2698 2693 v_4(:,nb_reaction_4) = b002(:) … … 2715 2710 nb_phot = nb_phot + 1 2716 2711 v_phot(:,nb_phot) = b004(:)*c(:,i_o2) 2717 2712 2718 2713 !--- b005: o(1d) + o3 -> o2 + o2 2719 2714 … … 2724 2719 nb_reaction_4 = nb_reaction_4 + 1 2725 2720 v_4(:,nb_reaction_4) = b005(:) 2726 2721 2727 2722 !--- b006: o(1d) + o3 -> o2 + o + o 2728 2723 … … 2733 2728 nb_reaction_4 = nb_reaction_4 + 1 2734 2729 v_4(:,nb_reaction_4) = b006(:) 2735 2730 2736 2731 !---------------------------------------------------------------------- 2737 2732 ! reactions des hox … … 2836 2831 ! jpl 2006 2837 2832 2833 ! do iz = 1,nz 2834 ! ak0 = 2.5*4.4E-32*(t(iz)/300.)**(-1.3) 2835 ! ak1 = 4.7E-11*(t(iz)/300.)**(-0.2) 2836 2837 ! rate = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1) 2838 ! xpo = 1./(1. + alog10((ak0*conc(iz))/ak1)**2) 2839 ! c011(iz) = rate*0.6**xpo 2840 ! end do 2841 2842 ! jpl 2019 2843 2838 2844 do iz = 1,nz 2839 ak0 = 2. 5*4.4E-32*(t(iz)/300.)**(-1.3)2840 ak1 = 4.7E-11*(t(iz)/300.)**(-0.2)2845 ak0 = 2.4*5.3e-32*(t(iz)/298.)**(-1.8) 2846 ak1 = 9.5e-11*(t(iz)/298.)**(0.4) 2841 2847 2842 2848 rate = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1) … … 2926 2932 ! jpl 2006 2927 2933 2928 d001(:) = 5.1e-12*exp(210./t(:)) 2934 ! d001(:) = 5.1e-12*exp(210./t(:)) 2935 2936 ! jpl 2019 2937 2938 ! For the sake of simplicity, it is assumed that the association 2939 ! yield (kf) gives the same product as the chemical activation yield 2940 ! (kca). Thus the only products are no + o2. There is no production 2941 ! of no3. 2942 2943 do iz = 1,nz 2944 2945 ! association 2946 2947 k0 = 2.5*3.4e-31*(298./t(iz))**(1.6) 2948 kinf = 2.3e-11*(298./t(iz))**(0.2) 2949 2950 kf = (kinf*k0*conc(iz)/(kinf + k0*conc(iz))) & 2951 *0.6**(1. + (log10(k0*conc(iz)/kinf))**2.)**(-1.0) 2952 2953 ! chemical acitvation 2954 2955 kint = 5.3e-12*exp(200./t(iz)) 2956 2957 kca = kint*(1. - kf/kinf) 2958 2959 ! total : association + chemical activation 2960 2961 d001(iz) = kf + kca 2962 2963 end do 2929 2964 2930 2965 nb_reaction_4 = nb_reaction_4 + 1 … … 2944 2979 ! jpl 2011 2945 2980 2946 d003(:) = 3.3e-12*exp(270./t(:)) 2981 ! d003(:) = 3.3e-12*exp(270./t(:)) 2982 2983 ! jpl 2019 2984 2985 d003(:) = 3.44e-12*exp(270./t(:)) 2947 2986 2948 2987 nb_reaction_4 = nb_reaction_4 + 1 … … 2963 3002 ! jpl 2011 2964 3003 2965 d005(:) = 1.5e-11*exp(-3600./t(:)) 3004 ! d005(:) = 1.5e-11*exp(-3600./t(:)) 3005 3006 ! jpl 2019 3007 3008 d005(:) = 3.3e-12*exp(-3150./t(:)) 2966 3009 2967 3010 nb_reaction_4 = nb_reaction_4 + 1 … … 2972 3015 ! jpl 2011 2973 3016 2974 d006(:) = 4.0e-10*exp(-340./t(:)) 3017 ! d006(:) = 4.0e-10*exp(-340./t(:)) 3018 3019 ! jpl 2019 3020 3021 d006(:) = 1.35e-10 2975 3022 2976 3023 nb_reaction_4 = nb_reaction_4 + 1 … … 2979 3026 !--- d007: n + o -> no 2980 3027 2981 d007(:) = 2.8e-17*(300./t(:))**0.53028 d007(:) = 1.9e-17*(300./t(:))**(0.5)*exp(1-(0.57/(t(:)**(0.5)))) 2982 3029 2983 3030 nb_reaction_4 = nb_reaction_4 + 1 … … 3046 3093 ! jpl 2015 3047 3094 3048 do iz = 1,nz3095 ! do iz = 1,nz 3049 3096 3050 3097 ! branch 1 : oh + co -> h + co2 3051 3098 3052 rate1 = 1.5e-13*(t(iz)/300.)**(0.0)3099 ! rate1 = 1.5e-13*(t(iz)/300.)**(0.0) 3053 3100 3054 3101 ! branch 2 : oh + co + m -> hoco + m 3055 3102 3056 ak0 = 5.9e-33*(t(iz)/300.)**(-1.0)3057 ak1 = 1.1e-12*(t(iz)/300.)**(1.3)3058 rate2 = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1)3059 xpo = 1./(1. + alog10((ak0*conc(iz))/ak1)**2)3060 3061 e001(iz) = rate1 + rate2*0.6**xpo3062 end do3103 ! ak0 = 5.9e-33*(t(iz)/300.)**(-1.0) 3104 ! ak1 = 1.1e-12*(t(iz)/300.)**(1.3) 3105 ! rate2 = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1) 3106 ! xpo = 1./(1. + alog10((ak0*conc(iz))/ak1)**2) 3107 3108 ! e001(iz) = rate1 + rate2*0.6**xpo 3109 ! end do 3063 3110 3064 3111 ! joshi et al., 2006 … … 3082 3129 ! end do 3083 3130 3131 ! jpl 2019 3132 3133 do iz = 1,nz 3134 3135 ! association 3136 3137 k0 = 2.5*6.9e-33*(398./t(iz))**(2.1) 3138 kinf = 1.1e-12*(298./t(iz))**(-1.3) 3139 3140 kf = (kinf*k0*conc(iz)/(kinf + k0*conc(iz))) & 3141 *0.6**(1. + (log10(k0*conc(iz)/kinf))**2.)**(-1.0) 3142 3143 ! chemical activation 3144 3145 kint = 1.85e-13*exp(-65/t(iz)) 3146 3147 kca = kint*(1. - kf/kinf) 3148 3149 ! total : association + chemical activation 3150 3151 e001(iz) = kf + kca 3152 3153 end do 3154 3084 3155 nb_reaction_4 = nb_reaction_4 + 1 3085 3156 v_4(:,nb_reaction_4) = e001(:) … … 3279 3350 !--- f018: clo + ho2 -> hocl + o2 3280 3351 3281 f018(:) = 2.7E-12*exp(220./t(:)) 3352 ! jpl 2019 3353 3354 f018(:) = 2.6E-12*exp(290./t(:)) 3282 3355 3283 3356 nb_reaction_4 = nb_reaction_4 + 1 … … 3515 3588 !--- g005: so + oh -> so2 + h 3516 3589 3517 ! jpl 201 53518 3519 g005(:) = 2. 7E-11*exp(335./t(:))3590 ! jpl 2019 3591 3592 g005(:) = 2.6E-11*exp(330./t(:)) 3520 3593 3521 3594 nb_reaction_4 = nb_reaction_4 + 1 … … 3557 3630 !--- g009: so2 + o + co2 -> so3 + co2 3558 3631 3559 ! jpl 20113560 3632 ! Naido 2005 3561 3633 3634 do iz = 1,nz 3635 ak0 = 5.*9.5*1.E-23*(t(iz)**(-3.0))*EXP(-2400./t(iz)) 3636 ak1 = 6.1*1.E-13*EXP(-850./t(iz)) 3637 rate = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1) 3638 xpo = 1./(1. + alog10((ak0*conc(iz))/ak1)**2) 3639 fc = 0.558*EXP(-t(iz)/316.)+0.442*EXP(-t(iz)/7442.) 3640 g009(iz) = rate*fc**xpo 3641 end do 3642 3643 ! jpl 2019 3644 3562 3645 ! do iz = 1,nz 3563 ! ak0 = 2.5*1.8 E-33*(t(iz)/300.)**(2.0)3564 ! ak1 = 4. 2E-14*(t(iz)/300.)**(1.8)3646 ! ak0 = 2.5*1.8e-33*(t(iz)/298)**(2.0) 3647 ! ak1 = 4.1e-14*(t(iz)/298)**(1.8) 3565 3648 3566 3649 ! rate = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1) 3567 ! xpo = 1./(1. + alog10((ak0*conc(iz))/ak1)**2)3650 ! xpo = 1./(1. + log10((ak0*conc(iz))/ak1)**2) 3568 3651 ! g009(iz) = rate*0.6**xpo 3569 ! g009(iz) = 0.0E+03570 3652 ! end do 3571 3572 do iz = 1,nz3573 ak0 = 5.*9.5*1.E-23*(t(iz)**(-3.0))*EXP(-2400./t(iz))3574 ak1 = 6.1*1.E-13*EXP(-850./t(iz))3575 rate = (ak0*conc(iz))/(1. + ak0*conc(iz)/ak1)3576 xpo = 1./(1. + alog10((ak0*conc(iz))/ak1)**2)3577 fc = 0.558*EXP(-t(iz)/316.)+0.442*EXP(-t(iz)/7442.)3578 g009(iz) = rate*fc**xpo3579 ! g009(iz) = 0.0E+03580 end do3581 3653 3582 3654 nb_reaction_4 = nb_reaction_4 + 1 -
trunk/LMDZ.VENUS/libf/phyvenus/photolysis_mod.F90
r2836 r2851 5 5 ! photolysis 6 6 7 integer, save :: nphot = 2 2! number of photolysis7 integer, save :: nphot = 23 ! number of photolysis 8 8 9 9 !$OMP THREADPRIVATE(nphot) 10 10 11 integer, parameter :: nabs = 19! number of absorbing gases11 integer, parameter :: nabs = 20 ! number of absorbing gases 12 12 13 13 ! spectral grid -
trunk/LMDZ.VENUS/libf/phyvenus/photolysis_online.F
r2836 r2851 3 3 subroutine photolysis_online(nlayer, nb_phot_max, 4 4 $ alt, press, temp, zmmean, 5 $ i_co2, i_co, i_o, i_o1d, i_o2, i_o3, 5 $ i_co2, i_co, i_o, i_o1d, i_o2, i_o3,i_h2, 6 6 $ i_oh, i_ho2, i_h2o2, i_h2o,i_h,i_hcl, 7 7 $ i_cl2, i_hocl, i_so2, i_so, i_so3, … … 23 23 $ i_cl2, i_hocl, i_so2, i_so, i_so3, i_clo, 24 24 $ i_ocs, i_cl, i_cocl2, i_h2so4, i_no2, 25 $ i_no, i_n2, i_n2d 25 $ i_no, i_n2, i_n2d, i_h2 26 26 27 27 real, dimension(nlayer), intent(in) :: press, temp, zmmean ! pressure (hpa)/temperature (k)/mean molecular mass (g.mol-1) … … 67 67 $ j_h2o, j_h2o2, j_ho2, j_h, j_hcl, j_cl2, j_hocl, j_so2, 68 68 $ j_so, j_so3, j_clo, j_ocs, j_cocl2, j_h2so4, j_no2, 69 $ j_no, j_n2 69 $ j_no, j_n2, j_h2 70 70 71 71 integer :: a_o2, a_co2, a_o3, a_h2o, a_h2o2, a_ho2, a_hcl, a_cl2, 72 72 $ a_hocl, a_so2, a_so, a_so3, a_clo, a_ocs, a_cocl2, 73 $ a_h2so4, a_no2, a_no, a_n2 73 $ a_h2so4, a_no2, a_no, a_n2, a_h2 74 74 integer :: i, ilay, iw, ialt 75 75 real :: deltaj … … 81 81 a_co2 = 2 ! co2 82 82 a_o3 = 3 ! o3 83 a_h2o = 4 ! h2o 84 a_h2o2 = 5 ! h2o2 85 a_ho2 = 6 ! ho2 86 a_hcl = 7 ! hcl 87 a_cl2 = 8 ! cl2 88 a_hocl = 9 ! hocl 89 a_so2 = 10 ! so2 90 a_so = 11 ! so 91 a_so3 = 12 ! so3 92 a_clo = 13 ! clo 93 a_ocs = 14 ! ocs 94 a_cocl2 = 15 ! cocl2 95 a_h2so4 = 16 ! h2so4 96 a_no2 = 17 ! no2 97 a_no = 18 ! no 98 a_n2 = 19 ! n2 83 a_h2 = 4 ! h2 84 a_h2o = 5 ! h2o 85 a_h2o2 = 6 ! h2o2 86 a_ho2 = 7 ! ho2 87 a_hcl = 8 ! hcl 88 a_cl2 = 9 ! cl2 89 a_hocl = 10 ! hocl 90 a_so2 = 11 ! so2 91 a_so = 12 ! so 92 a_so3 = 13 ! so3 93 a_clo = 14 ! clo 94 a_ocs = 15 ! ocs 95 a_cocl2 = 16 ! cocl2 96 a_h2so4 = 17 ! h2so4 97 a_no2 = 18 ! no2 98 a_no = 19 ! no 99 a_n2 = 20 ! n2 99 100 100 101 ! photodissociation rates numbering. … … 107 108 j_o3_o1d = 5 ! o3 + hv -> o2 + o(1d) 108 109 j_o3_o = 6 ! o3 + hv -> o2 + o 109 j_h2o = 7 ! h2o + hv -> h + oh 110 j_ho2 = 8 ! ho2 + hv -> oh + o 111 j_h2o2 = 9 ! h2o2 + hv -> oh + oh 112 j_hcl = 10 ! hcl + hv -> h + cl 113 j_cl2 = 11 ! cl2 + hv -> cl + cl 114 j_hocl = 12 ! hocl + hv -> oh + cl 115 j_so2 = 13 ! so2 + hv -> so + o 116 j_so = 14 ! so + hv -> s + o 117 j_so3 = 15 ! so3 + hv -> so2 + o 118 j_clo = 16 ! clo + hv -> cl + o 119 j_ocs = 17 ! ocs + hv -> co + s 120 j_cocl2 = 18 ! cocl2 + hv -> 2cl + co 121 j_h2so4 = 19 ! h2so4 + hv -> so3 + h2o 122 j_no2 = 20 ! no2 + hv -> no + o 123 j_no = 21 ! no + hv -> n + o 124 j_n2 = 22 ! n2 + hv -> n(2d) + n 110 j_h2 = 7 ! h2 + hv -> h + h 111 j_h2o = 8 ! h2o + hv -> h + oh 112 j_ho2 = 9 ! ho2 + hv -> oh + o 113 j_h2o2 = 10 ! h2o2 + hv -> oh + oh 114 j_hcl = 11 ! hcl + hv -> h + cl 115 j_cl2 = 12 ! cl2 + hv -> cl + cl 116 j_hocl = 13 ! hocl + hv -> oh + cl 117 j_so2 = 14 ! so2 + hv -> so + o 118 j_so = 15 ! so + hv -> s + o 119 j_so3 = 16 ! so3 + hv -> so2 + o 120 j_clo = 17 ! clo + hv -> cl + o 121 j_ocs = 18 ! ocs + hv -> co + s 122 j_cocl2 = 19 ! cocl2 + hv -> 2cl + co 123 j_h2so4 = 20 ! h2so4 + hv -> so3 + h2o 124 j_no2 = 21 ! no2 + hv -> no + o 125 j_no = 22 ! no + hv -> n + o 126 j_n2 = 23 ! n2 + hv -> n(2d) + n 125 127 126 128 ! j_hdo_od = ! hdo + hv -> od + h … … 174 176 do ilay = 1,nlayer 175 177 do iw = 1,nw-1 178 dtgas(ilay,iw,a_h2) = colinc(ilay)*rm(ilay,i_h2)*xsh2(iw) 176 179 dtgas(ilay,iw,a_h2o) = colinc(ilay)*rm(ilay,i_h2o)*xsh2o(iw) 177 180 dtgas(ilay,iw,a_ho2) = colinc(ilay)*rm(ilay,i_ho2)*xsho2(iw) … … 209 212 do ilay = 1,nlayer 210 213 do iw = 1,nw-1 214 sj(ilay,iw,j_h2) = xsh2(iw) ! h2 211 215 sj(ilay,iw,j_h2o) = xsh2o(iw) ! h2o 212 216 sj(ilay,iw,j_ho2) = xsho2(iw) ! ho2
Note: See TracChangeset
for help on using the changeset viewer.