bemio module documentation

The following sections contain autogenerated documentation from the docstrings within the bemio code and provide useful information for bemio users and developers alike. Examples describing how to use the various modules, classes, and functions are provided where possible.

bemio.io

wamit

This moduel provides functionality to read and interact with WAMIT simulation output data

class bemio.io.wamit.WamitOutput(out_file, density=1000.0, gravity=9.81, scale=False)[source]

Class to read and interact with WAMIT simulation data

Parameters:
out_file : str
Name of the wamit .out output file. In order to read scattering and Froud-Krylof forces the .3sc (scatterin) and .3fk (Froud-Krylof) coefficinet files must have the same base name as the .out file.
density : float, optional
Water density used to scale the hydrodynamic coefficient data
gravity : float, optional
Acceleration due to gravity used to scale the hydrodynamic coefficient dataaaa
scale : bool, optional
Boolean value to determine if the hydrodynamic data is scaled. See the bemio.data_structures.bem.scale function for more information
Examples

The user can create a WamitOutput data object directily as show below, or the bemio.io.wamit.read function can be used. The following example assumes there is a WAMIT output file named wamit.out.

>>> wamit_data = WamitOtuput(out_file=wamit.out)
bemio.io.wamit.read(out_file, density=1000.0, gravity=9.81, scale=False)[source]

Function to read WAMIT data into a data object of type(WamitOutput)

Parameters:
out_file : str
Name of the wamit .out output file. In order to read scattering and Froud-Krylof forces the .3sc (scatterin) and .3fk (Froud-Krylof) coefficinet files must have the same base name as the .out file.
density : float, optional
Water density used to scale the hydrodynamic coefficient data
gravity : float, optional
Acceleration due to gravity used to scale the hydrodynamic coefficient data
scale : bool, optional
Boolean value to determine if the hydrodynamic data is scaled. See the bemio.data_structures.bem.scale function for more information
Returns:
wamit_data
A WamitData object that contains the data from the WAMIT .out file specified
Examples:

The following example assumes there is a WAMIT output file named wamit.out

>>> wamit_data = read(out_file=wamit.out)

aqwa

class bemio.io.aqwa.AqwaOutput(hydro_file, list_file, scale=False)[source]

Class to read and interact with AQWA simulation data

Parameters:
hydro_file : str
Name of the AQWA .AH1 output file
list_file : str
Name of the AQWA .LIS output file
scale : bool, optional
Boolean value to determine if the hydrodynamic data is scaled. See the bemio.data_structures.bem.scale function for more information
Examples:

The following example assumes there are AQWA output files named aqwa.LIS and aqwa.LS1

>>> aqwa_data = AqwaOutput(hydro_file=aqwa.LS1, list_file=aqwa.LIS)
bemio.io.aqwa.read(hydro_file, list_file, scale=False)[source]

Function to read AQWA data into a data object of type(AqwaOutput)

Parameters:
hydro_file : str
Name of the AQWA .AH1 output file
list_file : str
Name of the AQWA .LIS output file
scale : bool, optional
Boolean value to determine if the hydrodynamic data is scaled. See the bemio.data_structures.bem.scale function for more information
Returns:
aqwa_data
A AqwaOutput object that contains hydrodynamic data
Examples:

The following example assumes there are AQWA output files named aqwa.LIS and aqwa.LS1

>>> aqwa_data = read(hydro_file=aqwa.LS1, list_file=aqwa.LIS)

nemoh

class bemio.io.nemoh.NemohOutput(sim_dir='./', cal_file='Nemoh.cal', results_dir='Results', mesh_dir='Mesh', scale=False)[source]

Class to read and interact with NEMOH simulation data

Parameters:
sim_dir : str, optional
Directory where NEMOH simulation results are located
cal_file : str, optional
Name of NEMOH .cal file
results_dir : float, optional
Name of the directory that contains the NEMOH results files
mesh_dir : str, optional
Name of the directory that contains the NEMOH mesh files
scale : bool, optional
Boolean value to determine if the hydrodynamic data is scaled. See the bemio.data_structures.bem.scale function for more information
Examples

The following example assumes that a NEMOH simulation was run and that there is data ./Results and ./Mesh directories. The Nemoh.cal file is assumed to be located at ./Nemoh.cal

>>> nemoh_data = NemohOtuput()
read_hydrostatics(file, body_num)[source]

Function to read NEMOH hydrostatic data

