예제 #1
0
    def __init__(
        self,
        z=None,  # grid (input)
        kappa=None,  # diffusivity profile (input)
        bs=0.025,  # surface buoyancy bound. cond (input)
        bbot=0.0,  # bottom buoyancy boundary condition (input)  
        bzbot=None,  # bottom strat. as alternative boundary condition (input) 
        b=0.0,  # Buoyancy profile (input, output)
        Area=None,  # Horizontal area (can be function of depth)
        N2min=1e-7  # Minimum strat. for conv adjustment
    ):
        r"""
    Parameters
    ----------

    z : ndarray
        Vertical depth levels of column grid. Units: m
    kappa : float, function, or ndarray
            Vertical diffusivity profile. Units: m\ :sup:`2`/s
    bs : float
         Surface level buoyancy boundary condition. Units: m/s\ :sup:`2`
    bbot : float; optional
           Bottom level buoyancy boundary condition. Units: m/s\ :sup:`2`
    bzbot : float; optional
            Bottom level buoyancy stratification. Can be used as an alternative to **bbot**. Units: s\ :sup:`-2`
    b : float, function, or ndarray
        Initial vertical buoyancy profile. Recalculated on model run. Units: m/s
    Area : float, function, or ndarray
           Horizontal area of basin. Units: m\ :sup:`2`
    N2min : float; optional
            Minimum stratification for convective adjustment. Units: s\ :sup:`-1`
    """

        # initialize grid:
        if isinstance(z, np.ndarray) and len(z) > 0:
            self.z = z
        else:
            raise TypeError('z needs to be numpy array providing grid levels')

        self.kappa = make_func(kappa, self.z, 'kappa')
        self.Area = make_func(Area, self.z, 'Area')

        self.bs = bs
        self.bbot = bbot
        self.bzbot = bzbot

        self.N2min = N2min

        self.b = make_array(b, self.z, 'b')

        if check_numpy_version():
            self.bz = np.gradient(self.b, z)
        else:
            self.bz = 0. * z  # notice that this is just for initialization of ode solver
예제 #2
0
  def test_update(self, psi):
    b1 = 10.0
    b2 = 50.0
    z = psi.z
    psi.update(b1, b2)

    b1_func = make_func(b1, z, 'b1')
    b2_func = make_func(b2, z, 'b2')

    assert (all(psi.b1(z) == b1_func(z)))
    assert (all(psi.b2(z) == b2_func(z)))
예제 #3
0
파일: test_psi_SO.py 프로젝트: pymoc/pymoc
  def test_update(self, psi_so):
    b = 10.0
    bs = 50.0
    z = psi_so.z
    y = psi_so.y
    psi_so.update(b, bs)

    b_func = make_func(b, z, 'b')
    bs_func = make_func(bs, y, 'bs')

    assert (all(psi_so.b(z) == b_func(z)))
    assert (all(psi_so.bs(y) == bs_func(y)))
예제 #4
0
  def test_psi_thermwind_init(self, psi_config):
    if not 'z' in psi_config or not isinstance(psi_config['z'], np.ndarray
                                               ) or not len(psi_config['z']):
      with pytest.raises(TypeError) as zinfo:
        Psi_Thermwind(**psi_config)
      assert (
          str(zinfo.value) == "z needs to be numpy array providing grid levels"
      )
      return

    psi = Psi_Thermwind(**psi_config)
    for k in ['z', 'f', 'sol_init', 'b1', 'b2']:
      assert hasattr(psi, k)

    psi_signature = funcsigs.signature(Psi_Thermwind)

    for k in ['f']:
      assert getattr(psi, k) == (
          psi_config[k] if k in psi_config and psi_config[k] else
          psi_signature.parameters[k].default
      )

    if not 'sol_init' in psi_config:
      testing.assert_array_equal(psi.sol_init, np.zeros((2, len(psi.z))))

    for k in ['b1', 'b2']:
      f = getattr(psi, k)
      ft = make_func(psi_config[k], psi_config['z'], k)
      for z in psi.z:
        assert f(z) == ft(z)
