Changeset 126 for trunk/mesoscale/LMD_LES_MARS/modif_mars
- Timestamp:
- May 23, 2011, 8:40:24 PM (14 years ago)
- Location:
- trunk/mesoscale/LMD_LES_MARS/modif_mars
- Files:
-
- 1 added
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/mesoscale/LMD_LES_MARS/modif_mars/Registry.EM
r96 r126 118 118 state real QH2O_ICE ikjftb scalar 1 - i01rhusdf=(bdy_interp:dt) "QH2O_ICE" "Water ice mixing ratio" "kg kg-1" 119 119 state real QDUST ikjftb scalar 1 - i01rhusdf=(bdy_interp:dt) "QDUST" "Dust mixing ratio" "kg kg-1" 120 state real qtrac1 ikjftb scalar 1 - i01rhusdf=(bdy_interp:dt) "qtrac1" "Decaying tracer 1" "kg kg-1" 120 121 #### 121 122 #### … … 1488 1489 package water mars==1 - moist:qv;scalar:qh2o,qh2o_ice 1489 1490 package dust mars==2 - moist:qv;scalar:qdust 1491 package radioac mars==20 - scalar:qtrac1 1490 1492 ##### MARS OPTIONS 1491 1493 ##### MARS OPTIONS -
trunk/mesoscale/LMD_LES_MARS/modif_mars/module_initialize_les.F
r92 r126 798 798 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n) 799 799 800 ! input therm data (element 0 is the ground so it's n+1 but n is 1000 anyway so...) 801 802 real r_therm(n),cp_therm(n),p_therm(n),rho_therm(n),t_therm(n) 803 800 804 ! diagnostics 801 805 … … 816 820 call read_sounding( p_surf, th_surf, qv_surf, & 817 821 h_input, th_input, qv_input, u_input, v_input,n, nl, debug ) 822 823 ! and the therm : 824 825 call read_therm(r_therm,cp_therm,p_therm,rho_therm,t_therm,n) 826 827 ! To use r/cp as defined above, one has to recompute teta from T (default MCD computes 828 ! teta for a variable r/cp) 829 830 do k=1,nl 831 th_input(k) = t_therm(k)*(p1000mb/p_therm(k))**(r/cp) 832 enddo 833 th_surf = t_therm(1)*(p1000mb/p_therm(1))**(r/cp) 834 ! ----- 818 835 819 836 if(dry) then … … 957 974 end subroutine read_sounding 958 975 976 subroutine read_therm(r,cp,p,rho,t,n) 977 implicit none 978 integer n 979 real r(n),cp(n),p(n),rho(n),t(n) 980 logical end_of_file 981 982 integer k 983 984 ! first element is the surface 985 986 open(unit=11,file='input_therm',form='formatted',status='old') 987 rewind(11) 988 end_of_file = .false. 989 k = 0 990 do while (.not. end_of_file) 991 992 read(11,*,end=101) r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1) 993 write(*,*) k, r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1) 994 k = k+1 995 go to 112 996 101 end_of_file = .true. 997 112 continue 998 enddo 999 1000 close(unit=11,status = 'keep') 1001 1002 end subroutine read_therm 1003 959 1004 END MODULE module_initialize_ideal
Note: See TracChangeset
for help on using the changeset viewer.