Line data Source code
1 : module ref_pres
2 : !--------------------------------------------------------------------------
3 : !
4 : ! Provides access to reference pressures for use by the physics
5 : ! parameterizations. The pressures are provided by the dynamical core
6 : ! since it determines the grid used by the physics.
7 : !
8 : ! Note that the init method for this module is called before the init
9 : ! method in physpkg; therefore, most physics modules can use these
10 : ! reference pressures during their init phases.
11 : !
12 : !--------------------------------------------------------------------------
13 :
14 : use shr_kind_mod, only: r8=>shr_kind_r8
15 : use ppgrid, only: pver, pverp
16 : use cam_history_support, only: add_vert_coord
17 : use cam_logfile, only: iulog
18 : use error_messages, only: alloc_err
19 :
20 : implicit none
21 : public
22 : save
23 :
24 : ! Reference pressures (Pa)
25 : real(r8), protected :: pref_edge(pverp) ! Layer edges
26 : real(r8), protected :: pref_mid(pver) ! Layer midpoints
27 : real(r8), protected :: pref_mid_norm(pver) ! Layer midpoints normalized by
28 : ! surface pressure ('eta' coordinate)
29 :
30 : real(r8), protected :: ptop_ref ! Top of model
31 : real(r8), protected :: psurf_ref ! reference pressure
32 :
33 : ! Number of top levels using pure pressure representation
34 : integer, protected :: num_pr_lev
35 :
36 : ! Pressure used to set troposphere cloud physics top (Pa)
37 : real(r8), protected :: trop_cloud_top_press = 0._r8
38 : ! Top level for troposphere cloud physics
39 : integer, protected :: trop_cloud_top_lev
40 :
41 : ! Pressure used to set MAM process top (Pa)
42 : real(r8), protected :: clim_modal_aero_top_press = 0._r8
43 : ! Top level for MAM processes that impact climate
44 : integer, protected :: clim_modal_aero_top_lev
45 :
46 : ! Molecular diffusion is calculated only if the model top is below this
47 : ! pressure (Pa).
48 : real(r8), protected :: do_molec_press = 0.1_r8
49 : ! Pressure used to set bottom of molecular diffusion region (Pa).
50 : real(r8), protected :: molec_diff_bot_press = 50._r8
51 : ! Flag for molecular diffusion, and molecular diffusion level index.
52 : logical, protected :: do_molec_diff = .false.
53 : integer, protected :: nbot_molec = 0
54 :
55 : ! Data for the trop_pref coordinate. It is the target of a pointer in a hist_coord_t
56 : ! object in the cam_history_support module. It is associated by the call to add_vert_coord.
57 : real(r8), private, allocatable, target :: trop_pref(:)
58 : real(r8), private, allocatable, target :: trop_prefi(:)
59 :
60 : !====================================================================================
61 : contains
62 : !====================================================================================
63 :
64 1536 : subroutine ref_pres_readnl(nlfile)
65 :
66 : use spmd_utils, only: masterproc
67 : use cam_abortutils, only: endrun
68 : use namelist_utils, only: find_group_name
69 : use units, only: getunit, freeunit
70 : use mpishorthand
71 :
72 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
73 :
74 : ! Local variables
75 : integer :: unitn, ierr
76 : character(len=*), parameter :: subname = 'ref_pres_readnl'
77 :
78 : namelist /ref_pres_nl/ trop_cloud_top_press, clim_modal_aero_top_press,&
79 : do_molec_press, molec_diff_bot_press
80 : !-----------------------------------------------------------------------------
81 :
82 1536 : if (masterproc) then
83 2 : unitn = getunit()
84 2 : open( unitn, file=trim(nlfile), status='old' )
85 2 : call find_group_name(unitn, 'ref_pres_nl', status=ierr)
86 2 : if (ierr == 0) then
87 2 : read(unitn, ref_pres_nl, iostat=ierr)
88 2 : if (ierr /= 0) then
89 0 : call endrun(subname // ':: ERROR reading namelist')
90 : end if
91 : end if
92 2 : close(unitn)
93 2 : call freeunit(unitn)
94 :
95 : ! Check that top for modal aerosols is not lower than
96 : ! top for clouds.
97 2 : if (clim_modal_aero_top_press > trop_cloud_top_press) &
98 : call endrun("ERROR: clim_modal_aero_top press must be less &
99 0 : &than or equal to trop_cloud_top_press.")
100 : end if
101 :
102 : #ifdef SPMD
103 : ! Broadcast namelist variables
104 1536 : call mpibcast(trop_cloud_top_press, 1 , mpir8, 0, mpicom)
105 1536 : call mpibcast(clim_modal_aero_top_press, 1 , mpir8, 0, mpicom)
106 1536 : call mpibcast(do_molec_press, 1 , mpir8, 0, mpicom)
107 1536 : call mpibcast(molec_diff_bot_press, 1 , mpir8, 0, mpicom)
108 : #endif
109 :
110 1536 : end subroutine ref_pres_readnl
111 :
112 : !====================================================================================
113 :
114 1536 : subroutine ref_pres_init(pref_edge_in, pref_mid_in, num_pr_lev_in)
115 :
116 : ! Initialize reference pressures
117 :
118 : ! arguments
119 : real(r8), intent(in) :: pref_edge_in(:) ! reference pressure at layer edges (Pa)
120 : real(r8), intent(in) :: pref_mid_in(:) ! reference pressure at layer midpoints (Pa)
121 : integer, intent(in) :: num_pr_lev_in ! number of top levels using pure pressure representation
122 :
123 : ! local variables
124 : integer :: nlev
125 : integer :: istat
126 : character(len=*), parameter :: sub = 'ref_pres_init'
127 : !---------------------------------------------------------------------------
128 :
129 43008 : pref_edge = pref_edge_in
130 41472 : pref_mid = pref_mid_in
131 1536 : num_pr_lev = num_pr_lev_in
132 :
133 1536 : ptop_ref = pref_edge(1)
134 1536 : psurf_ref = pref_edge(pverp)
135 :
136 41472 : pref_mid_norm = pref_mid/psurf_ref
137 :
138 : ! Find level corresponding to the top of troposphere clouds.
139 : trop_cloud_top_lev = press_lim_idx(trop_cloud_top_press, &
140 1536 : top=.true.)
141 :
142 : ! Find level corresponding to the top for MAM processes.
143 : clim_modal_aero_top_lev = press_lim_idx(clim_modal_aero_top_press, &
144 1536 : top=.true.)
145 :
146 : ! Find level corresponding to the molecular diffusion bottom.
147 1536 : do_molec_diff = (ptop_ref < do_molec_press)
148 1536 : if (do_molec_diff) then
149 : nbot_molec = press_lim_idx(molec_diff_bot_press, &
150 0 : top=.false.)
151 : end if
152 :
153 : ! Add vertical coordinates to history file for use with outputs that are only
154 : ! computed in the subdomain bounded by the top of troposphere clouds.
155 1536 : nlev = pver - trop_cloud_top_lev + 1
156 :
157 4608 : allocate(trop_pref(nlev), stat=istat)
158 1536 : call alloc_err(istat, sub, 'trop_pref', nlev)
159 43008 : trop_pref = pref_mid(trop_cloud_top_lev:)*0.01_r8 ! convert Pa to hPa
160 :
161 : call add_vert_coord('trop_pref', nlev, 'troposphere reference pressures', &
162 1536 : 'hPa', trop_pref, positive='down')
163 :
164 4608 : allocate(trop_prefi(nlev+1), stat=istat)
165 1536 : call alloc_err(istat, sub, 'trop_prefi', nlev+1)
166 44544 : trop_prefi = pref_edge(trop_cloud_top_lev:)*0.01_r8 ! convert Pa to hPa
167 :
168 : call add_vert_coord('trop_prefi', nlev+1, 'troposphere reference pressures (interfaces)', &
169 1536 : 'hPa', trop_prefi, positive='down')
170 :
171 1536 : end subroutine ref_pres_init
172 :
173 : !====================================================================================
174 :
175 : ! Convert pressure limiters to the appropriate level.
176 3072 : pure function press_lim_idx(p, top) result(k_lim)
177 : ! Pressure
178 : real(r8), intent(in) :: p
179 : ! Is this a top or bottom limit?
180 : logical, intent(in) :: top
181 : integer :: k_lim, k
182 :
183 3072 : if (top) then
184 3072 : k_lim = pver+1
185 3072 : do k = 1, pver
186 3072 : if (pref_mid(k) > p) then
187 : k_lim = k
188 : exit
189 : end if
190 : end do
191 : else
192 0 : k_lim = 0
193 0 : do k = pver, 1, -1
194 0 : if (pref_mid(k) < p) then
195 : k_lim = k
196 : exit
197 : end if
198 : end do
199 : end if
200 :
201 3072 : end function press_lim_idx
202 :
203 : !====================================================================================
204 :
205 : end module ref_pres
|