source: lmdz_wrf/trunk/WRFV3/phys/module_cam_physconst.F @ 1361

Last change on this file since 1361 was 1, checked in by lfita, 10 years ago
  • -- --- Opening of the WRF+LMDZ coupling repository --- -- -

WRF: version v3.3
LMDZ: version v1818

More details in:

File size: 9.9 KB
Line 
1#define WRF_PORT
2
3module physconst
4!------------------------------------------------------------------------
5! Based on physconst.F90 from CAM
6! Ported to WRF by William.Gustafson@pnl.gov, Nov. 2009
7! Updated to version from CESM 1.0.1 Nov. 2010
8!------------------------------------------------------------------------
9
10   ! Physical constants.  Use CCSM shared values whenever available.
11
12   use shr_kind_mod, only: r8 => shr_kind_r8
13   use shr_const_mod, only: shr_const_g,      shr_const_stebol, shr_const_tkfrz,  &
14                            shr_const_mwdair, shr_const_rdair,  shr_const_mwwv,   &
15                            shr_const_latice, shr_const_latvap, shr_const_cpdair, &
16                            shr_const_rhofw,  shr_const_cpwv,   shr_const_rgas,   &
17                            shr_const_karman, shr_const_pstd,   shr_const_rhodair,&
18                            shr_const_avogad, shr_const_boltz,  shr_const_cpfw,   &
19                            shr_const_rwv,    shr_const_zvir,   shr_const_pi,     &
20                            shr_const_rearth, shr_const_sday,   shr_const_cday,   &
21                            shr_const_spval         
22   implicit none
23   private
24   public  :: physconst_readnl
25   save
26   ! Constants based off share code or defined in physconst
27
28   real(r8), public, parameter :: avogad      = shr_const_avogad     ! Avogadro's number (molecules/kmole)
29   real(r8), public, parameter :: boltz       = shr_const_boltz      ! Boltzman's constant (J/K/molecule)
30   real(r8), public, parameter :: cday        = shr_const_cday       ! sec in calendar day ~ sec
31   real(r8), public, parameter :: cpair       = shr_const_cpdair     ! specific heat of dry air (J/K/kg)
32   real(r8), public, parameter :: cpliq       = shr_const_cpfw       ! specific heat of fresh h2o (J/K/kg)
33   real(r8), public, parameter :: karman      = shr_const_karman     ! Von Karman constant
34   real(r8), public, parameter :: latice      = shr_const_latice     ! Latent heat of fusion (J/kg)
35   real(r8), public, parameter :: latvap      = shr_const_latvap     ! Latent heat of vaporization (J/kg)
36   real(r8), public, parameter :: pi          = shr_const_pi         ! 3.14...
37   real(r8), public, parameter :: pstd        = shr_const_pstd       ! Standard pressure (Pascals)
38   real(r8), public, parameter :: r_universal = shr_const_rgas       ! Universal gas constant (J/K/kmol)
39   real(r8), public, parameter :: rhoh2o      = shr_const_rhofw      ! Density of liquid water (STP)
40   real(r8), public, parameter :: spval       = shr_const_spval      !special value
41   real(r8), public, parameter :: stebol      = shr_const_stebol     ! Stefan-Boltzmann's constant (W/m^2/K^4)
42   real(r8), public, parameter :: tmelt       = shr_const_tkfrz      ! Freezing point of water (K)
43
44   real(r8), public, parameter :: c0          = 2.99792458e8_r8      ! Speed of light in a vacuum (m/s)
45   real(r8), public, parameter :: planck      = 6.6260755e-34_r8     ! Planck's constant (J.s)
46
47   ! Molecular weights
48   real(r8), public, parameter :: mwco2       =  44._r8             ! molecular weight co2
49   real(r8), public, parameter :: mwn2o       =  44._r8             ! molecular weight n2o
50   real(r8), public, parameter :: mwch4       =  16._r8             ! molecular weight ch4
51   real(r8), public, parameter :: mwf11       = 136._r8             ! molecular weight cfc11
52   real(r8), public, parameter :: mwf12       = 120._r8             ! molecular weight cfc12
53   real(r8), public, parameter :: mwo3        =  48._r8             ! molecular weight O3
54   real(r8), public, parameter :: mwso2       =  64._r8
55   real(r8), public, parameter :: mwso4       =  96._r8
56   real(r8), public, parameter :: mwh2o2      =  34._r8
57   real(r8), public, parameter :: mwdms       =  62._r8
58
59
60   ! modifiable physical constants for aquaplanet
61
62   real(r8), public           :: gravit       = shr_const_g     ! gravitational acceleration (m/s**2)
63   real(r8), public           :: sday         = shr_const_sday  ! sec in siderial day ~ sec
64   real(r8), public           :: mwh2o        = shr_const_mwwv  ! molecular weight h2o
65   real(r8), public           :: cpwv         = shr_const_cpwv  ! specific heat of water vapor (J/K/kg)
66   real(r8), public           :: mwdry        = shr_const_mwdair! molecular weight dry air
67   real(r8), public           :: rearth       = shr_const_rearth! radius of earth (m)
68
69!---------------  Variables below here are derived from those above -----------------------
70
71   real(r8), public           :: rga          = 1._r8/shr_const_g                 ! reciprocal of gravit
72   real(r8), public           :: ra           = 1._r8/shr_const_rearth            ! reciprocal of earth radius
73   real(r8), public           :: omega        = 2.0_R8*shr_const_pi/shr_const_sday! earth rot ~ rad/sec
74   real(r8), public           :: rh2o         = shr_const_rgas/shr_const_mwwv     ! Water vapor gas constant ~ J/K/kg
75   real(r8), public           :: rair         = shr_const_rdair   ! Dry air gas constant     ~ J/K/kg
76   real(r8), public           :: epsilo       = shr_const_mwwv/shr_const_mwdair   ! ratio of h2o to dry air molecular weights
77   real(r8), public           :: zvir         = (shr_const_rwv/shr_const_rdair)-1.0_R8 ! (rh2o/rair) - 1
78   real(r8), public           :: cpvir        = (shr_const_cpwv/shr_const_cpdair)-1.0_R8 ! CPWV/CPDAIR - 1.0
79   real(r8), public           :: rhodair      = shr_const_pstd/(shr_const_rdair*shr_const_tkfrz)
80   real(r8), public           :: cappa        = (shr_const_rgas/shr_const_mwdair)/shr_const_cpdair  ! R/Cp
81   real(r8), public           :: ez           ! Coriolis expansion coeff -> omega/sqrt(0.375)   
82   real(r8), public           :: Cpd_on_Cpv   = shr_const_cpdair/shr_const_cpwv
83                         
84
85!================================================================================================
86contains
87!================================================================================================
88
89   ! Read namelist variables.
90   subroutine physconst_readnl(nlfile)
91#ifndef WRF_PORT
92      use namelist_utils,  only: find_group_name
93      use units,           only: getunit, freeunit
94      use mpishorthand
95      use spmd_utils,      only: masterproc
96      use abortutils,      only: endrun
97      use cam_logfile,     only: iulog
98#endif     
99      character(len=*), intent(in) :: nlfile  ! filepath for file containing namelist input
100#ifndef WRF_PORT
101      ! Local variables
102      integer :: unitn, ierr
103      character(len=*), parameter :: subname = 'physconst_readnl'
104      logical       newg, newsday, newmwh2o, newcpwv, newmwdry, newrearth
105
106      ! Physical constants needing to be reset (ie. for aqua planet experiments)
107      namelist /physconst_nl/  cpwv, gravit, mwdry, mwh2o, rearth, sday
108
109      !-----------------------------------------------------------------------------
110
111      if (masterproc) then
112         unitn = getunit()
113         open( unitn, file=trim(nlfile), status='old' )
114         call find_group_name(unitn, 'physconst_nl', status=ierr)
115         if (ierr == 0) then
116            read(unitn, physconst_nl, iostat=ierr)
117            if (ierr /= 0) then
118               call endrun(subname // ':: ERROR reading namelist')
119            end if
120         end if
121         close(unitn)
122         call freeunit(unitn)
123      end if
124
125#ifdef SPMD
126      ! Broadcast namelist variables
127      call mpibcast(cpwv,      1,                   mpir8,   0, mpicom)
128      call mpibcast(gravit,    1,                   mpir8,   0, mpicom)
129      call mpibcast(mwdry,     1,                   mpir8,   0, mpicom)
130      call mpibcast(mwh2o,     1,                   mpir8,   0, mpicom)
131      call mpibcast(rearth,    1,                   mpir8,   0, mpicom)
132      call mpibcast(sday,      1,                   mpir8,   0, mpicom)
133#endif
134
135
136     
137      newg     =  gravit .ne. shr_const_g
138      newsday  =  sday   .ne. shr_const_sday
139      newmwh2o =  mwh2o  .ne. shr_const_mwwv
140      newcpwv  =  cpwv   .ne. shr_const_cpwv
141      newmwdry =  mwdry  .ne. shr_const_mwdair
142      newrearth=  rearth .ne. shr_const_rearth
143     
144     
145     
146      if (newg .or. newsday .or. newmwh2o .or. newcpwv .or. newmwdry .or. newrearth) then
147         if (masterproc) then
148            write(iulog,*)'****************************************************************************'
149            write(iulog,*)'***    New Physical Constant Values set via namelist                     ***'
150            write(iulog,*)'***                                                                      ***'
151            write(iulog,*)'***    Physical Constant    Old Value                  New Value         ***'
152            if (newg)       write(iulog,*)'***       GRAVITY   ',shr_const_g,gravit,'***'
153            if (newsday)    write(iulog,*)'***       SDAY      ',shr_const_sday,sday,'***'
154            if (newmwh2o)   write(iulog,*)'***       MWH20     ',shr_const_mwwv,mwh2o,'***'
155            if (newcpwv)    write(iulog,*)'***       CPWV      ',shr_const_cpwv,cpwv,'***'
156            if (newmwdry)   write(iulog,*)'***       MWDRY     ',shr_const_mwdair,mwdry,'***'
157            if (newrearth)  write(iulog,*)'***       REARTH    ',shr_const_rearth,rearth,'***'
158            write(iulog,*)'****************************************************************************'
159         end if
160         rga         = 1._r8/gravit
161         ra          = 1._r8/rearth
162         omega       = 2.0_R8*pi/sday
163         cpvir       = cpwv/cpair - 1._r8
164         epsilo      = mwh2o/mwdry     
165         
166         !  rair and rh2o have to be defined before any of the variables that use them
167         
168         rair        = r_universal/mwdry
169         rh2o        = r_universal/mwh2o 
170         
171         cappa       = rair/cpair       
172         rhodair     = pstd/(rair*tmelt)
173         zvir        =  (rh2o/rair)-1.0_R8
174         ez          = omega / sqrt(0.375_r8)
175         Cpd_on_Cpv  = cpair/cpwv
176         
177      else     
178         ez          = omega / sqrt(0.375_r8)
179      end if
180#endif     
181    end subroutine physconst_readnl
182
183  end module physconst
Note: See TracBrowser for help on using the repository browser.