Parameters:
file : str
Name of the file containing the hydrostatic data
body_num : int
Number of the body corresponding to the Nemoh.cal input file
Returns:
The function does not directily return any variables, but calculates self.body[body_num].disp_vol (displace volume), self.body[body_num].wp_area (water plane area), and self.body[body_num].cb (center of gravity)
Examples:

This example assumes there is a file called Hydrostatics_1.dat that contains hydrodynamic data for body 1 in the Nemoh.cal file and that there is a NemohOutput object called nemoh_data already created.

>>> nemoh_data.read_hydrostatics(body_num=1,file='./Hydrostatics_1.dat')
read_kh(file, body_num)[source]

Function to read NEMOH linear spring stiffness data

Parameters:
file : str
Name of the file containing the linear spring stifness data
body_num : int
Number of the body corresponding to the Nemoh.cal input file
Returns:
The function does not directily return any variables, but calculates self.body[body_num].k (linear spring stiffness)
Examples:

This example assumes there is a file called KH_1.dat that contains linear spring stiffness data for body 1 in the Nemoh.cal file and that there is a NemohOutput object called nemoh_data already created.

>>> nemoh_data.read_kh(body_num=1, file='KH_1.dat')
bemio.io.nemoh.read(sim_dir='./', cal_file='Nemoh.cal', results_dir='Results', mesh_dir='Mesh', scale=False)[source]

Function to read NEMOH data into a data object of type(NemohOutput)

Parameters:
sim_dir : str, optional
Directory where NEMOH simulation results are located
cal_file : str, optional
Name of NEMOH .cal file
results_dir : float, optional
Name of the directory that contains the NEMOH results files
mesh_dir : str, optional
Name of the directory that contains the NEMOH mesh files
scale : bool, optional
Boolean value to determine if the hydrodynamic data is scaled. See the bemio.data_structures.bem.scale function for more information
Returns:
nemoh_data
A NemohOutput object that contains the hydrodynamic data
Examples

The following example assumes that a NEMOH simulation was run and that there is data ./Results and ./Mesh directories. The Nemoh.cal file is assumed to be located at ./Nemoh.cal

>>> nemoh_data = NemohOtuput()

output

bemio.io.output.write_hdf5(bemio_obj, out_file=None)[source]

Function that writes NEMOH, WAMIT, or NEMOH data to a standard human readable data format that uses the HDF5 format. This data can easily be input into various codes, such as MATLAB, Python, C++, etc. The data can easily be viewed using HDFVIEW.

Parameters:
data_object : {bemio.io.NemohOutput, bamio.io.WamitOutput, bemio.io.AqwaOutput}
A data object created using the bemio.io data readers
out_file : str, optional
The name of the output file. The file should have the .h5 file extension
Examples:

This example assumes there is a wamit output file called wamit.out that is read using the bemio.io.wamit.read function

>>> from bemio.io.wamit import read
>>> from bemio.io.output import write_hdf5
>>> wamit_data = read(out_file=wamit.out)
>>> write_hdf5(wamit_data)
Writing HDF5 data to ./wamit.h5

bemio.data_structures

bem

class bemio.data_structures.bem.HydrodynamicCoefficients[source]

bemio.io.bem internal calss to store hydrodynamic coefficients

class bemio.data_structures.bem.HydrodynamicData[source]

Class to store hydrodynamic data from NEMOH, WAMIT, and AQWA simulations

