source: trunk/LMDZ.TITAN/libf/phytitan/radcommon_h.F90 @ 2054

Last change on this file since 2054 was 2050, checked in by jvatant, 7 years ago

Major radiative transfer contribution : Add the 'corrk_recombin' option that allows to use
correlated-k for single species instead of pre-mix and enables more flexiblity for variable species.
-> Algorithm inspired from Lacis and Oinas 1991 and Amundsen et al 2016

+ Added 'recombin_corrk.F90' - Important improvements in 'sugas_corrk.90' and 'callcork.F90'

  • Must have the desired variable species as tracers -> TBD : Enable to force composition even if no tracers
  • To have decent CPU time recombining is not done on all gridpoints and wavelenghts but we calculate a gasi/v_recomb variable on the reference corrk-k T,P grid (only for T,P values who match the atmospheric conditions ) which is then processed as a standard pre-mix in optci/v routines, but updated every time tracers on the ref P grid have varied > 1%.


READ CAREFULY :

  • In case of 'corrk_recombin', the variable L_NREFVAR doesn't have the same meaning as before and doesn't stand for the different mixing ratios but the different species.
  • Input corr-k should be found in corrkdir within 'corrk_gcm_IR/VI_XXX.dat' and can contain a 'fixed' specie ( compulsory if you include self-broadening ) that MUST have been created with correct mixing ratios, or a variable specie for which mixing ratio MUST have been set to 1 ( no self-broadening then, assume it's a trace sepecie ) -> You can't neither have CIA of variable species included upstream in the corr-k
File size: 4.9 KB
Line 
1module radcommon_h
2      use radinc_h, only: L_NSPECTI, L_NSPECTV, L_NGAUSS, NTstar, NTstop
3      implicit none
4
5!----------------------------------------------------------------------C
6!
7!                             radcommon.h
8!
9!----------------------------------------------------------------------C
10!
11!  "Include" grid.h and radinc.h before this file in code that uses
12!  some or all of this common data set
13!
14!     WNOI       - Array of wavenumbers at the spectral interval
15!                  centers for the infrared.  Array is NSPECTI
16!                  elements long.
17!     DWNI       - Array of "delta wavenumber", i.e., the width,
18!                  in wavenumbers (cm^-1) of each IR spectral
19!                  interval.  NSPECTI elements long.
20!     WAVEI      - Array (NSPECTI elements long) of the wavelenght
21!                  (in microns) at the center of each IR spectral
22!                  interval.
23!     WNOV       - Array of wavenumbers at the spectral interval
24!                  center for the VISUAL.  Array is NSPECTV
25!                  elements long.
26!     DWNV       - Array of "delta wavenumber", i.e., the width,
27!                  in wavenumbers (cm^-1) of each VISUAL spectral
28!                  interval.  NSPECTV elements long.
29!     WAVEV      - Array (NSPECTV elements long) of the wavelenght
30!                  (in microns) at the center of each VISUAL spectral
31!                  interval.
32!     STELLARF   - Array (NSPECTV elements) of stellar flux (W/M^2) in
33!                  each spectral interval.  Values are for 1 AU,
34!                  scaled to the planetary distance elsewhere.
35!     TAURAY     - Array (NSPECTV elements) of the pressure-independent
36!                  part of Rayleigh scattering optical depth.
37!     FZEROI     - Fraction of zeros in the IR CO2 k-coefficients, for
38!                  each temperature, pressure, and spectral interval
39!     FZEROV     - Fraction of zeros in the VISUAL CO2 k-coefficients, for
40!                  each temperature, pressure, and spectral interval
41!
42
43      REAL*8 BWNI(L_NSPECTI+1), WNOI(L_NSPECTI), DWNI(L_NSPECTI), WAVEI(L_NSPECTI) !BWNI read by master in setspi
44      REAL*8 BWNV(L_NSPECTV+1), WNOV(L_NSPECTV), DWNV(L_NSPECTV), WAVEV(L_NSPECTV) !BWNV read by master in setspv
45      REAL*8 STELLARF(L_NSPECTV), TAURAY(L_NSPECTV)
46!$OMP THREADPRIVATE(WNOI,DWNI,WAVEI,&
47        !$OMP WNOV,DWNV,WAVEV,&
48        !$OMP STELLARF,TAURAY)
49
50      !! AS: introduced to avoid doing same computations again for continuum
51      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indi
52      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indv
53!$OMP THREADPRIVATE(indi,indv)
54
55      !!! ALLOCATABLE STUFF SO THAT DIMENSIONS ARE READ in *.dat FILES -- AS 12/2011 
56      REAL*8, DIMENSION(:,:,:,:,:), ALLOCATABLE :: gasi, gasv
57      REAL*8, DIMENSION(:), ALLOCATABLE :: PGASREF, TGASREF, PFGASREF, GWEIGHT
58 
59      ! For corrk recombining - JVO 18
60      REAL*8,  DIMENSION(:,:,:,:), ALLOCATABLE :: gasi_recomb, gasv_recomb
61      REAL*8,  DIMENSION(:,:),     ALLOCATABLE :: pqrold
62      REAL*8,  DIMENSION(:),       ALLOCATABLE :: w_cum
63      LOGICAL, DIMENSION(:,:),     ALLOCATABLE :: useptold
64      INTEGER, DIMENSION(:),       ALLOCATABLE :: permut_idx
65!$OMP THREADPRIVATE(gasi_recomb,gasv_recomb,pqrold,w_cum,useptold,permut_idx)
66     
67      LOGICAL, DIMENSION(:), ALLOCATABLE :: RADVAR_MASK ! for corrk recombin -> read by master in sugas_corrk
68      INTEGER, DIMENSION(:), ALLOCATABLE :: RADVAR_INDX ! for corrk recombin -> read by master in sugas_corrk
69     
70      real*8 FZEROI(L_NSPECTI)
71      real*8 FZEROV(L_NSPECTV)
72      real*8 pgasmin, pgasmax
73      real*8 tgasmin, tgasmax
74!$OMP THREADPRIVATE(gasi,gasv,&  !pgasref,tgasref,pfgasref, gweight read by master in sugas_corrk
75        !$OMP FZEROI,FZEROV)     !pgasmin,pgasmax,tgasmin,tgasmax read by master in sugas_corrk
76
77
78      REAL,SAVE :: tstellar ! Stellar brightness temperature (SW)
79
80      real*8,save :: planckir(L_NSPECTI,NTstop-NTstar+1)
81
82      real*8,save :: PTOP
83      real*8,save,allocatable :: TAUREF(:)
84
85      real*8,parameter :: UBARI = 0.5D0
86
87!$OMP THREADPRIVATE(QREFvis,QREFir,omegaREFvis,omegaREFir,&
88                !$OMP tstellar,planckir,PTOP)
89
90!     If the gas optical depth (top to the surface) is less than
91!     this value, we place that Gauss-point into the "zeros"
92!     channel.
93      real*8, parameter :: TLIMIT =  1.0D-30
94
95!     Factor to convert pressures from millibars to Pascals
96      real*8, parameter :: SCALEP = 1.00D+2
97
98      real*8, parameter :: sigma = 5.67032D-8
99      real*8, parameter :: grav = 6.672E-11
100
101      ! extinction of incoming sunlight (Saturn's rings, eclipses, etc...)
102      REAL, DIMENSION(:), ALLOCATABLE ,SAVE :: eclipse
103!$OMP THREADPRIVATE(eclipse)
104
105      ! Altitude-Latitude-dependent gravity
106      REAL, DIMENSION(:,:), ALLOCATABLE , SAVE :: gzlat ! This should be stored elsewhere ...
107      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: gzlat_ig
108      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: Cmk
109!$OMP THREADPRIVATE(gzlat,gzlat_ig,Cmk)
110
111end module radcommon_h
Note: See TracBrowser for help on using the repository browser.