source: trunk/LMDZ.GENERIC/libf/phystd/radcommon_h.F90 @ 2142

Last change on this file since 2142 was 2142, checked in by aslmd, 5 years ago

corrected a bug that prevented compiling: a removed variable was still mentionned in a OMP list

File size: 6.1 KB
Line 
1module radcommon_h
2      use radinc_h, only: L_NSPECTI, L_NSPECTV, NTstar, NTstop, &
3                          naerkind, nsizemax
4      implicit none
5
6!----------------------------------------------------------------------C
7!
8!                             radcommon.h
9!
10!----------------------------------------------------------------------C
11!
12!  "Include" grid.h and radinc.h before this file in code that uses
13!  some or all of this common data set
14!
15!     WNOI       - Array of wavenumbers at the spectral interval
16!                  centers for the infrared.  Array is NSPECTI
17!                  elements long.
18!     DWNI       - Array of "delta wavenumber", i.e., the width,
19!                  in wavenumbers (cm^-1) of each IR spectral
20!                  interval.  NSPECTI elements long.
21!     WAVEI      - Array (NSPECTI elements long) of the wavelenght
22!                  (in microns) at the center of each IR spectral
23!                  interval.
24!     WNOV       - Array of wavenumbers at the spectral interval
25!                  center for the VISUAL.  Array is NSPECTV
26!                  elements long.
27!     DWNV       - Array of "delta wavenumber", i.e., the width,
28!                  in wavenumbers (cm^-1) of each VISUAL spectral
29!                  interval.  NSPECTV elements long.
30!     WAVEV      - Array (NSPECTV elements long) of the wavelenght
31!                  (in microns) at the center of each VISUAL spectral
32!                  interval.
33!     STELLARF   - Array (NSPECTV elements) of stellar flux (W/M^2) in
34!                  each spectral interval.  Values are for 1 AU,
35!                  scaled to the planetary distance elsewhere.
36!     TAURAY     - Array (NSPECTV elements) of the pressure-independent
37!                  part of Rayleigh scattering optical depth.
38!     TAURAYVAR  - Array (NSPECTV elements) of the pressure-independent
39!                  part of Rayleigh scattering optical depth for the variable gas.
40!     FZEROI     - Fraction of zeros in the IR CO2 k-coefficients, for
41!                  each temperature, pressure, and spectral interval
42!     FZEROV     - Fraction of zeros in the VISUAL CO2 k-coefficients, for
43!                  each temperature, pressure, and spectral interval
44!
45!     AEROSOL RADIATIVE OPTICAL CONSTANTS
46!
47!   Shortwave
48!   ~~~~~~~~~
49!
50! For the "naerkind" kind of aerosol radiative properties :
51! QVISsQREF  :  Qext / Qext("longrefvis")
52! omegavis   :  single scattering albedo
53! gvis       :  assymetry factor
54!
55!   Longwave
56!   ~~~~~~~~
57!
58! For the "naerkind" kind of aerosol radiative properties :
59! QIRsQREF :  Qext / Qext("longrefvis")
60! omegaIR  :  mean single scattering albedo
61! gIR      :  mean assymetry factor
62
63      REAL*8 BWNI(L_NSPECTI+1), WNOI(L_NSPECTI), DWNI(L_NSPECTI), WAVEI(L_NSPECTI) !BWNI read by master in setspi
64      REAL*8 BWNV(L_NSPECTV+1), WNOV(L_NSPECTV), DWNV(L_NSPECTV), WAVEV(L_NSPECTV) !BWNV read by master in setspv
65      REAL*8 STELLARF(L_NSPECTV), TAURAY(L_NSPECTV), TAURAYVAR(L_NSPECTV)
66!$OMP THREADPRIVATE(WNOI,DWNI,WAVEI,&
67        !$OMP WNOV,DWNV,WAVEV,&
68        !$OMP STELLARF,TAURAY,TAURAYVAR)
69
70      REAL*8 blami(L_NSPECTI+1)
71      REAL*8 blamv(L_NSPECTV+1) ! these are needed by suaer.F90
72!$OMP THREADPRIVATE(blami,blamv)
73
74      !! AS: introduced to avoid doing same computations again for continuum
75      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indi
76      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indv
77!$OMP THREADPRIVATE(indi,indv)
78
79      !!! ALLOCATABLE STUFF SO THAT DIMENSIONS ARE READ in *.dat FILES -- AS 12/2011 
80      REAL*8, DIMENSION(:,:,:,:,:), ALLOCATABLE :: gasi, gasv
81      REAL*8, DIMENSION(:), ALLOCATABLE :: PGASREF, TGASREF, WREFVAR, PFGASREF, GWEIGHT
82      real*8 FZEROI(L_NSPECTI)
83      real*8 FZEROV(L_NSPECTV)
84      real*8 pgasmin, pgasmax
85      real*8 tgasmin, tgasmax
86!$OMP THREADPRIVATE(gasi,gasv,&  !wrefvar,pgasref,tgasref,pfgasref read by master in sugas_corrk
87        !$OMP FZEROI,FZEROV)     !pgasmin,pgasmax,tgasmin,tgasmax read by master in sugas_corrk
88
89      real QVISsQREF(L_NSPECTV,naerkind,nsizemax)
90      real omegavis(L_NSPECTV,naerkind,nsizemax)
91      real gvis(L_NSPECTV,naerkind,nsizemax)
92      real QIRsQREF(L_NSPECTI,naerkind,nsizemax)
93      real omegair(L_NSPECTI,naerkind,nsizemax)
94      real gir(L_NSPECTI,naerkind,nsizemax)
95!$OMP THREADPRIVATE(QVISsQREF,omegavis,gvis,QIRsQREF,omegair,gir)
96
97
98! Reference wavelengths used to compute reference optical depth (m)
99! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
100
101      REAL lamrefir(naerkind),lamrefvis(naerkind)
102
103! Actual number of grain size classes in each domain for a
104!   given aerosol:
105
106      INTEGER          :: nsize(naerkind,2)
107
108! Particle size axis (depend on the kind of aerosol and the
109!   radiation domain)
110
111      DOUBLE PRECISION :: radiustab(naerkind,2,nsizemax)
112!$OMP THREADPRIVATE(lamrefir,lamrefvis,radiustab) !nsize read by suaer_corrk
113
114! Extinction coefficient at reference wavelengths;
115!   These wavelengths are defined in aeroptproperties, and called
116!   longrefvis and longrefir.
117
118      REAL :: QREFvis(naerkind,nsizemax)
119      REAL :: QREFir(naerkind,nsizemax)
120!      REAL :: omegaREFvis(naerkind,nsizemax)
121      REAL :: omegaREFir(naerkind,nsizemax)
122
123      REAL,SAVE :: tstellar ! Stellar brightness temperature (SW)
124
125      real*8,save :: planckir(L_NSPECTI,NTstop-NTstar+1)
126
127      real*8,save :: PTOP
128
129      real*8,parameter :: UBARI = 0.5D0
130
131!$OMP THREADPRIVATE(QREFvis,QREFir,omegaREFir,&         ! gweight read by master in sugas_corrk
132                !$OMP tstellar,planckir,PTOP)
133
134!     If the gas optical depth (top to the surface) is less than
135!     this value, we place that Gauss-point into the "zeros"
136!     channel.
137      real*8, parameter :: TLIMIT =  1.0D-30
138
139!     Factor to convert pressures from millibars to Pascals
140      real*8, parameter :: SCALEP = 1.00D+2
141
142      real*8, parameter :: sigma = 5.67032D-8
143      real*8, parameter :: grav = 6.672E-11
144
145      real*8,save :: Cmk
146      real*8,save :: glat_ig
147!$OMP THREADPRIVATE(Cmk,glat_ig)
148
149      ! extinction of incoming sunlight (Saturn's rings, eclipses, etc...)
150      REAL, DIMENSION(:), ALLOCATABLE ,SAVE :: eclipse
151
152      !Latitude-dependent gravity
153      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: glat
154!$OMP THREADPRIVATE(glat,eclipse)
155
156end module radcommon_h
Note: See TracBrowser for help on using the repository browser.