polyadcirc.pyADCIRC package

Submodules

polyadcirc.pyADCIRC.basic module

This module contains a set of simple classes for use by

class polyadcirc.pyADCIRC.basic.MPI_for_no_mpi4py

Bases: object

A stand in class for when mpi4py is not installed

class polyadcirc.pyADCIRC.basic.comm_for_no_mpi4py

Bases: object

A set of stand ins for when mpi4py is not installed

Get_rank()
Return type:int
Returns:0
Get_size()
Return type:int
Returns:1
allgather(val)
Parameters:val (object) – object to allgather
Return type:object
Returns:val
allreduce(val1, op=None)
Parameters:
  • val1 (object) – object to allreduce
  • op – operation
Return type:

object

Returns:

val1

bcast(val, root=0)
Parameters:val (object) – object to broadcast
Return type:object
Returns:val
class polyadcirc.pyADCIRC.basic.location(x, y)

Bases: polyadcirc.pyADCIRC.basic.pickleable

Store x, y coordinates

class polyadcirc.pyADCIRC.basic.node(x, y, bathymetry)

Bases: polyadcirc.pyADCIRC.basic.location

Stores data specific to a single node

class polyadcirc.pyADCIRC.basic.pickleable

Bases: object

Class that makes objects easily pickleable via a dict

class polyadcirc.pyADCIRC.basic.time(dt, statim, rnday, dramp)

Bases: polyadcirc.pyADCIRC.basic.pickleable

Stores time data specific to ADCIRC model runs from the users and fort.15 input file

polyadcirc.pyADCIRC.convert_fort14_to_fort13 module

A set of methods for converting a fort.14 formatted file to a fort.13 formatted file.

exception polyadcirc.pyADCIRC.convert_fort14_to_fort13.Error

Bases: exceptions.Exception

Base class for exceptions in this module.

exception polyadcirc.pyADCIRC.convert_fort14_to_fort13.InputError(expr, msg)

Bases: polyadcirc.pyADCIRC.convert_fort14_to_fort13.Error

Exception raised for errors in the input.

polyadcirc.pyADCIRC.convert_fort14_to_fort13.convert(source, keep_flags=0, target='fort.13')

Reads in a fort5...14 file and produces a fort5...13 file. By default all remaining flagged nodes are removed. If the user desires to keep the flagged nodes, keep_flags should be set to 1.

Parameters:
  • source (string) – the fort5...14 file name, this is used to determine the nodal manningsn values
  • target (string) – the fort.13 files, this is used to determine the default manningsn value and relevant header information
  • keep_flags (int) –
    1. _noflags – creates a fort.13 formatted file and removes remaining flags for Griddata
    2. _flags – creates a fort.13 formatted file with remainings
      flags for Griddata intact
    3. _fillins – creates a fort.13 formatted file and fills in
      remiaing flags from Griddata with data from source if data is present
    4. _color – creates a fort.13 formatted file and color
      codes the data based on wheter or not is is a default value or not
polyadcirc.pyADCIRC.convert_fort14_to_fort13.convert_go(grid, folder_name=None, keep_flags=0)

See convert() where source is the final *.14 file produced by the bash script associated with grid in folder_name

Parameters:
  • gridgridInfo object
  • folder_name (string) – folder containing the fort5...14 file
  • keep_flags (int) –
    1. _noflags – creates a fort.13 formatted file and removes remaining flags for Griddata
    2. _flags – creates a fort.13 formatted file with remainings
      flags for Griddata intact
    3. _fillins – creates a fort.13 formatted file and fills in
      remiaing flags from Griddata with data from source if data is present
    4. _color – creates a fort.13 formatted file and color
      codes the data based on wheter or not is is a default value or not

polyadcirc.pyADCIRC.flag_fort14 module

