Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / convertor / nmm_interface_convertor / constants.F
bloba4b341cb83a215fd258fdb6d3fe70560bccf59e6
1 module constants
2 !$$$   module documentation block
3 !                .      .    .                                       .
4 ! module:    constants
5 !   prgmmr: treadon          org: np23                date: 2003-09-25
7 ! abstract:  This module contains the definition of various constants
8 !            used in the gsi code
10 ! program history log:
11 !   2003-09-25  treadon - original code
12 !   2004-03-02  treadon - allow global and regional constants to differ
13 !   2004-06-16  treadon - update documentation
14 !   2004-10-28  treadon - replace parameter tiny=1.e-12 with tiny_r_kind
15 !                         and tiny_single
16 !   2004-11-16  treadon - add huge_single, huge_r_kind parameters
17 !   2005-01-27  cucurull - add ione
18 !   2005-08-24  derber   - move cg_term to constants from qcmod
19 !   2006-03-07  treadon  - add rd_over_cp_mass
20 !   2006-05-18  treadon  - add huge_i_kind
21 !   2006-06-06       su  - add var-qc wgtlim, change value to 0.25 (ECMWF)
22 !   2006-07-28  derber   - add r1000
23 !   2007-03-20  rancic   - add r3600
25 ! Subroutines Included:
26 !   sub init_constants  - compute derived constants, set regional/global constants
28 ! Variable Definitions:
29 !   see below
31 ! attributes:
32 !   language: f90
33 !   machine:  ibm RS/6000 SP
35 !$$$
36   use kinds, only: r_single,r_kind,i_kind
37   implicit none
39 ! Declare constants
40   integer(i_kind) izero,ione
41   real(r_kind) rearth,grav,omega,rd,rv,cp,cv,cvap,cliq
42   real(r_kind) csol,hvap,hfus,psat,t0c,ttp,jcal,cp_mass,cg_term
43   real(r_kind) fv,deg2rad,rad2deg,pi,tiny_r_kind,huge_r_kind,huge_i_kind
44   real(r_kind) ozcon,rozcon,tpwcon,rd_over_g,rd_over_cp,g_over_rd
45   real(r_kind) amsua_clw_d1,amsua_clw_d2,constoz,zero,one,two,four
46   real(r_kind) one_tenth,quarter,three,five,rd_over_cp_mass
47   real(r_kind) rearth_equator,stndrd_atmos_ps,r1000
48   real(r_kind) r3600    
49   real(r_kind) semi_major_axis,semi_minor_axis,n_a,n_b
50   real(r_kind) eccentricity,grav_polar,grav_ratio
51   real(r_kind) grav_equator,earth_omega,grav_constant
52   real(r_kind) flattening,eccentricity_linear,somigliana
53   real(r_kind) dldt,dldti,hsub,psatk,tmix,xa,xai,xb,xbi
54   real(r_kind) eps,epsm1,omeps,wgtlim
55   real(r_kind) elocp,cpr,el2orc,cclimit,climit,epsq
56   real(r_kind) pcpeff0,pcpeff1,pcpeff2,pcpeff3,rcp,c0,delta
57   real(r_kind) h1000,factor1,factor2,rhcbot,rhctop,dx_max,dx_min,dx_inv
58   real(r_kind) h300,half,cmr,cws,ke2,row,rrow
59   real(r_single) zero_single,tiny_single,huge_single
62 ! Define constants common to global and regional applications
63 !           name     value                  description                     units
64 !           ----     -----                  -----------                     -----
65   parameter(rearth_equator= 6.37813662e6_r_kind) ! equatorial earth radius (m)
66   parameter(omega  = 7.2921e-5_r_kind)  !  angular velocity of earth       (1/s)
67   parameter(cp     = 1.0046e+3_r_kind)  !  specific heat of air @pressure  (J/kg/K)
68   parameter(cvap   = 1.8460e+3_r_kind)  !  specific heat of h2o vapor      (J/kg/K)
69   parameter(csol   = 2.1060e+3_r_kind)  !  specific heat of solid h2o (ice)(J/kg/K)
70   parameter(hvap   = 2.5000e+6_r_kind)  !  latent heat of h2o condensation (J/kg)
71   parameter(hfus   = 3.3358e+5_r_kind)  !  latent heat of h2o fusion       (J/kg)
72   parameter(psat   = 6.1078e+2_r_kind)  !  pressure at h2o triple point    (Pa)
73   parameter(t0c    = 2.7315e+2_r_kind)  !  temperature at zero celsius     (K)
74   parameter(ttp    = 2.7316e+2_r_kind)  !  temperature at h2o triple point (K)
75   parameter(jcal   = 4.1855e+0_r_kind)  !  joules per calorie              ()
76   parameter(stndrd_atmos_ps = 1013.25e2_r_kind) ! 1976 US standard atmosphere ps   (Pa)
78 ! Numeric constants
79   parameter(izero  = 0)
80   parameter(ione   = 1)
81   parameter(zero_single = 0.0_r_single)
82   parameter(zero   = 0.0_r_kind)
83   parameter(one_tenth  = 0.10_r_kind)
84   parameter(quarter= 0.25_r_kind)
85   parameter(one    = 1.0_r_kind)
86   parameter(two    = 2.0_r_kind)
87   parameter(three  = 3.0_r_kind)
88   parameter(four   = 4.0_r_kind)
89   parameter(five   = 5.0_r_kind)
90   parameter(r1000  = 1000.0_r_kind)
91   parameter(r3600  = 3600.0_r_kind)
93 ! Constants for gps refractivity
94   parameter(n_a=77.6_r_kind) !K/mb
95   parameter(n_b=3.73e+5_r_kind) !K^2/mb
97 ! Parameters below from WGS-84 model software inside GPS receivers.
98   parameter(semi_major_axis = 6378.1370e3_r_kind)    !                     (m)
99   parameter(semi_minor_axis = 6356.7523142e3_r_kind) !                     (m)
100   parameter(grav_polar = 9.8321849378_r_kind)        !                     (m/s2)
101   parameter(grav_equator = 9.7803253359_r_kind)      !                     (m/s2) 
102   parameter(earth_omega = 7.292115e-5_r_kind)        !                     (rad/s)
103   parameter(grav_constant = 3.986004418e14_r_kind)   !                     (m3/s2)
105 ! Derived geophysical constants
106   parameter(flattening = (semi_major_axis-semi_minor_axis)/semi_major_axis)!() 
107   parameter(somigliana = &
108        (semi_minor_axis/semi_major_axis) * (grav_polar/grav_equator) - one)!()
109   parameter(grav_ratio = (earth_omega*earth_omega * &
110        semi_major_axis*semi_major_axis * semi_minor_axis) / grav_constant) !()
112 ! Derived thermodynamic constants
113   parameter ( dldti = cvap-csol )
114   parameter ( hsub = hvap+hfus )
115   parameter ( psatk = psat*0.001_r_kind )
116   parameter ( tmix = ttp-20._r_kind )
117   parameter ( elocp = hvap/cp )
118   parameter ( rcp  = one/cp )
120 ! Constants used in GFS moist physics
121   parameter ( h300 = 300._r_kind )
122   parameter ( half = 0.5_r_kind )
123   parameter ( cclimit = 0.001_r_kind )
124   parameter ( climit = 1.e-20_r_kind)
125   parameter ( epsq = 2.e-12_r_kind )
126   parameter ( h1000 = 1000.0_r_kind)
127   parameter ( rhcbot=0.85_r_kind )
128   parameter ( rhctop=0.85_r_kind )
129   parameter ( dx_max=-8.8818363_r_kind )
130   parameter ( dx_min=-5.2574954_r_kind )
131   parameter ( dx_inv=one/(dx_max-dx_min) )
132   parameter ( c0=0.002_r_kind )
133   parameter ( delta=0.6077338_r_kind )
134   parameter ( pcpeff0=1.591_r_kind )
135   parameter ( pcpeff1=-0.639_r_kind )
136   parameter ( pcpeff2=0.0953_r_kind )
137   parameter ( pcpeff3=-0.00496_r_kind )
138   parameter ( cmr = one/0.0003_r_kind )
139   parameter ( cws = 0.025_r_kind )
140   parameter ( ke2 = 0.00002_r_kind )
141   parameter ( row = 1000._r_kind )
142   parameter ( rrow = one/row )
144 ! Constant used to process ozone
145   parameter ( constoz = 604229.0_r_kind)
147 ! Constants used in cloud liquid water correction for AMSU-A
148 ! brightness temperatures
149   parameter ( amsua_clw_d1 = 0.754_r_kind )
150   parameter ( amsua_clw_d2 = -2.265_r_kind )
152 ! Constants used for variational qc
153   parameter ( wgtlim = 0.25_r_kind)  ! Cutoff weight for concluding that obs has been
154                                      ! rejected by nonlinear qc. This limit is arbitrary
155                                      ! and DOES NOT affect nonlinear qc. It only affects
156                                      ! the printout which "counts" the number of obs that
157                                      ! "fail" nonlinear qc.  Observations counted as failing
158                                      ! nonlinear qc are still assimilated.  Their weight
159                                      ! relative to other observations is reduced. Changing
160                                      ! wgtlim does not alter the analysis, only
161                                      ! the nonlinear qc data "count"
163 contains
164   subroutine init_constants_derived
165 !$$$  subprogram documentation block
166 !                .      .    .                                       .
167 ! subprogram:    init_constants_derived          set derived constants
168 !     prgmmr:    treadon          org: np23           date: 2004-12-02
170 ! abstract:  This routine sets derived constants
172 ! program history log:
173 !   2004-12-02  treadon
174 !   2005-03-03  treadon - add implicit none
176 !   input argument list:
178 !   output argument list:
180 ! attributes:
181 !   language: f90
182 !   machine:  ibm rs/6000 sp
184 !$$$
185     implicit none
186     logical regional
187     real(r_kind) reradius,g,r_d,r_v,cliq_wrf
189 !   Trigonometric constants
190     pi      = acos(-one)
191     deg2rad = pi/180.0_r_kind
192     rad2deg = one/deg2rad
193     cg_term = (sqrt(two*pi))/two                  ! constant for variational qc
194     tiny_r_kind = tiny(zero)
195     huge_r_kind = huge(zero)
196     tiny_single = tiny(zero_single)
197     huge_single = huge(zero_single)
198     huge_i_kind = huge(izero)
200 !   Geophysical parameters used in conversion of geopotential to
201 !   geometric height
202     eccentricity_linear = sqrt(semi_major_axis**2 - semi_minor_axis**2)
203     eccentricity = eccentricity_linear / semi_major_axis
205     return
206   end subroutine init_constants_derived
208   subroutine init_constants(regional)
209 !$$$  subprogram documentation block
210 !                .      .    .                                       .
211 ! subprogram:    init_constants       set regional or global constants
212 !     prgmmr:    treadon          org: np23           date: 2004-03-02
214 ! abstract:  This routine sets constants specific to regional or global
215 !            applications of the gsi
217 ! program history log:
218 !   2004-03-02  treadon
219 !   2004-06-16  treadon, documentation
220 !   2004-10-28  treadon - use intrinsic TINY function to set value
221 !                         for smallest machine representable positive
222 !                         number
223 !   2004-12-03  treadon - move derived constants to init_constants_derived
224 !   2005-03-03  treadon - add implicit none
226 !   input argument list:
227 !     regional - if .true., set regional gsi constants;
228 !                otherwise (.false.), use global constants
230 !   output argument list:
232 ! attributes:
233 !   language: f90
234 !   machine:  ibm rs/6000 sp
236 !$$$
237     implicit none
238     logical regional
239     real(r_kind) reradius,g,r_d,r_v,cliq_wrf
241 !   Define regional constants here
242     if (regional) then
244 !      Name given to WRF constants
245        reradius = one/6370.e03_r_kind
246        g        = 9.81_r_kind
247        r_d      = 287.04_r_kind
248        r_v      = 461.6_r_kind
249        cliq_wrf = 4190.0_r_kind
250        cp_mass  = 1004.67_r_kind
252 !      Transfer WRF constants into unified GSI constants
253        rearth = one/reradius
254        grav   = g
255        rd     = r_d
256        rv     = r_v
257        cv     = cp-r_d
258        cliq   = cliq_wrf
259        rd_over_cp_mass = rd / cp_mass
261 !   Define global constants here
262     else
263        rearth = 6.3712e+6_r_kind
264        grav   = 9.80665e+0_r_kind
265        rd     = 2.8705e+2_r_kind
266        rv     = 4.6150e+2_r_kind
267        cv     = 7.1760e+2_r_kind
268        cliq   = 4.1855e+3_r_kind
269        cp_mass= zero
270        rd_over_cp_mass = zero
271     endif
274 !   Now define derived constants which depend on constants
275 !   which differ between global and regional applications.
277 !   Constants related to ozone assimilation
278     ozcon = grav*21.4e-9_r_kind
279     rozcon= one/ozcon
281 !   Constant used in vertical integral for precipitable water
282     tpwcon = 100.0_r_kind/grav
284 !   Derived atmospheric constants
285     fv         = rv/rd-one    ! used in virtual temperature equation 
286     dldt       = cvap-cliq
287     xa         = -(dldt/rv)
288     xai        = -(dldti/rv)
289     xb         = xa+hvap/(rv*ttp)
290     xbi        = xai+hsub/(rv*ttp)
291     eps        = rd/rv
292     epsm1      = rd/rv-one
293     omeps      = one-eps
294     factor1    = (cvap-cliq)/rv
295     factor2    = hvap/rv-factor1*t0c
296     cpr        = cp*rd
297     el2orc     = hvap*hvap/(rv*cp)
298     rd_over_g  = rd/grav
299     rd_over_cp = rd/cp
300     g_over_rd  = grav/rd
302     return
303   end subroutine init_constants
305 end module constants