Line data Source code
1 : module cam_snapshot
2 : !--------------------------------------------------------
3 : ! The purpose of this module is to handle taking the "snapshot" of CAM data.
4 : !
5 : ! This module writes out ALL the state, tend and pbuf fields. It also includes the cam_in and cam_out
6 : ! fields which are used within CAM
7 : !--------------------------------------------------------
8 :
9 : use shr_kind_mod, only: r8 => shr_kind_r8
10 : use cam_history, only: addfld, add_default, outfld
11 : use cam_history, only: cam_history_snapshot_deactivate, cam_history_snapshot_activate
12 : use cam_history_support, only: horiz_only
13 : use cam_abortutils, only: endrun
14 : use physics_buffer, only: physics_buffer_desc, pbuf_get_index, pbuf_get_field, pbuf_get_field_name
15 : use physics_types, only: physics_state, physics_tend, physics_ptend
16 : use camsrfexch, only: cam_out_t, cam_in_t
17 : use ppgrid, only: pcols, begchunk, endchunk
18 : use constituents, only: pcnst
19 : use phys_control, only: phys_getopts
20 : use cam_logfile, only: iulog
21 : use cam_snapshot_common, only: snapshot_type, cam_snapshot_deactivate, cam_snapshot_all_outfld, cam_snapshot_ptend_outfld
22 : use cam_snapshot_common, only: snapshot_type, cam_state_snapshot_init, cam_cnst_snapshot_init, cam_tend_snapshot_init
23 : use cam_snapshot_common, only: cam_ptend_snapshot_init, cam_in_snapshot_init, cam_out_snapshot_init
24 : use cam_snapshot_common, only: cam_pbuf_snapshot_init, snapshot_addfld
25 :
26 : implicit none
27 :
28 : private
29 :
30 : public :: cam_snapshot_init
31 : public :: cam_snapshot_all_outfld_tphysbc, cam_snapshot_all_outfld_tphysac
32 :
33 : private :: cam_tphysbc_snapshot_init, cam_tphysac_snapshot_init
34 :
35 : integer :: ntphysbc_var
36 : integer :: ntphysac_var
37 :
38 : integer :: cam_snapshot_before_num, cam_snapshot_after_num
39 :
40 : ! Note the maximum number of variables for each type
41 : type (snapshot_type) :: tphysbc_snapshot(30)
42 : type (snapshot_type) :: tphysac_snapshot(30)
43 :
44 : contains
45 :
46 0 : subroutine cam_snapshot_init(cam_in_arr, cam_out_arr, pbuf, index)
47 :
48 :
49 : !--------------------------------------------------------
50 : ! This subroutine does the addfld calls for ALL state, tend, ptend, and pbuf fields. It also includes the cam_in and cam_out
51 : ! elements which are used within CAM
52 : !--------------------------------------------------------
53 : type(cam_in_t), intent(in) :: cam_in_arr(begchunk:endchunk)
54 : type(cam_out_t), intent(in) :: cam_out_arr(begchunk:endchunk)
55 : type(physics_buffer_desc), pointer, intent(inout) :: pbuf(:,:)
56 : integer, intent(in) :: index
57 :
58 :
59 : call phys_getopts(cam_snapshot_before_num_out = cam_snapshot_before_num, &
60 1536 : cam_snapshot_after_num_out = cam_snapshot_after_num)
61 :
62 :
63 : ! Return if not turned on
64 1536 : if (cam_snapshot_before_num <= 0 .and. cam_snapshot_after_num <= 0) return ! No snapshot files are being requested
65 :
66 0 : call cam_state_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
67 0 : call cam_cnst_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
68 0 : call cam_tend_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
69 0 : call cam_ptend_snapshot_init(cam_snapshot_after_num)
70 0 : call cam_in_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num, cam_in_arr(index))
71 : call cam_out_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num, cam_out_arr(index))
72 0 : call cam_pbuf_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num, pbuf(:,index))
73 0 : call cam_tphysac_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
74 0 : call cam_tphysbc_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
75 :
76 : end subroutine cam_snapshot_init
77 :
78 0 : subroutine cam_snapshot_all_outfld_tphysbc(file_num, state, tend, cam_in, cam_out, pbuf, cmfmc, cmfcme, &
79 0 : zdu, rliq, rice, dlf, dlf2, rliq2, net_flx)
80 :
81 : use time_manager, only: is_first_step, is_first_restart_step
82 :
83 : !--------------------------------------------------------
84 : ! This subroutine does the outfld calls for ALL state, tend and pbuf fields for routines in tphysbc.
85 : ! It also includes the cam_in and cam_out elements which are used within CAM as well as variables which
86 : ! are local to tphysbc.
87 : !--------------------------------------------------------
88 :
89 : integer, intent(in) :: file_num
90 : type(physics_state), intent(in) :: state
91 : type(physics_tend), intent(in) :: tend
92 : type(cam_in_t), intent(in) :: cam_in
93 : type(cam_out_t), intent(in) :: cam_out
94 : type(physics_buffer_desc), pointer, intent(in) :: pbuf(:)
95 : real(r8), intent(in) :: cmfmc(:,:) ! convective mass flux
96 : real(r8), intent(in) :: cmfcme(:,:) ! cmf condensation - evaporation
97 : real(r8), intent(in) :: zdu(:,:) ! detraining mass flux from deep convection
98 : real(r8), intent(in) :: rliq(:) ! vertical integral of liquid not yet in q(ixcldliq)
99 : real(r8), intent(in) :: rice(:) ! vertical integral of ice not yet in q(ixcldice)
100 : real(r8), intent(in) :: dlf(:,:) ! local copy of DLFZM (copy so need to output)
101 : real(r8), intent(in) :: dlf2(:,:) ! Detraining cld H20 from shallow convections
102 : real(r8), intent(in) :: rliq2(:) ! vertical integral of liquid from shallow scheme
103 : real(r8), intent(in) :: net_flx(:)
104 :
105 : integer :: lchnk
106 :
107 : ! Return if the first timestep as not all fields may be filled in and this will cause a core dump
108 0 : if (is_first_step().or. is_first_restart_step()) return
109 :
110 : ! Return if not turned on
111 0 : if (cam_snapshot_before_num <= 0 .and. cam_snapshot_after_num <= 0) return ! No snapshot files are being requested
112 :
113 0 : lchnk = state%lchnk
114 :
115 0 : call outfld('tphysbc_cmfmc', cmfmc, pcols, lchnk)
116 0 : call outfld('tphysbc_cmfcme', cmfcme, pcols, lchnk)
117 0 : call outfld('tphysbc_zdu', zdu, pcols, lchnk)
118 0 : call outfld('tphysbc_rliq', rliq, pcols, lchnk)
119 0 : call outfld('tphysbc_rice', rice, pcols, lchnk)
120 0 : call outfld('tphysbc_dlf', dlf, pcols, lchnk)
121 0 : call outfld('tphysbc_dlf2', dlf2, pcols, lchnk)
122 0 : call outfld('tphysbc_rliq2', rliq2, pcols, lchnk)
123 0 : call outfld('tphysbc_net_flx', net_flx, pcols, lchnk)
124 :
125 0 : call cam_snapshot_all_outfld(file_num, state, tend, cam_in, cam_out, pbuf)
126 :
127 0 : end subroutine cam_snapshot_all_outfld_tphysbc
128 :
129 0 : subroutine cam_snapshot_all_outfld_tphysac(file_num, state, tend, cam_in, cam_out, pbuf, fh2o, surfric, obklen, flx_heat, &
130 0 : cmfmc, dlf, det_s, det_ice, net_flx)
131 :
132 0 : use time_manager, only: is_first_step
133 :
134 : !--------------------------------------------------------
135 : ! This subroutine does the outfld calls for ALL state, tend and pbuf fields for routines in tphysac.
136 : ! It also includes the cam_in and cam_out elements which are used within CAM as well as variables which
137 : ! are local to tphysac.
138 : !--------------------------------------------------------
139 :
140 : integer, intent(in) :: file_num
141 : type(physics_state), intent(in) :: state
142 : type(physics_tend), intent(in) :: tend
143 : type(cam_in_t), intent(in) :: cam_in
144 : type(cam_out_t), intent(in) :: cam_out
145 : type(physics_buffer_desc), pointer, intent(in) :: pbuf(:)
146 : real(r8), intent(in) :: fh2o(:) ! h2o flux to balance source from methane chemistry
147 : real(r8), intent(in) :: surfric(:) ! surface friction velocity
148 : real(r8), intent(in) :: obklen(:) ! Obukhov length
149 : real(r8), intent(in) :: flx_heat(:) ! Heat flux for check_energy_chng.
150 : real(r8), intent(in) :: cmfmc(:,:) ! convective mass flux
151 : real(r8), intent(in) :: dlf(:,:) ! local copy of DLFZM (copy so need to output)
152 : real(r8), intent(in) :: det_s(:) ! vertical integral of detrained static energy from ice
153 : real(r8), intent(in) :: det_ice(:) ! vertical integral of detrained ice
154 : real(r8), intent(in) :: net_flx(:)
155 :
156 : integer :: lchnk
157 :
158 : ! Return if the first timestep as not all fields may be filled in and this will cause a core dump
159 0 : if (is_first_step()) return
160 :
161 : ! Return if not turned on
162 0 : if (cam_snapshot_before_num <= 0 .and. cam_snapshot_after_num <= 0) return ! No snapshot files are being requested
163 :
164 0 : lchnk = state%lchnk
165 :
166 0 : call outfld('tphysac_fh2o', fh2o, pcols, lchnk)
167 0 : call outfld('tphysac_surfric', surfric, pcols, lchnk)
168 0 : call outfld('tphysac_obklen', obklen, pcols, lchnk)
169 0 : call outfld('tphysac_flx_heat', flx_heat, pcols, lchnk)
170 0 : call outfld('tphysac_cmfmc', cmfmc, pcols, lchnk)
171 0 : call outfld('tphysac_dlf', dlf, pcols, lchnk)
172 0 : call outfld('tphysac_det_s', det_s, pcols, lchnk)
173 0 : call outfld('tphysac_det_ice', det_ice, pcols, lchnk)
174 0 : call outfld('tphysac_net_flx', net_flx, pcols, lchnk)
175 :
176 0 : call cam_snapshot_all_outfld(file_num, state, tend, cam_in, cam_out, pbuf)
177 :
178 0 : end subroutine cam_snapshot_all_outfld_tphysac
179 :
180 0 : subroutine cam_tphysbc_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
181 :
182 : !--------------------------------------------------------
183 : ! This subroutine does the addfld calls for the misc tphysbc physics variables that are passed individually
184 : ! into physics packages
185 : !--------------------------------------------------------
186 :
187 : integer,intent(in) :: cam_snapshot_before_num, cam_snapshot_after_num
188 :
189 0 : ntphysbc_var = 0
190 :
191 : !--------------------------------------------------------
192 : ! Add the misc tphysbc variables to the output
193 : ! NOTE - flx_heat is added in tphysac
194 : !--------------------------------------------------------
195 :
196 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
197 0 : 'flx', 'tphysbc_flx_heat', 'unset', horiz_only)
198 :
199 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
200 0 : 'cmfmc', 'tphysbc_cmfmc', 'unset', 'lev')
201 :
202 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
203 0 : 'cmfcme', 'tphysbc_cmfcme', 'unset', 'lev')
204 :
205 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
206 0 : 'zdu', 'tphysbc_zdu', 'unset', 'lev')
207 :
208 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
209 0 : 'rliq', 'tphysbc_rliq', 'unset', horiz_only)
210 :
211 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
212 0 : 'rice', 'tphysbc_rice', 'unset', horiz_only)
213 :
214 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
215 0 : 'dlf', 'tphysbc_dlf', 'unset', 'lev')
216 :
217 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
218 0 : 'dlf2', 'tphysbc_dlf2', 'unset', 'lev')
219 :
220 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
221 0 : 'rliq2', 'tphysbc_rliq2', 'unset', horiz_only)
222 :
223 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
224 0 : 'net_flx', 'tphysbc_net_flx', 'unset', horiz_only)
225 :
226 :
227 0 : end subroutine cam_tphysbc_snapshot_init
228 :
229 0 : subroutine cam_tphysac_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
230 :
231 : !--------------------------------------------------------
232 : ! This subroutine does the addfld calls for the misc tphysac physics variables that are passed individually
233 : ! into physics packages
234 : !--------------------------------------------------------
235 :
236 : integer,intent(in) :: cam_snapshot_before_num, cam_snapshot_after_num
237 :
238 0 : ntphysac_var = 0
239 :
240 : !--------------------------------------------------------
241 : ! Add the misc tphysac variables to the output
242 : !--------------------------------------------------------
243 :
244 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
245 0 : 'fh2o', 'tphysac_fh2o', 'unset', horiz_only)
246 :
247 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
248 0 : 'surfric', 'tphysac_surfric', 'unset', horiz_only)
249 :
250 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
251 0 : 'obklen', 'tphysac_obklen', 'unset', horiz_only)
252 :
253 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
254 0 : 'flx', 'tphysac_flx_heat', 'unset', horiz_only)
255 :
256 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
257 0 : 'cmfmc', 'tphysac_cmfmc', 'unset', 'lev')
258 :
259 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
260 0 : 'zdu', 'tphysac_zdu', 'unset', 'lev')
261 :
262 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
263 0 : 'rliq', 'tphysac_rliq', 'unset', horiz_only)
264 :
265 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
266 0 : 'dlf', 'tphysac_dlf', 'unset', 'lev')
267 :
268 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
269 0 : 'dlf2', 'tphysac_dlf2', 'unset', 'lev')
270 :
271 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
272 0 : 'rliq2', 'tphysac_rliq2', 'unset', horiz_only)
273 :
274 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
275 0 : 'det_s', 'tphysac_det_s', 'unset', horiz_only)
276 :
277 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
278 0 : 'det_ice', 'tphysac_det_ice', 'unset', horiz_only)
279 :
280 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
281 0 : 'net_flx', 'tphysac_net_flx', 'unset', horiz_only)
282 :
283 :
284 0 : end subroutine cam_tphysac_snapshot_init
285 :
286 : end module cam_snapshot
|