Skip to content

Commit a8ec0d6

Browse files
Merge branch 'whannah/eam/update-zonal-mean-diagnostics' into next (PR #7438)
Update zonal mean diagnostic calculations This update was provided by Wandi Yu to fix the diagnostic zonal mean calculations, specifically the ones involving the product of two zonal mean perturbations (i.e. zm(u'v')). The previous approach was more explicit, but resulted in unphysical noise in the output. The updated method resolves the noise and matches what we expect from offline calculations. This also adds a timer around the call to phys_grid_ctem_diags() [BFB] * whannah/eam/update-zonal-mean-diagnostics: cosmetic fix cosmetic fix cosmetic fixes fix variable description fix units of zonal mean output fields add timer around phys_grid_ctem_diags update zonal mean diagnostic calculations
2 parents 5cc98cf + fc26370 commit a8ec0d6

File tree

2 files changed

+59
-48
lines changed

2 files changed

+59
-48
lines changed

components/eam/src/physics/cam/phys_grid_ctem.F90

Lines changed: 57 additions & 48 deletions
Original file line numberDiff line numberDiff line change
@@ -221,16 +221,15 @@ subroutine phys_grid_ctem_init
221221

222222
if (.not.do_tem_diags) return
223223

224-
call addfld ('PSzm',horiz_only, 'A','m s-1', 'Zonal-Mean surface pressure', gridname='ctem_zavg_phys' )
225-
call addfld ('Uzm', (/'lev'/), 'A','m s-1', 'Zonal-Mean zonal wind', gridname='ctem_zavg_phys' )
226-
call addfld ('Vzm', (/'lev'/), 'A','m s-1', 'Zonal-Mean meridional wind', gridname='ctem_zavg_phys' )
227-
call addfld ('Wzm', (/'lev'/), 'A','m s-1', 'Zonal-Mean vertical wind', gridname='ctem_zavg_phys' )
228-
call addfld ('THzm', (/'lev'/), 'A','K', 'Zonal-Mean potential temp', gridname='ctem_zavg_phys' )
229-
call addfld ('VTHzm',(/'lev'/), 'A','K m s-1','Meridional Heat Flux:', gridname='ctem_zavg_phys')
230-
call addfld ('WTHzm',(/'lev'/), 'A','K m s-1','Vertical Heat Flux:', gridname='ctem_zavg_phys')
231-
call addfld ('UVzm', (/'lev'/), 'A','m2 s-2', 'Meridional Flux of Zonal Momentum', gridname='ctem_zavg_phys')
232-
call addfld ('UWzm', (/'lev'/), 'A','m2 s-2', 'Vertical Flux of Zonal Momentum', gridname='ctem_zavg_phys')
233-
call addfld ('THphys',(/'lev'/), 'A', 'K', 'Potential temp', gridname='physgrid' )
224+
call addfld('PSzm', horiz_only, 'A','Pa', "Zonal-Mean Surface Pressure", gridname='ctem_zavg_phys' )
225+
call addfld('Uzm', (/'lev'/), 'A','m s-1', "Zonal-Mean Zonal Wind", gridname='ctem_zavg_phys' )
226+
call addfld('Vzm', (/'lev'/), 'A','m s-1', "Zonal-Mean Meridional Wind", gridname='ctem_zavg_phys' )
227+
call addfld('Wzm', (/'lev'/), 'A','Pa s-1', "Zonal-Mean Vertical Wind", gridname='ctem_zavg_phys' )
228+
call addfld('THzm', (/'lev'/), 'A','K', "Zonal-Mean Potential Temp", gridname='ctem_zavg_phys' )
229+
call addfld('VTHzm', (/'lev'/), 'A','K m s-1', "Zonal-Mean Meridional Heat Flux i.e. zonal_mean(v'theta')", gridname='ctem_zavg_phys')
230+
call addfld('WTHzm', (/'lev'/), 'A','K Pa s-1',"Zonal-Mean Vertical Heat Flux i.e. zonal_mean(w'theta')", gridname='ctem_zavg_phys')
231+
call addfld('UVzm', (/'lev'/), 'A','m2 s-2', "Zonal-Mean Meridional Flux of Zonal Momentum i.e. zonal_mean(u'v')", gridname='ctem_zavg_phys')
232+
call addfld('UWzm', (/'lev'/), 'A','Pa m s-2',"Zonal-Mean Vertical Flux of Zonal Momentum i.e. zonal_mean(u'w')", gridname='ctem_zavg_phys')
234233

235234
end subroutine phys_grid_ctem_init
236235

@@ -249,16 +248,20 @@ subroutine phys_grid_ctem_diags(phys_state)
249248
real(r8) :: v(pcols,pver,begchunk:endchunk)
250249
real(r8) :: w(pcols,pver,begchunk:endchunk)
251250

251+
real(r8) :: uv(pcols,pver,begchunk:endchunk)
252+
real(r8) :: uw(pcols,pver,begchunk:endchunk)
253+
real(r8) :: vth(pcols,pver,begchunk:endchunk)
254+
real(r8) :: wth(pcols,pver,begchunk:endchunk)
255+
252256
real(r8) :: pszm(pcols,begchunk:endchunk)
253257

254258
real(r8) :: uzm(pcols,pver,begchunk:endchunk)
255259
real(r8) :: vzm(pcols,pver,begchunk:endchunk)
256260
real(r8) :: wzm(pcols,pver,begchunk:endchunk)
257-
258-
real(r8) :: ud(pcols,pver,begchunk:endchunk)
259-
real(r8) :: vd(pcols,pver,begchunk:endchunk)
260-
real(r8) :: wd(pcols,pver,begchunk:endchunk)
261-
real(r8) :: thd(pcols,pver,begchunk:endchunk)
261+
real(r8) :: uvzm(pcols,pver,begchunk:endchunk)
262+
real(r8) :: uwzm(pcols,pver,begchunk:endchunk)
263+
real(r8) :: vthzm(pcols,pver,begchunk:endchunk)
264+
real(r8) :: wthzm(pcols,pver,begchunk:endchunk)
262265

263266
real(r8) :: uvp(pcols,pver,begchunk:endchunk)
264267
real(r8) :: uwp(pcols,pver,begchunk:endchunk)
@@ -296,32 +299,39 @@ subroutine phys_grid_ctem_diags(phys_state)
296299
u(:ncol,:,lchnk) = phys_state(lchnk)%u(:ncol,:)
297300
v(:ncol,:,lchnk) = phys_state(lchnk)%v(:ncol,:)
298301
end do
302+
! calculate the product of the variables
303+
do lchnk = begchunk,endchunk
304+
ncol = phys_state(lchnk)%ncol
305+
do k = 1,pver
306+
uv(:ncol,k,lchnk) = u(:ncol,k,lchnk) * v(:ncol,k,lchnk)
307+
uw(:ncol,k,lchnk) = u(:ncol,k,lchnk) * w(:ncol,k,lchnk)
308+
vth(:ncol,k,lchnk) = v(:ncol,k,lchnk) * theta(:ncol,k,lchnk)
309+
wth(:ncol,k,lchnk) = w(:ncol,k,lchnk) * theta(:ncol,k,lchnk)
310+
end do
311+
end do
299312

300313
! zonal means evaluated on the physics grid (3D) to be used in the deviations calculation below
301-
pszm(:,:) = zmean_fld_2D(ps(:,:))
302-
uzm(:,:,:) = zmean_fld_3D(u(:,:,:))
303-
vzm(:,:,:) = zmean_fld_3D(v(:,:,:))
304-
wzm(:,:,:) = zmean_fld_3D(w(:,:,:))
305-
thzm(:,:,:) = zmean_fld_3D(theta(:,:,:))
306-
307-
! diagnostic output
308-
do lchnk = begchunk, endchunk
309-
call outfld( 'THphys', theta(:,:,lchnk), pcols, lchnk)
310-
end do
314+
pszm(:,:) = zmean_fld_2D(ps(:,:))
315+
uzm(:,:,:) = zmean_fld_3D(u(:,:,:))
316+
vzm(:,:,:) = zmean_fld_3D(v(:,:,:))
317+
wzm(:,:,:) = zmean_fld_3D(w(:,:,:))
318+
thzm(:,:,:) = zmean_fld_3D(theta(:,:,:))
319+
uvzm(:,:,:) = zmean_fld_3D(uv(:,:,:))
320+
uwzm(:,:,:) = zmean_fld_3D(uw(:,:,:))
321+
vthzm(:,:,:) = zmean_fld_3D(vth(:,:,:))
322+
wthzm(:,:,:) = zmean_fld_3D(wth(:,:,:))
311323

312324
do lchnk = begchunk,endchunk
313325
ncol = phys_state(lchnk)%ncol
314326
do k = 1,pver
315-
! zonal deviations
316-
thd(:ncol,k,lchnk) = theta(:ncol,k,lchnk) - thzm(:ncol,k,lchnk)
317-
ud(:ncol,k,lchnk) = u(:ncol,k,lchnk) - uzm(:ncol,k,lchnk)
318-
vd(:ncol,k,lchnk) = v(:ncol,k,lchnk) - vzm(:ncol,k,lchnk)
319-
wd(:ncol,k,lchnk) = w(:ncol,k,lchnk) - wzm(:ncol,k,lchnk)
320-
! fluxes
321-
uvp(:ncol,k,lchnk) = ud(:ncol,k,lchnk) * vd(:ncol,k,lchnk)
322-
uwp(:ncol,k,lchnk) = ud(:ncol,k,lchnk) * wd(:ncol,k,lchnk)
323-
vthp(:ncol,k,lchnk) = vd(:ncol,k,lchnk) * thd(:ncol,k,lchnk)
324-
wthp(:ncol,k,lchnk) = wd(:ncol,k,lchnk) * thd(:ncol,k,lchnk)
327+
! zmean(u'v') = zmean(uv) - zmean(u)*zmean(v)
328+
! NOTE - the previous method involved explicitly calcalating the the perturbations from
329+
! the zonal mean (i.e. u - zmean(u) ), but this was found to be the source of excessive
330+
! noise in the diagnostic outputs. The new method yields much smoother results.
331+
uvp(:ncol,k,lchnk) = uvzm(:ncol,k,lchnk) - uzm(:ncol,k,lchnk) * vzm(:ncol,k,lchnk)
332+
uwp(:ncol,k,lchnk) = uwzm(:ncol,k,lchnk) - uzm(:ncol,k,lchnk) * wzm(:ncol,k,lchnk)
333+
vthp(:ncol,k,lchnk) = vthzm(:ncol,k,lchnk) - vzm(:ncol,k,lchnk) * thzm(:ncol,k,lchnk)
334+
wthp(:ncol,k,lchnk) = wthzm(:ncol,k,lchnk) - wzm(:ncol,k,lchnk) * thzm(:ncol,k,lchnk)
325335
end do
326336
end do
327337

@@ -331,9 +341,8 @@ subroutine phys_grid_ctem_diags(phys_state)
331341
call ZAobj%binAvg(vthp, vthza)
332342
call ZAobj%binAvg(wthp, wthza)
333343

334-
335-
if (any(abs(uvza)>1.e20_r8)) call endrun(prefix//'bad values in uvza')
336-
if (any(abs(uwza)>1.e20_r8)) call endrun(prefix//'bad values in uwza')
344+
if (any(abs(uvza) >1.e20_r8)) call endrun(prefix//'bad values in uvza')
345+
if (any(abs(uwza) >1.e20_r8)) call endrun(prefix//'bad values in uwza')
337346
if (any(abs(vthza)>1.e20_r8)) call endrun(prefix//'bad values in vthza')
338347
if (any(abs(wthza)>1.e20_r8)) call endrun(prefix//'bad values in wthza')
339348

@@ -344,20 +353,20 @@ subroutine phys_grid_ctem_diags(phys_state)
344353
call ZAobj%binAvg(thzm, thza)
345354

346355
if (any(abs(psza)>1.e20_r8)) call endrun(prefix//'bad values in psza')
347-
if (any(abs(uza)>1.e20_r8)) call endrun(prefix//'bad values in uza')
348-
if (any(abs(vza)>1.e20_r8)) call endrun(prefix//'bad values in vza')
349-
if (any(abs(wza)>1.e20_r8)) call endrun(prefix//'bad values in wza')
356+
if (any(abs(uza) >1.e20_r8)) call endrun(prefix//'bad values in uza')
357+
if (any(abs(vza) >1.e20_r8)) call endrun(prefix//'bad values in vza')
358+
if (any(abs(wza) >1.e20_r8)) call endrun(prefix//'bad values in wza')
350359
if (any(abs(thza)>1.e20_r8)) call endrun(prefix//'bad values in thza')
351360

352361
! diagnostic output
353362
do j = 1,nzalat
354-
call outfld('PSzm',psza(j),1,j)
355-
call outfld('Uzm',uza(j,:),1,j)
356-
call outfld('Vzm',vza(j,:),1,j)
357-
call outfld('Wzm',wza(j,:),1,j)
358-
call outfld('THzm',thza(j,:),1,j)
359-
call outfld('UVzm',uvza(j,:),1,j)
360-
call outfld('UWzm',uwza(j,:),1,j)
363+
call outfld('PSzm', psza(j), 1,j)
364+
call outfld('Uzm', uza(j,:), 1,j)
365+
call outfld('Vzm', vza(j,:), 1,j)
366+
call outfld('Wzm', wza(j,:), 1,j)
367+
call outfld('THzm', thza(j,:), 1,j)
368+
call outfld('UVzm', uvza(j,:), 1,j)
369+
call outfld('UWzm', uwza(j,:), 1,j)
361370
call outfld('VTHzm',vthza(j,:),1,j)
362371
call outfld('WTHzm',wthza(j,:),1,j)
363372
end do

components/eam/src/physics/cam/physpkg.F90

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -3235,7 +3235,9 @@ subroutine phys_timestep_init(phys_state, cam_out, pbuf2d)
32353235
if(Nudge_Model) call nudging_timestep_init(phys_state)
32363236

32373237
! Update Transformed Eularian Mean (TEM) diagnostics
3238+
call t_startf('phys_grid_ctem_diags')
32383239
call phys_grid_ctem_diags(phys_state)
3240+
call t_stopf('phys_grid_ctem_diags')
32393241

32403242
end subroutine phys_timestep_init
32413243

0 commit comments

Comments
 (0)