Ignore:
Timestamp:
Dec 3, 2021, 3:02:20 PM (3 years ago)
Author:
lrossi
Message:

MARS GCM

Utils: fixed a bug in concatnc. The program now properly carries over missing
values from the input files into the output file.

LR

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/util/concatnc.F90

    r2577 r2590  
    4848character (len=4) :: axis
    4949! axis:  "sol" or "ls" or "adls"
    50 integer :: nid,ierr,miss
     50integer :: nid,ierr,miss,validr
    5151! nid: [netcdf] file ID #
    5252! ierr: [netcdf] subroutine returned error code
     
    694694      ierr = NF_GET_VAR_REAL(nid,varid,var3d)
    695695      ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3d)
     696      validr=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)
    696697      miss=NF_GET_ATT_REAL(nid,varid,"missing_value",[missing])
    697       miss=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)
    698698
    699699      if (ierr.ne.NF_NOERR) then
     
    702702      endif
    703703
    704 ! In case there is a "valid_range" attribute
     704! In case there is a "missing_value" attribute
    705705      ! Write "valid_range" and "missing_value" attributes in output file
    706706      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)
    708708      endif
    709709
     
    14651465end subroutine sol2ls
    14661466!******************************************************************************
    1467 subroutine  missing_value(nout,nvarid,valid_range,missing)
     1467subroutine  missing_value(nout,nvarid,validr,miss,valid_range,missing)
    14681468!==============================================================================
    14691469! Purpose:
     
    14871487integer, intent(in) :: nvarid
    14881488! varid: [netcdf] variable ID #
     1489integer, intent(in) :: validr
     1490! validr : [netcdf] routines return code for "valid_range" attribute
     1491integer, intent(in) :: miss
     1492! miss : [netcdf] routines return code for "missing_value" attribute
    14891493real, dimension(2), intent(in) :: valid_range
    14901494! valid_range(2): [netcdf] "valid_range" attribute (min and max)
     
    15091513
    15101514! 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
     1515if (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
    15191525endif
    15201526
    15211527! 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
     1528if (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
    15301538endif
    15311539
Note: See TracChangeset for help on using the changeset viewer.