Changeset 2590
- Timestamp:
- Dec 3, 2021, 3:02:20 PM (3 years ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/README
r2589 r2590 3532 3532 => new flag: meteo_flux 3533 3533 => new input file in datadir: Meteo_flux_Plane_v2.dat 3534 3535 == 03/12/2021 == LR 3536 Utils: fixed a bug in concatnc. The program now properly carries over missing 3537 values from the input files into the output file. -
trunk/LMDZ.MARS/util/concatnc.F90
r2577 r2590 48 48 character (len=4) :: axis 49 49 ! axis: "sol" or "ls" or "adls" 50 integer :: nid,ierr,miss 50 integer :: nid,ierr,miss,validr 51 51 ! nid: [netcdf] file ID # 52 52 ! ierr: [netcdf] subroutine returned error code … … 694 694 ierr = NF_GET_VAR_REAL(nid,varid,var3d) 695 695 ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3d) 696 validr=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range) 696 697 miss=NF_GET_ATT_REAL(nid,varid,"missing_value",[missing]) 697 miss=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)698 698 699 699 if (ierr.ne.NF_NOERR) then … … 702 702 endif 703 703 704 ! In case there is a " valid_range" attribute704 ! In case there is a "missing_value" attribute 705 705 ! Write "valid_range" and "missing_value" attributes in output file 706 706 if (miss.eq.NF_NOERR) then 707 call missing_value(nout,varidout,valid _range,missing)707 call missing_value(nout,varidout,validr,miss,valid_range,missing) 708 708 endif 709 709 … … 1465 1465 end subroutine sol2ls 1466 1466 !****************************************************************************** 1467 subroutine missing_value(nout,nvarid,valid _range,missing)1467 subroutine missing_value(nout,nvarid,validr,miss,valid_range,missing) 1468 1468 !============================================================================== 1469 1469 ! Purpose: … … 1487 1487 integer, intent(in) :: nvarid 1488 1488 ! varid: [netcdf] variable ID # 1489 integer, intent(in) :: validr 1490 ! validr : [netcdf] routines return code for "valid_range" attribute 1491 integer, intent(in) :: miss 1492 ! miss : [netcdf] routines return code for "missing_value" attribute 1489 1493 real, dimension(2), intent(in) :: valid_range 1490 1494 ! valid_range(2): [netcdf] "valid_range" attribute (min and max) … … 1509 1513 1510 1514 ! Write valid_range() attribute 1511 ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range) 1512 1513 if (ierr.ne.NF_NOERR) then 1514 print*,'missing_value: valid_range attribution failed' 1515 print*, NF_STRERROR(ierr) 1516 !write(*,*) 'NF_NOERR', NF_NOERR 1517 !CALL abort 1518 stop 1515 if (validr.eq.NF_NOERR) then 1516 ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range) 1517 1518 if (ierr.ne.NF_NOERR) then 1519 print*,'missing_value: valid_range attribution failed' 1520 print*, NF_STRERROR(ierr) 1521 !write(*,*) 'NF_NOERR', NF_NOERR 1522 !CALL abort 1523 stop 1524 endif 1519 1525 endif 1520 1526 1521 1527 ! Write "missing_value" attribute 1522 ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,[missing]) 1523 1524 if (ierr.NE.NF_NOERR) then 1525 print*, 'missing_value: missing value attribution failed' 1526 print*, NF_STRERROR(ierr) 1527 ! WRITE(*,*) 'NF_NOERR', NF_NOERR 1528 ! CALL abort 1529 stop 1528 if (miss.eq.NF_NOERR) then 1529 ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,[missing]) 1530 1531 if (ierr.NE.NF_NOERR) then 1532 print*, 'missing_value: missing value attribution failed' 1533 print*, NF_STRERROR(ierr) 1534 ! WRITE(*,*) 'NF_NOERR', NF_NOERR 1535 ! CALL abort 1536 stop 1537 endif 1530 1538 endif 1531 1539
Note: See TracChangeset
for help on using the changeset viewer.