source: trunk/LMDZ.TITAN/libf/phytitan/radinc_h.F90 @ 3567

Last change on this file since 3567 was 2366, checked in by jvatant, 5 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: 3.9 KB
Line 
1module radinc_h
2
3  implicit none
4
5  include "bands.h"
6  include "scatterers.h"
7
8!======================================================================
9!
10!     RADINC.H
11!
12!     Includes for the radiation code; RADIATION LAYERS, LEVELS,
13!     number of spectral intervals. . .
14!
15!======================================================================
16
17!     RADIATION parameters
18
19!     In radiation code, layer 1 corresponds to the stratosphere.  Level
20!     1 is the top of the stratosphere.  The dummy layer is at the same
21!     temperature as the (vertically isothermal) stratosphere, and
22!     any time it is explicitly needed, the appropriate quantities will
23!     be dealt with (aka "top". . .)
24
25!     L_NLEVRAD corresponds to the surface - i.e., the GCM Level that
26!     is at the surface.  PLEV(L_NLEVRAD) = P(J,I)+PTROP,
27!     PLEV(2) = PTROP, PLEV(1) = ptop
28
29!     L_NLAYRAD is the number of radiation code layers
30!     L_NLEVRAD is the number of radiation code levels.  Level N is the
31!               top of layer N.
32!
33!     L_NSPECTI is the number of IR spectral intervals
34!     L_NSPECTV is the number of Visual(or Solar) spectral intervals
35!     L_NGAUSS  is the number of Gauss points for K-coefficients
36!               GAUSS POINT 17 (aka the last one) is the special case
37!
38!     L_NPREF   is the number of reference pressures that the
39!               k-coefficients are calculated on
40!     L_PINT    is the number of Lagrange interpolated reference
41!               pressures for the gas k-coefficients - now for a
42!               finer p-grid than before
43!     L_NTREF   is the number of reference temperatures for the
44!               k-coefficients
45!     L_TAUMAX  is the largest optical depth - larger ones are set
46!               to this value
47!
48!     L_REFVAR 
49!               EITHER (corrk_recombin = false)
50!      _          The number of different mixing ratio values for
51!     / \         the k-coefficients. Variable component of the mixture
52!    /   \        can in princple be anything: currently it's H2O.
53!   /  !  \
54!  /_______\    OR (corrk_recombin = true)
55!                 The TOTAL number of species to recombine when we
56!                 recombine corr-k instead of pre-mixed values.
57!----------------------------------------------------------------------
58
59      integer,save :: L_NLAYRAD  ! = nbp_lev ! set by ini_radinc_h
60      integer,save :: L_LEVELS   ! = 2*(nbp_lev-1)+3 ! set by ini_radinc_h
61      integer,save :: L_NLEVRAD  ! = nbp_lev+1 ! set by ini_radinc_h
62!$OMP THREADPRIVATE(L_NLAYRAD,L_LEVELS,L_NLEVRAD)
63
64      ! These are set in sugas_corrk
65      ! [uses allocatable arrays] -- AS 12/2011
66      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
67
68      integer, parameter :: L_NSPECTI = NBinfrared
69      integer, parameter :: L_NSPECTV = NBvisible
70
71      real,    parameter :: L_TAUMAX  = 35
72
73      ! For Planck function integration:
74      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
75      ! Integration boundary temperatures are NTstart/NTfac and Ntstop/NTfac
76      ! -- JVO 20 : Now read boundary T and integration dT as inputs in callphys.def
77      !             NTstart, Nstop and NTfac then set by ini_radinc_h
78      !             Smart user can adjust values depending he's running hot or cold atm
79      !             Default is wide range : 30K-1500K, with 0.1K step
80      !             ->  NTstart=300, Nstop=15000, NTfac=10
81      integer :: NTstart, NTstop
82      real*8  :: NTfac
83
84
85contains
86
87  subroutine ini_radinc_h(nbp_lev,tplanckmin,tplanckmax,dtplanck)
88  ! Initialize module variables
89  implicit none
90  integer,intent(in) :: nbp_lev
91  real*8, intent(in) :: tplanckmin
92  real*8, intent(in) :: tplanckmax
93  real*8, intent(in) :: dtplanck
94 
95  L_NLAYRAD = nbp_lev
96  L_LEVELS  = 2*(nbp_lev-1)+3
97  L_NLEVRAD = nbp_lev+1
98
99  NTfac   = 1.D0 / dtplanck
100  NTstart = int(tplanckmin * NTfac)
101  NTstop  = int(tplanckmax * NTfac)
102 
103  end subroutine
104
105end module radinc_h
Note: See TracBrowser for help on using the repository browser.