예제 #5
0
    def update(self, b1=None, b2=None):
        r"""
    Update the vertical buoyancy profiles from the southern basin and northern region.

    Parameters
    ----------

    b1 : float, function, or ndarray; optional
         Vertical buoyancy profile from the southern basin. Units: m/s\ :sup:`2`
    b2 : float, function, or ndarray; optional
         Vertical buoyancy profile from the northern basin, representing the
         deepwater formation region. Units: m/s\ :sup:`2`
    """
        # update buoyancy profiles
        if b1 is not None:
            self.b1 = make_func(b1, self.z, 'b1')
        if b2 is not None:
            self.b2 = make_func(b2, self.z, 'b2')
예제 #6
0
    def update(self, b=None, bs=None):
        r"""
    Update the vertical buoyancy profile and surface buoyancy, based on changes
    in the adjoining basin and/or in the surface boundary conditions.

    Parameters
    ----------

    b : float, function, or ndarray
        Vertical buoyancy profile from the adjoining basin, on the north
        side of the ACC. Units: m/s\ :sup:`2`
    bs : float, function, or ndarray
         Surface level buoyancy boundary condition. Can be a constant,
         or an array or function in y. Units: m/s\ :sup:`2`
    """

        if b is not None:
            self.b = make_func(b, self.z, 'b')
        if bs is not None:
            self.bs = make_func(bs, self.y, 'bs')
예제 #7
0
    def __init__(
            self,
            f=1.2e-4,  # Coriolis parameter (input)
            z=None,  # grid (input)
            sol_init=None,  # Initial conditions for ODE solver (input)
            b1=None,  # Buoyancy in the basin (input, output)
            b2=0.,  # Buoyancy in the deep water formation region (input, output)
    ):
        r"""
    Parameters
    ----------

    f : float
        Coriolis parameter. Units s\ :sup:`-1`
    z : ndarray
        Vertical depth levels of overturning grid. Units: m
    sol_init : ndarray; optional
               Initial guess at the solution to the thermal wind overturning streamfunction. Units: [...]
    b1 : float, function, or ndarray; optional
         Vertical buoyancy profile from the southern basin. Units: m/s\ :sup:`2`
    b2 : float, function, or ndarray; optional
         Vertical buoyancy profile from the northern basin, representing the
         deepwater formation region. Units: m/s\ :sup:`2`
    """

        self.f = f
        # initialize grid:
        if isinstance(z, np.ndarray):
            self.z = z
            nz = np.size(z)
        else:
            raise TypeError('z needs to be numpy array providing grid levels')

        self.b1 = make_func(b1, self.z, 'b1')
        self.b2 = make_func(b2, self.z, 'b2')

        # Set initial conditions for BVP solver
        if sol_init is None:
            self.sol_init = np.zeros((2, nz))
        else:
            self.sol_init = sol_init
예제 #8
0
    def calc_GM(self):
        r"""
    Compute the eddy (Gent & Mcwilliams) transport based on the meridionally
    averaged isopycnal slope.

    .. math::
      \begin{aligned}
      \Psi_{GM} &= K_{GM}\cdot s \\
      s\left(b\right) &\equiv \frac{z_B\left(b\right)}{L_y - y_{SO}\left(b\right)}
      \end{aligned}

    Where :math:`z_B\left(b\right)` is the depth of isopycnals of density class :math:`b`
    in the adjoining basin, and :math:`y_{SO}\left(b\right)` is the outcropping latitude
    of isopycnals of density class :math:`b` in the Southern Ocean, available via
    :meth:`pymoc.modules.Psi_SO.ys`.

    Returns
    -------

    Psi_GM : ndarray
             An array representing the meridional average of the eddy transport at
             each vertical level of the Southern Ocean model.

    """

        dy_atz = 0 * self.z
        eps = 0.1  # minimum dy (in meters) (to avoid div. by 0)
        for ii in range(0, np.size(self.z)):
            dy_atz[ii] = max(self.y[-1] - self.ys(self.b(self.z[ii])), eps)
        bottaper = self.calc_bottom_taper(self.Htaperbot, self.z)
        toptaper = self.calc_top_taper(self.Htapertop, self.z)
        if self.c is not None:
            temp = make_func(
                self.KGM * self.z / dy_atz * self.L * toptaper * bottaper,
                self.z, 'psiGM')
            N2 = self.calc_N2()

            def ode(z, y):
                return np.vstack((y[1], N2(z) / self.c**2. * (y[0] - temp(z))))

            #Solve the boundary value problem
            res = integrate.solve_bvp(ode, self.bc_GM, self.z,
                                      np.zeros((2, np.size(self.z))))
            # return solution interpolated onto original grid
            temp = res.sol(self.z)[0, :]
        else:
            temp = self.KGM * np.maximum(
                self.z / dy_atz, -self.smax) * self.L * toptaper * bottaper
        # limit Psi_GM to -Psi_Ek on isopycnals that don't outcrop:
        idx = dy_atz > self.y[-1] - self.y[0]
        temp[idx] = np.maximum(temp[idx], -self.Psi_Ek[idx] * 1e6)
        return temp
