source: trunk/LMDZ.GENERIC/libf/phystd/radinc_h.F90 @ 1477

Last change on this file since 1477 was 1315, checked in by milmd, 10 years ago

LMDZ.GENERIC. OpenMP directives added in generic physic. When running in pure OpenMP or hybrid OpenMP/MPI, may have some bugs with condense_cloud and wstats routines.

File size: 3.7 KB
Line 
1      module radinc_h
2
3      implicit none
4
5#include "dimensions.h"
6#include "bands.h"
7#include "scatterers.h"
8
9!======================================================================
10!
11!     RADINC.H
12!
13!     Includes for the radiation code; RADIATION LAYERS, LEVELS,
14!     number of spectral intervals. . .
15!
16!======================================================================
17
18!     RADIATION parameters
19
20!     In radiation code, layer 1 corresponds to the stratosphere.  Level
21!     1 is the top of the stratosphere.  The dummy layer is at the same
22!     temperature as the (vertically isothermal) stratosphere, and
23!     any time it is explicitly needed, the appropriate quantities will
24!     be dealt with (aka "top". . .)
25
26!     L_NLEVRAD corresponds to the surface - i.e., the GCM Level that
27!     is at the surface.  PLEV(L_NLEVRAD) = P(J,I)+PTROP,
28!     PLEV(2) = PTROP, PLEV(1) = ptop
29
30!     L_NLAYRAD is the number of radiation code layers
31!     L_NLEVRAD is the number of radiation code levels.  Level N is the
32!               top of layer N.
33!
34!     L_NSPECTI is the number of IR spectral intervals
35!     L_NSPECTV is the number of Visual(or Solar) spectral intervals
36!     L_NGAUSS  is the number of Gauss points for K-coefficients
37!               GAUSS POINT 17 (aka the last one) is the special case
38!
39!     L_NPREF   is the number of reference pressures that the
40!               k-coefficients are calculated on
41!     L_PINT    is the number of Lagrange interpolated reference
42!               pressures for the gas k-coefficients - now for a
43!               smaller p-grid than before
44!     L_NTREF   is the number of reference temperatures for the
45!               k-coefficients
46!     L_TAUMAX  is the largest optical depth - larger ones are set
47!               to this value
48!
49!     L_REFVAR  The number of different mixing ratio values for
50!               the k-coefficients. Variable component of the mixture
51!               can in princple be anything: currently it's H2O.
52!
53!     NAERKIND  The number of radiatively active aerosol types
54!
55!     NSIZEMAX  The maximum number of aerosol particle sizes
56!
57!----------------------------------------------------------------------
58
59      integer, parameter :: L_NLAYRAD  = llm
60      integer, parameter :: L_LEVELS   = 2*(llm-1)+3
61      integer, parameter :: L_NLEVRAD  = llm+1
62
63      ! These are set in sugas_corrk
64      ! [uses allocatable arrays] -- AS 12/2011
65      integer :: L_NPREF, L_NTREF, L_REFVAR, L_PINT   !L_NPREF, L_NTREF, L_REFVAR, L_PINT read by master in sugas_corrk
66
67      integer, parameter :: L_NGAUSS  = 17
68
69      integer, parameter :: L_NSPECTI = NBinfrared
70      integer, parameter :: L_NSPECTV = NBvisible
71
72!      integer, parameter :: NAERKIND  = 2 ! set in scatterers.h
73      real,    parameter :: L_TAUMAX  = 35
74
75      ! For Planck function integration:
76      ! equivalent temperatures are 1/NTfac of these values
77      integer, parameter :: NTstar = 500
78      integer, parameter :: NTstop = 15000 ! new default for all non hot Jupiter runs
79      real*8, parameter :: NTfac = 1.0D+1 
80      !integer, parameter :: NTstar = 1000
81      !integer, parameter :: NTstop = 25000
82      !real*8,parameter :: NTfac = 5.0D+1   
83      !integer, parameter :: NTstar = 2000
84      !integer, parameter :: NTstop = 50000
85      !real*8,parameter :: NTfac = 1.0D+2   
86
87      ! Maximum number of grain size classes for aerosol convolution:
88      ! This must correspond to size of largest dataset used for aerosol
89      ! optical properties in datagcm folder.
90      integer, parameter :: nsizemax = 60
91
92      character (len=100) :: corrkdir
93      save corrkdir
94!$OMP THREADPRIVATE(corrkdir)
95
96      character (len=100) :: banddir
97      save banddir
98!$OMP THREADPRIVATE(banddir)
99
100      end module radinc_h
Note: See TracBrowser for help on using the repository browser.