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

Last change on this file since 3529 was 2366, checked in by jvatant, 4 years ago

Titan GCM : Major maintenance catching up commits from the generic including :

  • r2356 and 2354 removing obsolete old dynamical core
  • various minor addition to physics and gestion of phys_state_var_mode, especially in dyn1d
  • adding MESOSCALE CPP keys around chemistry and microphysics (disabled in mesoscale for now)
File size: 4.8 KB
Line 
1module radcommon_h
2      use radinc_h, only: L_NSPECTI, L_NSPECTV, NTstart, 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, DIMENSION(:,:), ALLOCATABLE, SAVE :: planckir
81
82      real*8,save :: PTOP
83!$OMP THREADPRIVATE(tstellar,planckir,PTOP)
84
85      real*8,parameter :: UBARI = 0.5D0
86
87!     If the gas optical depth (top to the surface) is less than
88!     this value, we place that Gauss-point into the "zeros"
89!     channel.
90      real*8, parameter :: TLIMIT =  1.0D-30
91
92!     Factor to convert pressures from millibars to Pascals
93      real*8, parameter :: SCALEP = 1.00D+2
94
95      real*8, parameter :: sigma = 5.67032D-8
96      real*8, parameter :: grav = 6.672E-11
97
98      ! extinction of incoming sunlight (Saturn's rings, eclipses, etc...)
99      REAL, DIMENSION(:), ALLOCATABLE ,SAVE :: eclipse
100!$OMP THREADPRIVATE(eclipse)
101
102      ! Altitude-Latitude-dependent gravity
103      REAL, DIMENSION(:,:), ALLOCATABLE , SAVE :: gzlat ! This should be stored elsewhere ...
104      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: gzlat_ig
105      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: Cmk
106!$OMP THREADPRIVATE(gzlat,gzlat_ig,Cmk)
107
108end module radcommon_h
Note: See TracBrowser for help on using the repository browser.