Flags all of the nodes in a grid file accordining to the flagging scheme for Griddata_v1.32.F90

  1. Averaging Flag_value: -999.0 (x1), -950.0 (x2), -960.0 (x4), -970.0 (x8) This scheme evaluate the average to count the GIS database points which are inside of Grid scale. And this scheme can use different Grid scale size classified by flag_value.
  2. Highest point Flag_value: -888.0 This scheme pick up the highest point of GIS database inside of Grid_scale.
  3. Nearest point Flag_value: -777.0 This scheme pick up the nearest point of GIS database inside of Grid_scale.
polyadcirc.pyADCIRC.flag_fort14.flag_fort14(grid_file_name='fort.14', avg_scheme=2)

Modifiy grid_file_name so that all of the nodes are flagged appropriately for Griddata program and save to flagged_grid_file_name.

Parameters:
  • grid_file_name (string) – name of fort.14 formatted file
  • avg_scheme (int) – flag to choose which averaging scheme to use
Return type:

string

Returns:

flagged file path

polyadcirc.pyADCIRC.flag_fort14.flag_fort14_go(grid, avg_scheme=2)

Given a gridInfo object create a flagged version of the grid.file_name and save to grid.file_name

Parameters:
  • gridpolyadcirc.pyGriddata.gridInfo.gridObject
  • avg_scheme (int) – flag to choose which averaging scheme to use

See flag_fort14()

Return type:string
Returns:flagged file path

polyadcirc.pyADCIRC.fort13_management module

This module fort13_management handles the reading/writing of fort.13 formatted files.

polyadcirc.pyADCIRC.fort13_management.read_default(data, path=None, file_name='fort.13')

Read in default manningsn from a *.13 file and update data accordingly

Parameters:
Return type:

polyadcirc.run_framework.domain

Returns:

a reference to the polyadcirc.run_framework.domain object containing the default value

polyadcirc.pyADCIRC.fort13_management.read_manningsn(fid)

Reads in a node value line from fid

Parameters:fid (file) – the file object to be read from we might want to re write this as the fromstring method returns an array :rtype: numpy.ndarray
Returns:Returns an array([node, value])
polyadcirc.pyADCIRC.fort13_management.read_nodal_attr(data, path=None, file_name='fort.13', nums=None)

Load in nodal attributes from a *.13 file (only does Manning’s n for now) and return a dictonary (like a MATLAB struct) with these attributes).

Parameters:
Return type:

dict

Returns:

dictionary of Manning’s n values

polyadcirc.pyADCIRC.fort13_management.read_nodal_attr_dict(path=None, file_name='fort.13')

Load in nodal attributes from a *.13 file (only does Manning’s n for now) and return a dictonary (like a MATLAB struct) with these attributes).

Parameters:
  • path (string or None) – the directory containing the fort.13 to be read in
  • file_name (string) – the name of the fort.13 formatted file
Return type:

dict

Returns:

dictionary of Manning’s n

polyadcirc.pyADCIRC.fort13_management.read_node_num(path=None, file_name='fort.13')

Read in the number of nodes in the mesh from a *.13 file

Parameters:
  • path (string or None) – directory containing file_name
  • file_name (string) – file name
Return type:

int

Returns:

number of nodes in the mesh

polyadcirc.pyADCIRC.fort13_management.update_mann(data, path=None, default=None, file_name='fort.13')

Write out fort.13 to path with the attributes contained in Data.

Parameters:
  • data (numpy.ndarray or dict) – containing the nodal attribute information
  • path (string or None) – the directory to which the fort.13 file will be written
  • default (None or float) – default value
  • file_name (string) – the name of the fort.13 formatted file
polyadcirc.pyADCIRC.fort13_management.write_manningsn(fid, node, value)

Write out a formated node value line to fid

Parameters:
  • fid (file) – the file object to be written to
  • node (int) – the node number
  • value (float) – the nodal value
Return type:

string

Returns:

formatted string for a line in a fort.13 formatted file

polyadcirc.pyADCIRC.fort14_management module

This module, fort14_management handles the reading/writing of fort.14 formatted files.

polyadcirc.pyADCIRC.fort14_management.clean(grid_object, folder_name=None)

Removes all except for the most recent *.14 files in folder_name or the current working direcory

