API Reference

This section provides reference documentation for the Python API of the ADCIRC Subgrid Preprocessor. While the package is primarily designed for command-line use, the underlying Python classes can be used for custom workflows and advanced applications.

Note

This API is designed for internal use and may change in future versions. For stable interfaces, use the command-line tools described in other sections.

Core Classes

SubgridPreprocessor

The main processing class that orchestrates subgrid calculation.

class SubgridPreprocessor(input_file: InputFile, window_memory: int = 64)

Main class for processing subgrid data.

Parameters:
  • input_file (InputFile) – Configuration object containing processing parameters

  • window_memory (int, optional) – Memory limit for raster processing windows (MB)

Methods:

process() None

Execute the main subgrid processing workflow.

This method: - Loads input datasets (mesh, DEM, land cover) - Processes each mesh node to compute φ relationships - Calculates Manning’s roughness corrections - Stores results in internal data structures

write() None

Write processed subgrid data to NetCDF output file.

The output file contains: - φ factors as functions of water level - Average Manning’s coefficients - Various correction coefficients - Mesh coordinate information

Example Usage:

from AdcircSubgrid import SubgridPreprocessor, InputFile

# Create configuration
config = InputFile('config.yaml')

# Initialize processor
processor = SubgridPreprocessor(config, window_memory=128)

# Process and write results
processor.process()
processor.write()

InputFile

Configuration file parser and validator.

class InputFile(input_file: str)

Parser for YAML configuration files with schema validation.

Parameters:

input_file (str) – Path to YAML configuration file

Methods:

data() dict

Get validated configuration data.

Returns:

Dictionary containing configuration parameters

Return type:

dict

Example Usage:

from AdcircSubgrid import InputFile

# Load and validate configuration
config = InputFile('my_config.yaml')

# Access configuration data
config_data = config.data()
mesh_file = config_data['input']['adcirc_mesh']
output_file = config_data['output']['filename']

Mesh

ADCIRC mesh file reader and processor.

class Mesh(mesh_file: str)

Class for reading and processing ADCIRC mesh files (fort.14).

Parameters:

mesh_file (str) – Path to ADCIRC mesh file

Attributes:

nodes

Number of nodes in the mesh.

Type:

int

elements

Number of elements in the mesh.

Type:

int

x

Node x-coordinates.

Type:

numpy.ndarray

y

Node y-coordinates.

Type:

numpy.ndarray

depth

Node bathymetric depths.

Type:

numpy.ndarray

Example Usage:

from AdcircSubgrid import Mesh

# Load mesh
mesh = Mesh('fort.14')

# Access mesh properties
print(f"Mesh has {mesh.nodes} nodes and {mesh.elements} elements")

# Access coordinates and bathymetry
coordinates = list(zip(mesh.x, mesh.y))
depths = mesh.depth

Raster

Raster data reader and processor for DEM and land cover data.

class Raster(filename: str)

Generic raster data handler using GDAL.

Parameters:

filename (str) – Path to raster file

Methods:

get_values_at_points(x: np.ndarray, y: np.ndarray) np.ndarray

Extract raster values at specified coordinate points.

Parameters:
  • x (numpy.ndarray) – X-coordinates of query points

  • y (numpy.ndarray) – Y-coordinates of query points

Returns:

Raster values at query points

Return type:

numpy.ndarray

get_bbox() tuple

Get raster bounding box.

Returns:

Bounding box as (minx, miny, maxx, maxy)

Return type:

tuple

Example Usage:

from AdcircSubgrid import Raster
import numpy as np

# Load DEM
dem = Raster('elevation.tif')

# Sample elevations at specific points
x_points = np.array([-95.0, -94.5])
y_points = np.array([29.0, 29.5])
elevations = dem.get_values_at_points(x_points, y_points)

LookupTable

Manning’s roughness coefficient lookup table handler.

class LookupTable

Class for managing land cover to Manning’s roughness mappings.

Static Methods:

static ccap_lookup() np.ndarray

Generate default CCAP lookup table.

Returns:

Array mapping CCAP class codes to Manning’s coefficients

Return type:

numpy.ndarray

static from_file(filename: str) np.ndarray

Load custom lookup table from CSV file.

Parameters:

filename (str) – Path to CSV lookup file

Returns:

Array mapping class codes to Manning’s coefficients

Return type:

numpy.ndarray

Example Usage:

from AdcircSubgrid import LookupTable

# Use default CCAP lookup
ccap_table = LookupTable.ccap_lookup()

# Or load custom lookup
custom_table = LookupTable.from_file('custom_manning.csv')

# Get Manning's coefficient for land cover class
manning_n = ccap_table[class_code]

SubgridData

Container class for computed subgrid correction data.

class SubgridData(n_nodes: int, n_levels: int)

Storage container for subgrid calculations.

Parameters:
  • n_nodes (int) – Number of mesh nodes

  • n_levels (int) – Number of phi levels

Attributes:

phi

Phi factors for each node and water level.

Type:

numpy.ndarray

Shape:

(n_nodes, n_levels)

manning_avg

Average Manning’s coefficients for each node.

Type:

numpy.ndarray

Shape:

(n_nodes,)

phi_levels

Water surface elevation levels.

Type:

numpy.ndarray

Shape:

(n_levels,)

SubgridOutputFile

NetCDF output file writer for subgrid data.

class SubgridOutputFile(filename: str, mesh: Mesh, subgrid_data: SubgridData, input_config: dict)

Writer for subgrid NetCDF output files.