Attribuites:
rho : float
Water density
g : float
Acceleration due to gravity
wave_dir : np.array
Array of wave directions used to calculate hydrodynamic coefficients
num_bodies : int
Number of bodies in the BEM simulation
cg : np.array, shape = [3,1]
Center of gravity
cg : np.array, shape = [3,1]
Center of buoyancy
k : np.array, shape = [6,6]
Linear hydrostatic restoring stiffness
T : np.array
Wave periods considered in BEM simulation
w : np.array
Wave frequencies considered in BEM simulation
wp_area : float
Water plan area
buoy_force : float
Buoyancy force acting on the body in the equelibrium position
disp_volume : float
Water volume displaced by the body in the equelibrium position
water_depth : float
Water depth considered in the BEM simulation
body_num : int
Body number of the body in the BEM simulation
scaled : bool
Flag that indicates if the hydrodynamic coefficients have been scaled. See the bemio.data_structures.bem.scale for more information.
name : string
Name of the body in the BEM simulation
bem_code : string
Name of the BWM code used to generate the hydrodynamic coefficients.
bem_raw_data : string
BEM output file
am : HydrodynamicCoefficients
Object containing the added mass coefficients. Attribuites of the am object are:
am.all : np.array
Three dimensional array containing the added mass data. The array is comprised of added mass matricies for each frequency (w) of the BEM simulation. am.all[0,0,:] contains the added mass matrix corresponding to w[0], am.all[1,1,:] contains the added mass matrix corresponding to w[1], and so on.
am.inf : np.array
Two dimensional array containing the infiniate frequency added mass data.
am.zero : np.array
Two dimensional array containing the zero frequency added mass data.
ex : HydrodynamicCoefficients
Object containing excitation force coefficinets. Attribuites of the ex object are:
ex.mag : np.array
Magnitude of the excitation force for each frequency (w) of the BEM simulation. ex.re[0,:,:] contains the coefficients corresponding to w[0], ex.mag[1,:,:] contains the coefficients corresponding to w[1], and so on.
ex.phase : np.array
Phase angle of the excitation force magnitude. The data structure is organized the same way as ex.mag.
ex.re : np.array
Real component of the excitation force. The data structure is organized the same way as ex.mag.
ex.im : np.array
Imaginary component of the excitation force. The data structure is organized the same way as ex.mag.
ex.sc : HydrodynamicCoefficients
A data object that contains wave scattering force coefficnets. The orginization of the ex.sc coefficients is the same as the ex object.
ex.fk : HydrodynamicCoefficients
A data object that contains Froud-Krylof force coefficnets. The orginization of the ex.fk coefficients is the same as the ex object.
ex.rao : HydrodynamicCoefficients
Object that contains response amplitude operator information.
ex.ssy : HydrodynamicCoefficients
Object that contains WAMIT specific SSY data. See the WAMIT users manual for more information. (Note: This variable description needs to be improved)
Note:
This object is created by the bemio.io data readers and its purpose is to contain hydrodynamic coefficients from NEMOH, WAMIT, and AQWA data is a standard format.
calc_irf_excitation(t_end=100.0, n_t=1001, n_w=1001, w_min=None, w_max=None)[source]

Function to calculate the excitation impulse response function. For more information please see equation 2.12 in:

Dynamics Modeling and Loads Analysis of an Offshore Floating Wind Turbine, J. M Jonkman, 2007, NREL/TP-500-41958
Args:
t_end : float
Calculation range for the IRF. The IRF is calculated from -t_end to t_end
n_t : int
Number of timesteps in the IRF
n_w : int
Number of frequecy steps to use in the IRF calculation. If the frequency steps are different from the loaded BEM data, the excitation coefficinets are interpolated.
Returns:

No variables are directily returned by thi function. The following internal variables are calculated:

self.ex.irf.t : np.array
Timesteps at which the IRF is calculated
self.ex.irf.w : np.array
Frequency steps used to calculate the IRF
self.ex.irf.f : np.array
The excitation force IRF. This is a
3 dimensional np.array:
dim 1: number of degrees of freedom of the floating body, typically 6, dim 2: 1, dim 3: size(self.ex.irf.t)
Note:
When using this function it is important to look at the resulting IRF to insure that it is physical. The user may need to adjust the t_end, n_t, and n_w values in order to generate a realistic IRF.
Examples:

Once a HydrodynamicData data object is created (called ‘hydro_data’ here), the calc_irf_excitation function can be called:

>>> hydro_data.calc_irf_excitation()

The data can be plotted using matplotlib

>>> import matplotlib.pyplot as plt
>>> plt.plot(hydro_data.ex.irf.t,hydro_data.ex.irf.f)
>>> plt.show()
calc_irf_radiation(t_end=100.0, n_t=1001, n_w=1001, w_min=None, w_max=None)[source]

Function to calculate the wave radiation impulse response function. For more information please see Section 2.4.2 in:

Dynamics Modeling and Loads Analysis of an Offshore Floating Wind Turbine, J. M Jonkman, 2007, NREL/TP-500-41958, http://www.nrel.gov/docs/fy08osti/41958.pdf

and Section 13.6-7 in:

