Scheduled Downtime
On Tuesday 24 October 2023 @ 5pm MT the forums will be in read only mode in preparation for the downtime. On Wednesday 25 October 2023 @ 5am MT, this website will be down for maintenance and expected to return online later in the morning.
Normal Operations
The forums are back online with normal operations. If you notice any issues or errors related to the forums, please reach out to help@ucar.edu

How to compute baroclinic pressure anomaly at z=0 or similar

ezaron

Ed Zaron
New Member
Hello Experts and Enthusiasts:

The definition of baroclinic pressure using a z-coordinate could be written as,

p'(z) = \int_{-H}^z \rho(s) ds - (1/H) \int_{-H}^0 \int_{-H}^z \rho(s) ds dz.

This is the "baroclinic" pressure in the sense that it has zero vertical average, so it is orthogonal to the barotropic pressure force in a Boussinesq model which assumes incompressibility in the continuity equation.

My goal is to compute p'(z=0)/(rho0*g) in MOM6, which I would like to interpret as the sea-surface height of baroclinic waves.

I notice that pbce is computed in MOM_PressureForce_Montgomery.F90, which seems to be the rho-coordinate (or integration-by-parts) equivalent of the first term in the above expression (and pbce further divides this by the water depth).

I would like to compute either p'(z=0) or a related quantity that I could compare with sea level anomaly measurements, but I am a MOM6 newbie and unsure about a few details involved in computing the vertical average of pbce:
(1) Is it correct to regard pbce(i,j,k) as being defined at the bottom interface of level k?
(2) I see some expressions for summing layer thickness in h(i,j,k), and these loops range from js-2 to je+2 and is-2 to ie+2 (I assume these are including halos). Is it necessary for me to include the halos when working with pbce(i,j,k)? (I am not going to use p'(z=0) in the dynamics, I will just write it out to the diagnostics.)
(3) I am not confident about the use of GV%H_to_Z and GV%Z_to_H and when these are needed. Is there documentation for these factors?

This is my current effort at computing p'(z=0), which I refer to ask bsl(i,j):
CS%bsl(i,j) = 0.5*CS%pbce(i,j,1)*h(i,j,k)*GV%H_to_Z
do k=2,nz ; do j=js-2,je+2 ; do i=is-2,ie+2
CS%bsl(i,j) = CS%bsl(i,j) + 0.5*(CS%pbce(i,j,k-1) + CS%pbce(i,j,k))*h(i,j,k)*GV%H_to_Z
enddo ; enddo ; enddo
do j=js-2,je+2 ; do i=is-2,ie+2
CS%bsl(i,j) = CS%pbce(i,j,nz)*G%bathyT(i,j) - CS%bsl(i,j)
enddo ; enddo

I understand that my definition of the baroclinic pressure is derived from dynamics which involve a linear surface boundary condition at z=0, and I am unclear if there may be subtleties related to the nonlinear free surface in MOM6.

All the best,
Ed
 

awallcraft

Alan Wallcraft
New Member
A cleaner version, with the initial line inside the loops:

Code:
do j=js-2,je+2 ; do i=is-2,ie+2
  CS%bsl(i,j) = 0.5*CS%pbce(i,j,1)*h(i,j,k)*GV%H_to_Z
  do k=2,nz
    CS%bsl(i,j) = CS%bsl(i,j) + 0.5*(CS%pbce(i,j,k-1) + CS%pbce(i,j,k))*h(i,j,k)*GV%H_to_Z
  enddo
  CS%bsl(i,j) = CS%pbce(i,j,nz)*G%bathyT(i,j) - CS%bsl(i,j)
enddo ; enddo

I don't know if pbce is what you want.

Another place to look for inspiration is subroutine calculate_vertical_integrals in MOM_diagnostics.F90 for the calculation of pbo. This is probably where your diagnostic calculation should end up, and it uses do j=js,je ; do i=is,ie. Note that a non-zero p_surf would be from atmospheric pressure forcing (say).

Alan
 
Top