Parameters:
  • grid_objectgridInfo
  • folder_name (string) – folder to clean
Return type:

list

Returns:

list of fort.14 files in folder_name

polyadcirc.pyADCIRC.fort14_management.flag(grid_file_name='fort.14', avg_scheme=2)

Modifiy grid_file_name so that all of the nodes are flagged appropriately for Griddata program and save to flagged_grid_file_name.

Parameters:
  • grid_file_name (string) – name of fort.14 formatted file
  • avg_scheme (int) – averaging scheme flag

See flag_fort14()

polyadcirc.pyADCIRC.fort14_management.flag_go(grid, avg_scheme=2)

Given a gridInfo object create a flagged version of the grid.file_name[8:] and save to grid.file_name

Parameters:
  • gridgridInfo
  • avg_scheme (int) – averaging scheme flag

See flag_fort14()

Return type:string
Returns:flagged file name
polyadcirc.pyADCIRC.fort14_management.is_flagged(grid)
Parameters:gridgridInfo
Return type:bool
Returns:true if flagged_grid_file_name exists and false if it doesn’t exist
polyadcirc.pyADCIRC.fort14_management.read_spatial_grid(data, path=None, make_domain_map=False)

Reads in a fort.14 file in path and updates data

Parameters:
  • data (polyadcirc.run_framework.domain) – python object to save the fort.14 data to
  • path (string or None) – path to the``fort.14`` fortmatted file
  • make_domain_map (bool) – flag for whether or not to make a node to element map
Returns:

reference to data

polyadcirc.pyADCIRC.fort14_management.read_spatial_grid_header(data, path=None)

Reads in the the number of nodes and elements from the fort.14 file in path and updates data

Parameters:
Returns:

reference to data

polyadcirc.pyADCIRC.fort14_management.update(data, bathymetry=None, path=None, file_name='fort.14')

Write out bathymetry in data to fort.14 formated file, by updating path/file_name accordingly

Parameters:
  • data (polyadcirc.run_framework.domain) – python object to save the fort.14 data to
  • bathymetry (array or None) – if None then use the bathymetry in data
  • path (string or None) – path to the``fort.14`` fortmatted file
  • file_name (string) – file name

polyadcirc.pyADCIRC.fort15_management module

This module, fort15_management, handles the reading/writing of fort.15 formatted files.

polyadcirc.pyADCIRC.fort15_management.array_to_loc_list(station_array)

Convert a numpy.ndarray into a list of location.

Parameters:station_array (numpy.ndarray) – an array of shape (n, 2) where n is the number of station locations
Return type:list
Returns:list of location
polyadcirc.pyADCIRC.fort15_management.fake_stations(domain, num_stat)

Given a domain and a number of stations creates a grid of stations with num_stat stations in the x and y directions.

Parameters:
  • domain (domain) – physical domain used to define the minimum and maximum extent of the grid of stations to be created
  • num_stat (int or iterable of length 2) – number of stations in the x and y directions
Return type:

list

Returns:

list of location

polyadcirc.pyADCIRC.fort15_management.loc_list_to_array(station_list)

Convert a list of station locations to an array of x and y values.

Parameters:station_list (list) – list of location
Return type:numpy.ndarray of shape (n, 2)
Returns:array of station locations
polyadcirc.pyADCIRC.fort15_management.read_recording_data(data, path=None)

Reads in fort.15 and stores the following in data as a time object (DT, STATIM, RNDAY, DRAMP)

Calulcates and stores:

  • recording[key] = (meas_locs, total_obs, irtype)
  • stations[key] = list() of locations
Parameters:
  • data (domain) – object to store reference to time
  • path (string or None) – directory containing fort.15 file
Returns:

reference to data.recording and data.stations

Return type:

time

polyadcirc.pyADCIRC.fort15_management.set_ihot(ihot, path=None)

Set the hot start parameter to ihot

Parameters:ihot (int) – Hot start paramter, fort.## file to read hot start dtat from
polyadcirc.pyADCIRC.fort15_management.set_write_hot(nhstar, nhsinc, path=None)

