Changeset 2327 for trunk/LMDZ.MARS/util


Ignore:
Timestamp:
May 18, 2020, 6:25:32 PM (5 years ago)
Author:
abierjon
Message:

Mars GCM:
Some improvements in aeroptical.F90 :
1) add the possibility to input a diagfi_P.nc file (with pressure as altitude coordinates) ;
2) replace all "title" attributes into "long_name" attributes in accordance with ticket #48

AB

File:
1 edited

Legend:

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

    r2318 r2327  
    7474character(len=16) :: tmpvarname                  ! temporarily stores a variable name
    7575real,dimension(:,:,:,:,:),allocatable :: opa_aer ! Aerosols opacities [1/km]
    76 real :: missval                                  ! Value to put in outfile when the reff is out of bounds
    77 PARAMETER(missval=1E+20)
     76real,dimension(:),allocatable :: missval         ! Value to put in outfile when the reff is out of bounds
     77                                                 ! or when there is a mising value in input file
    7878
    7979
     
    179179allocate(mmr(naerkind,lonlen,latlen,altlen,timelen))
    180180allocate(reff(naerkind,lonlen,latlen,altlen,timelen))
     181
     182! Initialize missing value
     183allocate(missval(naerkind))
     184missval(:)=1.e+20
    181185
    182186! Initialize aerok to .true. for all aerosols
     
    201205  call status_check(ierr,error_text)
    202206  write(*,*) "Dust mass mixing ratio loaded from "//trim(gcmfile)
     207  ! Get missing value
     208  ierr=nf90_get_att(gcmfid,tmpvarid,"missing_value",missval(1))
     209  if (ierr.ne.nf90_noerr) then
     210    missval(1) = 1.e+20
     211  endif
    203212endif
    204213
     
    238247  call status_check(ierr,error_text)
    239248  write(*,*) "Water ice mass mixing ratio loaded from "//trim(gcmfile)
     249  ! Get missing value
     250  ierr=nf90_get_att(gcmfid,tmpvarid,"missing_value",missval(2))
     251  if (ierr.ne.nf90_noerr) then
     252    missval(2) = 1.e+20
     253  endif
    240254endif
    241255
     
    532546    ierr=nf90_put_att(outfid,tmpvaridout,"units","opacity/km")
    533547    ierr=nf90_put_att(outfid,tmpvaridout,"refwavelength",wvl_val)
    534     ierr=nf90_put_att(outfid,tmpvaridout,"missing_value",missval)
    535     write(*,*)"with missing value = ",missval
     548    ierr=nf90_put_att(outfid,tmpvaridout,"missing_value",missval(iaer))
     549    write(*,*)"with missing value = ",missval(iaer)
    536550   
    537551    ! End netcdf define mode
     
    557571
    558572           ! NB: this should especially handle cases when reff=0
    559            opa_aer(iaer,ilon,ilat,ialt,it)=missval
     573           opa_aer(iaer,ilon,ilat,ialt,it)=missval(iaer)
     574           
     575         else if (mmr(iaer,ilon,ilat,ialt,it).eq.missval(iaer)) then
     576           ! if there is a missing value in input file
     577           opa_aer(iaer,ilon,ilat,ialt,it)=missval(iaer)
     578           
    560579         else
    561580           if (radiusdyn(isize).eq.reff(iaer,ilon,ilat,ialt,it)) then
     
    10491068integer, intent(in) :: lonlen ! # of grid points along longitude
    10501069integer, intent(in) :: latlen ! # of grid points along latitude
    1051 integer, intent(in) :: altlen ! # of grid points along latitude
     1070integer, intent(in) :: altlen ! # of grid points along altitude
    10521071integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers
    10531072integer, intent(in) :: outfid ! NetCDF output file ID
     
    12041223
    12051224!==============================================================================
    1206 ! 2.2. Hybrid coordinates ap() , bp(), aps() and bps()
     1225! 2.1 Hybrid coordinates ap() , bp(), aps() and bps()
    12071226!==============================================================================
    12081227if(hybrid) then
     
    12161235  endif
    12171236  ! Write the attributes
    1218   ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid pressure at midlayers")
     1237  ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid pressure at midlayers")
    12191238  ierr=nf90_put_att(outfid,tmpvarid,"units"," ")
    12201239  ! End netcdf define mode
     
    12361255  endif
    12371256  ! Write the attributes
    1238   ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid sigma at midlayers")
     1257  ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid sigma at midlayers")
    12391258  ierr=nf90_put_att(outfid,tmpvarid,"units"," ")
    12401259  ! End netcdf define mode
     
    12581277    endif
    12591278    ! Write the attributes
    1260     ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid sigma at interlayers")
     1279    ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid sigma at interlayers")
    12611280    ierr=nf90_put_att(outfid,tmpvarid,"units"," ")
    12621281    ! End netcdf define mode
     
    12791298    endif
    12801299    ! Write the attributes
    1281     ierr=nf90_put_att(outfid,tmpvarid,"title","hybrid sigma at interlayers")
     1300    ierr=nf90_put_att(outfid,tmpvarid,"long_name","hybrid sigma at interlayers")
    12821301    ierr=nf90_put_att(outfid,tmpvarid,"units"," ")
    12831302    ! End netcdf define mode
     
    13011320  endif
    13021321  ! Write the attributes
    1303   ierr=nf90_put_att(outfid,tmpvarid,"title","sigma at midlayers")
     1322  ierr=nf90_put_att(outfid,tmpvarid,"long_name","sigma at midlayers")
    13041323  ierr=nf90_put_att(outfid,tmpvarid,"units"," ")
    13051324  ! End netcdf define mode
     
    13141333
    13151334!==============================================================================
    1316 ! 2.2. aire() and phisinit()
     1335! 2.2 aire() and phisinit()
    13171336!==============================================================================
    13181337
     
    13271346  endif
    13281347  ! Write the attributes
    1329   ierr=nf90_put_att(outfid,tmpvarid,"title","Mesh area")
     1348  ierr=nf90_put_att(outfid,tmpvarid,"long_name","Mesh area")
    13301349  ierr=nf90_put_att(outfid,tmpvarid,"units","m2")
    13311350  ! End netcdf define mode
     
    13391358endif ! of if (area)
    13401359
    1341 IF (phis) THEN
     1360if (phis) then
    13421361
    13431362  ! Switch to netcdf define mode
     
    13491368  endif
    13501369  ! Write the attributes
    1351   ierr=nf90_put_att(outfid,tmpvarid,"title","Ground level geopotential")
     1370  ierr=nf90_put_att(outfid,tmpvarid,"long_name","Ground level geopotential")
    13521371  ierr=nf90_put_att(outfid,tmpvarid,"units"," ")
    13531372  ! End netcdf define mode
     
    13601379  endif
    13611380
    1362 ENDIF ! of IF (phis)
     1381endif ! of if (phis)
    13631382
    13641383
Note: See TracChangeset for help on using the changeset viewer.