예제 #9
0
    def calc_N2(self):
        r"""
    Calculate the buouyancy (Brunt-Väisällä) frequency profile for the Southern Ocean

    Returns
    -------

    N2 : function
         A depth dependent function that returns the buoyancy frequency :math:`N^2` for a given depth :math:`z`.

    """

        dz = self.z[1:] - self.z[:-1]
        N2 = np.zeros(np.size(self.z))
        b = self.b(self.z)

        N2[1:-1] = (b[2:] - b[:-2]) / (dz[1:] + dz[:-1])
        N2[0] = (b[1] - b[0]) / dz[0]
        N2[-1] = (b[-1] - b[-2]) / dz[-1]

        return make_func(N2, self.z, 'N2')
예제 #10
0
    def solve_equi(self, wA):
        r"""
    Solve for the equilibrium buoyancy profile, given a specified vertical
    velocity profile, and pre-set surface and bottom boundary conditions, based
    on the system of equations defined by :meth:`pymoc.modules.Column.ode`.

    Parameters
    ----------

    wA : ndarray
         Area integrated velocity profile for the equilibrium solution. Units: m\ :sup:`3`/s

    """

        self.wA = make_func(wA, self.z, 'w')
        sol_init = np.zeros((2, np.size(self.z)))
        sol_init[0, :] = self.b
        sol_init[1, :] = self.bz
        res = integrate.solve_bvp(self.ode, self.bc, self.z, sol_init)
        # interpolate solution for b and db/dz onto original grid
        self.b = res.sol(self.z)[0, :]
        self.bz = res.sol(self.z)[1, :]
예제 #11
0
파일: test_psi_SO.py 프로젝트: pymoc/pymoc
  def test_psi_so_init(self, psi_so_config):
    if not 'z' in psi_so_config or not isinstance(
        psi_so_config['z'], np.ndarray
    ) or not len(psi_so_config['z']):
      with pytest.raises(TypeError) as zinfo:
        Psi_SO(**psi_so_config)
      assert (
          str(zinfo.value) == "z needs to be numpy array providing grid levels"
      )
      return
    if not 'y' in psi_so_config or not isinstance(
        psi_so_config['y'], np.ndarray
    ) or not len(psi_so_config['y']):
      with pytest.raises(TypeError) as yinfo:
        Psi_SO(**psi_so_config)
      assert (
          str(yinfo.value) ==
          "y needs to be numpy array providing horizontal grid (or boundaries) of ACC"
      )
      return

    # Implicit exceptions thrown if these params are missing, should change these to explicit checks
    for k in ['b', 'bs', 'tau']:
      if not k in psi_so_config:
        with pytest.raises(TypeError) as kinfo:
          Psi_SO(**psi_so_config)
        assert (
            str(kinfo.value) == "('" + k +
            "', 'needs to be either function, numpy array, or float')"
        )
        return

    psi_so = Psi_SO(**psi_so_config)
    for k in [
        'z', 'y', 'b', 'bs', 'tau', 'f', 'rho', 'L', 'KGM', 'c', 'bvp_with_Ek',
        'Hsill', 'HEk', 'Htapertop', 'Htaperbot', 'smax'
    ]:
      assert hasattr(psi_so, k)

    psi_so_signature = funcsigs.signature(Psi_SO)

    for k in [
        'f', 'rho', 'L', 'KGM', 'c', 'bvp_with_Ek', 'Hsill', 'HEk',
        'Htapertop', 'Htaperbot', 'smax'
    ]:
      assert getattr(psi_so, k) == (
          psi_so_config[k] if k in psi_so_config and psi_so_config[k] else
          psi_so_signature.parameters[k].default
      )

    for k in ['b']:
      f = getattr(psi_so, k)
      ft = make_func(psi_so_config[k], psi_so_config['z'], k)
      for z in psi_so.z:
        assert f(z) == ft(z)

    for k in ['bs', 'tau']:
      f = getattr(psi_so, k)
      ft = make_func(psi_so_config[k], psi_so_config['y'], k)
      for y in psi_so.y:
        assert f(y) == ft(y)
