def ddeta( var, grid, z=None, hcoord=None, scoord=None, hboundary="extend", hfill_value=np.nan, sboundary="extend", sfill_value=np.nan, attrs=None, ): """Calculate d/deta for a variable. Inputs ------ var: DataArray, ndarray Variable to operate on. grid: xgcm.grid Grid object associated with var. z: DataArray, ndarray, optional Depth [m]. If None, use z coordinate attached to var. hcoord: string, optional Name of horizontal grid to interpolate output to. Options are 'rho', 'psi', 'u', 'v'. scoord: string, optional Name of vertical grid to interpolate output to. Options are 's_rho', 's_w', 'rho', 'w'. hboundary: string, optional Passed to `grid` method calls; horizontal boundary selection for calculating horizontal derivative of var. This same value will be used for grid changes too. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. hfill_value: float, optional Passed to `grid` method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. sboundary: string, optional Passed to `grid` method calls; vertical boundary selection for calculating horizontal derivative of var. This same value will be used for vertical grid changes too. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. sfill_value: float, optional Passed to `grid` method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. attrs: dict, optional Dictionary of attributes to add to resultant arrays. Requires that q is DataArray. For example: `attrs={'name': 'varname', 'long_name': 'longvarname', 'units': 'units'}` Returns ------- DataArray or ndarray of dqdeta, the gradient of q in the eta-direction with attributes altered to reflect calculation. Notes ----- dqdeta = dqdy*dzdz - dqdz*dzdy Derivatives are taken in the ROMS curvilinear grid native eta-direction. These derivatives properly account for the fact that ROMS vertical coordinates are s coordinates and therefore can vary in time and space. This will alter the number of points in the eta and s dimensions. Example usage ------------- >>> xroms.ddeta(ds.salt, grid) """ var = xroms.hgrad( var, grid, which="eta", z=z, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, attrs=attrs, ) return var
def ertel( phi, u, v, f, grid, hcoord="rho", scoord="s_rho", hboundary="extend", hfill_value=None, sboundary="extend", sfill_value=None, ): """Calculate Ertel potential vorticity of phi. Inputs ------ phi: DataArray Conservative tracer. Usually this would be the buoyancy but could be another approximately conservative tracer. The buoyancy can be calculated as: >>> xroms.buoyancy(temp, salt, 0) and then input as `phi`. u: DataArray xi component of velocity [m/s] v: DataArray eta component of velocity [m/s] f: DataArray Coriolis parameter [1/s] grid: xgcm.grid Grid object associated with u, v hcoord: string, optional. Name of horizontal grid to interpolate output to. Options are 'rho', 'psi', 'u', 'v'. scoord: string, optional. Name of vertical grid to interpolate output to. Options are 's_rho', 's_w', 'rho', 'w'. hboundary: string, optional Passed to `grid` method calls; horizontal boundary selection for calculating horizontal derivatives of phi and for calculating relative vorticity. This same value will be used for all horizontal grid changes too. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. hfill_value: float, optional Passed to `grid` method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. sboundary: string, optional Passed to `grid` method calls; vertical boundary selection for calculating horizontal and vertical derivatives of phi, and for calculating relative vorticity. This same value will be used for all vertical grid changes too. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. sfill_value: float, optional Passed to `grid` method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. Returns ------- DataArray of the Ertel potential vorticity for the input tracer. Notes ----- epv = -v_z * phi_x + u_z * phi_y + (f + v_x - u_y) * phi_z This is not set up to accept different boundary choices for different variables. Example usage: >>> xroms.ertel(ds.dye_01, ds.u, ds.v, ds.f, ds.attrs['grid'], scoord='s_w'); """ assert isinstance(phi, xr.DataArray), "phi must be DataArray" assert isinstance(u, xr.DataArray), "u must be DataArray" assert isinstance(v, xr.DataArray), "v must be DataArray" assert isinstance(f, xr.DataArray), "f must be DataArray" phi_xi, phi_eta = xroms.hgrad( phi, grid, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) phi_xi = xroms.to_grid( phi_xi, grid, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) phi_eta = xroms.to_grid( phi_eta, grid, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) phi_z = xroms.ddz( phi, grid, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) # vertical shear (horizontal components of vorticity) u_z = xroms.ddz( u, grid, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) v_z = xroms.ddz( v, grid, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) # vertical component of vorticity vort = relative_vorticity( u, v, grid, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) vort = xroms.to_grid( vort, grid, hcoord=hcoord, scoord=scoord, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) # combine terms to get the ertel potential vorticity epv = -v_z * phi_xi + u_z * phi_eta + (f + vort) * phi_z attrs = { "name": "ertel", "long_name": "ertel potential vorticity", "units": "tracer/(m*s)", "grid": grid, } epv = xroms.to_grid( epv, grid, hcoord=hcoord, scoord=scoord, attrs=attrs, hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) return epv
def uv_geostrophic(zeta, f, grid, hboundary="extend", hfill_value=None, which="both"): """Calculate geostrophic velocities from zeta [m/s] Inputs ------ zeta: DataArray sea surface height [m] f: DataArray or ndarray Coriolis parameter [1/s] grid: xgcm.grid Grid object associated with zeta hboundary: string, optional Passed to `grid` method calls; horizontal boundary selection for moving f to rho grid. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. hfill_value: float, optional Passed to `grid` method calls; horizontal boundary selection for moving f to rho grid. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. which: string, optional Which components of geostrophic velocity to return. * 'both': return both components of hgrad * 'xi': return only xi-direction. * 'eta': return only eta-direction. Returns ------- DataArrays of components of geostrophic velocity calculated on their respective grids. Notes ----- vg = g * zeta_eta / (d eta * f) # on v grid ug = -g * zeta_xi / (d xi * f) # on u grid Translation to Python of Matlab copy of surf_geostr_vel of IRD Roms_Tools. Example usage ------------- >>> xroms.uv_geostrophic(ds.zeta, ds.f, grid) """ assert isinstance(zeta, xr.DataArray), "zeta must be DataArray" assert isinstance(f, xr.DataArray), "f must be DataArray" if which in ["both", "xi"]: # calculate derivatives of zeta dzetadxi = xroms.hgrad(zeta, grid, which="xi") # calculate geostrophic velocities ug = ( -g * dzetadxi / xroms.to_u(f, grid, hboundary=hboundary, hfill_value=hfill_value)) ug.attrs["name"] = "u_geo" ug.attrs["long_name"] = "geostrophic u velocity" ug.attrs["units"] = "m/s" # inherits grid from T ug.name = ug.attrs["name"] if which in ["both", "eta"]: # calculate derivatives of zeta dzetadeta = xroms.hgrad(zeta, grid, which="eta") # calculate geostrophic velocities vg = ( g * dzetadeta / xroms.to_v(f, grid, hboundary=hboundary, hfill_value=hfill_value)) vg.attrs["name"] = "v_geo" vg.attrs["long_name"] = "geostrophic v velocity" vg.attrs["units"] = "m/s" # inherits grid from T vg.name = vg.attrs["name"] if which == "both": return ug, vg elif which == "xi": return ug elif which == "eta": return vg else: print("nothing being returned from uv_geostrophic")
def relative_vorticity( u, v, grid, hboundary="extend", hfill_value=None, sboundary="extend", sfill_value=None, ): """Calculate the vertical component of the relative vorticity [1/s] Inputs ------ u: DataArray xi component of velocity [m/s] v: DataArray eta component of velocity [m/s] grid: xgcm.grid Grid object associated with u, v hboundary: string, optional Passed to `grid` method calls; horizontal boundary selection for calculating horizontal derivatives of u and v. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. hfill_value: float, optional Passed to `grid` method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. sboundary: string, optional Passed to `grid` method calls; vertical boundary selection for calculating horizontal derivatives of u and v. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. sfill_value: float, optional Passed to `grid` method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. Returns ------- DataArray of vertical component of relative vorticity psi/w grids. Notes ----- relative_vorticity = v_x - u_y Example usage ------------- >>> xroms.relative_vorticity(u, v, grid) """ assert isinstance(u, xr.DataArray), "u must be DataArray" assert isinstance(v, xr.DataArray), "v must be DataArray" dvdxi = xroms.hgrad( v, grid, which="xi", hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) dudeta = xroms.hgrad( u, grid, which="eta", hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) var = dvdxi - dudeta var.attrs["name"] = "vort" var.attrs["long_name"] = "vertical component of vorticity" var.attrs["units"] = "1/s" # inherits grid from T var.name = var.attrs["name"] return var
def M2( rho, grid, rho0=1025.0, hboundary="extend", hfill_value=None, sboundary="fill", sfill_value=np.nan, z=None, ): """Calculate the horizontal buoyancy gradient. Inputs ------ rho: DataArray Density [kg/m^3] grid: xgcm.grid Grid object associated with rho rho0: int, float, optional Reference density [kg/m^3]. hboundary: string, optional Passed to `grid` method calls; horizontal boundary selection for calculating horizontal derivatives of rho. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. hfill_value: float, optional Passed to `grid` method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. sboundary: string, optional Passed to `grid` method calls; vertical boundary selection for calculating horizontal derivatives of rho. From xgcm documentation: A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation. * 'fill': Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.) * 'extend': Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition. sfill_value: float, optional Passed to `grid` method calls; vertical boundary fill value associated with sboundary input. From xgcm documentation: The value to use in the boundary condition with `boundary='fill'`. z: DataArray, optional Depths [m] associated with rho. If None, use z coordinate attached to temperature. Returns ------- DataArray of the horizontal buoyancy gradient on rho/w grids. Output is `[T,Z,Y,X]`. Notes ----- M2 = g/rho0 * sqrt(d(rho)/dxi^2 + d(rho)deta^2) g=9.81 [m/s^2] Example usage ------------- >>> xroms.M2(rho, grid) """ assert isinstance(rho, xr.DataArray), "rho must be DataArray" # calculate spatial derivatives of density drhodxi, drhodeta = xroms.hgrad( rho, grid, which="both", hcoord="rho", hboundary=hboundary, hfill_value=hfill_value, sboundary=sboundary, sfill_value=sfill_value, ) # combine var = np.sqrt(drhodxi**2 + drhodeta**2) * g / rho0 var.attrs["name"] = "M2" var.attrs["long_name"] = "horizontal buoyancy gradient" var.attrs["units"] = "1/s" var.attrs["grid"] = grid var.name = var.attrs["name"] return var