XBeach Procedural Workflow¶
This notebook demonstrates the procedural approach to setting up an XBeach model using rompy-xbeach. Each component is defined step-by-step in Python code.
What is rompy-xbeach?¶
rompy-xbeach is a Python package for configuring and running XBeach nearshore wave models. It provides:
- Type-safe configuration using Pydantic models
- Data interfaces for generating forcing from various sources (NetCDF, CSV, etc.)
- Grid and bathymetry tools with interpolation and extension capabilities
- Integration with rompy for model run management
Workflow Overview¶
- Define the model grid - spatial extent and resolution
- Prepare bathymetry - interpolate source data onto the grid
- Set up forcing - wave boundaries, wind, and tide
- Configure physics - wave breaking, friction, etc.
- Generate workspace - create XBeach input files
For the declarative approach using YAML configuration, see the companion notebook.
from pathlib import Path
from datetime import timedelta
import shutil
import pandas as pd
import xarray as xr
import matplotlib.pyplot as plt
from cartopy import crs as ccrs
import warnings
warnings.filterwarnings("ignore")
Setup¶
# Data directory (using test data from rompy-xbeach)
DATA_DIR = Path("../../../../rompy-xbeach/tests/data")
# Output directory
OUT_DIR = Path("example-procedural-output")
if OUT_DIR.exists():
shutil.rmtree(OUT_DIR)
OUT_DIR.mkdir()
Step 1: Define Model Times¶
The TimeRange defines the simulation period. This controls:
- Start and end times for forcing data extraction
- Simulation duration (
tstopin XBeach)
from rompy.core.time import TimeRange
times = TimeRange(
start="2023-01-01T00:00",
end="2023-01-01T12:00",
interval="1h",
)
print(f"Simulation period: {times.start} to {times.end}")
print(f"Duration: {times.duration}")
Simulation period: 2023-01-01 00:00:00 to 2023-01-01 12:00:00 Duration: 12:00:00
Step 2: Define Model Grid¶
The RegularGrid defines the computational domain in XBeach coordinates.
The grid is defined by:
- Origin (
ori): Reference point in geographic coordinates - Rotation (
alfa): Grid rotation angle (degrees from North) - Resolution (
dx,dy): Cell sizes in meters - Extent (
nx,ny): Number of cells in x and y directions - CRS: Projected coordinate system for the grid
from rompy_xbeach.grid import RegularGrid
grid = RegularGrid(
ori=dict(x=115.594239, y=-32.641104, crs="EPSG:4326"),
alfa=347.0, # Rotation angle (degrees from North)
dx=10.0, # Cross-shore resolution (m)
dy=15.0, # Alongshore resolution (m)
nx=230, # Number of cross-shore cells
ny=220, # Number of alongshore cells
crs="EPSG:28350", # Projected CRS (MGA Zone 50)
)
# Visualize the grid
ax = grid.plot(scale="f")
ax.set_title("Model Grid")
Text(0.5, 1.0, 'Model Grid')
Step 3: Prepare Bathymetry¶
Bathymetry is prepared using XBeachBathy which handles:
- Loading source data (GeoTIFF, NetCDF, etc.)
- Interpolation onto the model grid
- Seaward extension to ensure uniform depth at the offshore boundary
- Lateral smoothing at the domain edges
3.1 Inspect Source Data¶
# Load and visualize the source bathymetry
bathyfile = DATA_DIR / "bathy.tif"
bathy_source = xr.open_dataset(bathyfile, engine="rasterio")
fig, ax = plt.subplots(figsize=(10, 8), subplot_kw=dict(projection=grid.projection))
bathy_source.band_data.isel(band=0).plot(
ax=ax,
transform=ccrs.PlateCarree(),
cmap="terrain",
vmin=-20, vmax=20,
cbar_kwargs=dict(label="Elevation (m)"),
)
grid.plot(ax=ax, grid_kwargs=dict(facecolor="none", edgecolor="red", alpha=0.8))
ax.set_title("Source Bathymetry with Model Grid")
Text(0.5, 1.0, 'Source Bathymetry with Model Grid')
3.2 Configure Bathymetry Processing¶
from rompy_xbeach.source import SourceGeotiff
from rompy_xbeach.interpolate import RegularGridInterpolator
from rompy_xbeach.data.bathy import XBeachBathy, SeawardExtensionLinear
# Source data
source = SourceGeotiff(filename=bathyfile)
# Interpolation method
interpolator = RegularGridInterpolator(
kwargs=dict(method="linear", fill_value=None)
)
# Seaward extension (extends grid offshore to reach target depth)
extension = SeawardExtensionLinear(
depth=25.0, # Target depth at offshore boundary (m)
slope=0.1, # Extension slope
)
# Create bathymetry object
bathy = XBeachBathy(
source=source,
posdwn=False, # Depths are negative (elevation positive up)
interpolator=interpolator,
extension=extension,
left=5, # Lateral smoothing cells (left edge)
right=5, # Lateral smoothing cells (right edge)
)
3.3 Generate and Visualize Bathymetry¶
# Generate bathymetry files
destdir = OUT_DIR / "bathy_check"
destdir.mkdir(exist_ok=True)
xfile, yfile, depfile, grid_extended = bathy.get(destdir=destdir, grid=grid)
# Load and visualize
dset = xr.Dataset.xbeach.from_xbeach(depfile, grid_extended)
dset.xbeach.plot_model_bathy(grid_extended, posdwn=False)
plt.title("Interpolated Model Bathymetry")
Text(0.5, 1.0, 'Interpolated Model Bathymetry')
Step 4: Set Up Forcing¶
XBeach requires forcing for:
- Wave boundary - offshore wave conditions
- Wind - surface wind stress
- Tide - water level variations
4.1 Wave Boundary¶
Wave boundaries can be generated from spectral data or integrated parameters.
Here we use BoundaryStationSpectraSwan to create SWAN-format spectral files.
from rompy_xbeach.source import SourceCRSWavespectra
from rompy_xbeach.data.boundary import BoundaryStationSpectraSwan
# Define spectral data source
wave_source = SourceCRSWavespectra(
uri=DATA_DIR / "ww3-spectra-20230101-short.nc",
reader="read_ww3",
)
# Create wave boundary
wave = BoundaryStationSpectraSwan(source=wave_source)
Visualize Wave Data¶
# Plot wave heights at first timestep
ds = wave_source.open().sel(time=times.start).squeeze()
fig, ax = plt.subplots(figsize=(10, 8), subplot_kw=dict(projection=grid.projection))
p = ax.scatter(
ds.lon, ds.lat,
s=100, c=ds.spec.hs(),
cmap="viridis", vmin=2.0, vmax=3.6,
transform=ccrs.PlateCarree(),
)
ax = grid.plot(
ax=ax,
scale="h",
grid_kwargs=dict(facecolor="black", edgecolor="red", alpha=1.0),
set_extent=False,
set_gridlines=False,
show_offshore=False,
show_origin=False,
)
plt.colorbar(p, label="Hs (m)")
ax.set_title(f"Wave Height at {times.start}")
ax.set_extent([114.3, 116, -33.7, -31.5], crs=ccrs.PlateCarree())
4.2 Wind Forcing¶
Wind forcing is extracted from station or gridded data at the grid centre or offshore boundary.
from rompy_xbeach.source import SourceCRSFile
from rompy_xbeach.data.wind import WindStation, WindVector
# Define wind data source
wind_source = SourceCRSFile(
uri=DATA_DIR / "smc-params-20230101.nc",
crs=4326,
)
# Create wind forcing (using vector components)
wind = WindStation(
source=wind_source,
coords=dict(s="seapoint"),
wind_vars=WindVector(u="uwnd", v="vwnd"),
)
Visualize Wind Data¶
ds = wind_source.open().sel(time=times.start).squeeze()
wspd = (ds.uwnd**2 + ds.vwnd**2)**0.5
fig, ax = plt.subplots(figsize=(10, 8), subplot_kw=dict(projection=grid.projection))
p = ax.scatter(
ds.longitude,
ds.latitude,
s=50,
c=wspd,
cmap="viridis",
transform=ccrs.PlateCarree(),
)
ax = grid.plot(
ax=ax,
scale="h",
grid_kwargs=dict(facecolor="black", edgecolor="red", alpha=1.0),
set_extent=False,
set_gridlines=False,
show_offshore=False,
show_origin=False,
)
plt.colorbar(p, label="Wind Speed (m/s)")
ax.set_title(f"Wind Speed at {times.start}")
ax.set_extent([114.3, 116, -33.7, -31.5], crs=ccrs.PlateCarree())
4.3 Tide Forcing¶
Tide forcing is generated from tidal constituents using oceantide.
from rompy_xbeach.source import SourceCRSOceantide
from rompy_xbeach.data.waterlevel import TideConsGrid
# Define tidal constituent source
tide_source = SourceCRSOceantide(
reader="read_otis_binary",
kwargs=dict(
gfile=DATA_DIR / "swaus_tide_cons/grid_m2s2n2k2k1o1p1q1mmmf",
hfile=DATA_DIR / "swaus_tide_cons/h_m2s2n2k2k1o1p1q1mmmf",
ufile=DATA_DIR / "swaus_tide_cons/u_m2s2n2k2k1o1p1q1mmmf",
),
crs=4326,
)
# Create tide forcing
tide = TideConsGrid(
source=tide_source,
coords=dict(x="lon", y="lat"),
)
Step 5: Configure the Model¶
The Config object brings together all components and defines XBeach parameters.
Parameters are organized into logical groups:
- physics - wave breaking, friction, numerics
- flow_boundary - boundary conditions
- output - output variables and timing
from rompy_xbeach.config import Config, DataInterface
config = Config(
# Grid and bathymetry
grid=grid,
bathy=bathy,
# Forcing data
input=DataInterface(
wave=wave,
wind=wind,
tide=tide,
),
# Flow boundary conditions
flow_boundary=dict(
front="abs_2d",
back="abs_2d",
left="neumann",
right="neumann",
),
# Tide boundary
tide_boundary=dict(
zs0=0.0,
paulrevere="land",
),
# Output configuration
output=dict(
outputformat="netcdf",
ncfilename="xboutput.nc",
tstart=0,
tintm=3600.0,
),
)
Step 6: Generate Workspace¶
The ModelRun object generates the complete XBeach workspace with all input files.
from rompy.model import ModelRun
modelrun = ModelRun(
run_id="procedural_example",
period=times,
output_dir=OUT_DIR,
config=config,
)
# Generate the workspace
rundir = modelrun()
2026-01-11 12:32:10 [INFO] rompy.model : ┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
2026-01-11 12:32:10 [INFO] rompy.model : ┃ MODEL RUN CONFIGURATION ┃
2026-01-11 12:32:10 [INFO] rompy.model : ┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
2026-01-11 12:32:10 [INFO] rompy.model : Run ID : procedural_example
2026-01-11 12:32:10 [INFO] rompy.model : Model Type : Config
2026-01-11 12:32:10 [INFO] rompy.model : Start Time : 2023-01-01T00:00:00
2026-01-11 12:32:10 [INFO] rompy.model : End Time : 2023-01-01T12:00:00
2026-01-11 12:32:10 [INFO] rompy.model : Duration : 12 hours
2026-01-11 12:32:10 [INFO] rompy.model : Time Interval : 1:00:00
2026-01-11 12:32:10 [INFO] rompy.model : Output Directory : example-procedural-output
2026-01-11 12:32:10 [INFO] rompy.model : ┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
2026-01-11 12:32:10 [INFO] rompy.model : ┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
2026-01-11 12:32:10 [INFO] rompy.model :
2026-01-11 12:32:10 [INFO] root : ┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
2026-01-11 12:32:10 [INFO] root : ┃ MODEL CONFIGURATION (Config) ┃
2026-01-11 12:32:10 [INFO] root : ┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
2026-01-11 12:32:10 [INFO] root :
2026-01-11 12:32:10 [INFO] rompy.model : Config:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: xbeach
2026-01-11 12:32:10 [INFO] rompy.model : template: /source/csiro/rompy-xbeach/src/rompy_xbeach/templates/base
2026-01-11 12:32:10 [INFO] rompy.model : checkout: main
2026-01-11 12:32:10 [INFO] rompy.model : grid:
2026-01-11 12:32:10 [INFO] rompy.model : grid_type: base
2026-01-11 12:32:10 [INFO] rompy.model : model_type: regular
2026-01-11 12:32:10 [INFO] rompy.model : ori:
2026-01-11 12:32:10 [INFO] rompy.model : x: 115.594239
2026-01-11 12:32:10 [INFO] rompy.model : y: -32.641104
2026-01-11 12:32:10 [INFO] rompy.model : crs: EPSG:4326
2026-01-11 12:32:10 [INFO] rompy.model : alfa: 347.0
2026-01-11 12:32:10 [INFO] rompy.model : dx: 10.0
2026-01-11 12:32:10 [INFO] rompy.model : dy: 15.0
2026-01-11 12:32:10 [INFO] rompy.model : nx: 230
2026-01-11 12:32:10 [INFO] rompy.model : ny: 220
2026-01-11 12:32:10 [INFO] rompy.model : crs: EPSG:28350
2026-01-11 12:32:10 [INFO] rompy.model : bathy:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: xbeach_bathy
2026-01-11 12:32:10 [INFO] rompy.model : id: data
2026-01-11 12:32:10 [INFO] rompy.model : source:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: geotiff
2026-01-11 12:32:10 [INFO] rompy.model : filename: ../../../../rompy-xbeach/tests/data/bathy.tif
2026-01-11 12:32:10 [INFO] rompy.model : band: 1
2026-01-11 12:32:10 [INFO] rompy.model : kwargs: {}
2026-01-11 12:32:10 [INFO] rompy.model : filter:
2026-01-11 12:32:10 [INFO] rompy.model : sort: {}
2026-01-11 12:32:10 [INFO] rompy.model : subset: {}
2026-01-11 12:32:10 [INFO] rompy.model : crop: {}
2026-01-11 12:32:10 [INFO] rompy.model : timenorm: {}
2026-01-11 12:32:10 [INFO] rompy.model : rename: {}
2026-01-11 12:32:10 [INFO] rompy.model : derived: {}
2026-01-11 12:32:10 [INFO] rompy.model : variables:
2026-01-11 12:32:10 [INFO] rompy.model : [0]: data
2026-01-11 12:32:10 [INFO] rompy.model : coords:
2026-01-11 12:32:10 [INFO] rompy.model : t: time
2026-01-11 12:32:10 [INFO] rompy.model : x: longitude
2026-01-11 12:32:10 [INFO] rompy.model : y: latitude
2026-01-11 12:32:10 [INFO] rompy.model : z: None
2026-01-11 12:32:10 [INFO] rompy.model : s: None
2026-01-11 12:32:10 [INFO] rompy.model : crop_data: True
2026-01-11 12:32:10 [INFO] rompy.model : buffer: 0.0
2026-01-11 12:32:10 [INFO] rompy.model : time_buffer:
2026-01-11 12:32:10 [INFO] rompy.model : [0]: 0
2026-01-11 12:32:10 [INFO] rompy.model : [1]: 0
2026-01-11 12:32:10 [INFO] rompy.model : interpolator:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: scipy_regular_grid
2026-01-11 12:32:10 [INFO] rompy.model : kwargs:
2026-01-11 12:32:10 [INFO] rompy.model : method: linear
2026-01-11 12:32:10 [INFO] rompy.model : fill_value: None
2026-01-11 12:32:10 [INFO] rompy.model : posdwn: False
2026-01-11 12:32:10 [INFO] rompy.model : left: 5
2026-01-11 12:32:10 [INFO] rompy.model : right: 5
2026-01-11 12:32:10 [INFO] rompy.model : extension:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: linear
2026-01-11 12:32:10 [INFO] rompy.model : depth: 25.0
2026-01-11 12:32:10 [INFO] rompy.model : slope: 0.1
2026-01-11 12:32:10 [INFO] rompy.model : interpolate_na: True
2026-01-11 12:32:10 [INFO] rompy.model : interpolate_na_kwargs: {}
2026-01-11 12:32:10 [INFO] rompy.model : input:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: data
2026-01-11 12:32:10 [INFO] rompy.model : wave: {}
2026-01-11 12:32:10 [INFO] rompy.model : wind:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: wind_station
2026-01-11 12:32:10 [INFO] rompy.model : id: wind
2026-01-11 12:32:10 [INFO] rompy.model : source:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: file
2026-01-11 12:32:10 [INFO] rompy.model : uri: ../../../../rompy-xbeach/tests/data/smc-params-20230101.nc
2026-01-11 12:32:10 [INFO] rompy.model : kwargs: {}
2026-01-11 12:32:10 [INFO] rompy.model : variable: None
2026-01-11 12:32:10 [INFO] rompy.model : crs: EPSG:4326
2026-01-11 12:32:10 [INFO] rompy.model : x_dim: x
2026-01-11 12:32:10 [INFO] rompy.model : y_dim: y
2026-01-11 12:32:10 [INFO] rompy.model : filter:
2026-01-11 12:32:10 [INFO] rompy.model : sort: {}
2026-01-11 12:32:10 [INFO] rompy.model : subset: {}
2026-01-11 12:32:10 [INFO] rompy.model : crop: {}
2026-01-11 12:32:10 [INFO] rompy.model : timenorm: {}
2026-01-11 12:32:10 [INFO] rompy.model : rename: {}
2026-01-11 12:32:10 [INFO] rompy.model : derived: {}
2026-01-11 12:32:10 [INFO] rompy.model : variables:
2026-01-11 12:32:10 [INFO] rompy.model : [0]: longitude
2026-01-11 12:32:10 [INFO] rompy.model : [1]: latitude
2026-01-11 12:32:10 [INFO] rompy.model : [2]: uwnd
2026-01-11 12:32:10 [INFO] rompy.model : [3]: vwnd
2026-01-11 12:32:10 [INFO] rompy.model : coords:
2026-01-11 12:32:10 [INFO] rompy.model : t: time
2026-01-11 12:32:10 [INFO] rompy.model : x: longitude
2026-01-11 12:32:10 [INFO] rompy.model : y: latitude
2026-01-11 12:32:10 [INFO] rompy.model : z: None
2026-01-11 12:32:10 [INFO] rompy.model : s: seapoint
2026-01-11 12:32:10 [INFO] rompy.model : crop_data: True
2026-01-11 12:32:10 [INFO] rompy.model : buffer: 0.0
2026-01-11 12:32:10 [INFO] rompy.model : time_buffer:
2026-01-11 12:32:10 [INFO] rompy.model : [0]: 1
2026-01-11 12:32:10 [INFO] rompy.model : [1]: 1
2026-01-11 12:32:10 [INFO] rompy.model : location: centre
2026-01-11 12:32:10 [INFO] rompy.model : sel_method: idw
2026-01-11 12:32:10 [INFO] rompy.model : sel_method_kwargs: {}
2026-01-11 12:32:10 [INFO] rompy.model : wind_vars:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: wind_vector
2026-01-11 12:32:10 [INFO] rompy.model : u: uwnd
2026-01-11 12:32:10 [INFO] rompy.model : v: vwnd
2026-01-11 12:32:10 [INFO] rompy.model : tide:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: tide_cons_grid
2026-01-11 12:32:10 [INFO] rompy.model : id: tide
2026-01-11 12:32:10 [INFO] rompy.model : source:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: oceantide
2026-01-11 12:32:10 [INFO] rompy.model : reader: read_otis_binary
2026-01-11 12:32:10 [INFO] rompy.model : kwargs:
2026-01-11 12:32:10 [INFO] rompy.model : gfile: ../../../../rompy-xbeach/tests/data/swaus_tide_cons/grid_m2s2n2k2k1o1p1q1mmmf
2026-01-11 12:32:10 [INFO] rompy.model : hfile: ../../../../rompy-xbeach/tests/data/swaus_tide_cons/h_m2s2n2k2k1o1p1q1mmmf
2026-01-11 12:32:10 [INFO] rompy.model : ufile: ../../../../rompy-xbeach/tests/data/swaus_tide_cons/u_m2s2n2k2k1o1p1q1mmmf
2026-01-11 12:32:10 [INFO] rompy.model : crs: EPSG:4326
2026-01-11 12:32:10 [INFO] rompy.model : x_dim: lon
2026-01-11 12:32:10 [INFO] rompy.model : y_dim: lat
2026-01-11 12:32:10 [INFO] rompy.model : filter:
2026-01-11 12:32:10 [INFO] rompy.model : sort: {}
2026-01-11 12:32:10 [INFO] rompy.model : subset: {}
2026-01-11 12:32:10 [INFO] rompy.model : crop: {}
2026-01-11 12:32:10 [INFO] rompy.model : timenorm: {}
2026-01-11 12:32:10 [INFO] rompy.model : rename: {}
2026-01-11 12:32:10 [INFO] rompy.model : derived: {}
2026-01-11 12:32:10 [INFO] rompy.model : variables:
2026-01-11 12:32:10 [INFO] rompy.model : [0]: h
2026-01-11 12:32:10 [INFO] rompy.model : coords:
2026-01-11 12:32:10 [INFO] rompy.model : t: time
2026-01-11 12:32:10 [INFO] rompy.model : x: lon
2026-01-11 12:32:10 [INFO] rompy.model : y: lat
2026-01-11 12:32:10 [INFO] rompy.model : z: None
2026-01-11 12:32:10 [INFO] rompy.model : s: None
2026-01-11 12:32:10 [INFO] rompy.model : crop_data: True
2026-01-11 12:32:10 [INFO] rompy.model : buffer: 0.0
2026-01-11 12:32:10 [INFO] rompy.model : time_buffer:
2026-01-11 12:32:10 [INFO] rompy.model : [0]: 1
2026-01-11 12:32:10 [INFO] rompy.model : [1]: 1
2026-01-11 12:32:10 [INFO] rompy.model : location: centre
2026-01-11 12:32:10 [INFO] rompy.model : sel_method: sel
2026-01-11 12:32:10 [INFO] rompy.model : sel_method_kwargs:
2026-01-11 12:32:10 [INFO] rompy.model : method: nearest
2026-01-11 12:32:10 [INFO] rompy.model : tideloc: 1
2026-01-11 12:32:10 [INFO] rompy.model : freq: 1h
2026-01-11 12:32:10 [INFO] rompy.model : physics:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: physics
2026-01-11 12:32:10 [INFO] rompy.model : wavemodel: None
2026-01-11 12:32:10 [INFO] rompy.model : advection: None
2026-01-11 12:32:10 [INFO] rompy.model : avalanching: None
2026-01-11 12:32:10 [INFO] rompy.model : bedfriction: None
2026-01-11 12:32:10 [INFO] rompy.model : cyclic: None
2026-01-11 12:32:10 [INFO] rompy.model : flow: None
2026-01-11 12:32:10 [INFO] rompy.model : gwflow: None
2026-01-11 12:32:10 [INFO] rompy.model : lwave: None
2026-01-11 12:32:10 [INFO] rompy.model : roller: None
2026-01-11 12:32:10 [INFO] rompy.model : setbathy: None
2026-01-11 12:32:10 [INFO] rompy.model : ships: None
2026-01-11 12:32:10 [INFO] rompy.model : single_dir: None
2026-01-11 12:32:10 [INFO] rompy.model : snells: None
2026-01-11 12:32:10 [INFO] rompy.model : swave: None
2026-01-11 12:32:10 [INFO] rompy.model : swrunup: None
2026-01-11 12:32:10 [INFO] rompy.model : vegetation: None
2026-01-11 12:32:10 [INFO] rompy.model : viscosity: None
2026-01-11 12:32:10 [INFO] rompy.model : wci: None
2026-01-11 12:32:10 [INFO] rompy.model : wind: None
2026-01-11 12:32:10 [INFO] rompy.model : flow_numerics: None
2026-01-11 12:32:10 [INFO] rompy.model : wave_numerics: None
2026-01-11 12:32:10 [INFO] rompy.model : constants: None
2026-01-11 12:32:10 [INFO] rompy.model : coriolis: None
2026-01-11 12:32:10 [INFO] rompy.model : sediment:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: sediment
2026-01-11 12:32:10 [INFO] rompy.model : sedtrans: None
2026-01-11 12:32:10 [INFO] rompy.model : numerics: None
2026-01-11 12:32:10 [INFO] rompy.model : morphology: None
2026-01-11 12:32:10 [INFO] rompy.model : bed_composition: None
2026-01-11 12:32:10 [INFO] rompy.model : bed_update: None
2026-01-11 12:32:10 [INFO] rompy.model : groundwater: None
2026-01-11 12:32:10 [INFO] rompy.model : q3d: None
2026-01-11 12:32:10 [INFO] rompy.model : mpi:
2026-01-11 12:32:10 [INFO] rompy.model : mmpi: None
2026-01-11 12:32:10 [INFO] rompy.model : mpiboundary: None
2026-01-11 12:32:10 [INFO] rompy.model : nmpi: None
2026-01-11 12:32:10 [INFO] rompy.model : output:
2026-01-11 12:32:10 [INFO] rompy.model : model_type: output
2026-01-11 12:32:10 [INFO] rompy.model : outputformat: netcdf
2026-01-11 12:32:10 [INFO] rompy.model : outputprecision: None
2026-01-11 12:32:10 [INFO] rompy.model : ncfilename: xboutput.nc
2026-01-11 12:32:10 [INFO] rompy.model : nrugdepth: None
2026-01-11 12:32:10 [INFO] rompy.model : rugdepth: None
2026-01-11 12:32:10 [INFO] rompy.model : timings: None
2026-01-11 12:32:10 [INFO] rompy.model : tstart: 0.0
2026-01-11 12:32:10 [INFO] rompy.model : tintc: None
2026-01-11 12:32:10 [INFO] rompy.model : tintg: None
2026-01-11 12:32:10 [INFO] rompy.model : tintm: 3600.0
2026-01-11 12:32:10 [INFO] rompy.model : tintp: None
2026-01-11 12:32:10 [INFO] rompy.model : tinth: None
2026-01-11 12:32:10 [INFO] rompy.model : writehotstart: None
2026-01-11 12:32:10 [INFO] rompy.model : tsglobal: None
2026-01-11 12:32:10 [INFO] rompy.model : tsmean: None
2026-01-11 12:32:10 [INFO] rompy.model : tspoint: None
2026-01-11 12:32:10 [INFO] rompy.model : projection: None
2026-01-11 12:32:10 [INFO] rompy.model : remdryoutput: None
2026-01-11 12:32:10 [INFO] rompy.model : rotate: None
2026-01-11 12:32:10 [INFO] rompy.model : flow_boundary:
2026-01-11 12:32:10 [INFO] rompy.model : front: abs_2d
2026-01-11 12:32:10 [INFO] rompy.model : back: abs_2d
2026-01-11 12:32:10 [INFO] rompy.model : left: neumann
2026-01-11 12:32:10 [INFO] rompy.model : right: neumann
2026-01-11 12:32:10 [INFO] rompy.model : lateralwave: None
2026-01-11 12:32:10 [INFO] rompy.model : nc: None
2026-01-11 12:32:10 [INFO] rompy.model : highcomp: None
2026-01-11 12:32:10 [INFO] rompy.model : tide_boundary:
2026-01-11 12:32:10 [INFO] rompy.model : tideloc: None
2026-01-11 12:32:10 [INFO] rompy.model : tidetype: None
2026-01-11 12:32:10 [INFO] rompy.model : zs0: 0.0
2026-01-11 12:32:10 [INFO] rompy.model : paulrevere: land
2026-01-11 12:32:10 [INFO] rompy.model : hotstart: None
2026-01-11 12:32:10 [INFO] rompy.model : tunits: None
2026-01-11 12:32:10 [INFO] rompy.model :
2026-01-11 12:32:10 [INFO] rompy.model : ┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
2026-01-11 12:32:10 [INFO] rompy.model : ┃ STARTING MODEL GENERATION ┃
2026-01-11 12:32:10 [INFO] rompy.model : ┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
2026-01-11 12:32:10 [INFO] rompy.model : Preparing input files in example-procedural-output
2026-01-11 12:32:10 [INFO] rompy.model : Processing model configuration...
2026-01-11 12:32:10 [INFO] rompy.model : Running configuration callable...
2026-01-11 12:32:10 [INFO] rompy_xbeach.config : Generating wave boundary data
2026-01-11 12:32:10 [INFO] rompy_xbeach.config : Generating wind forcing data
2026-01-11 12:32:10 [INFO] rompy_xbeach.config : Generating tide forcing data
2026-01-11 12:32:12 [INFO] rompy.model : Rendering model templates to example-procedural-output/procedural_example...
2026-01-11 12:32:12 [INFO] rompy.core.render : Template source: /source/csiro/rompy-xbeach/src/rompy_xbeach/templates/base
2026-01-11 12:32:12 [INFO] rompy.core.render : Output directory: example-procedural-output
2026-01-11 12:32:12 [INFO] rompy.core.render : Using template version: main
2026-01-11 12:32:12 [INFO] rompy.core.render : • Locating template repository...
2026-01-11 12:32:12 [INFO] rompy.core.render : Template repository located at: /source/csiro/rompy-xbeach/src/rompy_xbeach/templates/base
2026-01-11 12:32:12 [INFO] rompy.core.render : • Generating files from template...
2026-01-11 12:32:12 [INFO] rompy.core.render : • Rendering time: 0.02 seconds
2026-01-11 12:32:12 [INFO] rompy.core.render : • Total process time: 0.02 seconds
2026-01-11 12:32:12 [INFO] rompy.core.render : • Files created: 7
2026-01-11 12:32:12 [INFO] rompy.core.render : • Output location: /source/csiro/rompy-notebooks/notebooks/xbeach/as_python_script/example-procedural-output/procedural_example
2026-01-11 12:32:12 [INFO] rompy.model :
2026-01-11 12:32:12 [INFO] rompy.model : ┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
2026-01-11 12:32:12 [INFO] rompy.model : ┃ MODEL GENERATION COMPLETE ┃
2026-01-11 12:32:12 [INFO] rompy.model : ┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
2026-01-11 12:32:12 [INFO] rompy.model : Model files generated at: /source/csiro/rompy-notebooks/notebooks/xbeach/as_python_script/example-procedural-output/procedural_example
Step 7: Inspect Generated Files¶
# List generated files
modeldir = Path(modelrun.output_dir) / modelrun.run_id
print("Generated files:")
for f in sorted(modeldir.glob("*")):
print(f" {f.name}")
Generated files: bathy.txt params.txt swan-20230101T000000.txt tide-20230101T000000-20230101T120000.txt wind-20230101T000000-20230101T120000.txt xdata.txt ydata.txt
# View params.txt
params_file = modeldir / "params.txt"
print(params_file.read_text()[:2000]) # First 2000 chars
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%% XBeach parameter settings input file %%% %%% Date: 2026-01-10 23:32:10.253112+00:00 %%% User: rguedes-XPS-13-9350 %%% Template: /source/csiro/rompy-xbeach/src/rompy_xbeach/templates/base %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% tstop = 43200.0 tunits = seconds since 2023-01-01 00:00:00 wbctype = swan bcfile = swan-20230101T000000.txt windfile = wind-20230101T000000-20230101T120000.txt zs0file = tide-20230101T000000-20230101T120000.txt tideloc = 1 tidelen = 13 front = abs_2d back = abs_2d left = neumann right = neumann zs0 = 0.0 paulrevere = land posdwn = -1 nx = 241 ny = 229 dx = 10.0 dy = 15.0 xori = 368011.2066131959 yori = 6387580.505638544 alfa = 347.0 projection = +proj=utm +zone=50 +south +ellps=GRS80 +units=m +no_defs +type=crs depfile = bathy.txt outputformat = netcdf ncfilename = xboutput.nc tstart = 0.0 tintm = 3600.0
Inspect Forcing Files¶
def read_forcing(filepath, cols):
"""Read XBeach forcing file."""
df = pd.read_csv(filepath, header=None, sep=r"\s+", names=cols)
if "tsec" in cols:
df["time"] = [times.start + timedelta(seconds=s) for s in df.tsec]
df.set_index("time", inplace=True)
return df
Wind File¶
windfile = modeldir / config.params.get("windfile")
if windfile.exists():
df_wind = read_forcing(windfile, ["tsec", "wspd", "wdir"])
fig, (ax1, ax2) = plt.subplots(2, 1, figsize=(12, 5), sharex=True)
df_wind.wspd.plot(ax=ax1, marker="o")
ax1.set_ylabel("Wind Speed (m/s)")
ax1.grid(True, alpha=0.3)
df_wind.wdir.plot(ax=ax2, marker="o")
ax2.set_ylabel("Wind Direction (°)")
ax2.grid(True, alpha=0.3)
plt.suptitle("Generated Wind Forcing")
plt.tight_layout()
Tide File¶
tidefile = modeldir / config.params.get("zs0file")
if tidefile.exists():
df_tide = read_forcing(tidefile, ["tsec", "zs"])
fig, ax = plt.subplots(figsize=(12, 4))
df_tide.zs.plot(ax=ax, marker="o")
ax.set_ylabel("Water Level (m)")
ax.set_title("Generated Tide Forcing")
ax.grid(True, alpha=0.3)
Wave Boundary File¶
bcfile = modeldir / config.params.get("bcfile")
if bcfile.exists():
print(f"Wave boundary file: {bcfile.name}")
print(bcfile.read_text()[:500])
Wave boundary file: swan-20230101T000000.txt
SWAN 1 Swan standard spectral file
$ Created by wavespectra
$
LONLAT locations in spherical coordinates
1 number of locations
115.598402 -32.626715
AFREQ absolute frequencies in Hz
31 number of frequencies
0.03700
0.04070
0.04477
0.04925
0.05417
0.05959
0.06555
0.07210
0.07931
from wavespectra import read_swan
ds = read_swan(bcfile)
p = ds.squeeze().spec.plot()
Summary¶
This notebook demonstrated the procedural workflow for setting up an XBeach model:
- TimeRange - Define simulation period
- RegularGrid - Define computational domain
- XBeachBathy - Prepare bathymetry with interpolation and extension
- Data interfaces - Set up wave, wind, and tide forcing
- Config - Combine all components with physics parameters
- ModelRun - Generate the complete XBeach workspace
Next Steps¶
- See Example Declarative for YAML-based configuration
- See component tutorials in
components/for detailed parameter options - See data interface tutorials in
data-interface/for forcing options
Running XBeach¶
To run the model, execute XBeach in the generated workspace:
cd example-procedural-output/procedural_example
xbeach params.txt