WAMIT v7.0 Users Manual, http://www.wamit.com/manual.htm
Args:
t_end : float
Calculation range for the IRF. The IRF is calculated from 0 to t_end.
n_t : int
Number of timesteps in the IRF
n_w : int
Number of frequecy steps to use in the IRF calculation. If the frequency steps are different from the loaded BEM data, the radiation damping coefficinets are interpolated.
w_min : float, optional
Minimum frequency to use in IRF calculation. If no vlaue is given, this defaults to the minimum frequency from the BEM simulations
w_max : float, optional
Maximum frequency to use in IRF calculation. If no vlaue is given, this defaults to the maximum frequency from the BEM simulations
Returns:

No variables are directily returned by thi function. The following internal variables are calculated:

self.rd.irf.t : np.array
Timesteps at which the IRF is calculated
self.rd.irf.w : np.array
Frequency steps used to calculate the IRF
self.rd.irf.L : np.array
The wave radiation IRF. This is a 3 dimensional np.array with dim 1: Number of degrees of freedom for the given body, dim 2: Number of degrees of freedom in the entire floating system. For example, two 6 degree of freedom bodies will result in a dimension of 12, dim 3: size(self.ex.irf.t)
self.rd.irf.L : np.array
d/dt(self.rd.irf.L)
Note:
When using this function it is important to look at the resulting IRF to insure that it is physical. The user may need to adjust the t_end, n_t, and n_w values in order to generate a realistic IRF.
calc_ss_radiation(max_order=10, r2_thresh=0.95)[source]

Function to calculate the state space reailization of the wave radiation IRF.

Args:
max_order : int
Maximum order of the state space reailization fit
r2_thresh : float
The R^2 threshold used for the fit. THe value must be between 0 and 1
Returns:

No variables are directily returned by thi function. The following internal variables are calculated:

Ass :
time-invariant state matrix
Bss :
time-invariant input matrix
Css :
time-invariant output matrix
Dss :
time-invariant feedthrough matrix
k_ss_est :
Impusle response function as cacluated from state space approximation
status :
status of the realization, 0 - zero hydrodynamic oefficients, 1 - state space realization meets R2 thresholdm 2 - state space realization does not meet R2 threshold and at ss_max limit

Examples:

scale(scale=None)[source]

Function to scale the hydrodynamic coefficient.

Args:
scale (bool): Boolean operater to determin if hydrodynamic data should be scaled. If scale is True self.am (added mass) coeffcients are scaled by self.rho*self.g, self.ex (excitation) coefficinets are scaled by self.rho*self.g, and self.rd (radiation damping) coefficnets are scaled by self.rho*self.w. If scale is False self.am (added mass) coeffcients are scaled by 1./(self.rho*self.g), self.ex (excitation) coefficinets are scaled by 1./(self.rho*self.g), and self.rd (radiation damping) coefficnets are scaled by 1./(self.rho*self.w).
Returns:
No variables are directily returned by thi function. Hydrodynamic coefficnets are scaled as described above.
Note:
The bemio.io.nemoh, bemio.io.wamit, and bemio.io.aqwa functions read data and return it with the scale == False. If the user wishes to return data from these functions with scale == True, they shoud specify the scale argument when the data is read or should call the scale function after read. Regardless, for consistency, the data should be scaled before the calc_irf_radiation, calc_irf_excitation, or calc_ss_radiation functions are called.
Examples:

Once a HydrodynamicData data object is created (called ‘hydro_data’ here), the scalen function can be called:

>>> hydro_data.scale(scale=False)
class bemio.data_structures.bem.ImpulseResponseFunction[source]

bemio.io.bem internal calss to store impulse response function data

class bemio.data_structures.bem.Raw[source]

bemio.io.bem internal calss to store various data

class bemio.data_structures.bem.StateSpaceRealization[source]

bemio.io.bem internal calss to store state space realization data

bemio.data_structures.bem.generate_file_names(out_file)[source]

Internal bemio function to generate filenames needed by hydroData module

Parameters:
out_file : str
Name of hydrodynamic data file
Returns:
files : dictionary
A dictionary of file names used by bemio

wave_excitation

class bemio.data_structures.wave_excitation.ImpulseResponseFunction[source]

Internal bemio object to contain impulse response function (IRF) data

class bemio.data_structures.wave_excitation.WaveElevationTimeSeries[source]

Internal bemio object to contain wave elevation time series data

class bemio.data_structures.wave_excitation.WaveExcitationConvolution(irf, irf_t, eta, eta_t)[source]

Object for calculating wave excitation force time history using the convolution method