예제 #12
0
파일: test_psi_SO.py 프로젝트: pymoc/pymoc
  def test_calc_GM_no_bvp(self, psi_so):
    # Test constant isopycnal sloping
    dy_atz = 2001000.0
    GM = [psi_so.L * psi_so.KGM * z / dy_atz for z in psi_so.z]
    psi_so.Psi_Ek = psi_so.calc_Ekman()
    assert (
        all(
            np.around(GM, decimals=3) ==
            np.around(psi_so.calc_GM(), decimals=3)
        )
    )

    # Test minimum dy limit
    dy_atz = 0.1
    psi_so.b = make_func(np.linspace(0.3, 0.2, 81), psi_so.z, 'b')
    GM = [-psi_so.L * psi_so.KGM * psi_so.smax for z in psi_so.z]
    GM[-1] = 0
    assert (
        all(
            np.around(GM, decimals=3) ==
            np.around(psi_so.calc_GM(), decimals=3)
        )
    )

    # Test maximum isopycnal slope limit
    dy_atz = 2001000.0
    psi_so.smax = -0.01
    psi_so.b = make_func(np.linspace(0.03, -0.001, 81), psi_so.z, 'b')
    GM = [psi_so.L * psi_so.KGM * 0.01 for _ in psi_so.z]
    psi_so.Psi_Ek = psi_so.calc_Ekman()
    assert (
        all(
            np.around(GM, decimals=3) ==
            np.around(psi_so.calc_GM(), decimals=3)
        )
    )

    # Test quadratic bottom tapering of streamfunction
    psi_so.smax = 0.01
    psi_so.Htaperbot = 800.0
    dy_atz = 2001000.0
    GM = [psi_so.L * psi_so.KGM * z / dy_atz for z in psi_so.z]
    GM[0:17] = [(1.0 - (-3200.0 - psi_so.z[i])**2 / 6.4e5) * GM[i]
                for i in range(0, 17)]
    psi_so.Psi_Ek = psi_so.calc_Ekman()
    assert (
        all(
            np.around(GM, decimals=3) ==
            np.around(psi_so.calc_GM(), decimals=3)
        )
    )
    psi_so.Htaperbot = None

    # Test quadratic top tapering of streamfunction
    psi_so.smax = 0.01
    psi_so.Htapertop = 500.0
    dy_atz = 2001000.0
    GM = [psi_so.L * psi_so.KGM * z / dy_atz for z in psi_so.z]
    GM[-11:] = [(1.0 - (500.0 + psi_so.z[i])**2 / 2.5e5) * GM[i]
                for i in range(70, 81)]
    psi_so.Psi_Ek = psi_so.calc_Ekman()
    assert (
        all(
            np.around(GM, decimals=3) ==
            np.around(psi_so.calc_GM(), decimals=3)
        )
    )
    psi_so.Htapertop = None

    # Test ekman streamfunction limiting for non-outcropping isopycnals
    psi_so = Psi_SO(
        **{
            'z': np.asarray(np.linspace(-4000, 0, 81)),
            'y': np.asarray(np.linspace(0, 2.0e3, 51)),
            'b': np.linspace(0.03, -0.001, 81),
            'bs': np.linspace(0.05, 0.10, 51),
            'tau': 0.12
        }
    )
    dy_atz = 2001000.0
    GM = [-psi_so.L * psi_so.KGM * 0.01 for _ in psi_so.z]
    GM[-1] = 0
    psi_ek = np.asarray([gm + np.abs(gm / 2.0) for gm in GM])
    psi_so.Psi_Ek = psi_ek
    assert (
        all(
            np.around(-psi_ek * 1e6, decimals=3) ==
            np.around(psi_so.calc_GM(), decimals=3)
        )
    )
