source: trunk/LMDZ.VENUS/libf/phyvenus/radinc_h.F90 @ 3461

Last change on this file since 3461 was 2560, checked in by slebonnois, 3 years ago

SL: Implementation of SW computation based on generic model. Switch between this new SW module or old module that reads R. Haus tables implemented with a key (solarchoice)

File size: 4.2 KB
Line 
1module radinc_h
2
3  implicit none
4
5
6!-----------------------------------------------------------------------
7!   INCLUDE 'bands.h'
8!
9!   bands.h contains the visible & infrared bands in the model
10!
11!   NBinfrared = L_NSPECTI in the model
12!   NBvisible  = L_NSPECTV in the model
13!-----------------------------------------------------------------------
14
15      INTEGER, parameter :: NBinfrared=32
16      INTEGER, parameter :: NBvisible=36
17
18!-----------------------------------------------------------------------
19!   INCLUDE 'scatterers.h'
20!
21!======================================================================
22!
23!     RADINC.H
24!
25!     Includes for the radiation code; RADIATION LAYERS, LEVELS,
26!     number of spectral intervals. . .
27!
28!======================================================================
29
30!     RADIATION parameters
31
32!     In radiation code, layer 1 corresponds to the stratosphere.  Level
33!     1 is the top of the stratosphere.  The dummy layer is at the same
34!     temperature as the (vertically isothermal) stratosphere, and
35!     any time it is explicitly needed, the appropriate quantities will
36!     be dealt with (aka "top". . .)
37
38!     L_NLEVRAD corresponds to the surface - i.e., the GCM Level that
39!     is at the surface.  PLEV(L_NLEVRAD) = P(J,I)+PTROP,
40!     PLEV(2) = PTROP, PLEV(1) = ptop
41
42!     L_NLAYRAD is the number of radiation code layers
43!     L_NLEVRAD is the number of radiation code levels.  Level N is the
44!               top of layer N.
45!
46!     L_NSPECTI is the number of IR spectral intervals
47!     L_NSPECTV is the number of Visual(or Solar) spectral intervals
48!     L_NGAUSS  is the number of Gauss points for K-coefficients
49!               GAUSS POINT 17 (aka the last one) is the special case
50!
51!     L_NPREF   is the number of reference pressures that the
52!               k-coefficients are calculated on
53!     L_PINT    is the number of Lagrange interpolated reference
54!               pressures for the gas k-coefficients - now for a
55!               smaller p-grid than before
56!     L_NTREF   is the number of reference temperatures for the
57!               k-coefficients
58!     L_TAUMAX  is the largest optical depth - larger ones are set
59!               to this value
60!
61!     L_REFVAR  The number of different mixing ratio values for
62!               the k-coefficients. Variable component of the mixture
63!               can in princple be anything: currently it's H2O.
64!
65!     NAERKIND  The number of radiatively active aerosol types
66!
67!     NSIZEMAX  The maximum number of aerosol particle sizes
68!
69!----------------------------------------------------------------------
70
71      integer,save :: L_NLAYRAD  ! = nbp_lev ! set by ini_radinc_h
72      integer,save :: L_LEVELS   ! = 2*(nbp_lev-1)+3 ! set by ini_radinc_h
73      integer,save :: L_NLEVRAD  ! = nbp_lev+1 ! set by ini_radinc_h
74!$OMP THREADPRIVATE(L_NLAYRAD,L_LEVELS,L_NLEVRAD)
75
76      ! These are set in sugas_corrk
77      ! [uses allocatable arrays] -- AS 12/2011
78      integer :: L_NPREF, L_NTREF, L_REFVAR, L_PINT, L_NGAUSS  !L_NPREF, L_NTREF, L_REFVAR, L_PINT, L_NGAUSS read by master in sugas_corrk
79
80      integer, parameter :: L_NSPECTV = NBvisible
81
82      integer, parameter :: NAERKIND  = 5
83      real,    parameter :: L_TAUMAX  = 35
84
85      ! For Planck function integration:
86      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
87      ! Integration boundary temperatures are NTstart/NTfac and Ntstop/NTfac
88      ! -- JVO 20 : Now read boundary T and integration dT as inputs in callphys.def
89      !             NTstart, Nstop and NTfac then set by ini_radinc_h
90      !             Smart user can adjust values depending he's running hot or cold atm
91      !             Default is wide range : 30K-1500K, with 0.1K step
92      !             ->  NTstart=300, Nstop=15000, NTfac=10
93      integer :: NTstart, NTstop
94      real*8  :: NTfac
95
96      ! Maximum number of grain size classes for aerosol convolution:
97      ! This must correspond to size of largest dataset used for aerosol
98      ! optical properties in datagcm folder.
99      integer, parameter :: nsizemax = 60
100
101contains
102
103  subroutine ini_radinc_h(nbp_lev)
104  ! Initialize module variables
105  implicit none
106  integer,intent(in) :: nbp_lev
107 
108  L_NLAYRAD = nbp_lev
109  L_LEVELS  = 2*(nbp_lev-1)+3
110  L_NLEVRAD = nbp_lev+1
111
112  NTfac   = 10
113  NTstart = 300
114  NTstop  = 15000
115 
116  end subroutine
117
118end module radinc_h
Note: See TracBrowser for help on using the repository browser.