Changeset 1776 in lmdz_wrf
- Timestamp:
- Feb 15, 2018, 7:18:41 PM (7 years ago)
- Location:
- trunk/tools
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/tools/diag_tools.py
r1773 r1776 910 910 return zmla, zmladims, zmlavdims 911 911 912 913 def Forcompute_zwind(ua, va, zsl, uas, vas, hgt, sina, cosa, zval, dimns, dimvns): 914 """ Function to compute the wind at a given height following the power law method 915 Forcompute_zwind(ua, va, zsl, uas, vas, hgt, sina, cosa, zval, dimns, dimvns) 916 [ua]= x-component of unstaggered 3D wind (assuming [[t],z,y,x]) [ms-1] 917 [va]= y-component of unstaggered 3D wind (assuming [[t],z,y,x]) [ms-1] 918 [zsl]= height above sea level [m] 919 [uas]= x-component of unstaggered 10 m wind (assuming [[t],z,y,x]) [ms-1] 920 [vas]= y-component of unstaggered 10 m wind (assuming [[t],z,y,x]) [ms-1] 921 [hgt]= topographical height (assuming [m] 922 [sina]= local sine of map rotation [1.] 923 [cosa]= local cosine of map rotation [1.] 924 [zval]= desired height for winds [m] 925 [dimns]= list of the name of the dimensions of [ua] 926 [dimvns]= list of the name of the variables with the values of the 927 dimensions of [ua] 928 """ 929 fname = 'Forcompute_zwind' 930 931 vardims = dimns[:] 932 varvdims = dimvns[:] 933 934 if len(ua.shape) == 4: 935 var1= np.zeros((ua.shape[0],ua.shape[2],ua.shape[3]), dtype=np.float) 936 var2= np.zeros((ua.shape[0],ua.shape[2],ua.shape[3]), dtype=np.float) 937 938 dx = ua.shape[3] 939 dy = ua.shape[2] 940 dz = ua.shape[1] 941 dt = ua.shape[0] 942 vardims.pop(1) 943 varvdims.pop(1) 944 945 pvar1, pvar2= fdin.module_fordiagnostics.compute_zwind4d(ua=ua.transpose(), \ 946 va=va[:].transpose(), zsl=zsl[:].transpose(), uas=uas.transpose(), \ 947 vas=vas.transpose(), hgt=hgt.transpose(), sina=sina.transpose(), \ 948 cosa=cosa.transpose(), zextrap=zval, d1=dx, d2=dy, d3=dz, d4=dt) 949 var1 = pvar1.transpose() 950 var2 = pvar2.transpose() 951 else: 952 print errormsg 953 print ' ' + fname + ': rank', len(ua.shape), 'not ready !!' 954 print ' it only computes 4D [t,z,y,x] rank values' 955 quit(-1) 956 957 return var1, var2, vardims, varvdims 912 958 913 959 def compute_OMEGAw(omega, p, t, dimns, dimvns): -
trunk/tools/diagnostics.inf
r1773 r1776 34 34 ua, WRFua, U@V@SINALPHA@COSALPHA 35 35 va, WRFva, U@V@SINALPHA@COSALPHA 36 uavaz, WRFzwind, U@V@WRFgeop@U10@V10@HGT@SINALPHA@COSALPHA@z=100. 36 37 wa, OMEGAw, vitw@pres@temp 37 38 wds, TSwds, u@v … … 45 46 zg, WRFght, PH@PHB 46 47 zmla, WRFzmlagen, T@QVAPOR@WRFgeop@HGT 47 -
trunk/tools/diagnostics.py
r1773 r1776 81 81 'WRFmrso', 'WRFp', \ 82 82 'WRFpsl_ptarget', 'WRFrvors', 'WRFslw', 'ws', 'wds', 'wss', 'WRFheight', \ 83 'WRFheightrel', 'WRFua', 'WRFva' ]83 'WRFheightrel', 'WRFua', 'WRFva', 'WRzwind'] 84 84 85 85 methods = ['accum', 'deaccum'] … … 344 344 gen.searchInlist(NONcheckingvars, depv) and \ 345 345 not gen.searchInlist(methods, depv) and not depvars[0] == 'deaccum' \ 346 and not depvars[0] == 'accum' :346 and not depvars[0] == 'accum' and not depv[0:2] == 'z=': 347 347 print errormsg 348 348 print ' ' + main + ": file '" + opts.ncfile + \ … … 1123 1123 ncvar.insert_variable(ncobj, 'zmla', diagout, diagoutd, diagoutvd, newnc) 1124 1124 1125 # WRFzwind wind extrapolation at a given hieght computation from WRF U, V, WRFgeop, 1126 # U10, V10, HGT, SINALPHA, COSALPHA, z=[zval] 1127 elif diagn == 'WRFzwind': 1128 var0 = ncobj.variables[depvars[0]][:] 1129 var1 = ncobj.variables[depvars[1]][:] 1130 dimz = var0.shape[1] 1131 var2 = WRFgeop[:,1:dimz+1,:,:]/9.8 1132 var3 = ncobj.variables[depvars[3]][:] 1133 var4 = ncobj.variables[depvars[4]][:] 1134 var5 = ncobj.variables[depvars[5]][0,:,:] 1135 var6 = ncobj.variables[depvars[6]][0,:,:] 1136 var7 = ncobj.variables[depvars[7]][0,:,:] 1137 var8 = np.float(depvars[8].split('=')[1]) 1138 1139 # un-staggering 3D winds 1140 unstgdims = [var0.shape[0], var0.shape[1], var0.shape[2], var0.shape[3]-1] 1141 va = np.zeros(tuple(unstgdims), dtype=np.float) 1142 unvar0 = np.zeros(tuple(unstgdims), dtype=np.float) 1143 unvar1 = np.zeros(tuple(unstgdims), dtype=np.float) 1144 unvar0 = 0.5*(var0[:,:,:,0:var0.shape[3]-1] + var0[:,:,:,1:var0.shape[3]]) 1145 unvar1 = 0.5*(var1[:,:,0:var1.shape[2]-1,:] + var1[:,:,1:var1.shape[2],:]) 1146 1147 diagout1, diagout2, diagoutd, diagoutvd = diag.Forcompute_zwind(unvar0, \ 1148 unvar1, var2, var3, var4, var5, var6, var7, var8, dnames, dvnames) 1149 1150 # Removing the nonChecking variable-dimensions from the initial list 1151 varsadd = [] 1152 for nonvd in NONchkvardims: 1153 if gen.searchInlist(dvnames,nonvd): diagoutvd.remove(nonvd) 1154 varsadd.append(nonvd) 1155 1156 ncvar.insert_variable(ncobj, 'uaz', diagout1, diagoutd, diagoutvd, newnc) 1157 ncvar.insert_variable(ncobj, 'vaz', diagout2, diagoutd, diagoutvd, newnc) 1158 1125 1159 else: 1126 1160 print errormsg -
trunk/tools/module_ForDiagnostics.f90
r1773 r1776 28 28 ! compute_zint4D: Subroutine to vertically integrate a 4D variable in any vertical coordinates 29 29 ! compute_zmla_generic4D: Subroutine to compute pbl-height following a generic method 30 ! compute_zwind4D: Subroutine to compute extrapolate the wind at a given height following the 'power law' methodology 30 31 31 32 !!! … … 642 643 END SUBROUTINE compute_zmla_generic4D 643 644 645 SUBROUTINE compute_zwind4D(ua, va, zsl, uas, vas, hgt, sina, cosa, zextrap, uaz, vaz, d1, d2, d3, d4) 646 ! Subroutine to compute extrapolate the wind at a given height following the 'power law' methodology 647 648 IMPLICIT NONE 649 650 INTEGER, INTENT(in) :: d1, d2, d3, d4 651 REAL(r_k), DIMENSION(d1,d2,d3,d4), INTENT(in) :: ua, va, zsl 652 REAL(r_k), DIMENSION(d1,d2,d4), INTENT(in) :: uas, vas 653 REAL(r_k), DIMENSION(d1,d2), INTENT(in) :: hgt, sina, cosa 654 REAL(r_k), INTENT(in) :: zextrap 655 REAL(r_k), DIMENSION(d1,d2,d4), INTENT(out) :: uaz, vaz 656 657 ! Local 658 INTEGER :: i, j, it 659 660 !!!!!!! Variables 661 ! tpot: potential air temperature [K] 662 ! qratio: water vapour mixing ratio [kgkg-1] 663 ! z: height above sea level [m] 664 ! hgt: terrain height [m] 665 ! sina, cosa: local sine and cosine of map rotation [1.] 666 ! zmla3D: boundary layer height from surface [m] 667 668 fname = 'compute_zwind4D' 669 670 DO i=1, d1 671 DO j=1, d2 672 DO it=1, d4 673 CALL var_zwind(d3, ua(i,j,:,it), va(i,j,:,it), zsl(i,j,:,it), uas(i,j,it), vas(i,j,it), & 674 hgt(i,j), sina(i,j), cosa(i,j), zextrap, uaz(i,j,it), vaz(i,j,it)) 675 END DO 676 END DO 677 END DO 678 679 RETURN 680 681 END SUBROUTINE compute_zwind4D 644 682 645 683 END MODULE module_ForDiagnostics -
trunk/tools/module_ForDiagnosticsVars.f90
r1773 r1776 29 29 ! var_clt: total cloudiness [0,1] 30 30 ! var_zmla_generic: Subroutine to compute pbl-height following a generic method 31 ! var_zwind: Subroutine to extrapolate the wind at a given height following the 'power law' methodology 31 32 ! VirtualTemp1D: Function for 1D calculation of virtual temperaure 32 33 ! VirtualTemperature: WRF's AFWA method to compute virtual temperature … … 1026 1027 END SUBROUTINE var_zmla_generic 1027 1028 1029 SUBROUTINE var_zwind(d1, u, v, z, u10, v10, sa, ca, topo, newz, unewz, vnewz) 1030 ! Subroutine to extrapolate the wind at a given height following the 'power law' methodology 1031 ! wss[newz] = wss[z1]*(newz/z1)**alpha 1032 ! alpha = (ln(wss[z2])-ln(wss[z1]))/(ln(z2)-ln(z1)) 1033 ! AFTER: Phd Thesis: 1034 ! Benedicte Jourdier. Ressource eolienne en France metropolitaine : methodes dâevaluation du 1035 ! potentiel, variabilite et tendances. Climatologie. Ecole Doctorale Polytechnique, 2015. French 1036 ! 1037 IMPLICIT NONE 1038 1039 INTEGER, INTENT(in) :: d1 1040 REAL(r_k), DIMENSION(d1), INTENT(in) :: u,v,z 1041 REAL(r_k), INTENT(in) :: u10, v10, topo, sa, ca, newz 1042 REAL(r_k), INTENT(out) :: unewz, vnewz 1043 1044 ! Local 1045 INTEGER :: inear 1046 REAL(r_k) :: zaground 1047 REAL(r_k), DIMENSION(2) :: v1, v2, zz, alpha, uvnewz 1048 1049 !!!!!!! Variables 1050 ! u,v: vertical wind components [ms-1] 1051 ! z: height above surface [m] 1052 ! u10,v10: 10-m wind components [ms-1] 1053 ! topo: topographical height [m] 1054 ! sa, ca: local sine and cosine of map rotation [1.] 1055 ! newz: desired height above grpund of extrapolation 1056 ! unewz,vnewz: Wind compoonents at the given height [ms-1] 1057 1058 fname = 'var_zwind' 1059 1060 PRINT *,' ilev zaground newz z[ilev+1] z[ilev+2] _______' 1061 IF (z(1) < newz ) THEN 1062 DO inear = 1,d1-2 1063 zaground = z(inear+2) 1064 PRINT *, inear, z(inear), newz, z(inear+1), z(inear+2) 1065 IF ( zaground >= newz) EXIT 1066 END DO 1067 ELSE 1068 PRINT *, 1, z(1), newz, z(2), z(3), ' z(1) > newz' 1069 inear = d1 - 2 1070 END IF 1071 1072 IF (inear == d1-2) THEN 1073 ! No vertical pair of levels is below newz, using 10m wind as first value and the first level as the second 1074 v1(1) = u10 1075 v1(2) = v10 1076 v2(1) = u(1) 1077 v2(2) = v(1) 1078 zz(1) = 10. 1079 zz(2) = z(1) 1080 ELSE 1081 v1(1) = u(inear) 1082 v1(2) = v(inear) 1083 v2(1) = u(inear+1) 1084 v2(2) = v(inear+1) 1085 zz(1) = z(inear) - topo 1086 zz(2) = z(inear+1) - topo 1087 END IF 1088 1089 ! Computing for each component 1090 alpha = (LOG(ABS(v2))-LOG(ABS(v1)))/(LOG(zz(2))-LOG(zz(1))) 1091 PRINT *,' Computing with v1:', v1, ' ms-1 v2:', v2, ' ms-1' 1092 PRINT *,' z1:', zz(1), 'm z2:', zz(2), ' m' 1093 PRINT *,' alhpa u:', alpha(1), ' alpha 2:', alpha(2) 1094 1095 uvnewz = v1*(newz/zz(1))**alpha 1096 ! Earth-rotation 1097 unewz = uvnewz(1)*ca - uvnewz(2)*sa 1098 vnewz = uvnewz(1)*sa + uvnewz(2)*ca 1099 1100 PRINT *,' result vz:', uvnewz 1101 1102 !STOP 1103 1104 RETURN 1105 1106 END SUBROUTINE var_zwind 1107 1028 1108 END MODULE module_ForDiagnosticsVars -
trunk/tools/variables_values.dat
r1762 r1776 462 462 U10, uas, eastward_wind, -30., 30., eastward|10m|wind, ms-1, seismic 463 463 Vent zonal 10m, uas, eastward_wind, -30., 30., eastward|10m|wind, ms-1, seismic 464 uaz, uaz, eastward_wind, -30., 30., eastward|wind|at|z|value, ms-1, seismic 464 465 va, va, northward_wind, -30., 30., northward|wind, ms-1, seismic 465 466 vitv, va, northward_wind, -30., 30., northward|wind, ms-1, seismic … … 472 473 V10, vas, northward_wind, -30., 30., northward|10m|wind, ms-1, seismic 473 474 Vent meridien 10m, vas, northward_wind, -30., 30., northward|10m|wind, ms-1, seismic 475 vaz, vaz, northward_wind, -30., 30., northward|wind|at|z|value, ms-1, seismic 474 476 veget, veget, vegetation_type, 1, 13, type|of|plant|functional|type, 1, Reds 475 477 wakedeltaq, wakedeltaq, wake_delta_vapor, -0.003, 0.003, wake|delta|mixing|ratio, kgkg-1s-1, seismic
Note: See TracChangeset
for help on using the changeset viewer.