예제 #13
0
    def __init__(
            self,
            z=None,  # vertical grid (array, in)
            y=None,  # horizontal grid (array, in)  
            b=None,  # buoyancy profile at northern end of ACC (function, array or float, in)
            bs=None,  # surface buoyancy (function, array or float, in)
            tau=None,  # surface wind stress (function, array or float, in)  
            f=1.2e-4,  # Coriolis parameter (in)
            rho=1030,  # Density of sea water (in)
            L=1e7,  # Zonal length of the ACC (in)
            KGM=1e3,  # GM coefficient (in)
            c=None,  # phase speed for F2010 BVP smoother of GM streamfunction
            bvp_with_Ek=False,  # if true, apply boundary condition that Psi_GM=-Psi_EK in F2010 BVP smoother 
            Hsill=None,  # height (in m above ocean floor) of the "sill", where Psi_Ek is tapered
            HEk=None,  # depth of surface Ekman layer
            Htapertop=None,  # A quadratic tapering of the GM streamfunction at the surface
            Htaperbot=None,  # A quadratic tapering of the GM streamfunction at the bottom
            smax=0.01,  # maximum slope for clipping of GM streamfunction
    ):
        r"""
    Parameters
    ----------

    z : ndarray
        Vertical depth levels of overturning grid. Units: m
    y : ndarray
        Meridional overturning grid. Units: m
    b : float, function, or ndarray
        Vertical buoyancy profile from the adjoining basin, on the north
        side of the ACC. Units: m/s\ :sup:`2`
    bs : float, function, or ndarray
         Surface level buoyancy boundary condition. Can be a constant,
         or an array or function in y. Units: m/s\ :sup:`2`
    tau : float, function, or ndarray
          Surface wind stress. Can be a constant, or an array or function
          in y. Units: N/m\ :sup:`2`
    f : float
        Coriolis parameter. Units s\ :sup:`-1`
    rho : float
          Density of sea water for Boussinesq approximation. Units: kg/m\ :sup:`3`
    L : float
        Zonal length of the modeled ACC. Units: m
    KGM : float
          Gent & McWilliams (GM) eddy diffusivity coefficient. Units: 
    c : float
        Phase speed cutoff for smoothing when solving the GM boundary value problem. Units: m/s 
    bvp_with_Ek : logical
            Whether to enforce the boundary condition that Psi_GM=-Psi_Ek at the ocean
            surface and bottom when solving the boundary value problem for the GM streamfunction.
    Hsill : float
            Height above the bottom at which the Ekman streamfunction is tapered. Units: m
    Hek : float
          Depth of the surface Ekman layer. Units: m
    Htapertop : float
                Height of the quadratic surface tapering layer for the GM streamfunction. Units: m
    Htaperbot : float
                Height of the quadratic bottom tapering layer for the GM streamfunction. Units: m
    smax : float
           Maximum slope of the GM streamfunction, above which Psi_GM is clipped. Units: m\ :sup:`-1`
    """

        # initialize grid:
        if isinstance(z, np.ndarray):
            self.z = z
        else:
            raise TypeError('z needs to be numpy array providing grid levels')

        if isinstance(y, np.ndarray):
            self.y = y
        else:
            raise TypeError(
                'y needs to be numpy array providing horizontal grid (or boundaries) of ACC'
            )

        self.b = make_func(b, self.z, 'b')
        self.bs = make_func(bs, self.y, 'bs')
        self.tau = make_func(tau, self.y, 'tau')
        self.f = f
        self.rho = rho
        self.L = L
        self.KGM = KGM
        self.c = c
        self.bvp_with_Ek = bvp_with_Ek
        self.Hsill = Hsill
        self.HEk = HEk
        self.Htapertop = Htapertop
        self.Htaperbot = Htaperbot
        self.smax = smax
예제 #14
0
 def make_func(self, myst, name, xin): # Seems unecessary to define a method that already exists identically as a function, no?
   return make_func(myst, xin, name)
예제 #15
0
 def make_func(self, myst, name, xin):
     return make_func(myst, xin, name)