Changeset 1665 for trunk/LMDZ.VENUS
- Timestamp:
- Feb 14, 2017, 2:02:02 PM (8 years ago)
- Location:
- trunk/LMDZ.VENUS/libf/phyvenus
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F
r1661 r1665 143 143 144 144 REAL flxmw(klon,klev) 145 146 ! A VIRER POUR DYNAMICO 145 147 REAL,INTENT(IN) :: plevmoy(klev+1) ! planet-averaged mean pressure (Pa) at interfaces 146 148 REAL,INTENT(IN) :: tmoy(klev) ! planet-averaged mean temperature (Pa) at mid-layers … … 180 182 181 183 c Pour calcul GW drag oro et nonoro: CALCUL de N2: 182 real zdzlev(klon,klev) 183 real ztlev(klon,klev),zpklev(klon,klev) 184 real ztetalay(klon,klev),ztetalev(klon,klev) 185 real zdtetalev(klon,klev) 184 real zdtlev(klon,klev),zdzlev(klon,klev) 185 real ztlev(klon,klev) 186 186 real zn2(klon,klev) ! BV^2 at plev 187 187 … … 540 540 do j=1,klev 541 541 rnew(ig,j)=R 542 cpnew(ig,j)=cpdet(t moy(j))542 cpnew(ig,j)=cpdet(t(ig,j)) 543 543 mmean(ig,j)=RMD 544 544 akknew(ig,j)=1.e-4 545 rho(ig,j)=pplay(ig,j)*mmean(ig,j)*1e-3/(rnew(ig,j)*t(ig,j)) 545 546 enddo 546 547 c stop … … 664 665 nmicro=0 665 666 if (ok_cloud .AND. (cl_scheme.eq.1)) nmicro=2 666 if (ok_cloud .AND. (cl_scheme.eq.2)) nmicro= 8667 if (ok_cloud .AND. (cl_scheme.eq.2)) nmicro=12 667 668 668 669 c CAS 1D POUR MICROPHYS Aurelien … … 1479 1480 c 1480 1481 if (ok_orodr.or.ok_gw_nonoro) then 1481 c CALCUL DE N2 1482 1483 c CALCUL DE N2 1484 c UTILISE LA RELATION ENTRE N2 ET STABILITE 1485 c N2 = RG/T (dT/dz+RG/cp(T)) 1486 c ET DONC EN N'UTILISE QUE LA TEMPERATURE, PAS teta. 1487 1482 1488 do i=1,klon 1483 1489 do k=2,klev 1484 1490 ztlev(i,k) = (t_seri(i,k)+t_seri(i,k-1))/2. 1485 zpklev(i,k) = sqrt(ppk(i,k)*ppk(i,k-1))1486 1491 enddo 1487 1492 enddo 1488 call t2tpot(klon*klev,ztlev, ztetalev,zpklev)1489 call t2tpot(klon*klev,t_seri,ztetalay,ppk)1490 1493 do i=1,klon 1491 1494 do k=2,klev 1492 zdtetalev(i,k) = ztetalay(i,k)-ztetalay(i,k-1) 1493 zdzlev(i,k) = (zphi(i,k)-zphi(i,k-1))/RG 1494 zn2(i,k) = RG*zdtetalev(i,k)/(ztetalev(i,k)*zdzlev(i,k)) 1495 ztlev(i,k) = (t_seri(i,k)+t_seri(i,k-1))/2. 1496 zdtlev(i,k) = t_seri(i,k)-t_seri(i,k-1) 1497 zdzlev(i,k) = (zphi(i,k)-zphi(i,k-1))/RG 1498 zn2(i,k) = RG/ztlev(i,k) * ( zdtlev(i,k)/zdzlev(i,k) 1499 . + RG/cpdet(ztlev(i,k)) ) 1495 1500 zn2(i,k) = max(zn2(i,k),1.e-12) ! securite 1496 1501 enddo … … 1796 1801 ! <field id="..." /> in context_lmdz_physics.xml to be correctly used) 1797 1802 1798 CALL send_xios_field("phis",phis) 1803 ! 2D fields 1804 1805 CALL send_xios_field("phis",pphis) 1799 1806 cell_area_out(:)=cell_area(:) 1800 1807 if (is_north_pole_phy) cell_area_out(1)=cell_area(1)/nbp_lon … … 1803 1810 CALL send_xios_field("tsol",ftsol) 1804 1811 CALL send_xios_field("psol",paprs(:,1)) 1805 CALL send_xios_field("ue",ue) 1806 c VENUS: regardee a l envers!!!!!!!!!!!!!!! 1807 CALL send_xios_field("ve",-1.*ve) 1812 CALL send_xios_field("cdragh",cdragh) 1813 CALL send_xios_field("cdragm",cdragm) 1814 1815 CALL send_xios_field("tops",topsw) 1816 CALL send_xios_field("topl",toplw) 1817 CALL send_xios_field("sols",solsw) 1818 CALL send_xios_field("soll",sollw) 1819 1820 ! 3D fields 1808 1821 1809 1822 CALL send_xios_field("temp",t_seri) 1823 CALL send_xios_field("pres",pplay) 1824 CALL send_xios_field("geop",zphi) 1810 1825 CALL send_xios_field("vitu",u_seri) 1811 1826 c VENUS: regardee a l envers!!!!!!!!!!!!!!! 1812 1827 CALL send_xios_field("vitv",-1.*v_seri) 1828 CALL send_xios_field("vitw",omega) 1829 CALL send_xios_field("Kz",ycoefh) 1830 CALL send_xios_field("mmean",mmean) 1831 CALL send_xios_field("rho",rho) 1832 1833 CALL send_xios_field("dudyn",d_u_dyn) 1834 CALL send_xios_field("duvdf",d_u_vdf) 1835 c VENUS: regardee a l envers!!!!!!!!!!!!!!! 1836 CALL send_xios_field("dvvdf",-1.*d_v_vdf) 1837 CALL send_xios_field("duajs",d_u_ajs) 1838 CALL send_xios_field("dugwo",d_u_oro) 1839 CALL send_xios_field("dugwno",d_u_hin) 1840 CALL send_xios_field("dumolvis",d_u_molvis) 1841 c VENUS: regardee a l envers!!!!!!!!!!!!!!! 1842 CALL send_xios_field("dvmolvis",-1.*d_v_molvis) 1843 CALL send_xios_field("dtdyn",d_t_dyn) 1844 CALL send_xios_field("dtphy",d_t) 1845 CALL send_xios_field("dtvdf",d_t_vdf) 1846 CALL send_xios_field("dtajs",d_t_ajs) 1847 CALL send_xios_field("dtswr",dtsw) 1848 CALL send_xios_field("dtswrNLTE",d_t_nirco2) 1849 CALL send_xios_field("dtswrLTE",heat) 1850 CALL send_xios_field("dtlwr",dtlw) 1851 CALL send_xios_field("dtlwrNLTE",d_t_nlte) 1852 CALL send_xios_field("dtlwrLTE",-1.*cool) 1853 CALL send_xios_field("dteuv",d_t_euv) 1854 CALL send_xios_field("dtcond",d_t_conduc) 1855 CALL send_xios_field("dtec",d_t_ec) 1856 1857 CALL send_xios_field("SWnet",swnet) 1858 CALL send_xios_field("LWnet",lwnet) 1859 CALL send_xios_field("fluxvdf",fluxt) 1860 CALL send_xios_field("fluxdyn",flux_dyn) 1861 CALL send_xios_field("fluxajs",flux_ajs) 1862 CALL send_xios_field("fluxec",flux_ec) 1863 1864 c plusieurs traceurs !!!outputs in [vmr] 1865 IF (iflag_trac.eq.1) THEN 1866 DO iq=1,nqmax 1867 CALL send_xios_field(tname(iq),qx(:,:,iq)*mmean(:,:)/M_tr(iq)) 1868 ENDDO 1869 ENDIF 1870 1871 IF (callthermos .and. ok_chem) THEN 1872 CALL send_xios_field("d_qmoldifCO2",d_q_moldif(:,:,i_co2)) 1873 CALL send_xios_field("d_qmoldifO3p",d_q_moldif(:,:,i_o)) 1874 CALL send_xios_field("d_qmoldifN2",d_q_moldif(:,:,i_n2)) 1875 ENDIF 1813 1876 1814 1877 #endif -
trunk/LMDZ.VENUS/libf/phyvenus/phytrac_chimie.F
r1661 r1665 95 95 c============================================================= 96 96 97 IF (reinit_trac ) THEN97 IF (reinit_trac.and.ok_chem) THEN 98 98 PRINT*,'REINIT MIXING RATIO TRACEURS' 99 99 … … 119 119 trac(:,:,:)=1.0E-30 120 120 121 c !!! AVEC NUAGE !!! 122 trac(:,1:22,i_ocs)=3.E-6 123 trac(:,:,i_hcl)=0.4E-6 124 trac(:,1:22,i_so2)=10.E-6 125 trac(:,1:22,i_h2o)=30.0E-6 121 c !!! Verif traceurs !!! 122 if ( (i_ocs.ne.0).and.(i_co.ne.0).and.(i_hcl.ne.0) 123 & .and.(i_so2.ne.0).and.(i_h2o.ne.0).and.(i_n2.ne.0) 124 & .and.(i_co2.ne.0) ) then 125 trac(:,1:22,i_ocs)=3.E-6 126 trac(:,1:22,i_co)=25.E-6 127 trac(:,:,i_hcl)=0.4E-6 128 trac(:,1:22,i_so2)=9.E-6 129 trac(:,1:22,i_h2o)=30.0E-6 126 130 127 131 c remettre tous les traceurs du start => trac(:,:,:)=trac_sav(:,:,:) … … 129 133 c N2 n est pas encore une espece chimique du modele chimique 130 134 c traceur passif pour la chimie-transport 131 trac(:,:,i_n2)=0.35d-1135 trac(:,:,i_n2)=0.35d-1 132 136 133 137 !!!! GG: Initialization CO2 = 1 - qtot 134 138 !! It assures that vmr_tot = 1 135 139 c On a donc le CO2 qui est le restant d atmosphere Venus 136 140 trac_sum(:,:)=0.0 137 141 c SEULEMENT LES GAZ 138 142 DO iq=2,nqmax-nmicro … … 141 145 142 146 trac(:,:,i_co2)= 1-trac_sum(:,:) 147 else 148 write(*,*) "Réinitialisation des traceurs avec chimie " 149 write(*,*) "IL manque un traceur pour les options choisies" 150 stop 151 endif ! verif traceurs 143 152 144 153 c=============================================================
Note: See TracChangeset
for help on using the changeset viewer.