Changeset 829 for trunk/LMDZ.MARS/util/zrecast.F90
- Timestamp:
- Nov 6, 2012, 12:02:45 PM (12 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/util/zrecast.F90
r788 r829 1 1 module planet_const 2 ! planetary constants (set via init_planet_const routine) 3 implicit none 4 real :: a0 ! Mean planetary radius (a0=3396.E3 for Mars) 5 real :: g0 ! gravity at a0 (Mars: g0=3.7257964 ; Lemoine et al. 2001) 6 real :: Rmean ! reduced gaz constant (Mars: R=191) 7 end module planet_const 2 8 3 9 program zrecast … … 51 57 ! EM 03/2011 : Add possibility to have outputs as distance to center 52 58 ! of planet 59 ! EM 11/2012 : Adapted so it can be used on "generic" model outputs; planet 60 ! constants (radius, R, etc.) are now read from file 53 61 implicit none 54 62 … … 154 162 endif 155 163 164 ! load planet constants (radius, gravity, ...) 165 166 call init_planet_const(infid) 167 156 168 !=============================================================================== 157 169 ! 1.2 Get # and names of variables in input file … … 1743 1755 end program 1744 1756 1757 !=============================================================================== 1758 1759 subroutine init_planet_const(infid) 1760 ! initialize planetary constants using the "controle" array in the file 1761 ! if "cointrole" array not found in file, look for it in "diagfi.nc" 1762 use planet_const 1763 use netcdf 1764 implicit none 1765 1766 integer,intent(in) :: infid ! input file ID 1767 1768 ! local variables 1769 character(len=100) :: varname 1770 integer :: varid ! variable ID 1771 integer :: status 1772 integer :: infid2 ! for an alternate input file 1773 real :: controle(100) ! to store "controle" array 1774 1775 ! look for "controle" in input file 1776 infid2=infid ! initialization 1777 varname="controle" 1778 status=nf90_inq_varid(infid,trim(varname),varid) 1779 if (status.ne.nf90_noerr) then 1780 write(*,*) "init_planet_const: Failed to find ",trim(varname) 1781 write(*,*) " looking for it in file diagfi.nc" 1782 status=nf90_open("diagfi.nc",NF90_NOWRITE,infid2) 1783 if (status.ne.nf90_noerr) then 1784 write(*,*) " no diafi.nc file ... looking for diagfi1.nc" 1785 status=nf90_open("diagfi1.nc",NF90_NOWRITE,infid2) 1786 if (status.ne.nf90_noerr) then 1787 write(*,*) "might as well stop here..." 1788 stop 1789 endif 1790 endif 1791 status=nf90_inq_varid(infid2,trim(varname),varid) 1792 if (status.ne.nf90_noerr) then 1793 write(*,*) " Failed to find ",trim(varname)," in file!!" 1794 stop 1795 endif 1796 write(*,*) "OK, found ",trim(varname) 1797 endif 1798 status=nf90_get_var(infid2,varid,controle) 1799 if (status.ne.nf90_noerr) then 1800 write(*,*) "init_planet_const: Failed to load ",trim(varname) 1801 stop 1802 endif 1803 1804 ! now that we have "controle"; extract relevent informations 1805 a0=controle(5) ! = radius of planet (m) 1806 g0=controle(7) ! = gravity (m.s-2) at a0 1807 Rmean=1000.*8.314511/controle(8) ! controle(8)=mugaz = molar mass (g.mol-1) of atm. 1808 1809 end subroutine init_planet_const 1745 1810 1746 1811 !=============================================================================== … … 1752 1817 ! surface altitudes" corresponding to GCM atmospheric levels 1753 1818 !============================================================================== 1819 use planet_const, only : a0,g0 1754 1820 implicit none 1755 1821 !============================================================================== … … 1785 1851 1786 1852 ! Parameters needed to integrate hydrostatic equation: 1787 real,parameter :: g0=3.72579641853 !real,parameter :: g0=3.7257964 1788 1854 !g0: exact mean gravity at radius=3396.km (Lemoine et al. 2001) 1789 real,parameter :: a0=3396.E31855 !real,parameter :: a0=3396.E3 1790 1856 !a0: 'mean' Mars radius=3396.km 1791 1857 real :: gz … … 1871 1937 !=============================================================================== 1872 1938 1873 !#include"build_gcm_alt.F90"1874 1939 subroutine build_gcm_za(lonlength,latlength,altlength,timelength, & 1875 1940 phis,ps,press,temp,rho,za_gcm) … … 1878 1943 ! altitudes" corresponding to GCM atmospheric levels 1879 1944 !============================================================================== 1945 use planet_const, only : a0,g0 1880 1946 implicit none 1881 1947 !============================================================================== … … 1911 1977 1912 1978 ! Parameters needed to integrate hydrostatic equation: 1913 real,parameter :: g0=3.72579641979 !real,parameter :: g0=3.7257964 1914 1980 !g0: exact mean gravity at radius=3396.km (Lemoine et al. 2001) 1915 real,parameter :: a0=3396.E31981 !real,parameter :: a0=3396.E3 1916 1982 !a0: 'mean' Mars radius=3396.km 1917 1983 real :: gz … … 2062 2128 ! surface altitudes" corresponding to GCM atmospheric levels 2063 2129 !============================================================================== 2130 use planet_const, only : a0,g0,Rmean 2064 2131 implicit none 2065 2132 !============================================================================== … … 2086 2153 !=============================================================================== 2087 2154 real,dimension(:),allocatable :: sigma ! GCM sigma levels 2088 real,parameter :: R=191 ! molecular gas constant2155 !real,parameter :: R=191 ! molecular gas constant 2089 2156 !real,dimension(:,:,:,:),allocatable :: zlocal ! local above surface altitude 2090 2157 real :: Tmean ! "mean" value of temperature for a given level … … 2092 2159 2093 2160 ! Parameters needed to integrate hydrostatic equation: 2094 real,parameter :: g0=3.72579642161 !real,parameter :: g0=3.7257964 2095 2162 !g0: exact mean gravity at radius=3396.km (Lemoine et al. 2001) 2096 real,parameter :: a0=3396.E32163 !real,parameter :: a0=3396.E3 2097 2164 !a0: 'mean' Mars radius=3396.km 2098 2165 real :: gz … … 2114 2181 sigma(1)=press(iloop,jloop,1,tloop)/ps(iloop,jloop,tloop) 2115 2182 zs_gcm(iloop,jloop,1,tloop)=-log(sigma(1))* & 2116 R*temp(iloop,jloop,1,tloop)/g02183 Rmean*temp(iloop,jloop,1,tloop)/g0 2117 2184 ! za_gcm(iloop,jloop,1,tloop)=zlocal(iloop,jloop,1,tloop)+ & 2118 2185 ! phis(iloop,jloop)/g0 … … 2139 2206 ! compute zs_gcm(iloop,jloop,kloop,tloop) 2140 2207 zs_gcm(iloop,jloop,kloop,tloop)=zs_gcm(iloop,jloop,kloop-1,tloop)- & 2141 log(sigma(kloop)/sigma(kloop-1))*R *Tmean/gz2208 log(sigma(kloop)/sigma(kloop-1))*Rmean*Tmean/gz 2142 2209 2143 2210 ! compute za_gcm(kloop) … … 2157 2224 !=============================================================================== 2158 2225 2159 !#include"crude_gcm_alt.F90"2160 2226 subroutine crude_gcm_za(lonlength,latlength,altlength,timelength, & 2161 2227 phis,ps,press,temp,za_gcm) … … 2164 2230 ! altitudes" corresponding to GCM atmospheric levels 2165 2231 !============================================================================== 2232 use planet_const, only : a0,g0,Rmean 2166 2233 implicit none 2167 2234 !============================================================================== … … 2188 2255 !=============================================================================== 2189 2256 real,dimension(:),allocatable :: sigma ! GCM sigma levels 2190 real,parameter :: R=191 ! molecular gas constant2257 !real,parameter :: R=191 ! molecular gas constant 2191 2258 real,dimension(:,:,:,:),allocatable :: zlocal ! local above surface altitude 2192 2259 real :: Tmean ! "mean" value of temperature for a given level … … 2194 2261 2195 2262 ! Parameters needed to integrate hydrostatic equation: 2196 real,parameter :: g0=3.72579642263 !real,parameter :: g0=3.7257964 2197 2264 !g0: exact mean gravity at radius=3396.km (Lemoine et al. 2001) 2198 real,parameter :: a0=3396.E32265 !real,parameter :: a0=3396.E3 2199 2266 !a0: 'mean' Mars radius=3396.km 2200 2267 real :: gz … … 2216 2283 sigma(1)=press(iloop,jloop,1,tloop)/ps(iloop,jloop,tloop) 2217 2284 zlocal(iloop,jloop,1,tloop)=-log(sigma(1))* & 2218 R*temp(iloop,jloop,1,tloop)/g02285 Rmean*temp(iloop,jloop,1,tloop)/g0 2219 2286 za_gcm(iloop,jloop,1,tloop)=zlocal(iloop,jloop,1,tloop)+ & 2220 2287 phis(iloop,jloop)/g0 … … 2239 2306 ! compute zlocal(kloop) 2240 2307 zlocal(iloop,jloop,kloop,tloop)=zlocal(iloop,jloop,kloop-1,tloop)- & 2241 log(sigma(kloop)/sigma(kloop-1))*R *Tmean/gz2308 log(sigma(kloop)/sigma(kloop-1))*Rmean*Tmean/gz 2242 2309 2243 2310 ! compute za_gcm(kloop) … … 2332 2399 !=============================================================================== 2333 2400 2334 !#include"za_coord_interp.F90"2335 2401 subroutine z_coord_interp(lonlen,latlen,altlen,tlen,newaltlen, & 2336 2402 missing,z_gcm,gcmdata,flag,z_new,newdata) … … 2462 2528 ! if flag=0, and extrapolated (exponentially) if flag=1. 2463 2529 !============================================================================== 2530 use planet_const, only : a0,g0,Rmean 2464 2531 implicit none 2465 2532 !============================================================================== … … 2503 2570 real :: x,y ! temporary variables 2504 2571 integer :: iloop,jloop,kloop,tloop 2505 real,parameter :: g0=3.72579642572 !real,parameter :: g0=3.7257964 2506 2573 !g0: exact mean gravity at radius=3396.km (Lemoine et al. 2001) 2507 real,parameter :: a0=3396.E32574 !real,parameter :: a0=3396.E3 2508 2575 !a0: 'mean' Mars radius=3396.km 2509 real,parameter :: Rmean=191 ! molecular gas constant2576 !real,parameter :: Rmean=191 ! molecular gas constant 2510 2577 real :: gz 2511 2578 ! gz: gravitational acceleration at a given (above areoid) altitude
Note: See TracChangeset
for help on using the changeset viewer.