Ignore:
Timestamp:
Jan 19, 2017, 2:46:53 PM (8 years ago)
Author:
jvatant
Message:

Modifications to custom radiative transfer to Titan
+ Enables an altitude dependant gfrac for CIA computations

-> many radical changes in su_gases and co ..
-> read vertical CH4 profile with call_profilgases
-> Now you need a 'profile.def' that I will add in the deftank

+ Added interpolate CIA routines for CH4
+ Added temporary mean aerosol profile opacity routine (disr_haze)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.TITAN/libf/phytitan/su_gases.F90

    r1647 r1648  
    1 subroutine su_gases
     1subroutine su_gases(nlayer,tracer)
    22
    33  use gases_h
     
    55  implicit none
    66
    7   integer igas, ierr, count
    8 
     7  integer,intent(in) :: nlayer
     8  logical,intent(in) :: tracer
     9 
     10  integer igas, ierr
     11 
    912  !==================================================================
    1013  !     
     
    1720  !     R. Wordsworth (2011)
    1821  !     Allocatable arrays by A. Spiga (2011)
    19   !     
     22  !     Titan's version : J. Vatant d'Ollone (2017)
     23  !              * force igas_X labels and def gnom for N2,CH4,H2
     24  !              * get rid of variable species ( enrichment dimension will be set in corrk routines )
     25  !
    2026  !==================================================================
    2127
    2228!$OMP MASTER
    23   ! load gas names from file 'gases.def'
     29
     30
     31  ngasmx   = 3 ! N2, CH4, H2
     32
     33
     34  ! load reference level and reference molar fractions from file 'gases.def'
    2435  open(90,file='gases.def',status='old',form='formatted',iostat=ierr)
     36 
    2537  if (ierr.eq.0) then
    2638     write(*,*) "sugases.F90: reading file gases.def"
    27      read(90,*)
    28      read(90,*,iostat=ierr) ngasmx
     39     read(90,*) ! header
     40
     41     ! We allocate gfrac and we set gas molar fractions for the reference level only.
     42     ! This will be useful for the cpp_mu and rayleigh routines
     43     ! Other gas molar fractions are now set in routine callprofilgases
     44     
     45     write(*,*) 'sugases.F90: allocating and reading gas molar fractions from reference level in gases.def...'
     46     
     47     if(.not.allocated(gfrac)) allocate(gfrac(ngasmx,nlayer))
     48     
     49     read(90,*,iostat=ierr) nivref     
    2950     if (ierr.ne.0) then
    30         write(*,*) "sugases.F90: error reading number of gases"
     51        write(*,*) "sugases.F90: error reading reference level"
    3152        write(*,*) "   (first line of gases.def) "
    3253        call abort
    3354     endif
    34 
    35      print*,ngasmx, " gases found in gases.def. Allocating names and molar fractions..."
    36 
    37      if (.not.allocated(gnom)) allocate(gnom(ngasmx))
     55     
     56     print*, "layer", nivref, "is reference level found in gases.def ..."
     57     
    3858     do igas=1,ngasmx
    39         read(90,*,iostat=ierr) gnom(igas)
     59        read(90,*,iostat=ierr) gfrac(igas,nivref)
    4060        if (ierr.ne.0) then
    41            write(*,*) 'sugases.F90: error reading gas names in gases.def...'
     61           write(*,*) 'sugases.F90: error reading reference gas molar fractions in gases.def... aborting'
    4262           call abort
    4363        endif
    4464     enddo                  !of do igas=1,ngasmx
    4565
    46      vgas=0
    47      if(.not.allocated(gfrac)) allocate(gfrac(ngasmx))
    48      do igas=1,ngasmx
    49         read(90,*,iostat=ierr) gfrac(igas)
    50         if (ierr.ne.0) then
    51            write(*,*) 'sugases.F90: error reading gas molar fractions in gases.def...'
    52            call abort
    53         endif
    5466
    55         ! find variable gas (if any)
    56         if(gfrac(igas).eq.-1.0)then
    57            if(vgas.eq.0)then
    58               vgas=igas
    59            else
    60               print*,'You seem to be choosing two variable gases'
    61               print*,'Check that gases.def is correct'
    62               call abort
    63            endif
    64         endif
    65 
    66      enddo                  !of do igas=1,ngasmx
     67     ! We force gnom = (N2, CH4, H2) and igas_X for Titan
     68     
     69     igas_N2  = 1
     70     igas_CH4 = 2
     71     igas_H2  = 3
    6772
    6873
    69      ! assign the 'igas_X' labels
    70      count=0
    71      do igas=1,ngasmx
    72         if (trim(gnom(igas)).eq."H2_" .or. trim(gnom(igas)).eq."H2") then
    73            igas_H2=igas
    74            count=count+1
    75         elseif (trim(gnom(igas)).eq."N2_" .or. trim(gnom(igas)).eq."N2") then
    76            igas_N2=igas
    77            count=count+1
    78         elseif (trim(gnom(igas)).eq."CH4") then
    79            igas_CH4=igas
    80            count=count+1
    81         elseif (trim(gnom(igas)).eq."C2H6") then
    82            igas_C2H6=igas
    83            count=count+1
    84         elseif (trim(gnom(igas)).eq."C2H2") then
    85            igas_C2H2=igas
    86            count=count+1
    87         endif
    88      enddo
    89 
    90      if(count.ne.ngasmx)then
    91         print*,'Mismatch between ngas and number of recognised gases in sugas_corrk.F90.'
    92         print*,'Either we haven`t managed to assign all the gases, or there are duplicates.'
    93         print*,'Please try again.'
    94      endif
    95 
     74     if (.not.allocated(gnom)) allocate(gnom(ngasmx))
     75     gnom(igas_N2)  = "N2_"
     76     gnom(igas_CH4) = "CH4"
     77     gnom(igas_H2)  = "H2_"
     78     
     79     
    9680  else
    9781     write(*,*) 'Cannot find required file "gases.def"'
    9882     call abort
    9983  endif
     84 
    10085  close(90)
    10186!$OMP END MASTER
Note: See TracChangeset for help on using the changeset viewer.