Set the hot start file generation paramters to nhstar and nhsinc

Parameters:
  • nhstar (int) – Type of hot start file to write
  • nhsinc (int) – Frequency in timesteps to write data out to the hotstart file
polyadcirc.pyADCIRC.fort15_management.subdomain(fulldomain_path, subdomain_path)

Create a modifed fort.15 for the subdomain at subdomain_path.

  • Reduce RNDAY by 0.05%
  • Set NBFR to 0
  • Remove periodic forcing boundary frequencies
  • Remove recording stations outside of the subdomain grid and update the number of stations accordingly
Parameters:
  • flag (int) – type of subdomain 0 - ellipse, 1 - circle
  • fulldomain_path (string) – fulldomain dir containing fort.15 file
  • subdomain_path (string) – subdomain dir containing fort.15 file
polyadcirc.pyADCIRC.fort15_management.trim_locations(flag, subdomain_path, locs)

Remove locations outside of the subdomain from locs

Parameters:
  • flag (int) – type of subdomain 0 - ellipse, 1 - circle
  • subdomain_path (string) – subdomain dir containing fort.15 file
  • locs (list) – list of location objects
Return type:

list

Returns:

list of locations inside the subdomain

polyadcirc.pyADCIRC.fort15_management.trim_locations_circle(subdomain_path, locs)

Remove locations outside of the circular subdomain from locs

Parameters:
  • subdomain_path (string) – subdomain dir containing fort.15 file
  • locs (list) – list of location objects
Return type:

list

Returns:

list of locations inside the subdomain

polyadcirc.pyADCIRC.fort15_management.trim_locations_ellipse(subdomain_path, locs)

Remove locations outside of the elliptical subdomain from locs

Parameters:
  • subdomain_path (string) – subdomain dir containing fort.15 file
  • locs (list) – list of location objects
Return type:

list

Returns:

list of locations inside the subdomain

polyadcirc.pyADCIRC.fort1920_management module

This module is for the manipulation and creation of fort.19 and fort.20 files.

polyadcirc.pyADCIRC.fort1920_management.sin_wave(t_start, t_finish, amplitude, nnodes, time, periods=0.5, shift=0, timinc=None)

Creates data for a sine wave for forcing over the simulation time evert timinc.

Parameters:
  • t_start (float) – starting time of sine shaped wave in days
  • t_finish (float) – ending time of the sine shaped wave in days
  • amplitude (float) – amplitude of the sine shaped wave
  • nnodes (int) – number of nodes for this BC
  • time (time) – container for information from the fort.15
  • periods (float) – number of periods to include in the wave
  • shift (float) – number of periods to shift the wave
  • timinc (int) – time increment (secs) between consecutive sets of data
Return type:

tuple of (numpy.ndarray, numpy.ndarray, int)

Returns:

(times, values at times, timinc)

polyadcirc.pyADCIRC.fort1920_management.step_wave(t_start, t_finish, amplitude, nnodes, time, timinc=None)

Creates data for a square wave for forcing over the simulation time evert timinc.

Parameters:
  • t_start (float) – starting time of sine shaped wave in days
  • t_finish (float) – ending time of the step shaped wave in days
  • amplitude (float) – amplitude of the step shaped wave
  • nnodes (int) – number of nodes for this BC
  • time (time) – container for information from the fort.15
  • timinc (int) – time increment (secs) between consecutive sets of data
Return type:

tuple of (numpy.ndarray, numpy.ndarray, int)

Returns:

(times, values at times, timinc)

polyadcirc.pyADCIRC.fort1920_management.write_fort19(etiminc, esbin, file_name=None)

Write out a fort.19 formatted file to file_name

neta = total number of elevation specified boundary nodes

Parameters:
  • etimic – time increment (secs) between consecutive sets of elevation specified bounadry condition values
  • esbin (numpy.ndarray of shape (k, neta)) – elevation at specfied elevation nodes k

See also

ADCIRC Non-periodic Elevation Boundary Condition File