Parameters:
irf : np.array
Wave excitation force impulse response function.
irf_t : np.array
Time series corresponding to irf
eta : np.array
Wave elevation time series
eta_t : np.array
Time series corresponding to eta
Attribuites:
self.irf : ImpulseResponseFunction
Object containing excitation force IRF information
self.wave_elevation : WaveElevationTimeSeries
Object containing wave elevation time series data
self.excitation_force : WaveExcitationForce
Object containing wave excitation force data
class bemio.data_structures.wave_excitation.WaveExcitationForce[source]

Internal bemio object to contain wave excitation force data

bemio.data_structures.wave_excitation.convolution(irf, irf_t, eta, eta_t, dt=None)[source]

Function to calculate wave excitation force using the convolution method

Patrameters:
irf : np.array
Wave excitation force impulse response function.
irf_t : np.array
Time series corresponding to irf
eta : np.array
Wave elevation time series
eta_t : np.array
Time series corresponding to eta
dt : float, optional
Time step for calculating the
Returns:
excitation_force : WaveExcitationConvolution
This function returns a WaveExcitationConvolution object with the wave exciting force and other information. See the WaveExcitationConvolution for more information.
Example:

The following example assumes that variables irf, irf_t, eta, and eta_t of type type(np.array) exist in the workspace. The contents of these variables are described above.

Calculate excitation force using the convolution method

>>> ex = convolution(irf=irf, irf_t=irf_t, eta=eta, eta_t=eta_t)

Plot the data

>>> plt.figure()
>>> plt.plot(ex.excitation_force.t,ex.excitation_force.f)

bemio.mesh_utilities

mesh

This module serves provides utilities to work with the following mesh types:
  • WAMTI
  • Nemoh
  • VTK Polydata
  • STL files
The functionality provided includes:
  • Ability to read mesh the mesh formats listed above
  • Ability to convert between mesh formats listed above
  • Utilities to calculate volume, surface area, linear spring stiffness, and other related mesh parameters.
class bemio.mesh_utilities.mesh.PanelMesh(file_name)[source]

Class to store mesh data. All mesh data is currently read and stored as quad elements. Tri elements are supported, but are stored as quad elements with a repeated point.

Parameters:
file_name : str
Name of mesh file to read. Currently WAMIT (.gdf), Stereolithography (.stl), VTK PolyDATA (.vtp), and NEMOH (.dat) mesh formats are supported
Attribuites:
files : dictionary
Dictionary containg input and output file names
orig_type : str
Mesh type of input file
points : list
List of points that define the mesh. points[n] = [x coord, y coord, z coord].
faces : list
List of points that define connectivity for each face. face[n] = [point 1 index , point 2 index, point 3 index, point 4 index], where ‘point 1-4’ are integres that correspond to the point index in the points attribuite.
center_of_gravity : np.array
Center of gravity of floating body
center_of_buoyancy : np.array
Center of buoyancy
volume_vtk : np.array
Mesh volume determined using VTK
volume_x, y, and z : np.array
Mesh volume determined using internal bemio calculations
surface_area_vtk : float
Surface area determined using VTK
surface_area : float
Surface area determined using internal bemio calculations
normals : np.array
Cell normals. This arrays is size [faces.shape[0], 3]. normals[n] = [x, y, z] is a vector that defines the normal vector for face n
cell_surface_area : np.array
Cell survace area. This array is size [faces.shape[0], 3]. cell_surface_area[n] is the surface aras of face[n].
centroid : np.array
Cell centroid. This array is size [faces.shape[0], 3]. cell_surface_area[n] is centroid of face[n].
hydrostatic_stiffness : np.array
The linear hydrostatic stiffness matrix of the mesh assuming the water surface is at z=0
bounds : dictionary
The bounds of the mesh. bounds[‘min’] and bounds[‘max’] are the minimum and maximum mesh dimensions, respectively
calculate_center_of_gravity_vtk()[source]

Function to calculate the center of gravity

Note

The VTK Pytnon bindings must be installed to use this function

Examples:

This example assumes that a mesh has been read by bemio and mesh data is contained in a PanelMesh object called mesh

>>> mesh.calculate_center_of_gravity_vtk()
cell_surface_area

Getter for cell_surface_area

Calculated

center_of_buoyancy

Getter for the center_of_buoyancy variable.

Calculated as defined in Section 3.1 of the WAMIT v7.0 users manual.

hydrostatic_stiffness

