Changeset 1073 for trunk/LMDZ.MARS/util/zrecast.F90
- Timestamp:
- Oct 17, 2013, 4:19:26 PM (11 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/util/zrecast.F90
r909 r1073 60 60 ! constants (radius, R, etc.) are now read from file 61 61 ! TN 01/2013 : Adapted for large output files with at least 2 variables > 2 GiB 62 ! TN 10/2013 : Read and write controle field, if available 62 63 ! 63 64 implicit none … … 80 81 integer nbvar4dinfile ! # of 4D (lon,lat,alt,time) variables in input file 81 82 integer outfid ! NetCDF output file ID 82 integer lon_dimid,lat_dimid,alt_dimid,time_dimid ! NetCDF dimension IDs83 integer lon_varid,lat_varid,alt_varid,time_varid 83 integer lon_dimid,lat_dimid,alt_dimid,time_dimid,ctl_dimid ! NetCDF dimension IDs 84 integer lon_varid,lat_varid,alt_varid,time_varid,ctl_varid 84 85 integer gcm_layers_dimid ! NetCDF dimension ID for # of layers in GCM 85 86 integer sigma_varid,aps_varid,bps_varid … … 103 104 integer latlength ! # of grid points along latitude 104 105 integer altlength ! # of grid point along altitude (of input datasets) 106 real,dimension(:),allocatable :: ctl ! controle 107 integer ctllength ! # of grid points along controle 105 108 real,dimension(:),allocatable :: time ! time 106 109 integer timelength ! # of points along time … … 439 442 if (ierr.ne.NF_NOERR) then 440 443 stop "Error: Failed to get altitude length" 444 endif 445 endif 446 447 ! controle 448 ierr=NF_INQ_DIMID(infid,"index",tmpdimid) 449 if (ierr.ne.NF_NOERR) then 450 write(*,*) "Could not get controle dimension ID" 451 ctllength = 0 452 else 453 ierr=NF_INQ_VARID(infid,"controle",tmpvarid) 454 if (ierr.ne.NF_NOERR) then 455 stop "Error: Failed to get controle ID" 456 else 457 ierr=NF_INQ_DIMLEN(infid,tmpdimid,ctllength) 458 if (ierr.ne.NF_NOERR) then 459 stop "Error: Failed to get controle length" 460 else 461 allocate(ctl(ctllength),stat=ierr) 462 if (ierr.ne.0) then 463 write(*,*) "Error: Failed to allocate ctl(ctllength)" 464 write(*,*) " ctllength=",ctllength 465 stop 466 endif 467 ierr=NF_GET_VAR_REAL(infid,tmpvarid,ctl) 468 if (ierr.ne.NF_NOERR) then 469 stop "Error: Failed reading controle" 470 endif 471 endif 441 472 endif 442 473 endif … … 1030 1061 endif 1031 1062 1063 ! controle 1064 if (ctllength .ne. 0) then 1065 ierr=NF_DEF_DIM(outfid,"index",ctllength,ctl_dimid) 1066 if (ierr.ne.NF_NOERR) then 1067 write(*,*) "Error: Could not define controle dimension" 1068 write(*,*) NF_STRERROR(ierr) 1069 stop 1070 endif 1071 endif 1072 1032 1073 ! GCM layers (for sigma or aps and bps) 1033 1074 ierr=NF_DEF_DIM(outfid,"GCM_layers",altlength,gcm_layers_dimid) … … 1162 1203 endif 1163 1204 endif ! of if (have_sigma) 1205 1206 ! controle 1207 if (ctllength .ne. 0) then 1208 ierr=NF_DEF_VAR(outfid,"controle",NF_REAL,1,ctl_dimid,ctl_varid) 1209 if (ierr.ne.NF_NOERR) then 1210 write(*,*) "Error: Could not define controle variable" 1211 write(*,*) NF_STRERROR(ierr) 1212 stop 1213 endif 1214 1215 ! controle attributes 1216 text='Control parameters' 1217 ierr=NF_PUT_ATT_TEXT(outfid,ctl_varid,'title',len_trim(text),text) 1218 if (ierr.ne.NF_NOERR) then 1219 stop "Error: Problem writing title for controle" 1220 endif 1221 endif 1164 1222 1165 1223 ! GCM_layers … … 1509 1567 endif 1510 1568 1569 ! Write controle 1570 if (ctllength .ne. 0) then 1571 ierr=NF_PUT_VAR_REAL(outfid,ctl_varid,ctl) 1572 if (ierr.ne.NF_NOERR) then 1573 write(*,*) "Error: Could not write controle data to output file" 1574 write(*,*) NF_STRERROR(ierr) 1575 stop 1576 endif 1577 endif 1578 1511 1579 ! Write time 1512 1580 ierr=NF_PUT_VARA_REAL(outfid,time_varid,1,timelength,time) … … 1771 1839 subroutine init_planet_const(infid) 1772 1840 ! initialize planetary constants using the "controle" array in the file 1773 ! if "co introle" array not found in file, look for it in "diagfi.nc"1841 ! if "controle" array not found in file, look for it in "diagfi.nc" 1774 1842 use planet_const 1775 1843 use netcdf
Note: See TracChangeset
for help on using the changeset viewer.