Changeset 2131


Ignore:
Timestamp:
Apr 29, 2019, 2:52:42 PM (6 years ago)
Author:
jvatant
Message:

+ Add diagnostics of optical thickness, in 1D, if 'diagdtau' key is activated, it

outputs dtaui/v(altitude) in diagfi.nc for every narrowband (could be done with one var
but would require to be able to have writediag in 5D

--JVO

Location:
trunk
Files:
10 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/README

    r2128 r2131  
    14771477- new version of the thermal plume model (new entrainment and detrainment, removed alimentation)
    14781478- removal of useless files thermcell_alp, thermcell_dv2, thermcell_alim and thermcell_dry.F90
     1479
     1480== 29/04/2019 == JVO
     1481+ Add diagnostics of optical thickness, in 1D, if 'diagdtau' key is activated, it
     1482 outputs dtaui/v(altitude) in diagfi.nc for every narrowband (could be done with one var
     1483 but would require to be able to have writediag in 5D)
  • trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90

    r2062 r2131  
    2020      logical,save :: specOLR
    2121      logical,save :: kastprof
    22 !$OMP THREADPRIVATE(enertest,nonideal,meanOLR,kastprof)
     22      logical,save :: diagdtau
     23!$OMP THREADPRIVATE(enertest,nonideal,meanOLR,kastprof,diagdtau)
    2324      logical,save :: newtonian
    2425      logical,save :: check_cpp_match
  • trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90

    r2078 r2131  
    327327     write(*,*)" specOLR = ",specOLR
    328328
     329     write(*,*)"Output diagnostic optical thickness?"
     330     diagdtau=.false.
     331     call getin_p("diagdtau",diagdtau)
     332     write(*,*)" diagdtau = ",diagdtau
     333     ! sanity check
     334     if (diagdtau.and.(ngrid.ne.1)) then
     335       print*,"Diagnostic optical thickness can be output in 1D only !"
     336       print*,"Start again with diagdtau=.false."
     337       stop
     338     endif
     339
    329340     write(*,*)"Operate in kastprof mode?"
    330341     kastprof=.false.
  • trunk/LMDZ.GENERIC/libf/phystd/optci.F90

    r2032 r2131  
    1111  use radinc_h, only: L_LEVELS, L_NLAYRAD, L_NSPECTI, L_NGAUSS, &
    1212                      L_NLEVRAD, L_REFVAR, naerkind
    13   use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi,glat_ig
     13  use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi,glat_ig,gweight
    1414  use gases_h, only: gfrac, ngasmx, igas_N2, igas_He, igas_H2O, igas_H2
    1515  use comcstfi_mod, only: g, r, mugaz
    16   use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple
     16  use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,diagdtau
    1717  implicit none
    1818
     
    8484  !real*8 rho !! see test below
    8585
     86  ! variables for diags of optical depth
     87  real *8 INT_DTAU(L_NLAYRAD,L_NSPECTI)
     88  character*2 str2
     89
    8690  integer igas, jgas
    8791
     
    449453!  print*,DTAUI
    450454!  call abort
     455
     456
     457  !  Optical thickness for 1D diagnostics (added by JVO)
     458  if (diagdtau) then ! diagtau can be true only if 1D
     459    do l=1,L_NLAYRAD
     460       do nw=1,L_NSPECTI
     461        INT_DTAU(L,NW) = 0.0d+0
     462          DO NG=1,L_NGAUSS
     463             INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtaui(L,nw,ng)*gweight(NG)
     464          enddo
     465       enddo
     466    enddo
     467    do nw=1,L_NSPECTI
     468      write(str2,'(i2.2)') nw
     469      call writediagfi(1,'dtaui'//str2,'Layer optical thickness in IR band '//str2,'',1,int_dtau(L_NLAYRAD:1:-1,nw))
     470    enddo
     471  endif
     472 
    451473 
    452474
  • trunk/LMDZ.GENERIC/libf/phystd/optcv.F90

    r2032 r2131  
    1010
    1111  use radinc_h, only: L_NLAYRAD, L_NLEVRAD, L_LEVELS, L_NSPECTV, L_NGAUSS, L_REFVAR, NAERKIND
    12   use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig
     12  use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig,gweight
    1313  use gases_h, only: gfrac, ngasmx, igas_H2, igas_H2O, igas_He, igas_N2
    1414  use comcstfi_mod, only: g, r, mugaz
    15   use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,callgasvis
     15  use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,callgasvis,diagdtau
    1616
    1717  implicit none
     
    9191  ! variables for k in units m^-1
    9292  real*8 dz(L_LEVELS)
     93
     94  ! variables for diags of optical depth
     95  real*8 INT_DTAU(L_NLAYRAD,L_NSPECTV)
     96  character*2 str2
    9397
    9498  integer igas, jgas
     
    375379
    376380
     381  !  Optical thickness for 1D diagnostics (added by JVO)
     382  if (diagdtau) then ! diagtau can be true only if 1D
     383    do l=1,L_NLAYRAD
     384       do nw=1,L_NSPECTV
     385        INT_DTAU(L,NW) = 0.0d+0
     386          DO NG=1,L_NGAUSS
     387             INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtauv(L,nw,ng)*gweight(NG)
     388          enddo
     389       enddo
     390    enddo
     391    do nw=1,L_NSPECTV
     392      write(str2,'(i2.2)') nw
     393      call writediagfi(1,'dtauv'//str2,'Layer optical thickness in VI band '//str2,'',1,int_dtau(L_NLAYRAD:1:-1,nw))
     394    enddo
     395  endif
     396
     397
    377398end subroutine optcv
    378399
  • trunk/LMDZ.TITAN/README

    r2109 r2131  
    15571557+ Altitude of the last level at 1e7m from physics was certainly source of divergence
    15581558+ Sanity check for negative is moved from within mm_microphysic to the end of calmufi avoiding rounding pbs
     1559
     1560== 29/04/2019 == JVO
     1561+ Add diagnostics of optical thickness, in 1D, if 'diagdtau' key is activated, it
     1562 outputs dtaui/v(altitude) in diagfi.nc for every narrowband (could be done with one var
     1563 but would require to be able to have writediag in 5D)
  • trunk/LMDZ.TITAN/libf/phytitan/callkeys_mod.F90

    r2050 r2131  
    2727      logical,save :: meanOLR
    2828      logical,save :: specOLR
    29 !$OMP THREADPRIVATE(enertest,nonideal,meanOLR)
     29      logical,save :: diagdtau
     30!$OMP THREADPRIVATE(enertest,nonideal,meanOLR,specOLR,diagdtau)
    3031      logical,save :: newtonian
    3132      logical,save :: check_cpp_match
  • trunk/LMDZ.TITAN/libf/phytitan/inifis_mod.F90

    r2078 r2131  
    336336     call getin_p("specOLR",specOLR)
    337337     write(*,*)" specOLR = ",specOLR
     338
     339     write(*,*)"Output diagnostic optical thickness?"
     340     diagdtau=.false.
     341     call getin_p("diagdtau",diagdtau)
     342     write(*,*)" diagdtau = ",diagdtau
     343     ! sanity check
     344     if (diagdtau.and.(ngrid.ne.1)) then
     345       print*,"Diagnostic optical thickness can be output in 1D only !"
     346       print*,"Start again with diagdtau=.false."
     347       stop
     348     endif
     349
    338350
    339351     write(*,*)"Uniform absorption in radiative transfer?"
  • trunk/LMDZ.TITAN/libf/phytitan/optci.F90

    r2050 r2131  
    77  use gases_h
    88  use comcstfi_mod, only: r
    9   use callkeys_mod, only: continuum,graybody,corrk_recombin,      &
     9  use callkeys_mod, only: continuum,graybody,corrk_recombin,diagdtau,      &
    1010                          callclouds,callmufi,seashaze,uncoupl_optic_haze
    1111  use tracer_h, only : nmicro,nice
    12   use MMP_OPTICS
    1312
    1413  implicit none
     
    6160  real*8 ASF_T(L_LEVELS,L_NSPECTI)
    6261  real*8 INT_DTAU(L_NLAYRAD,L_NSPECTI)
    63   real*8 K_HAZE(L_NLAYRAD,L_NSPECTI)
    6462
    6563  CHARACTER*2  str2
     
    412410
    413411
    414 !  Titan's outputs (J.V.O, 2016)===============================================
    415 !      do l=1,L_NLAYRAD
    416 !         do nw=1,L_NSPECTI
    417 !          INT_DTAU(L,NW) = 0.0d+0
    418 !            DO NG=1,L_NGAUSS
    419 !               INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtaui(L,nw,ng)*gweight(NG)
    420 !            enddo
    421 !         enddo
    422 !      enddo
    423 
    424 !       do nw=1,L_NSPECTI
    425 !          write(str2,'(i2.2)') nw
    426 !         call writediagfi(1,'kgi'//str2,'Gaz extinction coefficient IR band '//str2,'m-1',1,int_dtau(L_NLAYRAD:1:-1,nw)/dz_lay(L_NLAYRAD:1:-1))
    427 !          call writediagfi(1,'khi'//str2,'Haze extinction coefficient IR band '//str2,'m-1',1,k_haze(L_NLAYRAD:1:-1,nw)/dz_lay(L_NLAYRAD:1:-1))       
    428 !       enddo 
    429  
    430 ! ============================================================================== 
     412  !  Optical thickness for 1D diagnostics (added by JVO)
     413  if (diagdtau) then ! diagtau can be true only if 1D
     414    do l=1,L_NLAYRAD
     415       do nw=1,L_NSPECTI
     416        INT_DTAU(L,NW) = 0.0d+0
     417          DO NG=1,L_NGAUSS
     418             INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtaui(L,nw,ng)*gweight(NG)
     419          enddo
     420       enddo
     421    enddo
     422    do nw=1,L_NSPECTI
     423      write(str2,'(i2.2)') nw
     424      call writediagfi(1,'dtaui'//str2,'Layer optical thickness in IR band '//str2,'',1,int_dtau(L_NLAYRAD:1:-1,nw))
     425    enddo
     426  endif
     427 
    431428
    432429  if(firstcall) firstcall = .false.
  • trunk/LMDZ.TITAN/libf/phytitan/optcv.F90

    r2095 r2131  
    77  use gases_h
    88  use comcstfi_mod, only: r
    9   use callkeys_mod, only: continuum,graybody,callgasvis,corrk_recombin,  &
     9  use callkeys_mod, only: continuum,graybody,callgasvis,corrk_recombin,diagdtau,  &
    1010                          callclouds,callmufi,seashaze,uncoupl_optic_haze
    1111  use tracer_h, only: nmicro,nice
     
    6969  real*8 ASF_T(L_LEVELS,L_NSPECTI)
    7070  real*8 INT_DTAU(L_NLAYRAD,L_NSPECTI)
    71   real*8 K_HAZE(L_NLAYRAD,L_NSPECTI)
    7271 
    7372  CHARACTER*2  str2
     
    391390
    392391
    393 !  Titan's outputs (JVO, 2016)===============================================
    394 !      do l=1,L_NLAYRAD
    395 !         do nw=1,L_NSPECTV
    396 !          INT_DTAU(L,NW) = 0.0d+0
    397 !            DO NG=1,L_NGAUSS
    398 !               INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtauv(L,nw,ng)*gweight(NG)
    399 !            enddo
    400 !         enddo
    401 !      enddo
    402 
    403 !       do nw=1,L_NSPECTV
    404 !          write(str2,'(i2.2)') nw
    405 !         call writediagfi(1,'kgv'//str2,'Gaz extinction coefficient VI band '//str2,'m-1',1,int_dtau(L_NLAYRAD:1:-1,nw)/dz_lay(L_NLAYRAD:1:-1))
    406 !          call writediagfi(1,'khv'//str2,'Haze extinction coefficient VI band '//str2,'m-1',1,k_haze(L_NLAYRAD:1:-1,nw)/dz_lay(L_NLAYRAD:1:-1))       
    407 !       enddo 
    408 
    409 ! ============================================================================== 
     392  !  Optical thickness for 1D diagnostics (added by JVO)
     393  if (diagdtau) then ! diagtau can be true only if 1D
     394    do l=1,L_NLAYRAD
     395       do nw=1,L_NSPECTV
     396        INT_DTAU(L,NW) = 0.0d+0
     397          DO NG=1,L_NGAUSS
     398             INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtauv(L,nw,ng)*gweight(NG)
     399          enddo
     400       enddo
     401    enddo
     402    do nw=1,L_NSPECTV
     403      write(str2,'(i2.2)') nw
     404      call writediagfi(1,'dtauv'//str2,'Layer optical thickness in VI band '//str2,'',1,int_dtau(L_NLAYRAD:1:-1,nw))
     405    enddo
     406  endif
     407
    410408
    411409  if(firstcall) firstcall = .false.
Note: See TracChangeset for help on using the changeset viewer.