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, flx_heat, cmfmc, cmfcme, &
79 0 : zdu, rliq, rice, dlf, dlf2, rliq2, det_s, det_ice, 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) :: flx_heat(:) ! Heat flux for check_energy_chng.
96 : real(r8), intent(in) :: cmfmc(:,:) ! convective mass flux
97 : real(r8), intent(in) :: cmfcme(:,:) ! cmf condensation - evaporation
98 : real(r8), intent(in) :: zdu(:,:) ! detraining mass flux from deep convection
99 : real(r8), intent(in) :: rliq(:) ! vertical integral of liquid not yet in q(ixcldliq)
100 : real(r8), intent(in) :: rice(:) ! vertical integral of ice not yet in q(ixcldice)
101 : real(r8), intent(in) :: dlf(:,:) ! local copy of DLFZM (copy so need to output)
102 : real(r8), intent(in) :: dlf2(:,:) ! Detraining cld H20 from shallow convections
103 : real(r8), intent(in) :: rliq2(:) ! vertical integral of liquid from shallow scheme
104 : real(r8), intent(in) :: det_s(:) ! vertical integral of detrained static energy from ice
105 : real(r8), intent(in) :: det_ice(:) ! vertical integral of detrained ice
106 : real(r8), intent(in) :: net_flx(:)
107 :
108 : integer :: lchnk
109 :
110 : ! Return if the first timestep as not all fields may be filled in and this will cause a core dump
111 0 : if (is_first_step().or. is_first_restart_step()) return
112 :
113 : ! Return if not turned on
114 0 : if (cam_snapshot_before_num <= 0 .and. cam_snapshot_after_num <= 0) return ! No snapshot files are being requested
115 :
116 0 : lchnk = state%lchnk
117 :
118 0 : call outfld('tphysbc_flx_heat', flx_heat, pcols, lchnk)
119 0 : call outfld('tphysbc_cmfmc', cmfmc, pcols, lchnk)
120 0 : call outfld('tphysbc_cmfcme', cmfcme, pcols, lchnk)
121 0 : call outfld('tphysbc_zdu', zdu, pcols, lchnk)
122 0 : call outfld('tphysbc_rliq', rliq, pcols, lchnk)
123 0 : call outfld('tphysbc_rice', rice, pcols, lchnk)
124 0 : call outfld('tphysbc_dlf', dlf, pcols, lchnk)
125 0 : call outfld('tphysbc_dlf2', dlf2, pcols, lchnk)
126 0 : call outfld('tphysbc_rliq2', rliq2, pcols, lchnk)
127 0 : call outfld('tphysbc_det_s', det_s, pcols, lchnk)
128 0 : call outfld('tphysbc_det_ice', det_ice, pcols, lchnk)
129 0 : call outfld('tphysbc_net_flx', net_flx, pcols, lchnk)
130 :
131 0 : call cam_snapshot_all_outfld(file_num, state, tend, cam_in, cam_out, pbuf)
132 :
133 0 : end subroutine cam_snapshot_all_outfld_tphysbc
134 :
135 0 : subroutine cam_snapshot_all_outfld_tphysac(file_num, state, tend, cam_in, cam_out, pbuf, fh2o, surfric, obklen, flx_heat)
136 :
137 0 : use time_manager, only: is_first_step
138 :
139 : !--------------------------------------------------------
140 : ! This subroutine does the outfld calls for ALL state, tend and pbuf fields for routines in tphysac.
141 : ! It also includes the cam_in and cam_out elements which are used within CAM as well as variables which
142 : ! are local to tphysac.
143 : !--------------------------------------------------------
144 :
145 : integer, intent(in) :: file_num
146 : type(physics_state), intent(in) :: state
147 : type(physics_tend), intent(in) :: tend
148 : type(cam_in_t), intent(in) :: cam_in
149 : type(cam_out_t), intent(in) :: cam_out
150 : type(physics_buffer_desc), pointer, intent(in) :: pbuf(:)
151 : real(r8), intent(in) :: fh2o(:) ! h2o flux to balance source from methane chemistry
152 : real(r8), intent(in) :: surfric(:) ! surface friction velocity
153 : real(r8), intent(in) :: obklen(:) ! Obukhov length
154 : real(r8), intent(in) :: flx_heat(:) ! Heat flux for check_energy_chng.
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 :
171 0 : call cam_snapshot_all_outfld(file_num, state, tend, cam_in, cam_out, pbuf)
172 :
173 0 : end subroutine cam_snapshot_all_outfld_tphysac
174 :
175 0 : subroutine cam_tphysbc_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
176 :
177 : !--------------------------------------------------------
178 : ! This subroutine does the addfld calls for the misc tphysbc physics variables that are passed individually
179 : ! into physics packages
180 : !--------------------------------------------------------
181 :
182 : integer,intent(in) :: cam_snapshot_before_num, cam_snapshot_after_num
183 :
184 0 : ntphysbc_var = 0
185 :
186 : !--------------------------------------------------------
187 : ! Add the misc tphysbc variables to the output
188 : ! NOTE - flx_heat is added in tphysac
189 : !--------------------------------------------------------
190 :
191 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
192 0 : 'flx', 'tphysbc_flx_heat', 'unset', horiz_only)
193 :
194 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
195 0 : 'cmfmc', 'tphysbc_cmfmc', 'unset', 'lev')
196 :
197 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
198 0 : 'cmfcme', 'tphysbc_cmfcme', 'unset', 'lev')
199 :
200 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
201 0 : 'zdu', 'tphysbc_zdu', 'unset', 'lev')
202 :
203 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
204 0 : 'rliq', 'tphysbc_rliq', 'unset', horiz_only)
205 :
206 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
207 0 : 'rice', 'tphysbc_rice', 'unset', horiz_only)
208 :
209 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
210 0 : 'dlf', 'tphysbc_dlf', 'unset', 'lev')
211 :
212 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
213 0 : 'dlf2', 'tphysbc_dlf2', 'unset', 'lev')
214 :
215 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
216 0 : 'rliq2', 'tphysbc_rliq2', 'unset', horiz_only)
217 :
218 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
219 0 : 'det_s', 'tphysbc_det_s', 'unset', horiz_only)
220 :
221 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
222 0 : 'det_ice', 'tphysbc_det_ice', 'unset', horiz_only)
223 :
224 : call snapshot_addfld( ntphysbc_var, tphysbc_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
225 0 : 'net_flx', 'tphysbc_net_flx', 'unset', horiz_only)
226 :
227 :
228 0 : end subroutine cam_tphysbc_snapshot_init
229 :
230 0 : subroutine cam_tphysac_snapshot_init(cam_snapshot_before_num, cam_snapshot_after_num)
231 :
232 : !--------------------------------------------------------
233 : ! This subroutine does the addfld calls for the misc tphysac physics variables that are passed individually
234 : ! into physics packages
235 : !--------------------------------------------------------
236 :
237 : integer,intent(in) :: cam_snapshot_before_num, cam_snapshot_after_num
238 :
239 0 : ntphysac_var = 0
240 :
241 : !--------------------------------------------------------
242 : ! Add the misc tphysac variables to the output
243 : !--------------------------------------------------------
244 :
245 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
246 0 : 'fh2o', 'tphysac_fh2o', 'unset', horiz_only)
247 :
248 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
249 0 : 'surfric', 'tphysac_surfric', 'unset', horiz_only)
250 :
251 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
252 0 : 'obklen', 'tphysac_obklen', 'unset', horiz_only)
253 :
254 : call snapshot_addfld( ntphysac_var, tphysac_snapshot, cam_snapshot_before_num, cam_snapshot_after_num, &
255 0 : 'flx', 'tphysac_flx_heat', 'unset', horiz_only)
256 :
257 0 : end subroutine cam_tphysac_snapshot_init
258 :
259 : end module cam_snapshot
|