Parameters:
  • filename (str) – Output filename

  • mesh (Mesh) – Mesh object with coordinate information

  • subgrid_data (SubgridData) – Computed subgrid corrections

  • input_config (dict) – Configuration dictionary

Methods:

write() None

Write subgrid data to NetCDF file.

Creates a CF-compliant NetCDF file with: - Node coordinates and connectivity - Phi factors and water levels - Manning’s coefficients - Metadata and provenance information

Utility Functions

Calculation Levels

Functions for determining water level distributions.

calculate_levels_linear(z_min: float, z_max: float, n_levels: int) np.ndarray

Generate linearly spaced water levels.

Parameters:
  • z_min (float) – Minimum water level

  • z_max (float) – Maximum water level

  • n_levels (int) – Number of levels

Returns:

Array of water levels

Return type:

numpy.ndarray

calculate_levels_histogram(elevations: np.ndarray, n_levels: int) np.ndarray

Generate water levels based on elevation histogram.

Parameters:
  • elevations (numpy.ndarray) – Array of elevation values

  • n_levels (int) – Number of levels

Returns:

Array of water levels based on elevation quantiles

Return type:

numpy.ndarray

JIT-Compiled Helpers

Numba-accelerated functions for performance-critical calculations.

nan_mean_jit(array: np.ndarray) float

Compute mean of array ignoring NaN values (JIT-compiled).

Parameters:

array (numpy.ndarray) – Input array

Returns:

Mean value excluding NaN

Return type:

float

nan_sum_jit(array: np.ndarray) float

Compute sum of array ignoring NaN values (JIT-compiled).

Parameters:

array (numpy.ndarray) – Input array

Returns:

Sum excluding NaN

Return type:

float

Visualization Functions

The visualization functionality is primarily accessed through command-line tools, but some functions are available for programmatic use:

plot_mesh(filename: str, variable: str, water_level: float, show: bool = False, output_filename: str = None, **kwargs)

Create mesh-level visualization of subgrid data.

Parameters:
  • filename (str) – Subgrid NetCDF filename

  • variable (str) – Variable name to plot

  • water_level (float) – Water level for visualization

  • show (bool, optional) – Display plot interactively

  • output_filename (str, optional) – Save plot to file

node_plot(filename: str, node: int, basis: str = 'phi', show: bool = False, save: str = None, index_base: int = 0)

Create node-specific visualization of phi relationships.

Parameters:
  • filename (str) – Subgrid NetCDF filename

  • node (int) – Node number to analyze

  • basis (str, optional) – Plot basis (‘phi’ or ‘wse’)

  • show (bool, optional) – Display plot interactively

  • save (str, optional) – Save plot to filename

  • index_base (int, optional) – Node numbering base (0 or 1)

Example Custom Workflow

Here’s an example of using the API for a custom processing workflow:

import numpy as np
from AdcircSubgrid import (
    InputFile, SubgridPreprocessor, Mesh, Raster,
    LookupTable, SubgridData, SubgridOutputFile
)

def custom_subgrid_workflow(config_file, custom_processing=False):
    """
    Example custom workflow using the Python API.
    """

    # Load configuration
    config = InputFile(config_file)
    config_data = config.data()

    if custom_processing:
        # Custom processing example

        # Load mesh manually
        mesh = Mesh(config_data['input']['adcirc_mesh'])

        # Load DEM
        dem = Raster(config_data['input']['dem'])

        # Load land cover
        landcover = Raster(config_data['input']['land_cover'])

        # Get lookup table
        if config_data['input']['manning_lookup'] == 'ccap':
            lookup = LookupTable.ccap_lookup()
        else:
            lookup = LookupTable.from_file(config_data['input']['manning_lookup'])

        # Initialize data container
        n_levels = config_data['options']['n_phi_levels']
        subgrid_data = SubgridData(mesh.nodes, n_levels)

        # Custom processing loop (simplified example)
        for i in range(mesh.nodes):
            # Extract elevations around node
            # (This is a simplified example - actual implementation is more complex)
            x, y = mesh.x[i], mesh.y[i]

            # Sample DEM in neighborhood of node
            # ... custom sampling logic ...

            # Compute phi relationships
            # ... custom phi calculation ...

            # Store results
            # subgrid_data.phi[i, :] = computed_phi
            # subgrid_data.manning_avg[i] = computed_manning

        # Write results
        output_writer = SubgridOutputFile(
            config_data['output']['filename'],
            mesh,
            subgrid_data,
            config_data
        )
        output_writer.write()

    else:
        # Standard processing
        processor = SubgridPreprocessor(config)
        processor.process()
        processor.write()

# Usage
custom_subgrid_workflow('config.yaml', custom_processing=False)

Error Handling

The API includes several custom exceptions:

exception SubgridError

Base exception class for subgrid processing errors.

exception ConfigurationError

Raised when configuration file validation fails.

exception MeshError

Raised when mesh file reading or processing fails.

exception RasterError

Raised when raster data processing fails.

Example error handling:

from AdcircSubgrid import SubgridPreprocessor, InputFile, ConfigurationError

try:
    config = InputFile('config.yaml')
    processor = SubgridPreprocessor(config)
    processor.process()
    processor.write()
except ConfigurationError as e:
    print(f"Configuration error: {e}")
except Exception as e:
    print(f"Processing error: {e}")

This API reference provides the foundation for advanced users who need to customize or extend the subgrid preprocessing functionality beyond what’s available through the command-line interface.