Line data Source code
1 : module camsrfexch
2 :
3 : !-----------------------------------------------------------------------
4 : ! Module to handle data that is exchanged between the CAM atmosphere
5 : ! model and the surface models (land, sea-ice, and ocean).
6 : !-----------------------------------------------------------------------
7 :
8 : use shr_kind_mod, only: r8 => shr_kind_r8, r4 => shr_kind_r4
9 : use constituents, only: pcnst
10 : use ppgrid, only: pcols, begchunk, endchunk
11 : use phys_grid, only: get_ncols_p, phys_grid_initialized
12 : use infnan, only: posinf, assignment(=)
13 : use cam_abortutils, only: endrun
14 : use cam_logfile, only: iulog
15 : use srf_field_check, only: active_Sl_ram1, active_Sl_fv, active_Sl_soilw, &
16 : active_Fall_flxdst1, active_Fall_flxvoc, active_Fall_flxfire, &
17 : active_Faxa_nhx, active_Faxa_noy
18 :
19 :
20 :
21 : implicit none
22 : private
23 :
24 : ! Public interfaces
25 : public atm2hub_alloc ! Atmosphere to surface data allocation method
26 : public hub2atm_alloc ! Merged hub surface to atmosphere data allocation method
27 : public atm2hub_deallocate
28 : public hub2atm_deallocate
29 : public cam_export
30 :
31 : ! Public data types
32 : public cam_out_t ! Data from atmosphere
33 : public cam_in_t ! Merged surface data
34 :
35 : !---------------------------------------------------------------------------
36 : ! This is the data that is sent from the atmosphere to the surface models
37 : !---------------------------------------------------------------------------
38 :
39 : type cam_out_t
40 : integer :: lchnk ! chunk index
41 : integer :: ncol ! number of columns in chunk
42 : real(r8) :: tbot(pcols) ! bot level temperature
43 : real(r8) :: zbot(pcols) ! bot level height above surface
44 : real(r8) :: topo(pcols) ! surface topographic height (m)
45 : real(r8) :: ubot(pcols) ! bot level u wind
46 : real(r8) :: vbot(pcols) ! bot level v wind
47 : real(r8) :: qbot(pcols,pcnst) ! bot level specific humidity
48 : real(r8) :: pbot(pcols) ! bot level pressure
49 : real(r8) :: rho(pcols) ! bot level density
50 : real(r8) :: netsw(pcols) !
51 : real(r8) :: flwds(pcols) !
52 : real(r8) :: precsc(pcols) !
53 : real(r8) :: precsl(pcols) !
54 : real(r8) :: precc(pcols) !
55 : real(r8) :: precl(pcols) !
56 : real(r8) :: soll(pcols) !
57 : real(r8) :: sols(pcols) !
58 : real(r8) :: solld(pcols) !
59 : real(r8) :: solsd(pcols) !
60 : real(r8) :: thbot(pcols) !
61 : real(r8) :: co2prog(pcols) ! prognostic co2
62 : real(r8) :: co2diag(pcols) ! diagnostic co2
63 : real(r8) :: ozone(pcols) ! surface ozone concentration (mole/mole)
64 : real(r8) :: lightning_flash_freq(pcols) ! cloud-to-ground lightning flash frequency (/min)
65 : real(r8) :: psl(pcols)
66 : real(r8) :: bcphiwet(pcols) ! wet deposition of hydrophilic black carbon
67 : real(r8) :: bcphidry(pcols) ! dry deposition of hydrophilic black carbon
68 : real(r8) :: bcphodry(pcols) ! dry deposition of hydrophobic black carbon
69 : real(r8) :: ocphiwet(pcols) ! wet deposition of hydrophilic organic carbon
70 : real(r8) :: ocphidry(pcols) ! dry deposition of hydrophilic organic carbon
71 : real(r8) :: ocphodry(pcols) ! dry deposition of hydrophobic organic carbon
72 : real(r8) :: dstwet1(pcols) ! wet deposition of dust (bin1)
73 : real(r8) :: dstdry1(pcols) ! dry deposition of dust (bin1)
74 : real(r8) :: dstwet2(pcols) ! wet deposition of dust (bin2)
75 : real(r8) :: dstdry2(pcols) ! dry deposition of dust (bin2)
76 : real(r8) :: dstwet3(pcols) ! wet deposition of dust (bin3)
77 : real(r8) :: dstdry3(pcols) ! dry deposition of dust (bin3)
78 : real(r8) :: dstwet4(pcols) ! wet deposition of dust (bin4)
79 : real(r8) :: dstdry4(pcols) ! dry deposition of dust (bin4)
80 : real(r8), pointer, dimension(:) :: nhx_nitrogen_flx ! nitrogen deposition fluxes (kgN/m2/s)
81 : real(r8), pointer, dimension(:) :: noy_nitrogen_flx ! nitrogen deposition fluxes (kgN/m2/s)
82 : end type cam_out_t
83 :
84 : !---------------------------------------------------------------------------
85 : ! This is the merged state of sea-ice, land and ocean surface parameterizations
86 : !---------------------------------------------------------------------------
87 :
88 : type cam_in_t
89 : integer :: lchnk ! chunk index
90 : integer :: ncol ! number of active columns
91 : real(r8) :: asdir(pcols) ! albedo: shortwave, direct
92 : real(r8) :: asdif(pcols) ! albedo: shortwave, diffuse
93 : real(r8) :: aldir(pcols) ! albedo: longwave, direct
94 : real(r8) :: aldif(pcols) ! albedo: longwave, diffuse
95 : real(r8) :: lwup(pcols) ! longwave up radiative flux
96 : real(r8) :: lhf(pcols) ! latent heat flux
97 : real(r8) :: shf(pcols) ! sensible heat flux
98 : real(r8) :: wsx(pcols) ! surface u-stress (N)
99 : real(r8) :: wsy(pcols) ! surface v-stress (N)
100 : real(r8) :: tref(pcols) ! ref height surface air temp
101 : real(r8) :: qref(pcols) ! ref height specific humidity
102 : real(r8) :: u10(pcols) ! 10m wind speed
103 : real(r8) :: ugustOut(pcols) ! gustiness added
104 : real(r8) :: u10withGusts(pcols) ! 10m wind speed with gusts added
105 : real(r8) :: ts(pcols) ! merged surface temp
106 : real(r8) :: sst(pcols) ! sea surface temp
107 : real(r8) :: snowhland(pcols) ! snow depth (liquid water equivalent) over land
108 : real(r8) :: snowhice(pcols) ! snow depth over ice
109 : real(r8) :: fco2_lnd(pcols) ! co2 flux from lnd
110 : real(r8) :: fco2_ocn(pcols) ! co2 flux from ocn
111 : real(r8) :: fdms(pcols) ! dms flux
112 : real(r8) :: landfrac(pcols) ! land area fraction
113 : real(r8) :: icefrac(pcols) ! sea-ice areal fraction
114 : real(r8) :: ocnfrac(pcols) ! ocean areal fraction
115 : real(r8) :: cflx(pcols,pcnst) ! constituent flux (emissions)
116 : real(r8) :: ustar(pcols) ! atm/ocn saved version of ustar
117 : real(r8) :: re(pcols) ! atm/ocn saved version of re
118 : real(r8) :: ssq(pcols) ! atm/ocn saved version of ssq
119 : real(r8), pointer, dimension(:) :: ram1 !aerodynamical resistance (s/m) (pcols)
120 : real(r8), pointer, dimension(:) :: fv !friction velocity (m/s) (pcols)
121 : real(r8), pointer, dimension(:) :: soilw !volumetric soil water (m3/m3)
122 : real(r8), pointer, dimension(:,:) :: depvel ! deposition velocities
123 : real(r8), pointer, dimension(:,:) :: dstflx ! dust fluxes
124 : real(r8), pointer, dimension(:,:) :: meganflx ! MEGAN fluxes
125 : real(r8), pointer, dimension(:,:) :: fireflx ! wild fire emissions
126 : real(r8), pointer, dimension(:) :: fireztop ! wild fire emissions vert distribution top
127 : end type cam_in_t
128 :
129 : !===============================================================================
130 : CONTAINS
131 : !===============================================================================
132 :
133 1536 : subroutine hub2atm_alloc( cam_in )
134 :
135 : ! Allocate space for the surface to atmosphere data type. And initialize
136 : ! the values.
137 :
138 : use shr_drydep_mod, only: n_drydep
139 : use shr_megan_mod, only: shr_megan_mechcomps_n
140 : use shr_fire_emis_mod,only: shr_fire_emis_mechcomps_n
141 :
142 : ! ARGUMENTS:
143 : type(cam_in_t), pointer :: cam_in(:) ! Merged surface state
144 :
145 : ! LOCAL VARIABLES:
146 : integer :: c ! chunk index
147 : integer :: ierror ! Error code
148 : character(len=*), parameter :: sub = 'hub2atm_alloc'
149 : !-----------------------------------------------------------------------
150 :
151 1536 : if ( .not. phys_grid_initialized() ) call endrun(sub//": phys_grid not called yet")
152 4608 : allocate (cam_in(begchunk:endchunk), stat=ierror)
153 1536 : if ( ierror /= 0 )then
154 0 : write(iulog,*) sub//': Allocation error: ', ierror
155 0 : call endrun(sub//': allocation error')
156 : end if
157 :
158 7728 : do c = begchunk,endchunk
159 6192 : nullify(cam_in(c)%ram1)
160 6192 : nullify(cam_in(c)%fv)
161 6192 : nullify(cam_in(c)%soilw)
162 6192 : nullify(cam_in(c)%depvel)
163 6192 : nullify(cam_in(c)%dstflx)
164 6192 : nullify(cam_in(c)%meganflx)
165 6192 : nullify(cam_in(c)%fireflx)
166 7728 : nullify(cam_in(c)%fireztop)
167 : enddo
168 7728 : do c = begchunk,endchunk
169 6192 : if (active_Sl_ram1) then
170 6192 : allocate (cam_in(c)%ram1(pcols), stat=ierror)
171 6192 : if ( ierror /= 0 ) call endrun(sub//': allocation error ram1')
172 : endif
173 6192 : if (active_Sl_fv) then
174 6192 : allocate (cam_in(c)%fv(pcols), stat=ierror)
175 6192 : if ( ierror /= 0 ) call endrun(sub//': allocation error fv')
176 : endif
177 6192 : if (active_Sl_soilw) then
178 0 : allocate (cam_in(c)%soilw(pcols), stat=ierror)
179 0 : if ( ierror /= 0 ) call endrun(sub//': allocation error soilw')
180 : end if
181 6192 : if (active_Fall_flxdst1) then
182 : ! Assume 4 bins from surface model ....
183 6192 : allocate (cam_in(c)%dstflx(pcols,4), stat=ierror)
184 6192 : if ( ierror /= 0 ) call endrun(sub//': allocation error dstflx')
185 : endif
186 7728 : if (active_Fall_flxvoc .and. shr_megan_mechcomps_n>0) then
187 0 : allocate (cam_in(c)%meganflx(pcols,shr_megan_mechcomps_n), stat=ierror)
188 0 : if ( ierror /= 0 ) call endrun(sub//': allocation error meganflx')
189 : endif
190 : end do
191 :
192 1536 : if (n_drydep>0) then
193 0 : do c = begchunk,endchunk
194 0 : allocate (cam_in(c)%depvel(pcols,n_drydep), stat=ierror)
195 0 : if ( ierror /= 0 ) call endrun(sub//': allocation error depvel')
196 : end do
197 : endif
198 :
199 1536 : if (active_Fall_flxfire .and. shr_fire_emis_mechcomps_n>0) then
200 0 : do c = begchunk,endchunk
201 0 : allocate(cam_in(c)%fireflx(pcols,shr_fire_emis_mechcomps_n), stat=ierror)
202 0 : if ( ierror /= 0 ) call endrun(sub//': allocation error fireflx')
203 0 : allocate(cam_in(c)%fireztop(pcols), stat=ierror)
204 0 : if ( ierror /= 0 ) call endrun(sub//': allocation error fireztop')
205 : enddo
206 : endif
207 :
208 7728 : do c = begchunk,endchunk
209 6192 : cam_in(c)%lchnk = c
210 6192 : cam_in(c)%ncol = get_ncols_p(c)
211 105264 : cam_in(c)%asdir (:) = 0._r8
212 105264 : cam_in(c)%asdif (:) = 0._r8
213 105264 : cam_in(c)%aldir (:) = 0._r8
214 105264 : cam_in(c)%aldif (:) = 0._r8
215 105264 : cam_in(c)%lwup (:) = 0._r8
216 105264 : cam_in(c)%lhf (:) = 0._r8
217 105264 : cam_in(c)%shf (:) = 0._r8
218 105264 : cam_in(c)%wsx (:) = 0._r8
219 105264 : cam_in(c)%wsy (:) = 0._r8
220 105264 : cam_in(c)%tref (:) = 0._r8
221 105264 : cam_in(c)%qref (:) = 0._r8
222 105264 : cam_in(c)%u10 (:) = 0._r8
223 105264 : cam_in(c)%ugustOut (:) = 0._r8
224 105264 : cam_in(c)%u10withGusts (:) = 0._r8
225 105264 : cam_in(c)%ts (:) = 0._r8
226 105264 : cam_in(c)%sst (:) = 0._r8
227 105264 : cam_in(c)%snowhland(:) = 0._r8
228 105264 : cam_in(c)%snowhice (:) = 0._r8
229 105264 : cam_in(c)%fco2_lnd (:) = 0._r8
230 105264 : cam_in(c)%fco2_ocn (:) = 0._r8
231 105264 : cam_in(c)%fdms (:) = 0._r8
232 6192 : cam_in(c)%landfrac (:) = posinf
233 6192 : cam_in(c)%icefrac (:) = posinf
234 6192 : cam_in(c)%ocnfrac (:) = posinf
235 :
236 6192 : if (associated(cam_in(c)%ram1)) &
237 105264 : cam_in(c)%ram1 (:) = 0.1_r8
238 6192 : if (associated(cam_in(c)%fv)) &
239 105264 : cam_in(c)%fv (:) = 0.1_r8
240 6192 : if (associated(cam_in(c)%soilw)) &
241 0 : cam_in(c)%soilw (:) = 0.0_r8
242 6192 : if (associated(cam_in(c)%dstflx)) &
243 427248 : cam_in(c)%dstflx(:,:) = 0.0_r8
244 6192 : if (associated(cam_in(c)%meganflx)) &
245 0 : cam_in(c)%meganflx(:,:) = 0.0_r8
246 :
247 321984 : cam_in(c)%cflx (:,:) = 0._r8
248 105264 : cam_in(c)%ustar (:) = 0._r8
249 105264 : cam_in(c)%re (:) = 0._r8
250 105264 : cam_in(c)%ssq (:) = 0._r8
251 6192 : if (n_drydep>0) then
252 0 : cam_in(c)%depvel (:,:) = 0._r8
253 : endif
254 7728 : if (active_Fall_flxfire .and. shr_fire_emis_mechcomps_n>0) then
255 0 : cam_in(c)%fireflx(:,:) = 0._r8
256 0 : cam_in(c)%fireztop(:) = 0._r8
257 : endif
258 : end do
259 :
260 1536 : end subroutine hub2atm_alloc
261 :
262 : !===============================================================================
263 :
264 1536 : subroutine atm2hub_alloc( cam_out )
265 :
266 : ! Allocate space for the atmosphere to surface data type. And initialize
267 : ! the values.
268 :
269 : ! ARGUMENTS:
270 : type(cam_out_t), pointer :: cam_out(:) ! Atmosphere to surface input
271 :
272 : ! LOCAL VARIABLES:
273 : integer :: c ! chunk index
274 : integer :: ierror ! Error code
275 : character(len=*), parameter :: sub = 'atm2hub_alloc'
276 : !-----------------------------------------------------------------------
277 :
278 1536 : if (.not. phys_grid_initialized()) call endrun(sub//": phys_grid not called yet")
279 4608 : allocate (cam_out(begchunk:endchunk), stat=ierror)
280 1536 : if ( ierror /= 0 )then
281 0 : write(iulog,*) sub//': Allocation error: ', ierror
282 0 : call endrun(sub//': allocation error: cam_out')
283 : end if
284 :
285 7728 : do c = begchunk,endchunk
286 6192 : cam_out(c)%lchnk = c
287 6192 : cam_out(c)%ncol = get_ncols_p(c)
288 105264 : cam_out(c)%tbot(:) = 0._r8
289 105264 : cam_out(c)%zbot(:) = 0._r8
290 105264 : cam_out(c)%topo(:) = 0._r8
291 105264 : cam_out(c)%ubot(:) = 0._r8
292 105264 : cam_out(c)%vbot(:) = 0._r8
293 321984 : cam_out(c)%qbot(:,:) = 0._r8
294 105264 : cam_out(c)%pbot(:) = 0._r8
295 105264 : cam_out(c)%rho(:) = 0._r8
296 105264 : cam_out(c)%netsw(:) = 0._r8
297 105264 : cam_out(c)%flwds(:) = 0._r8
298 105264 : cam_out(c)%precsc(:) = 0._r8
299 105264 : cam_out(c)%precsl(:) = 0._r8
300 105264 : cam_out(c)%precc(:) = 0._r8
301 105264 : cam_out(c)%precl(:) = 0._r8
302 105264 : cam_out(c)%soll(:) = 0._r8
303 105264 : cam_out(c)%sols(:) = 0._r8
304 105264 : cam_out(c)%solld(:) = 0._r8
305 105264 : cam_out(c)%solsd(:) = 0._r8
306 105264 : cam_out(c)%thbot(:) = 0._r8
307 105264 : cam_out(c)%co2prog(:) = 0._r8
308 105264 : cam_out(c)%co2diag(:) = 0._r8
309 105264 : cam_out(c)%ozone(:) = 0._r8
310 105264 : cam_out(c)%lightning_flash_freq(:) = 0._r8
311 105264 : cam_out(c)%psl(:) = 0._r8
312 105264 : cam_out(c)%bcphidry(:) = 0._r8
313 105264 : cam_out(c)%bcphodry(:) = 0._r8
314 105264 : cam_out(c)%bcphiwet(:) = 0._r8
315 105264 : cam_out(c)%ocphidry(:) = 0._r8
316 105264 : cam_out(c)%ocphodry(:) = 0._r8
317 105264 : cam_out(c)%ocphiwet(:) = 0._r8
318 105264 : cam_out(c)%dstdry1(:) = 0._r8
319 105264 : cam_out(c)%dstwet1(:) = 0._r8
320 105264 : cam_out(c)%dstdry2(:) = 0._r8
321 105264 : cam_out(c)%dstwet2(:) = 0._r8
322 105264 : cam_out(c)%dstdry3(:) = 0._r8
323 105264 : cam_out(c)%dstwet3(:) = 0._r8
324 105264 : cam_out(c)%dstdry4(:) = 0._r8
325 105264 : cam_out(c)%dstwet4(:) = 0._r8
326 :
327 6192 : nullify(cam_out(c)%nhx_nitrogen_flx)
328 6192 : allocate (cam_out(c)%nhx_nitrogen_flx(pcols), stat=ierror)
329 6192 : if ( ierror /= 0 ) call endrun(sub//': allocation error nhx_nitrogen_flx')
330 105264 : cam_out(c)%nhx_nitrogen_flx(:) = 0._r8
331 :
332 6192 : nullify(cam_out(c)%noy_nitrogen_flx)
333 6192 : allocate (cam_out(c)%noy_nitrogen_flx(pcols), stat=ierror)
334 6192 : if ( ierror /= 0 ) call endrun(sub//': allocation error noy_nitrogen_flx')
335 106800 : cam_out(c)%noy_nitrogen_flx(:) = 0._r8
336 : end do
337 :
338 1536 : end subroutine atm2hub_alloc
339 :
340 : !===============================================================================
341 :
342 1536 : subroutine atm2hub_deallocate(cam_out)
343 :
344 : type(cam_out_t), pointer :: cam_out(:) ! Atmosphere to surface input
345 : !-----------------------------------------------------------------------
346 :
347 1536 : if(associated(cam_out)) then
348 1536 : deallocate(cam_out)
349 : end if
350 1536 : nullify(cam_out)
351 :
352 1536 : end subroutine atm2hub_deallocate
353 :
354 : !===============================================================================
355 :
356 1536 : subroutine hub2atm_deallocate(cam_in)
357 :
358 : type(cam_in_t), pointer :: cam_in(:) ! Atmosphere to surface input
359 :
360 : integer :: c
361 : !-----------------------------------------------------------------------
362 :
363 1536 : if(associated(cam_in)) then
364 7728 : do c=begchunk,endchunk
365 6192 : if(associated(cam_in(c)%ram1)) then
366 6192 : deallocate(cam_in(c)%ram1)
367 6192 : nullify(cam_in(c)%ram1)
368 : end if
369 6192 : if(associated(cam_in(c)%fv)) then
370 6192 : deallocate(cam_in(c)%fv)
371 6192 : nullify(cam_in(c)%fv)
372 : end if
373 6192 : if(associated(cam_in(c)%soilw)) then
374 0 : deallocate(cam_in(c)%soilw)
375 0 : nullify(cam_in(c)%soilw)
376 : end if
377 6192 : if(associated(cam_in(c)%dstflx)) then
378 6192 : deallocate(cam_in(c)%dstflx)
379 6192 : nullify(cam_in(c)%dstflx)
380 : end if
381 6192 : if(associated(cam_in(c)%meganflx)) then
382 0 : deallocate(cam_in(c)%meganflx)
383 0 : nullify(cam_in(c)%meganflx)
384 : end if
385 7728 : if(associated(cam_in(c)%depvel)) then
386 0 : deallocate(cam_in(c)%depvel)
387 0 : nullify(cam_in(c)%depvel)
388 : end if
389 :
390 : enddo
391 :
392 1536 : deallocate(cam_in)
393 : end if
394 1536 : nullify(cam_in)
395 :
396 1536 : end subroutine hub2atm_deallocate
397 :
398 :
399 : !======================================================================
400 :
401 1495368 : subroutine cam_export(state,cam_out,pbuf)
402 :
403 : ! Transfer atmospheric fields into necessary surface data structures
404 :
405 : use physics_types, only: physics_state
406 : use ppgrid, only: pver
407 : use cam_history, only: outfld
408 : use chem_surfvals, only: chem_surfvals_get
409 : use co2_cycle, only: co2_transport, c_i
410 : use physconst, only: rair, mwdry, mwco2, gravit, mwo3
411 : use constituents, only: pcnst
412 : use physics_buffer, only: pbuf_get_index, pbuf_get_field, physics_buffer_desc
413 : use rad_constituents, only: rad_cnst_get_gas
414 : use cam_control_mod, only: simple_phys
415 :
416 : implicit none
417 :
418 : ! Input arguments
419 : type(physics_state), intent(in) :: state
420 : type (cam_out_t), intent(inout) :: cam_out
421 : type(physics_buffer_desc), pointer :: pbuf(:)
422 :
423 : ! Local variables
424 :
425 : integer :: i ! Longitude index
426 : integer :: m ! constituent index
427 : integer :: lchnk ! Chunk index
428 : integer :: ncol
429 : integer :: psl_idx
430 : integer :: prec_dp_idx, snow_dp_idx, prec_sh_idx, snow_sh_idx
431 : integer :: prec_sed_idx,snow_sed_idx,prec_pcw_idx,snow_pcw_idx
432 : integer :: srf_ozone_idx, lightning_idx
433 :
434 1495368 : real(r8), pointer :: psl(:)
435 :
436 1495368 : real(r8), pointer :: prec_dp(:) ! total precipitation from ZM convection
437 1495368 : real(r8), pointer :: snow_dp(:) ! snow from ZM convection
438 1495368 : real(r8), pointer :: prec_sh(:) ! total precipitation from Hack convection
439 1495368 : real(r8), pointer :: snow_sh(:) ! snow from Hack convection
440 1495368 : real(r8), pointer :: prec_sed(:) ! total precipitation from ZM convection
441 1495368 : real(r8), pointer :: snow_sed(:) ! snow from ZM convection
442 1495368 : real(r8), pointer :: prec_pcw(:) ! total precipitation from Hack convection
443 1495368 : real(r8), pointer :: snow_pcw(:) ! snow from Hack convection
444 1495368 : real(r8), pointer :: o3_ptr(:,:), srf_o3_ptr(:)
445 1495368 : real(r8), pointer :: lightning_ptr(:)
446 : !-----------------------------------------------------------------------
447 :
448 1495368 : lchnk = state%lchnk
449 1495368 : ncol = state%ncol
450 :
451 2990736 : psl_idx = pbuf_get_index('PSL')
452 1495368 : call pbuf_get_field(pbuf, psl_idx, psl)
453 :
454 1495368 : prec_dp_idx = pbuf_get_index('PREC_DP', errcode=i)
455 1495368 : snow_dp_idx = pbuf_get_index('SNOW_DP', errcode=i)
456 1495368 : prec_sh_idx = pbuf_get_index('PREC_SH', errcode=i)
457 1495368 : snow_sh_idx = pbuf_get_index('SNOW_SH', errcode=i)
458 1495368 : prec_sed_idx = pbuf_get_index('PREC_SED', errcode=i)
459 1495368 : snow_sed_idx = pbuf_get_index('SNOW_SED', errcode=i)
460 1495368 : prec_pcw_idx = pbuf_get_index('PREC_PCW', errcode=i)
461 1495368 : snow_pcw_idx = pbuf_get_index('SNOW_PCW', errcode=i)
462 1495368 : srf_ozone_idx = pbuf_get_index('SRFOZONE', errcode=i)
463 1495368 : lightning_idx = pbuf_get_index('LGHT_FLASH_FREQ', errcode=i)
464 :
465 1495368 : if (prec_dp_idx > 0) then
466 1495368 : call pbuf_get_field(pbuf, prec_dp_idx, prec_dp)
467 : end if
468 1495368 : if (snow_dp_idx > 0) then
469 1495368 : call pbuf_get_field(pbuf, snow_dp_idx, snow_dp)
470 : end if
471 1495368 : if (prec_sh_idx > 0) then
472 1495368 : call pbuf_get_field(pbuf, prec_sh_idx, prec_sh)
473 : end if
474 1495368 : if (snow_sh_idx > 0) then
475 1495368 : call pbuf_get_field(pbuf, snow_sh_idx, snow_sh)
476 : end if
477 1495368 : if (prec_sed_idx > 0) then
478 1495368 : call pbuf_get_field(pbuf, prec_sed_idx, prec_sed)
479 : end if
480 1495368 : if (snow_sed_idx > 0) then
481 1495368 : call pbuf_get_field(pbuf, snow_sed_idx, snow_sed)
482 : end if
483 1495368 : if (prec_pcw_idx > 0) then
484 1495368 : call pbuf_get_field(pbuf, prec_pcw_idx, prec_pcw)
485 : end if
486 1495368 : if (snow_pcw_idx > 0) then
487 1495368 : call pbuf_get_field(pbuf, snow_pcw_idx, snow_pcw)
488 : end if
489 :
490 24969168 : do i=1,ncol
491 23473800 : cam_out%tbot(i) = state%t(i,pver)
492 23473800 : cam_out%thbot(i) = state%t(i,pver) * state%exner(i,pver)
493 23473800 : cam_out%zbot(i) = state%zm(i,pver)
494 23473800 : cam_out%topo(i) = state%phis(i) / gravit
495 23473800 : cam_out%ubot(i) = state%u(i,pver)
496 23473800 : cam_out%vbot(i) = state%v(i,pver)
497 23473800 : cam_out%pbot(i) = state%pmid(i,pver)
498 23473800 : cam_out%psl(i) = psl(i)
499 24969168 : cam_out%rho(i) = cam_out%pbot(i)/(rair*cam_out%tbot(i))
500 : end do
501 5981472 : do m = 1, pcnst
502 76402872 : do i = 1, ncol
503 74907504 : cam_out%qbot(i,m) = state%q(i,pver,m)
504 : end do
505 : end do
506 :
507 24969168 : cam_out%co2diag(:ncol) = chem_surfvals_get('CO2VMR') * 1.0e+6_r8
508 1495368 : if (co2_transport()) then
509 0 : do i=1,ncol
510 0 : cam_out%co2prog(i) = state%q(i,pver,c_i(4)) * 1.0e+6_r8 *mwdry/mwco2
511 : end do
512 : end if
513 :
514 : ! get bottom layer ozone concentrations to export to surface models
515 1495368 : if (srf_ozone_idx > 0) then
516 0 : call pbuf_get_field(pbuf, srf_ozone_idx, srf_o3_ptr)
517 0 : cam_out%ozone(:ncol) = srf_o3_ptr(:ncol)
518 1495368 : else if (.not.simple_phys) then
519 1495368 : call rad_cnst_get_gas(0, 'O3', state, pbuf, o3_ptr)
520 24969168 : cam_out%ozone(:ncol) = o3_ptr(:ncol,pver) * mwdry/mwo3 ! mole/mole
521 : endif
522 :
523 : ! get cloud to ground lightning flash freq (/min) to export to surface models
524 1495368 : if (lightning_idx>0) then
525 1495368 : call pbuf_get_field(pbuf, lightning_idx, lightning_ptr)
526 24969168 : cam_out%lightning_flash_freq(:ncol) = lightning_ptr(:ncol)
527 : end if
528 :
529 : !
530 : ! Precipation and snow rates from shallow convection, deep convection and stratiform processes.
531 : ! Compute total convective and stratiform precipitation and snow rates
532 : !
533 24969168 : do i=1,ncol
534 23473800 : cam_out%precc (i) = 0._r8
535 23473800 : cam_out%precl (i) = 0._r8
536 23473800 : cam_out%precsc(i) = 0._r8
537 23473800 : cam_out%precsl(i) = 0._r8
538 23473800 : if (prec_dp_idx > 0) then
539 23473800 : cam_out%precc (i) = cam_out%precc (i) + prec_dp(i)
540 : end if
541 23473800 : if (prec_sh_idx > 0) then
542 23473800 : cam_out%precc (i) = cam_out%precc (i) + prec_sh(i)
543 : end if
544 23473800 : if (prec_sed_idx > 0) then
545 23473800 : cam_out%precl (i) = cam_out%precl (i) + prec_sed(i)
546 : end if
547 23473800 : if (prec_pcw_idx > 0) then
548 23473800 : cam_out%precl (i) = cam_out%precl (i) + prec_pcw(i)
549 : end if
550 23473800 : if (snow_dp_idx > 0) then
551 23473800 : cam_out%precsc(i) = cam_out%precsc(i) + snow_dp(i)
552 : end if
553 23473800 : if (snow_sh_idx > 0) then
554 23473800 : cam_out%precsc(i) = cam_out%precsc(i) + snow_sh(i)
555 : end if
556 23473800 : if (snow_sed_idx > 0) then
557 23473800 : cam_out%precsl(i) = cam_out%precsl(i) + snow_sed(i)
558 : end if
559 23473800 : if (snow_pcw_idx > 0) then
560 23473800 : cam_out%precsl(i) = cam_out%precsl(i) + snow_pcw(i)
561 : end if
562 :
563 : ! jrm These checks should not be necessary if they exist in the parameterizations
564 23473800 : if (cam_out%precc(i) .lt.0._r8) cam_out%precc(i)=0._r8
565 23473800 : if (cam_out%precl(i) .lt.0._r8) cam_out%precl(i)=0._r8
566 23473800 : if (cam_out%precsc(i).lt.0._r8) cam_out%precsc(i)=0._r8
567 23473800 : if (cam_out%precsl(i).lt.0._r8) cam_out%precsl(i)=0._r8
568 23473800 : if (cam_out%precsc(i).gt.cam_out%precc(i)) cam_out%precsc(i)=cam_out%precc(i)
569 24969168 : if (cam_out%precsl(i).gt.cam_out%precl(i)) cam_out%precsl(i)=cam_out%precl(i)
570 :
571 : end do
572 :
573 2990736 : end subroutine cam_export
574 :
575 1495368 : end module camsrfexch
|