Getter for the hydrostatic_stiffness variable.

Calculated as defined in Section 3.1 of the WAMIT v7.0 users manual.

open()[source]

Function to open a VTK PolyData object in the default viewer of your operating system.

Note

This function is only available for OSX and Linux systems and and requires you have a program installed that has the ability to open VTK PolyData (.vtp) files.

Example:

This example assumes that a mesh has been read by bemio and mesh data is contained in a PanelMesh object called mesh

>>> mesh.open()
scale(scale_vect)[source]

Function used to scale mesh objects in the x, y, and z directions.

Parameters:
scale_vect : list
A list that contains the x, y, and z scale factors for the mesh
Examples:

This example assumes that a mesh has been read by bemio and mesh data is contained in a PanelMesh object called mesh

Here is how to scale a mesh by a factor of 2 in the x direction and .5 in the y direction:

>>> mesh.scale(scale_vect=[2, 0.5, 1])
translate(translation_vect, translate_cog=True)[source]

Function used to translate mesh obvjects in the x, y, and z directions.

Parameters:
translation_vect : list
A list that contains the desired x, y, and z translation for the mesh
Examples:

This example assumes that a mesh has been read by bemio and mesh data is contained in a PanelMesh object called mesh

Here is how to translate a mesh by 2 in the x direction and .5 in the y direction:

>>> mesh.translate(scale_vect=[2, 0.5, 0])
view(color=[0.5, 1, 0.5], opacity=1.0, save_png=False, camera_pos=[50, 50, 50], interact=True)[source]

Function to view the mesh using the VTK library

Parameters:
color : list, optional
VTK color specification for the mesh
opackty : float, optional
VTK opacity for the mesh. Must be between 0. and 1.
save_png : bool
Boolean operater that determines if a .png image of the mesh is saved.
interact : bool, optional
Boolean operater that determines if the user can interact with the geometry (e.g. zoom and rotate) after it is displayed
camera_pos : list, optional
Camera position
Examples:

This example assumes that a mesh has been read by bemio and mesh data is contained in a PanelMesh object called mesh

>>> mesh.view()

The mesh view window must be closed in order to return command to the Python shell

write(mesh_format='VTP')[source]

Function to write NEMOH, WAMIT, or VTK PolyData formats.

Parameters:
mesh_format : string {‘VTP’, ‘WAMIT’, ‘NEMOH’}
Variable that specifies the mesh format to write.
Examples:

This example assumes that a mesh has been read by bemio and mesh data is contained in a PanelMesh object called mesh

Here is how a WAMIT mesh would be written >>> mesh.write(mesh_format=’WAMTI’)

bemio.mesh_utilities.mesh.collapse_to_plane(mesh_obj, plane_ind=2, plane_loc=-1e-05, cut_dir=1.0)[source]

Function to collapse points to a given plane

Note

This function is not yet implemented

bemio.mesh_utilities.mesh.cut_mesh(mesh_obj, plane_ind=2, plane_loc=-1e-05, cut_dir=1.0)[source]

Function to remove cells on one side of plane

Note

This function is still early in the stages of development and needs to be improved and made more robust.

Parameters:
mesh_obj : PanelMesh
Mesh object to cut
plane_ind : int, optional
Index of plane along which to cut the mesh, 0 == x, 1 == y, 2 == z
plane_loc : float
Location of the mesh cut
cut_dir : int, {1, -1}
Direction for the mesh cut
Returns:
cut_mesh : PanelMesh
Panel mesh object that has been cut as specified
Examples:
None avaiable to data
bemio.mesh_utilities.mesh.read(file_name)[source]

Function to read surface mesh files. Currently VTK PolyData (.vtk), WAMIT (.gdf), NEMOH (.dat), and Stereolithography (.stl) mesh formats are supported

Parameters:
file_name : str
Name of the mesh file
Returns:
mesh_data : PanelMesh
A PanelMesh object that contains the mesh data
Exmaples:

This example assumes that a VTK PlolyData mesh named mesh.vtp exists in the current working directory

>>> mesh = read('mesh.vtp')

The mesh can then be converted to another format using the write function. In this case a wamit mesh is created.

>>> mesh.write(mesh_format='WAMIT')

If the VTK python bindings are installed the mesh can be viewed using the following command:

>>> mesh.view()

If you would are using OSX or Linux and have Paraview installed you can view the file using the follwing command:

>>> mesh.open()