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

Last change on this file since 3693 was 3693, checked in by emillour, 3 months ago

Generic PCM:
Minor corrections:

  • about OpenMP in rad_common_h.F90 : (unclosed bracket in the THREADPRIVATE statement)
  • in interpolate_continuum.F90 : "filename" character string size must be specified as '*' (i.e. variable) and not a fixed number

Took the opportunity to also clean-up interpolate_continuum.F90 by
making it a module, adding some intent() statements, specifying
why saved variables are not threadprivate, get rid of useless "return"
statement at the end of routines, etc.
EM

File size: 6.7 KB
Line 
1module radcommon_h
2      use radinc_h, only: L_NSPECTI, L_NSPECTV, NTstart, NTstop, &
3                          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!     FZEROI     - Fraction of zeros in the IR CO2 k-coefficients, for
37!                  each temperature, pressure, and spectral interval
38!     FZEROV     - Fraction of zeros in the VISUAL CO2 k-coefficients, for
39!                  each temperature, pressure, and spectral interval
40!
41!     AEROSOL RADIATIVE OPTICAL CONSTANTS
42!
43!   Shortwave
44!   ~~~~~~~~~
45!
46! For the "naerkind" kind of aerosol radiative properties :
47! QVISsQREF  :  Qext / Qext("longrefvis")
48! omegavis   :  single scattering albedo
49! gvis       :  assymetry factor
50!
51!   Longwave
52!   ~~~~~~~~
53!
54! For the "naerkind" kind of aerosol radiative properties :
55! QIRsQREF :  Qext / Qext("longrefir")
56! omegaIR  :  mean single scattering albedo
57! gIR      :  mean assymetry factor
58!
59!
60! Note - QIRsQREF in the martian model was scaled to longrefvis,
61! here it is scaled to longrefir, which is actually a dummy parameter,
62! as we do not output scaled aerosol opacity ...
63!
64
65      REAL*8 BWNI(L_NSPECTI+1), WNOI(L_NSPECTI), DWNI(L_NSPECTI), WAVEI(L_NSPECTI) !BWNI read by master in setspi
66      REAL*8 BWNV(L_NSPECTV+1), WNOV(L_NSPECTV), DWNV(L_NSPECTV), WAVEV(L_NSPECTV) !BWNV read by master in setspv
67      REAL*8 STELLARF(L_NSPECTV)
68!$OMP THREADPRIVATE(WNOI,DWNI,WAVEI,&
69        !$OMP WNOV,DWNV,WAVEV,&
70        !$OMP STELLARF)
71
72      REAL*8 blami(L_NSPECTI+1)
73      REAL*8 blamv(L_NSPECTV+1) ! these are needed by suaer.F90
74!$OMP THREADPRIVATE(blami,blamv)
75
76      !! AS: introduced to avoid doing same computations again for continuum
77      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indi
78      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: indv
79!$OMP THREADPRIVATE(indi,indv)
80
81      !!! ALLOCATABLE STUFF SO THAT DIMENSIONS ARE READ in *.dat FILES -- AS 12/2011 
82      REAL*8, DIMENSION(:,:,:,:,:), ALLOCATABLE :: gasi, gasv
83      REAL*8, DIMENSION(:), ALLOCATABLE :: PGASREF, TGASREF, WREFVAR, PFGASREF, GWEIGHT
84      real*8 FZEROI(L_NSPECTI)
85      real*8 FZEROV(L_NSPECTV)
86      real*8 pgasmin, pgasmax
87      real*8 tgasmin, tgasmax
88!$OMP THREADPRIVATE(gasi,gasv,&  !wrefvar,pgasref,tgasref,pfgasref read by master in sugas_corrk
89        !$OMP FZEROI,FZEROV)     !pgasmin,pgasmax,tgasmin,tgasmax read by master in sugas_corrk
90
91      real,allocatable,save :: QVISsQREF(:,:,:)
92      real,allocatable,save :: omegavis(:,:,:)
93      real,allocatable,save :: gvis(:,:,:)
94      real,allocatable,save :: QIRsQREF(:,:,:)
95      real,allocatable,save :: omegair(:,:,:)
96      real,allocatable,save :: gir(:,:,:)
97!$OMP THREADPRIVATE(QVISsQREF,omegavis,gvis,QIRsQREF,omegair,gir)
98
99
100! Reference wavelengths used to compute reference optical depth (m)
101! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
102
103      real,allocatable,save :: lamrefir(:),lamrefvis(:)
104!$OMP THREADPRIVATE(lamrefir,lamrefvis)
105! Actual number of grain size classes in each domain for a
106!   given aerosol:
107
108      integer,allocatable,save :: nsize(:,:)
109!$OMP THREADPRIVATE(nsize) ! nsize filled by suaer_corrk
110
111! Particle size axis (depend on the kind of aerosol and the
112!   radiation domain)
113
114      double precision,allocatable,save :: radiustab(:,:,:)
115!$OMP THREADPRIVATE(radiustab)
116
117! Extinction coefficient at reference wavelengths;
118!   These wavelengths are defined in aeroptproperties, and called
119!   longrefvis and longrefir.
120
121      real,allocatable,save :: QREFvis(:,:)
122      real,allocatable,save :: QREFir(:,:)
123!      REAL :: omegaREFvis(naerkind,nsizemax)
124      real,allocatable,save :: omegaREFir(:,:)
125
126      REAL,SAVE :: tstellar ! Stellar brightness temperature (SW)
127
128      REAL*8, DIMENSION(:,:), ALLOCATABLE, SAVE :: planckir
129
130      real*8,save :: PTOP
131
132      real*8,parameter :: UBARI = 0.5D0
133
134!$OMP THREADPRIVATE(QREFvis,QREFir,omegaREFir,&         ! gweight read by master in sugas_corrk
135                !$OMP tstellar,planckir,PTOP)
136
137!     If the gas optical depth (top to the surface) is less than
138!     this value, we place that Gauss-point into the "zeros"
139!     channel.
140      real*8, parameter :: TLIMIT =  1.0D-30
141
142!     Factor to convert pressures from millibars to Pascals
143      real*8, parameter :: SCALEP = 1.00D+2
144
145      real*8, parameter :: sigma = 5.67032D-8
146      real*8, parameter :: grav = 6.672E-11
147
148      real*8,save :: Cmk
149      real*8,save :: glat_ig
150!$OMP THREADPRIVATE(Cmk,glat_ig)
151
152      ! extinction of incoming sunlight (Saturn's rings, eclipses, etc...)
153      REAL, DIMENSION(:), ALLOCATABLE ,SAVE :: eclipse
154
155      !Latitude-dependent gravity
156      REAL, DIMENSION(:), ALLOCATABLE , SAVE :: glat
157!$OMP THREADPRIVATE(glat,eclipse)
158
159contains
160     
161subroutine ini_radcommon_h(naerkind)
162  ! Initialize module variables
163  implicit none
164  integer,intent(in) :: naerkind
165
166  allocate(QVISsQREF(L_NSPECTV,naerkind,nsizemax))
167  allocate(omegavis(L_NSPECTV,naerkind,nsizemax))
168  allocate(gvis(L_NSPECTV,naerkind,nsizemax))
169  allocate(QIRsQREF(L_NSPECTI,naerkind,nsizemax))
170  allocate(omegair(L_NSPECTI,naerkind,nsizemax))
171  allocate(gir(L_NSPECTI,naerkind,nsizemax))
172 
173  allocate(lamrefir(naerkind))
174  allocate(lamrefvis(naerkind))
175  allocate(nsize(naerkind,2))
176  allocate(radiustab(naerkind,2,nsizemax))
177 
178  allocate(QREFvis(naerkind,nsizemax))
179  allocate(QREFir(naerkind,nsizemax))
180  allocate(omegaREFir(naerkind,nsizemax))
181 
182end subroutine ini_radcommon_h
183
184end module radcommon_h
Note: See TracBrowser for help on using the repository browser.