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

Last change on this file since 1966 was 1947, checked in by jvatant, 7 years ago

Enables altitude-depending gravity field g(z) (glat->gzlat) in physics
+ Can be dangerous ( disagreement with dyn) but important (compulsive !) to have correct altitudes in the chemistry
+ Can be activated with eff_gz flag in callphys.def
-- JVO

File size: 4.3 KB
RevLine 
[1529]1module radcommon_h
[1788]2      use radinc_h, only: L_NSPECTI, L_NSPECTV, L_NGAUSS, NTstar, NTstop
[135]3      implicit none
4
5!----------------------------------------------------------------------C
6!
7!                             radcommon.h
[1722]8!
[135]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
[1315]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
[1648]45      REAL*8 STELLARF(L_NSPECTV), TAURAY(L_NSPECTV)
[1315]46!$OMP THREADPRIVATE(WNOI,DWNI,WAVEI,&
47        !$OMP WNOV,DWNV,WAVEV,&
[1648]48        !$OMP STELLARF,TAURAY)
[135]49
[873]50      !! AS: introduced to avoid doing same computations again for continuum
[878]51      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indi
52      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indv
[1315]53!$OMP THREADPRIVATE(indi,indv)
[873]54
[470]55      !!! ALLOCATABLE STUFF SO THAT DIMENSIONS ARE READ in *.dat FILES -- AS 12/2011 
56      REAL*8, DIMENSION(:,:,:,:,:), ALLOCATABLE :: gasi, gasv
[1648]57      REAL*8, DIMENSION(:), ALLOCATABLE :: PGASREF, TGASREF, PFGASREF
[135]58      real*8 FZEROI(L_NSPECTI)
59      real*8 FZEROV(L_NSPECTV)
60      real*8 pgasmin, pgasmax
61      real*8 tgasmin, tgasmax
[1648]62!$OMP THREADPRIVATE(gasi,gasv,&  !pgasref,tgasref,pfgasref read by master in sugas_corrk
[1315]63        !$OMP FZEROI,FZEROV)     !pgasmin,pgasmax,tgasmin,tgasmax read by master in sugas_corrk
[135]64
65
[1529]66      REAL,SAVE :: tstellar ! Stellar brightness temperature (SW)
[135]67
[1529]68      real*8,save :: planckir(L_NSPECTI,NTstop-NTstar+1)
[135]69
[1529]70      real*8,save :: PTOP
71      real*8,save,allocatable :: TAUREF(:)
[135]72
[1529]73      real*8,parameter :: UBARI = 0.5D0
[135]74
[1529]75      real*8,save :: gweight(L_NGAUSS)
[1315]76!$OMP THREADPRIVATE(QREFvis,QREFir,omegaREFvis,omegaREFir,&     ! gweight read by master in sugas_corrk
[1722]77                !$OMP tstellar,planckir,PTOP)
[135]78
79!     If the gas optical depth (top to the surface) is less than
80!     this value, we place that Gauss-point into the "zeros"
81!     channel.
82      real*8, parameter :: TLIMIT =  1.0D-30
83
84!     Factor to convert pressures from millibars to Pascals
85      real*8, parameter :: SCALEP = 1.00D+2
86
[959]87      real*8, parameter :: sigma = 5.67032D-8
[1194]88      real*8, parameter :: grav = 6.672E-11
[135]89
[1133]90      ! extinction of incoming sunlight (Saturn's rings, eclipses, etc...)
[1529]91      REAL, DIMENSION(:), ALLOCATABLE ,SAVE :: eclipse
[1947]92!$OMP THREADPRIVATE(eclipse)
[135]93
[1947]94      ! Altitude-Latitude-dependent gravity
95      REAL, DIMENSION(:,:), ALLOCATABLE , SAVE :: gzlat ! This should be stored elsewhere ...
96      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: gzlat_ig
97      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: Cmk
98!$OMP THREADPRIVATE(gzlat,gzlat_ig,Cmk)
[1194]99
[1529]100end module radcommon_h
Note: See TracBrowser for help on using the repository browser.