Changeset 543 for trunk


Ignore:
Timestamp:
Feb 27, 2012, 10:01:25 AM (13 years ago)
Author:
aslmd
Message:

LMDZ.GENERIC: Temperature grid for Planck calculations can now be refined through the parameter NTfac in radinc_h

Location:
trunk/LMDZ.GENERIC
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/README

    r540 r543  
    595595- kcm1d no longer consistent with new code, but I haven't updated as I'm still working on it
    596596  and as far as I know noone else uses it. If that changes let me know.
     597
     598== 27/02/2012 == AS
     599- Temperature grid for Planck calculations can now be refined through the parameter NTfac in radinc_h.
     600  Default is NTfac = 1.0D-1, i.e. Delta T = 0.1 K
  • trunk/LMDZ.GENERIC/libf/phystd/gfluxi.F

    r253 r543  
    6262 
    6363
    64       IF (NLL .GT. NLP) STOP 'PARAMETER NL TOO SMALL IN GFLUXV'
     64      IF (NLL .GT. NLP) STOP 'PARAMETER NL TOO SMALL IN GFLUXI'
    6565
    6666      NLAYER = L_NLAYRAD
     
    7272        !NT2   = TLEV(2*L+2)*10.0D0-499
    7373        !NT    = TLEV(2*L)*10.0D0-499
    74         NT    = int(TLEV(2*L)*10.0D0)   - NTstar+1
    75         NT2   = int(TLEV(2*L+2)*10.0D0) - NTstar+1
     74        NT    = int(TLEV(2*L)*NTfac)   - NTstar+1
     75        NT2   = int(TLEV(2*L+2)*NTfac) - NTstar+1
    7676
    7777        B1(L) = (PLANCKIR(NW,NT2)-PLANCKIR(NW,NT))/DTAU(L)
     
    8787      !NT    = TLEV(2*L+1)*10.0D0-499
    8888      !NT2   = TLEV(2*L)*10.0D0-499
    89       NT    = int(TLEV(2*L+1)*10.0D0) - NTstar+1
    90       NT2   = int(TLEV(2*L)*10.0D0)   - NTstar+1
     89      NT    = int(TLEV(2*L+1)*NTfac) - NTstar+1
     90      NT2   = int(TLEV(2*L)*NTfac)   - NTstar+1
    9191      B1(L) = (PLANCKIR(NW,NT)-PLANCKIR(NW,NT2))/DTAU(L)
    9292      B0(L) = PLANCKIR(NW,NT2)
  • trunk/LMDZ.GENERIC/libf/phystd/radinc_h.F90

    r537 r543  
    7373
    7474      ! For Planck function integration:
    75       ! equivalent temperatures are 1/10 of these values
     75      ! equivalent temperatures are 1/NTfac of these values
    7676      integer, parameter :: NTstar = 500
    7777      integer, parameter :: NTstop = 9000 ! new default for all non hot Jupiter runs
     78      real*8, parameter :: NTfac = 1.0D+1 
     79      !integer, parameter :: NTstar = 1000
     80      !integer, parameter :: NTstop = 25000
     81      !real*8,parameter :: NTfac = 5.0D+1   
     82      !integer, parameter :: NTstar = 2000
     83      !integer, parameter :: NTstop = 50000
     84      !real*8,parameter :: NTfac = 1.0D+2   
    7885
    7986      ! Maximum number of grain size classes for aerosol convolution:
  • trunk/LMDZ.GENERIC/libf/phystd/setspi.F90

    r374 r543  
    2222!==================================================================
    2323
    24       use radinc_h,    only: L_NSPECTI,corrkdir,banddir,NTstar,NTstop
     24      use radinc_h,    only: L_NSPECTI,corrkdir,banddir,NTstar,NTstop,NTfac
    2525      use radcommon_h, only: BWNI,BLAMI,WNOI,DWNI,WAVEI,planckir,sigma
    2626      use datafile_mod, only: datadir
     
    139139      print*,''
    140140      print*,'setspi: Current Planck integration range:'
    141       print*,'T = ',dble(NTstar)/1.0D+1, ' to ',dble(NTstop)/1.0D+1,' K.'
     141      print*,'T = ',dble(NTstar)/NTfac, ' to ',dble(NTstop)/NTfac,' K.'
    142142
    143143      do NW=1,L_NSPECTI
     
    147147         bma = (b-a)/2.0
    148148         do nt=NTstar,NTstop
    149             T   = dble(NT)/1.0D+1
     149            T   = dble(NT)/NTfac
    150150            ans = 0.0D0
    151151
     
    164164         do nt=NTstar,NTstop
    165165            plancksum=0.0
    166             T=dble(NT)/1.0D+1
     166            T=dble(NT)/NTfac
    167167       
    168168            do NW=1,L_NSPECTI
     
    174174               planckir(NW,nt-NTstar+1)=     &
    175175                  planckir(NW,nt-NTstar+1)*  &
    176                           sigma*(dble(nt)/1.0D+1)**4/plancksum
     176                          sigma*(dble(nt)/NTfac)**4/plancksum
    177177            end do
    178178         end do
     
    188188         print*,'setspi: At lower limit:'
    189189         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
    190          print*,'actual sig*T^4   = ',sigma*(dble(nt)/1.0D+1)**4,' W m^-2'
     190         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
    191191         
    192192         ! check energy conservation at upper temperature boundary
     
    198198         print*,'setspi: At upper limit:'
    199199         print*,'in model sig*T^4 = ',plancksum,' W m^-2'
    200          print*,'actual sig*T^4   = ',sigma*(dble(nt)/1.0D+1)**4,' W m^-2'
     200         print*,'actual sig*T^4   = ',sigma*(dble(nt)/NTfac)**4,' W m^-2'
    201201         print*,''
    202202      endif
  • trunk/LMDZ.GENERIC/libf/phystd/sfluxi.F

    r526 r543  
    6767      TSURF = TLEV(L_LEVELS)
    6868
    69       NTS   = int(TSURF*10.0D0)-NTstar+1
    70       NTT   = int(TTOP *10.0D0)-NTstar+1
     69      NTS   = int(TSURF*NTfac)-NTstar+1
     70      NTT   = int(TTOP *NTfac)-NTstar+1
    7171!      NTS   = TSURF*10.0D0-499
    7272!      NTT   = TTOP *10.0D0-499
Note: See TracChangeset for help on using the changeset viewer.