Changeset 2528 for trunk/LMDZ.MARS/util/zrecast.F90
- Timestamp:
- May 27, 2021, 9:41:55 AM (3 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/util/zrecast.F90
r2432 r2528 62 62 ! TN 10/2013 : Read and write controle field, if available 63 63 ! FF 10/2020 : Force interpolation in log space for variable starting by "rho" 64 ! FF 05/2021 : Force interpolation in log space for variable starting by "num_" 65 ! FF 05/2021 : "phisinit" added to the putput file 66 ! FF 05/2021 : Conversion to pressure coordinate possible even if temperature 67 ! not available 64 68 ! 65 69 implicit none … … 84 88 integer outfid ! NetCDF output file ID 85 89 integer lon_dimid,lat_dimid,alt_dimid,time_dimid,ctl_dimid ! NetCDF dimension IDs 86 integer lon_varid,lat_varid,alt_varid,time_varid,ctl_varid 90 integer lon_varid,lat_varid,alt_varid,time_varid,ctl_varid,phis_varid 87 91 integer gcm_layers_dimid ! NetCDF dimension ID for # of layers in GCM 88 92 integer sigma_varid,aps_varid,bps_varid … … 94 98 integer,dimension(4) :: datashape ! shape of 4D datasets 95 99 integer,dimension(3) :: surfdatashape ! shape of 3D (surface+time) datasets 100 integer,dimension(2) :: phisdatashape ! shape of 2D (surface) datasets 96 101 97 102 real :: miss_val= 1.E+20 ! special "missing value" to specify missing data … … 140 145 real,dimension(:),allocatable :: zradius ! distance to center of planet 141 146 logical :: have_rho ! Flag: true if density 'rho' is available 147 logical :: have_temp ! Flag: true if temperature is available in some way 142 148 logical :: have_sigma ! Flag: true if sigma levels are known (false if hybrid 143 149 ! coordinates are used) … … 937 943 ierr=NF_INQ_VARID(infid,"teta",tmpvarid) 938 944 if (ierr.ne.NF_NOERR) then 939 stop "Error: Failed to get t or teta ID" 945 write(*,*) "Error: Failed to get temp, t or teta ID" 946 have_temp=.false. 947 if (ztype.ne.1) stop 940 948 endif 941 949 ierr=NF_GET_VAR_REAL(infid,tmpvarid,teta) 942 950 if (ierr.ne.NF_NOERR) then 943 stop "Error: Failed reading atmospheric temperature" 951 write(*,*) "Failed reading atmospheric temperature" 952 if (ztype.ne.1) stop 944 953 endif 945 954 … … 953 962 enddo 954 963 enddo 964 have_temp=.true. 955 965 956 966 endif … … 960 970 stop "Error: Failed reading atmospheric temperature" 961 971 endif 972 have_temp=.true. 962 973 endif 963 974 … … 1358 1369 endif 1359 1370 1371 ! Add phisinit in outputfile 1372 phisdatashape(1)=lon_dimid 1373 phisdatashape(2)=lat_dimid 1374 ierr=NF_DEF_VAR(outfid,"phisinit",NF_REAL,2,phisdatashape,phis_varid) 1375 text='Geopotential at the surface' 1376 ierr=NF_PUT_ATT_TEXT(outfid,phis_varid,'long_name',len_trim(text),text) 1377 if (ierr.ne.NF_NOERR) then 1378 stop "Error: Problem writing long_name for phisinit" 1379 endif 1380 text='m2.s2' 1381 ierr=NF_PUT_ATT_TEXT(outfid,phis_varid,'units',len_trim(text),text) 1382 if (ierr.ne.NF_NOERR) then 1383 stop "Error: Problem writing units for geopotential phisinit" 1384 endif 1385 1386 1387 1360 1388 ! 3.3.2 Define 3D variables (ie: surface+time variables) 1361 1389 … … 1386 1414 1387 1415 ! add pressure or zareoid 1388 if ( ztype.eq.1) then ! pressure vertical coordinate1416 if ((ztype.eq.1).and.(have_temp)) then ! pressure vertical coordinate 1389 1417 ! zareoid dataset 1390 1418 ierr=NF_DEF_VAR(outfid,"zareoid",NF_REAL,4,datashape,za_varid) … … 1671 1699 endif 1672 1700 1701 ! write phisinit 1702 ierr=NF_PUT_VAR_REAL(outfid,phis_varid,phisinit) 1703 if (ierr.ne.NF_NOERR) then 1704 write(*,*) "Error: Could not write phisinit to output file" 1705 write(*,*) NF_STRERROR(ierr) 1706 stop 1707 endif 1708 1673 1709 !=============================================================================== 1674 1710 ! 3.5 Write 3D variables … … 1810 1846 ! interpolate dataset onto new grid 1811 1847 ! check if variable is "rho" (to set flag for interpolation below) 1812 if ( var(i)(1:3).eq.'rho') then1848 if ((var(i)(1:3).eq.'rho').or.(var(i)(1:4).eq.'num_')) then 1813 1849 j=1 1814 1850 write(*,*) trim(var(i)), ' is interpolated in log scale' … … 1866 1902 ! interpolate dataset onto new grid 1867 1903 ! check if variable is "rho" (to set flag for interpolation below) 1868 if ( var(i)(1:3).eq.'rho') then1904 if ((var(i)(1:3).eq.'rho').or.(var(i)(1:4).eq.'num_')) then 1869 1905 j=1 1870 1906 write(*,*) trim(var(i)), ' is interpolated in log scale' … … 1922 1958 ! interpolate dataset onto new grid 1923 1959 ! check if variable is "rho" (to set flag for interpolation below) 1924 if ( var(i)(1:3).eq.'rho') then1960 if ((var(i)(1:3).eq.'rho').or.(var(i)(1:4).eq.'num_')) then 1925 1961 j=1 1926 1962 write(*,*) trim(var(i)), ' is interpolated in log scale' … … 2839 2875 else if (z_new(kloop).gt.z_gcm(iloop,jloop,altlen,tloop)) then 2840 2876 ! z_new(kloop) is above highest GCM level 2841 newdata(iloop,jloop,kloop,tloop)=gcmdata(iloop,jloop,altlen,tloop) 2877 ! newdata(iloop,jloop,kloop,tloop)=gcmdata(iloop,jloop,altlen,tloop) 2878 newdata(iloop,jloop,kloop,tloop)=missing 2842 2879 else ! z_new(kloop) is within range 2843 2880 x=z_new(kloop)
Note: See TracChangeset
for help on using the changeset viewer.