polyadcirc.pyADCIRC.fort1920_management.write_fort20(ftiminc, qnin, file_name=None)

Write out a fort.20 formatted file to file_name

nflbn = total number of normal flow specified boundary nodes

Parameters:
  • ftiminc (numpy.ndarray of shape (nflbn,)) – time increment (secs) between consecutive sets of normal flow specified bounadry condition values
  • qnin (numpy.ndarray of shape (k, nflbn)) – normal flow/unit width at specified nomarl flow node k

See also

ADCIRC Non-periodic Normal Flow Boundary Condition File

polyadcirc.pyADCIRC.output module

This module provides methods for retrieving data from ASCII ADCIRC formatted timeseries and non-timeseries data files and returning that data as numpy arrays.

polyadcirc.pyADCIRC.output.get_data_nts(kk, path, data, nts_data, file_names=['tinun.63'])

Retrieves data from a nontimeseries formatted files in path and adds data to nts_data

Parameters:
  • kk (int) – run number
  • path (string) – RF_directory_* path
  • datadomain
  • nts_data (dict) – reference to dict() to store data to
  • file_names (list) – list of ADCIRC output files to retrieve data from
polyadcirc.pyADCIRC.output.get_data_ts(kk, path, ts_data, time_obs, file_names=['fort.61'], timesteps=None, ihot=None)

Retrieves data from a timeseries formatted files in path and adds data to ts_data

Parameters:
  • kk (int) – run number
  • path (string) – RF_directory_* path
  • ts_data (dict) – reference to dict() to store data to
  • time_obs (dict) – reference to dict() to store time data to
  • file_names (list) – list of ADCIRC output files to
  • ihot (int) – hotstart flag (0, 67, 68)
  • timesteps (int) – number of timesteps to read
polyadcirc.pyADCIRC.output.get_nts_sr(path, data, file_name)

Retrieves data from a nontimeseries formatted file in path and adds data to nts_data

Parameters:
  • path (string) – RF_directory_* path
  • datadomain
  • file_name (string) – ADCIRC output file to retrieve data from
Return type:

numpy.ndarray

Returns:

array of dimensions (data.node_num,)

polyadcirc.pyADCIRC.output.get_ts_sr(path, file_name, get_time=False, timesteps=None, ihot=None)

Retrieves data from a timeseries formatted file in path and adds data to ts_data

Parameters:
  • path (string) – RF_directory_* path
  • file_name (string) – ADCIRC output file to retrieve data from
  • bool – flag for whether or not to record times of recordings
  • timesteps (int) – number of timesteps to read
  • ihot (int) – hotstart flag (0, 67, 68)
Return type:

numpy.ndarray

Returns:

array of dimensions (data.node_num,)

polyadcirc.pyADCIRC.plotADCIRC module

A set of functions for plotting data from runSet

polyadcirc.pyADCIRC.plotADCIRC.add_2d_axes_labels(fig=None, ics=1)

Add 2D axes labels for either cartesian or polar coordinates

Parameters:
  • fig (figure) – matplotlib.pyplot.figure object
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar)
polyadcirc.pyADCIRC.plotADCIRC.basis_functions(domain, bv_array, path=None, save=True, show=False, ics=1, ext='.png', cmap=<matplotlib.colors.LinearSegmentedColormap object>)

Given a bv_array containing basis functions, plot them

