mpas_tools.viz.paraview_extractor.extract_vtk
- mpas_tools.viz.paraview_extractor.extract_vtk(filename_pattern, variable_list='all', dimension_list=None, mesh_filename=None, blocking=10000, output_32bit=False, combine=False, append=False, out_dir='vtk_files', xtime='xtime', lonlat=False, time=None, ignore_time=False, topo_dim=None, topo_cell_index=None, include_mesh_vars=False, fc_region_mask=None, temp_dir='./culled_region', use_progress_bar=True)[source]
Extract fields from a time series of NetCDF files as VTK files for plotting in ParaView.
To extract fields across multiple files, passing in a regular expression for the filename pattern, for example
filename_pattern="hist.comp.*.nc"
.By default, time-independent fields on cells are written to a file
vtk_files/staticFieldsOnCells.vtp
and time-dependent fields on cells are written to
vtk_files/timeDependentFieldsOnCells.pvd vtk_files/time_series/timeDependentFieldsOnCells.0.vtp vtk_files/time_series/timeDependentFieldsOnCells.1.vtp ...
and similarly for edges and vertices. Time-independent fields can be included in each time step of the time-dependent fields for with
combine=True
. This allows time-dependent and -independent fields to be combined in filters within ParaView at the expense of considerable additional storage space.Indices for extra dimensions can either be supplied at runtime at a prompt (if
dimension_list=None
) or via a list of strings with the dimensions and associated indices. For each extra dimension, you can specify nothing for the indices (an empty string, meaning skip any fields with this dimension), a single index, a comma-separated list of indices, or a range of indices (separated by 1 or 2 colons). For example,dimension_list=['maxEdges=', 'nVertLeves=0:10:2', 'nParticles=0,2,4,6,8']
will ignore any fields with dimension
maxEdges
, extract every other layer from the first 10 vertical levels (each into its own field) and extract the 5 specified particles.An index array can also be specified in this way (and these can be mixed with integer indices in a comma-separated list but not in a colon-separated range):
dimension_list=['nVertLeves=0,maxLevelCell']
will extract fields from the first vertical level and the vertical level with index given by
maxLevelCell
.The extractor includes optional support for extracting geometry appropriate for displaying variables at the depth of a topographic feature (typically the top or bottom of the domain) for MPAS components with a spatially variable top or bottom index (e.g.
maxLevelCell
in MPAS-Ocean). This is accomplished with arguments such as:topo_dim='nVertLevels', topo_cell_index='maxLevelCell'
Fields on cells are sampled at the topographic index and the geometry includes polygons corresponding to edges so that vertical faces between adjacent cells can be displayed. Fields are extracted as normal except that they are sampled as point data rather than cell data, allowing computations in ParaView to display the topography. A mask field is also included indicating which parts of edge polygons correspond to the boundary of the domain (
boundaryMask == 1
) and which parts of cell and edge polygons are interior (boundaryMask == 0
). Together, this can be used to plot topography by using a calculator filter like the following:coords*(1.0 + 100.0/mag(coords)*((1 - boundaryMask)*(-bottomDepth) + 10.0*boundaryMask))
If this is entered into a Calculator Filter in ParaView with the “coordinate result” box checked, the result will to display the MPAS-Ocean topography, exaggerated by a factor of 100, with a value equivalent to 10 m along boundary points of edge polygons (a “water-tight” surface).
- Parameters:
filename_pattern (str) – MPAS Filename pattern
variable_list (list of str, optional) – List of variables to extract (‘all’ for all variables, ‘allOnCells’ for all variables on cells, etc.)”
dimension_list (list of str, optional) – A list of dimensions and associated indices
mesh_filename (str, optional) – MPAS Mesh filename. By default, the first file matching
filename_pattern
will be usedblocking (int, optional) – Size of blocks when reading MPAS file
output_32bit (bool, optional) – Whether the vtk files will be written using 32bit floats
combine (bool, optional) – Whether time-independent fields are written to each file along with time-dependent fields
append (bool, optional) – Whether only vtp files that do not already exist are written out
out_dir (str, optional) – The output directory
xtime (str, optional") – The name of the xtime variable or ‘none’ to extract Time dim without xtime
lonlat (bool, optional) – Whether the resulting points are in lon-lat space, not Cartesian
time (str, optional) – Indices for the time dimension
ignore_time (bool, optional) – Whether to ignore the Time dimension if it exists for files with a Time dimension but no xtime variable (e.g. mesh file)
topo_dim (str, optional) – Dimension and range for topography dimension
topo_cell_index (str, optional) – Index array indicating the bottom of the domain (default is the topo_dim-1 for all cells)
include_mesh_vars (bool, optional) – Whether to include mesh variables as well as time-series variables in the extraction
fc_region_mask (geometric_features.FeatureCollection, optional) – A feature collection used to define a mask. The MPAS data is culled to lie within the mask before conversion to VTK proceeds
temp_dir (str, optional) – If fc_region_mask is supplied, a temporary directory where the culled mesh and time series files are stored
use_progress_bar (bool, optional) – Whether to display progress bars (problematic in logging to a file)