|
| 1 | +import xarray as xr |
| 2 | +import numpy as np |
| 3 | +import scipy.sparse |
| 4 | +import scipy.sparse.linalg |
| 5 | + |
| 6 | +import logging |
| 7 | +import sys |
| 8 | +from mpas_tools.ocean.depth import compute_zmid |
| 9 | + |
| 10 | + |
| 11 | +def compute_barotropic_streamfunction(ds_mesh, ds, logger=None, |
| 12 | + min_depth=-5., max_depth=1.e4, |
| 13 | + prefix='timeMonthly_avg_', |
| 14 | + time_index=0): |
| 15 | + """ |
| 16 | + Compute barotropic streamfunction. Returns BSF in Sv on vertices. |
| 17 | +
|
| 18 | + Parameters |
| 19 | + ---------- |
| 20 | + ds_mesh : ``xarray.Dataset`` |
| 21 | + A dataset containing MPAS mesh variables |
| 22 | +
|
| 23 | + ds : ``xarray.Dataset`` |
| 24 | + A dataset containing MPAS output variables ``normalVelocity`` and |
| 25 | + ``layerThickness`` (possibly with a ``prefix``) |
| 26 | +
|
| 27 | + logger : ``logging.Logger``, optional |
| 28 | + A logger for the output if not stdout |
| 29 | +
|
| 30 | + min_depth : float, optional |
| 31 | + The minimum depth (positive down) to compute transport over |
| 32 | +
|
| 33 | + max_depth : float, optional |
| 34 | + The maximum depth (positive down) to compute transport over |
| 35 | +
|
| 36 | + prefix : str, optional |
| 37 | + The prefix on the ``normalVelocity`` and ``layerThickness`` variables |
| 38 | +
|
| 39 | + time_index : int, optional |
| 40 | + The time at which to index ``ds`` (if it has ``Time`` as a dimension) |
| 41 | + """ |
| 42 | + |
| 43 | + useStdout = logger is None |
| 44 | + if useStdout: |
| 45 | + logger = logging.getLogger() |
| 46 | + logger.addHandler(logging.StreamHandler(sys.stdout)) |
| 47 | + logger.setLevel(logging.INFO) |
| 48 | + |
| 49 | + inner_edges, transport = _compute_transport( |
| 50 | + ds_mesh, ds, min_depth=min_depth, max_depth=max_depth, prefix=prefix, |
| 51 | + time_index=time_index) |
| 52 | + logger.info('transport computed.') |
| 53 | + |
| 54 | + nvertices = ds_mesh.sizes['nVertices'] |
| 55 | + |
| 56 | + cells_on_vertex = ds_mesh.cellsOnVertex - 1 |
| 57 | + vertices_on_edge = ds_mesh.verticesOnEdge - 1 |
| 58 | + is_boundary_cov = cells_on_vertex == -1 |
| 59 | + boundary_vertices = np.logical_or(is_boundary_cov.isel(vertexDegree=0), |
| 60 | + is_boundary_cov.isel(vertexDegree=1)) |
| 61 | + boundary_vertices = np.logical_or(boundary_vertices, |
| 62 | + is_boundary_cov.isel(vertexDegree=2)) |
| 63 | + |
| 64 | + # convert from boolean mask to indices |
| 65 | + boundary_vertices = np.flatnonzero(boundary_vertices.values) |
| 66 | + |
| 67 | + n_boundary_vertices = len(boundary_vertices) |
| 68 | + n_inner_edges = len(inner_edges) |
| 69 | + |
| 70 | + indices = np.zeros((2, 2 * n_inner_edges + n_boundary_vertices), dtype=int) |
| 71 | + data = np.zeros(2 * n_inner_edges + n_boundary_vertices, dtype=float) |
| 72 | + |
| 73 | + # The difference between the streamfunction at vertices on an inner |
| 74 | + # edge should be equal to the transport |
| 75 | + v0 = vertices_on_edge.isel(nEdges=inner_edges, TWO=0).values |
| 76 | + v1 = vertices_on_edge.isel(nEdges=inner_edges, TWO=1).values |
| 77 | + |
| 78 | + ind = np.arange(n_inner_edges) |
| 79 | + indices[0, 2 * ind] = ind |
| 80 | + indices[1, 2 * ind] = v1 |
| 81 | + data[2 * ind] = 1. |
| 82 | + |
| 83 | + indices[0, 2 * ind + 1] = ind |
| 84 | + indices[1, 2 * ind + 1] = v0 |
| 85 | + data[2 * ind + 1] = -1. |
| 86 | + |
| 87 | + # the streamfunction should be zero at all boundary vertices |
| 88 | + ind = np.arange(n_boundary_vertices) |
| 89 | + indices[0, 2 * n_inner_edges + ind] = n_inner_edges + ind |
| 90 | + indices[1, 2 * n_inner_edges + ind] = boundary_vertices |
| 91 | + data[2 * n_inner_edges + ind] = 1. |
| 92 | + |
| 93 | + rhs = np.zeros(n_inner_edges + n_boundary_vertices, dtype=float) |
| 94 | + |
| 95 | + # convert to Sv |
| 96 | + ind = np.arange(n_inner_edges) |
| 97 | + rhs[ind] = 1e-6 * transport |
| 98 | + |
| 99 | + ind = np.arange(n_boundary_vertices) |
| 100 | + rhs[n_inner_edges + ind] = 0. |
| 101 | + |
| 102 | + matrix = scipy.sparse.csr_matrix( |
| 103 | + (data, indices), |
| 104 | + shape=(n_inner_edges + n_boundary_vertices, nvertices)) |
| 105 | + |
| 106 | + solution = scipy.sparse.linalg.lsqr(matrix, rhs) |
| 107 | + bsf_vertex = xr.DataArray(-solution[0], |
| 108 | + dims=('nVertices',)) |
| 109 | + |
| 110 | + return bsf_vertex |
| 111 | + |
| 112 | +def _compute_transport(ds_mesh, ds, min_depth, max_depth, prefix, |
| 113 | + time_index): |
| 114 | + |
| 115 | + cells_on_edge = ds_mesh.cellsOnEdge - 1 |
| 116 | + inner_edges = np.logical_and(cells_on_edge.isel(TWO=0) >= 0, |
| 117 | + cells_on_edge.isel(TWO=1) >= 0) |
| 118 | + |
| 119 | + if 'Time' in ds.dims: |
| 120 | + ds = ds.isel(Time=time_index) |
| 121 | + |
| 122 | + # convert from boolean mask to indices |
| 123 | + inner_edges = np.flatnonzero(inner_edges.values) |
| 124 | + |
| 125 | + cell0 = cells_on_edge.isel(nEdges=inner_edges, TWO=0) |
| 126 | + cell1 = cells_on_edge.isel(nEdges=inner_edges, TWO=1) |
| 127 | + |
| 128 | + normal_velocity = \ |
| 129 | + ds[f'{prefix}normalVelocity'].isel(nEdges=inner_edges) |
| 130 | + layer_thickness = ds[f'{prefix}layerThickness'] |
| 131 | + layer_thickness_edge = 0.5 * (layer_thickness.isel(nCells=cell0) + |
| 132 | + layer_thickness.isel(nCells=cell1)) |
| 133 | + |
| 134 | + n_vert_levels = ds.sizes['nVertLevels'] |
| 135 | + |
| 136 | + vert_index = xr.DataArray.from_dict( |
| 137 | + {'dims': ('nVertLevels',), 'data': np.arange(n_vert_levels)}) |
| 138 | + mask_bottom = (vert_index < ds_mesh.maxLevelCell).T |
| 139 | + mask_bottom_edge = 0.5 * (mask_bottom.isel(nCells=cell0) + |
| 140 | + mask_bottom.isel(nCells=cell1)) |
| 141 | + |
| 142 | + if 'zMid' not in ds.keys(): |
| 143 | + z_mid = compute_zmid(ds_mesh.bottomDepth, ds_mesh.maxLevelCell, |
| 144 | + ds_mesh.layerThickness) |
| 145 | + else: |
| 146 | + z_mid = ds.zMid |
| 147 | + z_mid_edge = 0.5 * (z_mid.isel(nCells=cell0) + |
| 148 | + z_mid.isel(nCells=cell1)) |
| 149 | + |
| 150 | + mask = np.logical_and(np.logical_and(z_mid_edge >= -max_depth, |
| 151 | + z_mid_edge <= -min_depth), |
| 152 | + mask_bottom_edge) |
| 153 | + normal_velocity = normal_velocity.where(mask) |
| 154 | + layer_thickness_edge = layer_thickness_edge.where(mask) |
| 155 | + transport = ds_mesh.dvEdge[inner_edges] * \ |
| 156 | + (layer_thickness_edge * normal_velocity).sum(dim='nVertLevels') |
| 157 | + |
| 158 | + return inner_edges, transport |
0 commit comments