Parameters:
  • domaindomain
  • bv_array (:classnp.array`) – array of basis vectors based on land classification
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.bathymetry(domain, path=None, save=True, show=False, mesh=False, contour=False, ics=1, ext='.png', title=False)

Given a domain, plot the bathymetry

Parameters:
  • domaindomain
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • mesh (bool) – flag for whether or not to show mesh
  • contour (bool) – use tripcolor() or tricontour()
  • ics (int) – coordinate system (1 cartisian, 2 polar)
  • ext (string) – file extesion
polyadcirc.pyADCIRC.plotADCIRC.colorbar(mappable=None)

Add a colorbar to the current figure/plot/subfigure/axes

Parameters:mappable – See matplotlib
polyadcirc.pyADCIRC.plotADCIRC.field(domain, z, title, clim=None, path=None, save=True, show=False, ics=1, ext='.png', cmap=<matplotlib.colors.LinearSegmentedColormap object>)

Given a domain, plot the nodal value z

Parameters:
  • domaindomain
  • znumpy.ndarray
  • title (string) – plot title
  • climnumpy.clim
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar coords)
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.get_Triangulation(domain, path=None, save=True, show=False, ics=1, ext='.png', title=False)
Parameters:
  • domaindomain
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – coordinate system (1 cartisian, 2 polar)
  • ext (string) – file extesion
Return type:

matplotlib.tri.Triangulation

Returns:

triangulation of domain

polyadcirc.pyADCIRC.plotADCIRC.load_fig_gen_cmap()

Todo

Color map to match FigureGen

polyadcirc.pyADCIRC.plotADCIRC.mean_field(domain, points, bv_dict, path=None, save=True, show=False, ics=1, ext='.png')

Given a bv_dict a set of random points, plot the r_fields generated in random_manningsn.runSet.run_points

Parameters:
  • domaindomain
  • points (numpy.ndarray) – weights for points at which the random domain was sampled
  • bv_dict (list of dict) – list of basis vectors based on land classification
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar coords)
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.nts_contour(nts_data, domain, keys=None, points=None, path=None, save=True, show=False, ics=1, ext='.png')

Plot the non timeseries data in nts_data as a set of 2D contour plots. To plot only a subset of this data list which ADCIRC Output types to plot in keys, and which runs in points.

Note

This only applies to 1D nts data (irtype = 1)

Parameters:
  • nts_data (dict) – nts_data from runSet
  • domaindomain
  • keys (list) – list of types of ADCIRC output data to plot
  • points (list) – list of runs to plot
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.nts_line_data(nts_data, keys=None, path=None, save=True, show=False, ext='.png')

Plot the non timeseries data in data with points as the x axis. To plot only a subset of this data list which ADCIRC Output types to plot in keys.

Note

This only applies to 1D nts data (irtype = 1)

Parameters:
  • nts_data (dict) – nts_data from runSet
  • keys (list) – list of types of ADCIRC output data to plot
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.nts_pcolor(nts_data, domain, keys=None, points=None, path=None, save=True, show=False, ics=1, ext='.png')

Plot the non timeseries data in nts_data as a set of 2D color plots. To plot only a subset of this data list which ADCIRC Output types to plot in keys, and which runs in points.

Note

This only applies to 1D nts data (irtype = 1)

Parameters:
  • nts_data (:class:dict) – nts_data from :class:~polyadcirc.run_framework.random_manningsn.runSet
  • domaindomain
  • keys (list) – list of types of ADCIRC output data to plot
  • points (list) – list of runs to plot
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.random_fields(domain, points, bv_dict, path=None, save=True, show=False, ics=1, ext='.png', cmap=<matplotlib.colors.LinearSegmentedColormap object>)

Given a bv_dict a set of random points, plot the r_fields generated in run_points()

Parameters:
  • domaindomain
  • points (numpy.ndarray) – weights for points at which the random domain was sampled
  • bv_dict (list of dict) – list of basis vectors based on land classification
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar coords)
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.save_show(full_name, save, show, ext)

Save or show the current figure

Parameters:
  • full_name (string) – path to save the figure
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.station_data(ts_data, time_obs, keys=None, stations=None, path=None, save=True, show=False, ext='.png')

Plot the station data in ts_data with time_obs as the x axis. To plot only a subset of this data list which ADCIRC Output types to plot in keys.

Parameters:
  • ts_data (dict) – ts_data from runSet
  • time_obs (dict) – time_obs from runSet
  • keys (list) – list of types of ADCIRC output data to plot
  • stations (dict) – dictonary of lists of stations to plot
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.station_locations(domain, path=None, bathy=False, save=True, show=False, ics=1, ext='.png', station_markers=None)

Given a domain, plot the observation stations

Parameters:
  • domaindomain
  • path (string or None) – directory to store plots
  • bathy (bool) – flat for whether or not to show bathymetry
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – coordinate system (1 cartisian, 2 polar)
  • ext (string) – file extesion
polyadcirc.pyADCIRC.plotADCIRC.ts_pcolor(ts_data, time_obs, domain, keys=None, points=None, path=None, save=True, show=False, ics=1, ext='.png')

Plot the timeseries data in ts_data as a series of 2D color plots. To plot only a subset of this data list which ADCIRC Output types to plot in keys and which runs in points.

Note

This only applies to 1D nts data (irtype = 1)

Parameters:
  • ts_data (dict) – ts_data from runSet
  • time_obs (:class:dict) – time_obs from runSet
  • domaindomain
  • keys (list) – list of types of ADCIRC output data to plot
  • points (list) – list of runs or points to plot
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar)
  • ext (string) – file extension
polyadcirc.pyADCIRC.plotADCIRC.ts_quiver(ts_data, time_obs, domain, keys=None, points=None, path=None, save=True, show=False, ics=1, ext='.png')

Plot the timeseries data in ts_data as a series of 2D quiver plots. To plot only a subset of this data list which ADCIRC Output types to plot in keys and which runs in points.

Note

This only applies to 2D nts data (irtype = 2)

Parameters:
  • ts_data (dict) – ts_data from runSet
  • time_obs (dict) – time_obs from runSet
  • domaindomain
  • keys (list) – list of types of ADCIRC output data to plot
  • points (list) – list of runs or points to plot
  • path (string or None) – directory to store plots
  • save (bool) – flag for whether or not to save plots
  • show (bool) – flag for whether or not to show plots
  • ics (int) – polar coordinate option (1 = cart coords, 2 = polar
  • ext (string) – file extension

Todo

maybe add plt.streamplot, but this requires a regularly spaced grid and is therefore more expenseive to do

polyadcirc.pyADCIRC.post_management module

This module controls the automatic writing of in.post* files

polyadcirc.pyADCIRC.post_management.post_script_ALL(path)

Write out a bash scrip to run adcpost with in.postALL and and save it to path

Parameters:path (string) – folder to save postALL.sh to
polyadcirc.pyADCIRC.post_management.post_script_n(path, n)

Write out a bash scrip to run adcpost with in.postn and save it to path

Parameters:path (string) – folder to save postn.sh to
polyadcirc.pyADCIRC.post_management.post_script_sub(path)

Write out a bash scrip to run adcpost with in.postsub and and save it to path

Parameters:path (string) – folder to save postsub.sh to
polyadcirc.pyADCIRC.post_management.write_ALL(path, hotfiles=False)

Write out a in.postALL file and save it to path

Parameters:
  • path (string) – folder to save in.postALL to
  • hotfiles (bool) – flag whether or not to post process hot start files
polyadcirc.pyADCIRC.post_management.write_multi(path, nums, hotfiles=False)

Write out a in.postmulti file and save it to path

Parameters:
  • path (string) – folder to save in.postn to
  • nums (list) – list of ADCPOST file type codes
  • hotfiles (bool) – flag whether or not to post process hot start files
polyadcirc.pyADCIRC.post_management.write_n(path, n, hotfiles=False)

Write out a in.postn file and save it to path

Parameters:
  • path (string) – folder to save in.postn to
  • n (int) – ADCPOST file type code
  • hotfiles (bool) – flag whether or not to post process hot start files
polyadcirc.pyADCIRC.post_management.write_sub(path, hotfiles=False)

Write out a in.postsub file and save it to path

Parameters:
  • path (string) – folder to save in.postsub to
  • hotfiles (bool) – flag whether or not to post process hot start files

polyadcirc.pyADCIRC.prep_management module

This module controls the automatic writing of in.prep* files

polyadcirc.pyADCIRC.prep_management.prep_script_12(path)

Write out a bash scrip to run adcprep with in.prep1 and in.prep2 and save it to path

Parameters:path (string) – folder to save prep12.sh to
polyadcirc.pyADCIRC.prep_management.prep_script_n(path, n)

Write out a bash scrip to run adcprep with in.prepn and save it to path

Parameters:path (string) – folder to save prepn.sh to
polyadcirc.pyADCIRC.prep_management.write_1(path, nprocs=12, nfile='fort.14')

Write out a in.prep1 file for nfile and save it to path

Parameters:
  • path (string) – folder to save in.prep1 to
  • nprocs (int) – number of processors a single PADCIRC run will be run on
  • nfile (string) – name of fort.14 formatted file
polyadcirc.pyADCIRC.prep_management.write_2(path, nprocs=12)

Write out a in.prep2 file for nfile and save it to path

Parameters:
  • path (string) – folder to save in.prep2 to
  • nprocs (int) – number of processors a single PADCIRC run will be run on
polyadcirc.pyADCIRC.prep_management.write_5(path, nprocs=12, nfile='fort.13')

Write out a in.prep5 file for nfile and save it to path

Parameters:
  • path (string) – folder to save in.prep5 to
  • nprocs (int) – number of processors a single PADCIRC run will be run on
  • nfile (string) – name of fort.13 formatted file
polyadcirc.pyADCIRC.prep_management.write_n(path, n, nprocs=12, nfile=None)

Write out a in.prepn file for nfile and save it to path

Parameters:
  • path (string) – folder to save in.prepn to
  • n (int) – number of the in.prepn file
  • nprocs (int) – number of processors a single PADCIRC run will be run on
  • nfile (string) – name of fort.nn formatted file

polyadcirc.pyADCIRC.volume module

This module contains methods used in calculating the volume of water present in an ADCIRC simulation.

Todo

Some of these routines could be parallelized.

polyadcirc.pyADCIRC.volume.element_volume(domain, element, elevation)

Calculates the volume of water contained an element with a given sea surface elevation at each of the nodes.

Parameters:
  • domaindomain
  • element (array_like) – list of nodes defining an element
  • elevation – eta, sea surface height (NOT WATER COLUMN HEIGHT)
Return type:

double

Returns:

volume

polyadcirc.pyADCIRC.volume.side(domain, element, side_num, elevation)

Calculates dot(x, n*A) where x is the barycenter, n is the normal vector to the surface defined by z and the element verticies, and A is the area of the surface defined by z and the element vertices.

Parameters:
  • domaindomain
  • element (array_like) – list of nodes defining an element
  • elevation – eta, sea surface height (NOT WATER COLUMN HEIGHT)
Return type:

double

Returns:

dot(x, n*A)

polyadcirc.pyADCIRC.volume.sub_volume(domain, elevation, elements)

Calculates the total volume of water contained in an ADCIRC simulation for a given set of elements with sea surface height given by elevation.

Parameters:
  • domaindomain
  • elevation – eta, sea surface height (NOT WATER COLUMN HEIGHT)
  • elements – list of element numbers to calcuate volumes for
Return type:

tuple

Returns:

total volume, element-wise volume

polyadcirc.pyADCIRC.volume.total_volume(domain, elevation)

Calculates the total volume of water contained in an ADCIRC simulation with sea surface height given by elevation.

Parameters:
  • domaindomain
  • elevation – eta, sea surface height (NOT WATER COLUMN HEIGHT)
Return type:

tuple

Returns:

total volume, element-wise volume

polyadcirc.pyADCIRC.volume.triangle(domain, element, z, norm_dir=1.0)

Calculates dot(x, n*A) where x is the barycenter, n is the normal vector to the surface defined by z and the element verticies, and A is the area of the surface defined by z and the element vertices.

Parameters:
  • domaindomain
  • element (array_like) – list of nodes defining an element
  • z (numpy.ndarray) – z-coordinate relative to the geiod, z = eta OR z = -h
  • norm_dir (double) – 1.0 up, -1.0 down, direction of the normal vector
Return type:

double

Returns:

dot(x, n*A)

Module contents

This subpackage contains