Skip to content

Commit bbf52c3

Browse files
Apply the initial changes
1 parent d2a3aa7 commit bbf52c3

3 files changed

Lines changed: 105 additions & 56 deletions

File tree

uxarray/core/zonal.py

Lines changed: 2 additions & 26 deletions
Original file line numberDiff line numberDiff line change
@@ -3,10 +3,6 @@
33

44

55
from uxarray.grid.integrate import _zonal_face_weights, _zonal_face_weights_robust
6-
from uxarray.grid.utils import (
7-
_get_cartesian_face_edge_nodes,
8-
_get_lonlat_rad_face_edge_nodes,
9-
)
106

117

128

@@ -22,28 +18,8 @@ def _compute_non_conservative_zonal_mean(uxda, latitudes, use_robust_weights=Fal
2218
# Create a NumPy array for storing results
2319
result = np.zeros(shape, dtype=uxda.dtype)
2420

25-
26-
faces_edges_cartesian = _get_cartesian_face_edge_nodes(
27-
uxgrid.face_node_connectivity.values,
28-
uxgrid.n_face,
29-
uxgrid.n_max_face_edges,
30-
uxgrid.node_x.values,
31-
uxgrid.node_y.values,
32-
uxgrid.node_z.values,
33-
)
34-
35-
36-
faces_edges_lonlat_rad = _get_lonlat_rad_face_edge_nodes(
37-
uxgrid.face_node_connectivity.values,
38-
uxgrid.n_face,
39-
uxgrid.n_max_face_edges,
40-
uxgrid.node_lon.values,
41-
uxgrid.node_lat.values,
42-
)
43-
44-
# Check if the bounds are available
45-
if "bounds" not in uxgrid._ds:
46-
uxgrid._populate_bounds()
21+
bounds = uxgrid.bounds.values
22+
faces_edges_cartesian = uxgrid.face_edges_cartesian.values
4723

4824
for i, lat in enumerate(latitudes):
4925
face_indices = uxda.uxgrid.get_faces_at_constant_latitude(lat)

uxarray/grid/geometry.py

Lines changed: 88 additions & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -1346,13 +1346,91 @@ def compute_temp_latlon_array(
13461346
return temp_latlon_array
13471347

13481348

1349+
def _populate_face_edges_cartesian(grid, return_array=False):
1350+
faces_edges_cartesian = _get_cartesian_face_edge_nodes(
1351+
grid.face_node_connectivity.values,
1352+
grid.n_face,
1353+
grid.n_max_face_edges,
1354+
grid.node_x.values,
1355+
grid.node_y.values,
1356+
grid.node_z.values,
1357+
)
1358+
1359+
faces_edges_cartesian_xarray = xr.DataArray(
1360+
faces_edges_cartesian,
1361+
dims=["n_face", "n_max_face_edges", "two", "three"],
1362+
attrs={
1363+
"cf_role": "face_edges_cartesian",
1364+
"_FillValue": INT_FILL_VALUE,
1365+
"long_name": "Provide the Cartesian coordinates of the edge for each face",
1366+
"start_index": INT_DTYPE(0)
1367+
},
1368+
)
1369+
1370+
if return_array:
1371+
return faces_edges_cartesian_xarray
1372+
else:
1373+
grid._ds["face_edges_cartesian"] = faces_edges_cartesian_xarray
1374+
1375+
1376+
def _populate_faces_edges_cartesian(grid, return_array=False):
1377+
faces_edges_cartesian = _get_cartesian_face_edge_nodes(
1378+
grid.face_node_connectivity.values,
1379+
grid.n_face,
1380+
grid.n_max_face_edges,
1381+
grid.node_x.values,
1382+
grid.node_y.values,
1383+
grid.node_z.values,
1384+
)
1385+
1386+
faces_edges_cartesian_xarray = xr.DataArray(
1387+
faces_edges_cartesian,
1388+
dims=["n_face", "n_max_face_edges", "two", "three"],
1389+
attrs={
1390+
"cf_role": "faces_edges_cartesian",
1391+
"_FillValue": INT_FILL_VALUE,
1392+
"long_name": "Provide the Cartesian coordinates of the edge for each face",
1393+
"start_index": INT_DTYPE(0)
1394+
},
1395+
)
1396+
1397+
if return_array:
1398+
return faces_edges_cartesian_xarray
1399+
else:
1400+
grid._ds["faces_edges_cartesian"] = faces_edges_cartesian_xarray
1401+
1402+
1403+
def _populate_faces_edges_spherical(grid, return_array=False):
1404+
faces_edges_lonlat_rad = _get_lonlat_rad_face_edge_nodes(
1405+
grid.face_node_connectivity.values,
1406+
grid.n_face,
1407+
grid.n_max_face_edges,
1408+
grid.node_lon.values,
1409+
grid.node_lat.values,
1410+
)
1411+
1412+
faces_edges_lonlat_rad_xarray = xr.DataArray(
1413+
faces_edges_lonlat_rad,
1414+
dims=["n_face", "n_max_face_edges", "two", "two"],
1415+
attrs={
1416+
"cf_role": "faces_edges_spherical",
1417+
"_FillValue": INT_FILL_VALUE,
1418+
"long_name": "Provide the Spherical coordinates (lon lat) in radians of the edge for each face",
1419+
"start_index": INT_DTYPE(0)
1420+
},
1421+
)
1422+
1423+
if return_array:
1424+
return faces_edges_lonlat_rad_xarray
1425+
else:
1426+
grid._ds["faces_edges_spherical"] = faces_edges_lonlat_rad_xarray
1427+
1428+
13491429
def _populate_bounds(
1350-
grid,
1351-
faces_edges_cartesian: np.ndarray = None,
1352-
faces_edges_lonlat_rad: np.ndarray = None,
1353-
is_latlonface: bool = False,
1354-
is_face_GCA_list=None,
1355-
return_array=False
1430+
grid,
1431+
is_latlonface: bool = False,
1432+
is_face_GCA_list=None,
1433+
return_array=False
13561434
):
13571435
"""Populates the bounds of the grid based on the geometry of its faces,
13581436
taking into account special conditions such as faces crossing the
@@ -1430,26 +1508,11 @@ def _populate_bounds(
14301508
# TODO: Is it possible to have a more flexible normalization? (Avoid duplicated normalizations)
14311509
grid.normalize_cartesian_coordinates()
14321510

1433-
# Prepare data for Numba functions
1434-
# TODO: remove the duplicate call for this function in both zonal average and bounds
1435-
if faces_edges_cartesian is None:
1436-
faces_edges_cartesian = _get_cartesian_face_edge_nodes(
1437-
grid.face_node_connectivity.values,
1438-
grid.n_face,
1439-
grid.n_max_face_edges,
1440-
grid.node_x.values,
1441-
grid.node_y.values,
1442-
grid.node_z.values,
1443-
)
1511+
# If the face_edges_cartesian and face_edges_lonlat_rad didn't exist,
1512+
# The following call will automatically populate them
1513+
faces_edges_cartesian = grid.faces_edges_cartesian.values
14441514

1445-
if faces_edges_lonlat_rad is None:
1446-
faces_edges_lonlat_rad = _get_lonlat_rad_face_edge_nodes(
1447-
grid.face_node_connectivity.values,
1448-
grid.n_face,
1449-
grid.n_max_face_edges,
1450-
grid.node_lon.values,
1451-
grid.node_lat.values,
1452-
)
1515+
faces_edges_lonlat_rad = grid.faces_edges_spherical.values
14531516

14541517
n_nodes_per_face = grid.n_nodes_per_face.values
14551518

uxarray/grid/grid.py

Lines changed: 15 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -71,6 +71,8 @@
7171
_grid_to_matplotlib_polycollection,
7272
_grid_to_matplotlib_linecollection,
7373
_populate_bounds,
74+
_populate_faces_edges_cartesian,
75+
_populate_faces_edges_spherical,
7476
_construct_boundary_edge_indices,
7577
compute_temp_latlon_array,
7678
_find_faces,
@@ -1444,20 +1446,24 @@ def face_areas(self, value):
14441446
self._ds["face_areas"] = value
14451447

14461448
@property
1447-
def face_edges_cartesian(self):
1449+
def faces_edges_cartesian(self):
14481450
"""Cartesian Coordinates for each Face Edge.
14491451
14501452
Dimensions ``(n_face, n_max_face_edges, 2, 3)``
14511453
"""
1452-
pass
1454+
if "faces_edges_cartesian" not in self._ds:
1455+
_populate_faces_edges_cartesian(self)
1456+
return self._ds["faces_edges_cartesian"]
14531457

14541458
@property
1455-
def face_edges_spherical(self):
1456-
"""Latitude Longitude Bounds for each Face in radians.
1459+
def faces_edges_spherical(self):
1460+
"""Latitude Longitude Coordinates for each Face in radians.
14571461
14581462
Dimensions ``(n_face, n_max_face_edges, 2, 2)``
14591463
"""
1460-
pass
1464+
if "faces_edges_spherical" not in self._ds:
1465+
_populate_faces_edges_spherical(self)
1466+
return self._ds["faces_edges_spherical"]
14611467

14621468
@property
14631469
def bounds(self):
@@ -1472,6 +1478,10 @@ def bounds(self):
14721478
"This initial execution will be significantly longer.",
14731479
RuntimeWarning,
14741480
)
1481+
if "faces_edges_cartesian" not in self._ds:
1482+
_populate_faces_edges_cartesian(self)
1483+
if "faces_edges_spherical" not in self._ds:
1484+
_populate_faces_edges_spherical(self)
14751485
_populate_bounds(self)
14761486
return self._ds["bounds"]
14771487

0 commit comments

Comments
 (0)