Changeset 1378 for trunk/MESOSCALE/LMD_MM_MARS/SRC
- Timestamp:
- Feb 18, 2015, 3:36:22 PM (10 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/MESOSCALE/LMD_MM_MARS/SRC/LES/modif_mars/module_initialize_les.F
r1236 r1378 115 115 REAL :: lon_input, lat_input, alt_input, tsurf_input 116 116 ! for mode 3 117 REAL, DIMENSION(nl_max) :: profdustq,profdustn ,dust_p_level117 REAL, DIMENSION(nl_max) :: profdustq,profdustn 118 118 !!MARS 119 119 … … 757 757 ENDDO 758 758 759 !!!!! MARS INITIALIZE MODE 3 FOR LES 760 !! RESTART=config_flags%restart -> this routine is not used in restart 761 762 IF ( ( config_flags%mars == 3 ) ) then 763 764 write (*,*) 'MARS MODE 3 INITIALIZATION, READING INPUT_DUST (module_initialize_les.F)' 765 ! load a profile of dust (same for all points) 766 call read_dust(profdustq,profdustn,dust_p_level,nl_in) 767 p_level = grid%znu(1)*(pd_surf - grid%p_top) + grid%p_top 768 IF (dust_p_level(1) .lt. p_level) then !input profile needs rescaling to avoid a plateau. This happens when you use different sources to initialize the pressure of the LES and the input_dust. Typicaly: you use the MCD for input_therm and different runs for input_dust, that dont have the same surface pressure ! This trick is ok because we dont want to initialize with a very precise profile of dust, just a realistic one. 769 dust_p_level = dust_p_level * p_level/dust_p_level(1) 770 ENDIF 771 DO k=1,kte-1 759 !!!!! MARS 760 761 ! interpolate water vapor 762 if ( ( config_flags%mars == 1 ) & 763 .OR. ( config_flags%mars == 11 ) & 764 .OR. ( config_flags%mars == 12 ) ) then 765 print *, '**** INTERPOLATE HV **** RANK 2 in SCALAR' 766 DO k=1,kte-1 772 767 p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top 773 print *, '**** INTERPOLATE DUSTQ **** RANK 2 in SCALAR' 774 scalar(its:ite,k,jts:jte,2) = interp_0( profdustq, dust_p_level, p_level, nl_in ) 775 print *, '**** INTERPOLATE DUSTN **** RANK 3 in SCALAR' 776 scalar(its:ite,k,jts:jte,3) = interp_0( profdustn, dust_p_level, p_level, nl_in ) 777 ENDDO 778 ENDIF 768 scalar(its:ite,k,jts:jte,2) = interp_0( qv, pd_in, p_level, nl_in ) 769 scalar(its:ite,k,jts:jte,3) = 0. 770 !! water ice is set to 0 (was put into water vapor when building prof from MCD) 771 ENDDO 772 print *, "WATER VAPOR",scalar(its,:,jts,2) 773 endif 774 775 ! interpolate qdust 776 if ( ( config_flags%mars == 11 ) & 777 .OR. ( config_flags%mars == 12 ) ) then 778 call read_dust(profdustq,profdustn,nl_in) 779 print *, '**** INTERPOLATE DUSTQ **** RANK 4 in SCALAR' 780 print *, '**** INTERPOLATE DUSTN **** RANK 5 in SCALAR' 781 DO k=1,kte-1 782 p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top 783 scalar(its:ite,k,jts:jte,4) = interp_0( profdustq, pd_in, p_level, nl_in ) 784 scalar(its:ite,k,jts:jte,5) = interp_0( profdustn, pd_in, p_level, nl_in ) 785 ENDDO 786 print *, "DUST Q", scalar(its,:,jts,4) 787 print *, "DUST N", scalar(its,:,jts,5) 788 endif 789 790 if ( config_flags%mars == 12 ) then 791 scalar(its:ite,1:kte-1,jts:jte,6) = 0. 792 scalar(its:ite,1:kte-1,jts:jte,7) = 0. 793 endif 779 794 780 795 !!!!! MARS … … 866 881 ! ----- 867 882 868 if(dry) then869 do k=1,nl870 qv_input(k) = 0.871 enddo872 endif883 !if(dry) then 884 ! do k=1,nl 885 ! qv_input(k) = 0. 886 ! enddo 887 !endif 873 888 874 889 if(debug) write(6,*) ' number of input levels = ',nl … … 884 899 885 900 do k=1,nl 901 !!!!!!!!!!!!!! MARS 902 !! from mol/mol to kg/kg 903 qv_input(k) = qv_input(k)*18./mwdry 886 904 qv_input(k) = 0.001*qv_input(k) 887 905 enddo … … 889 907 p_surf = 100.*p_surf ! convert to pascals 890 908 qvf = 1. + rvovrd*qv_input(1) 909 !!MARS 910 qvf = 1. 891 911 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm)) 892 912 pi_surf = (p_surf/p1000mb)**(r/cp) … … 904 924 qvf = 1. + rvovrd*qv_input(1) 905 925 qvf1 = 1. + qv_input(1) 926 !!MARS 927 qvf = 1. 928 qvf1 = 1. 906 929 rho_input(1) = rho_surf 907 930 dz = h_input(1) … … 921 944 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k))) 922 945 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here 923 946 !!MARS 947 qvf = 1. 948 qvf1 = 1. 949 !!MARS 924 950 do it=1,10 925 951 pm_input(k) = pm_input(k-1) & … … 1035 1061 1036 1062 !!cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc 1037 subroutine read_dust(pdustq,pdustn, dust_p_level,n)1063 subroutine read_dust(pdustq,pdustn,n) 1038 1064 implicit none 1039 1065 integer n 1040 real pdustq(n),pdustn(n) ,dust_p_level(n)1066 real pdustq(n),pdustn(n) 1041 1067 logical end_of_file 1042 1068 … … 1051 1077 do while (.not. end_of_file) 1052 1078 1053 read(11,*,end=102) pdustq(k+1),pdustn(k+1) , dust_p_level(k+1)1054 write(*,*) k, pdustq(k+1), pdustn(k+1), dust_p_level(k+1)1079 read(11,*,end=102) pdustq(k+1),pdustn(k+1) 1080 write(*,*) k,pdustq(k+1),pdustn(k+1) 1055 1081 k = k+1 1056 1082 go to 113 … … 1064 1090 !!cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc 1065 1091 1066 1067 1092 END MODULE module_initialize_ideal
Note: See TracChangeset
for help on using the changeset viewer.