edg_acoustics
Documentation about edg_acoustics (in init.py file).
Submodules
Package Contents
Classes
Acoustics simulation data structure for running a DG acoustics simulation. |
|
Mesh data structure generated from common mesh file formats or geometry files of GMSH. |
|
Abstract base class for boundary conditions. |
|
Setup absorptive boundary condition of a DG acoustics simulation for a specific scenario. |
|
Setup initial condition of a DG acoustics simulation for a specific scenario. |
|
Setup a monopole source for a specific scenario. |
|
abstract base class for fluxes. |
|
Calculation of upwind fluxes. |
|
Base class for time integrators. |
|
Class for time integrator of Taylor-series time integration scheme. |
|
Postprocessor for monopole source simulation results. |
Attributes
- class edg_acoustics.AcousticsSimulation(rho0, c0, Nx, mesh, BC_list, node_tolerance=NODETOL)[source]
Acoustics simulation data structure for running a DG acoustics simulation.
AcousticsSimulationcontains the domain discretization and sets up the DG finite element discretisation for the solution to the acoustic wave propagation.- Parameters:
rho0 (float) – the density of the medium in which the acoustic wave propagates.
c0 (float) – the speed of sound in the medium in which the acoustic wave propagates.
Nx (int) – the polynomial degree of the approximating DG finite element space used to solve the acoustic wave propagation problem.
mesh (edg_acoustics.Mesh) – the mesh object containing the mesh information for the domain discretisation.
BC_list (dict[str, int]) – a dictionary containing the definition of the boundary conditions that are present in the mesh. BC_list must contain the same keys as
edg_acoustics.Mesh.BC_triangles, i.e., all boundary conditions in the mesh must have an associated boundary condition definition.node_tolerance (float) – the tolerance used to check if a node belongs to a facet or not. <default>:
edg_acoustics.acoustics_simulation.NODETOL
- Raises:
ValueError – If BC_list[‘my_label’] is not present in the mesh, an error is raised. If a label is present in the mesh but not in BC_list, an error is raised.
- BC_list
a dictionary containing the definition of the boundary conditions that are present in the mesh. BC_list must contain the same keys as
edg_acoustics.Mesh.BC_triangles, i.e., all boundary conditions in the mesh must have an associated boundary condition definition.
- BCnode
List of boundary map nodes, each element being a dictionary with keys (values) [‘label’ (int), ‘map’ (numpy.ndarray), ‘vmap’ (numpy.ndarray)].
- Dr
the reference element differentiation matrices containing the discrete representation of \(\frac{\partial}{\partial r}\), in the rst reference element coordinate system.
- Type:
- Ds
the reference element differentiation matrices containing the discrete representation of \(\frac{\partial}{\partial s}\), in the rst reference element coordinate system.
- Type:
- Dt
the reference element differentiation matrices containing the discrete representation of \(\frac{\partial}{\partial t}\), in the rst reference element coordinate system.
- Type:
- dim
the geometric dimension of the space where the acoustic problem is solved. Always set to 3.
- Type:
- dtscale
the time step scale based on the mesh size measure, is set to the minimum diameter of the inscribed spheres in all elements.
- Type:
- Fmask
[4, Nfp]array containing indices of nodes per surface of the tetrahedron.- Type:
- Fscale
[4*Nfp, N_tets]ratio of surface to volume Jacobian of facial node.- Type:
- J
[Np, N_tets]The determinant of the Jacobian matrix for the coordinate transformation, at the collocation nodes.- Type:
- lift
[Np, 4*Nfp]an array containing the product of inverse of the mass matrix (3D) with the face-mass matrices (2D).- Type:
- mesh
the mesh object containing the mesh information for the domain discretisation.
- Type:
- Nx
the polynomial degree of the approximating DG finite element space used to solve the acoustic wave propagation problem.
- Type:
- node_tolerance
the tolerance used to check if a node belongs to a facet or not. <default>:
edg_acoustics.acoustics_simulation.NODETOL.- Type:
- rst
the reference element coordinates \((r, s, t)\) of the collocation points. Physical coordinates
xyzare obtained by mapping for each element therstcoordinates of the reference element into the physical domain. rst[0] contains the r-coordinates, rst[1] contains the s-coordinates, rst[2] contains the t-coordinates.- Type:
- rst_xyz
[3, 3, Np, N_tets]The derivative of the local coordinates \(R = (r, s, t)\) with respect to the physical coordinates \(X = (x, y, z)\), i.e., \(\frac{\partial R}{\partial X}\), at the collocation nodes. Specifically:rst_xyz[0, 0]: rx (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(r\) with respect to the physical coordinates \(x\), i.e., \(\frac{\partial r}{\partial x}\), at the collocation nodes.rst_xyz[1, 0]: sx (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(s\) with respect to the physical coordinates \(x\), i.e., \(\frac{\partial s}{\partial x}\), at the collocation nodes.rst_xyz[2, 0]: tx (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(t\) with respect to the physical coordinates \(x\), i.e., \(\frac{\partial t}{\partial x}\), at the collocation nodes.rst_xyz[0, 1]: ry (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(r\) with respect to the physical coordinates \(y\), i.e., \(\frac{\partial r}{\partial y}\), at the collocation nodes.rst_xyz[1, 1]: sy (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(s\) with respect to the physical coordinates \(y\), i.e., \(\frac{\partial s}{\partial y}\), at the collocation nodes.rst_xyz[2, 1]: ty (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(t\) with respect to the physical coordinates \(y\), i.e., \(\frac{\partial t}{\partial y}\), at the collocation nodes.rst_xyz[0, 2]: rz (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(r\) with respect to the physical coordinates \(z\), i.e., \(\frac{\partial r}{\partial z}\), at the collocation nodes.rst_xyz[1, 2]: sz (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(s\) with respect to the physical coordinates \(z\), i.e., \(\frac{\partial s}{\partial z}\), at the collocation nodes.rst_xyz[2, 2]: tz (numpy.ndarray):
[Np, N_tets]The derivative of the local coordinates \(t\) with respect to the physical coordinates \(z\), i.e., \(\frac{\partial t}{\partial z}\), at the collocation nodes.
- Type:
- sJ
[4*Nfp, N_tets]The determinant of the surface Jacobian matrix at each of the collocation nodes, for each of the 4 faces of theN_tetselements.- Type:
- V
[Np, Np]vandermonde matrix of the orthonormal basis functions on the reference simplex element. Polynomial basis can exactly represent polynomials up to degreeNx. Consider the set ofNp3D nodes, with the coordinates of each node \(i\) equal to \((r_{i}, s_{i}, t_{i})\), inrst, and the set of \(m\) orthonormal basis functions, the vandermonde matrix will be \(V_{i,j} = f_{j}(r_{i}, s_{i}, t_{i})\).- Type:
- xyz
[3, Np, N_tets]the physical space coordinates \((x, y, z)\) of the collocation points of each element of the mesh. xyz[0] contains the x-coordinates, xyz[1] contains the y-coordinates, xyz[2] contains the z-coordinates.- Type:
- n_xyz
[3, 4*Np, N_tets]the outwards normals \(\vec{n}= (n_x, n_y, n_z)\) at each collocation point on the element faces. Specifically:n_xyz[0, :] (numpy.ndarray): nx
[4*Nfp, N_tets]The \(x\)-component of the outward normal \(\vec{n}\) at each of theNfpnodes on each of the 4 facets of each of theN_tetselements.n_xyz[1, :] (numpy.ndarray): ny
[4*Nfp, N_tets]The \(y\)-component of the outward normal \(\vec{n}\) at each of theNfpnodes on each of the 4 facets of each of theN_tetselements.n_xyz[2, :] (numpy.ndarray): nz
[4*Nfp, N_tets]The \(z\)-component of the outward normal \(\vec{n}\) at each of theNfpnodes on each of the 4 facets of each of theN_tetselements.
- Type:
- vmapM
[4*Nfp*N_tets, 1]an array containing the global indices for interior values.- Type:
- vmapP
[4*Nfp*N_tets, 1]an array containing the global indices for exterior values.- Type:
- init_local_system()[source]
Call the methods to initialise the local system and compute all the attributes of the AcousticsSimulation object.
- static compute_Np(Nx)[source]
Computes the number of collocation nodes for basis of polynomial degree
Nx.
- static compute_Nfp(Nx)[source]
Computes the number of collocation nodes lying on a face of the elements for basis of polynomial degree
Nx.
- static compute_Nx_from_Np(Np)[source]
Computes the polynomial degree
Nxof basis from the number of collocation pointsNp.
- static check_BC_list(BC_list, mesh)[source]
Check if BC_list is compatible with mesh. Given a mesh with a set of boundary conditions specified in mesh.BC_triangles, check if the list of boundary conditions specification, BC_list, is compatible. By compatible we mean that all boundary conditions (keys) in BC_list exist in mesh.BC_labels, and vice-versa.
- static compute_collocation_nodes(EToV, vertices, Nx, dim=3)[source]
Compute reference element \((r, s, t)\) coordinates of collocation points
rstand the physical domainxyzcoordinates for each element.- Parameters:
EToV (numpy.ndarray) – See
edg_acoustics.Mesh.EToV.vertices (numpy.ndarray) – see
edg_acoustics.Mesh.vertices.
- Returns:
- static compute_van_der_monde_matrix(Nx, rst, dim=3)[source]
- Compute vandermonde matrix of the orthonormal basis functions on the reference simplex element. Polynomial basis
can exactly represent polynomials up to degree
Nx. Consider the set ofNp3D nodes, with the coordinates of each node \(i\) equal to \((r_{i}, s_{i}, t_{i})\), inrst, and the set of \(m\) orthonormal basis functions, the vandermonde matrix will be \(V_{i,j} = f_{j}(r_{i}, s_{i}, t_{i})\).
- Parameters:
- Returns:
V (numpy.ndarray) – see
V.- Return type:
- static compute_derivative_matrix(Nx, rst, dim=3)[source]
Compute the derivative matrix on the reference element.
- static compute_Fmask(rst, node_tol)[source]
Find all the :attr:
Nfpface nodes that lie on each surface.- Parameters:
rst (numpy.ndarray) – see
rst.node_tol (float) – see
node_tolerance.
- Returns:
Fmask (numpy.ndarray) – see
Fmask.
- static compute_lift(V, rst, Fmask)[source]
Compute the lift matrix.
- Parameters:
V (numpy.ndarray) – see
V.rst (numpy.ndarray) – see
rst.Fmask (numpy.ndarray) – see
Fmask.
- Returns:
lift (numpy.ndarray) – see
lift.
- static geometric_factors_3d(xyz, Dr, Ds, Dt)[source]
Compute the metric elements for the local mappings of the elements.
- Parameters:
xyz (numpy.ndarray) – see
xyz.Dr (numpy.ndarray) – see
Dr.Ds (numpy.ndarray) – see
Ds.Dt (numpy.ndarray) – see
Dt.
- Returns:
- static normals_3d(xyz, rst_xyz, J, Fmask)[source]
Compute outward pointing normals at element’s faces as well as surface Jacobians.
- Parameters:
xyz (numpy.ndarray) – see
xyz.rst_xyz (numpy.ndarray) – see
rst_xyz.J (numpy.ndarray) – see
J.Fmask (numpy.ndarray) – see
Fmask.
- Returns:
- static build_maps_3d(xyz, EToE, EToF, Fmask, node_tol)[source]
Find connectivity for nodes given per surface in all elements
- Parameters:
xyz (numpy.ndarray) – see
xyz.EToE (numpy.ndarray) – see
edg_acoustics.Mesh.EToE.EToF (numpy.ndarray) – see
edg_acoustics.Mesh.EToF.Fmask (numpy.ndarray) – see
Fmask.node_tol (float) – see
node_tolerance.
- Returns:
- static ismember_col(a, b)[source]
find the indices of the columns of a that are in b
- Parameters:
a (numpy.ndarray) – matrix a to be checked
b (numpy.ndarray) – matrix b to be checked against
- Returns:
indices (numpy.ndarray) – boolean indices of the columns of a that are in b
- static build_BCmaps_3d(BC_list, EToV, vmapM, BC_triangles, Nx)[source]
Build specialized nodal maps for various types of boundary conditions,specified in BC_list
- Parameters:
EToV (numpy.ndarray) – see
edg_acoustics.Mesh.EToV.vmapM (numpy.ndarray) – see
vmapM.BC_triangles (dict[str, numpy.ndarray]) – see
edg_acoustics.mesh.Mesh.BC_triangles.
- Returns:
BCnode (list[dict]) – see
BCnode.
- static diameter_3d(Fscale)[source]
Compute the minimum diameter of the inscribed spheres in all elements.
- Parameters:
Fscale (numpy.ndarray) – see
Fscale.- Returns:
diameter (float) – minimum diameter of the inscribed spheres in all elements.
- grad_3d(U, axis)[source]
Compute partial derivative dU/dx, dU/dy, dU/dz, or gradient dU/dx + dU/dy + dU/dz
- Parameters:
U (numpy.ndarray) –
[Np, N_tets]the acoustic variables that needs to be differentiated.axis (str) – the axis to be differentiated w.r.t, e.g. ‘x’, ‘y’, ‘z’, ‘xyz’
- Returns:
dUdx (numpy.ndarray) –
[Np, N_tets]derivatives \(\frac{\partial U}{\partial x}\) at every nodal point, if axis is ‘x’.dUdy (numpy.ndarray) –
[Np, N_tets]derivatives \(\frac{\partial U}{\partial y}\) at every nodal point, if axis is ‘y’.dUdz (numpy.ndarray) –
[Np, N_tets]derivatives \(\frac{\partial U}{\partial z}\) at every nodal point, if axis is ‘z’.Tuple of gradient (numpy.ndarray) –
[Np, N_tets]gradient (\(\frac{\partial U}{\partial x}, \frac{\partial U}{\partial y}, \frac{\partial U}{\partial z}\)), if axis is ‘xyz’.
- static locate_simplex(node_coordinates, EToV, rec, methodLocate='scipy')[source]
Locate the simplices containing the sample points.
- Parameters:
node_coordinates (numpy.ndarray) – see
edg_acoustics.Mesh.vertices.EToV (numpy.ndarray) – see
edg_acoustics.Mesh.EToV.rec (numpy.ndarray) – An (N_rec x 3) array containing the (x, y, z) coordinates of N_rec microphone locations.
methodLocate (str) – search method to locate the simplices containing the sample points. Available methods are ‘scipy’ and ‘brute_force’. brutal force approach, adopted from https://stackoverflow.com/questions/25179693/how-to-check-whether-the-point-is-in-the-tetrahedron-or-not/60745339#60745339
- Returns:
nodeindex (numpy.ndarray) – indices of simplices containing the N_rec microphone points
- sample3D(methodLocate='scipy')[source]
Compute interpolation weights required to interpolate the nodal data to the sample (i.e., microphone location)
- Parameters:
methodLocate (str) – search method to locate the simplices containing the sample points. Available methods are ‘scipy’ and ‘brute_force’. brutal force approach, adopted from https://stackoverflow.com/questions/25179693/how-to-check-whether-the-point-is-in-the-tetrahedron-or-not/60745339#60745339
- Returns:
sampleWeight (numpy.ndarray) –
[N_rec, Np]interpolation weights required to interpolate the nodal data to the sample (i.e., microphone location).nodeindex (numpy.ndarray) –
[N_rec, ]index of simplices that contatin the sample (microphone) points.
- init_IC(IC)[source]
setup the initial condition and save it to the
AcousticsSimulationclass.- Parameters:
IC (edg_acoustics.InitialCondition) – the initial condition object.
- init_BC(BC)[source]
load the boundary condition and save it to the
AcousticsSimulationclass.- Parameters:
BC (edg_acoustics.BoundaryCondition) – the boundary condition object.
- init_rec(rec, methodLocate='scipy')[source]
load the receiver locations and save it to the
AcousticsSimulationclass. Then compute the interpolation weights required to interpolate the nodal data to the sample (i.e., microphone location).- Parameters:
rec (numpy.ndarray) – An
[3, N_rec]array containing the (x, y, z) coordinates of N_rec microphone locations.methodLocate (str) – search method to locate the simplices containing the sample points. Available methods are ‘scipy’ and ‘brute_force’. brutal force approach, adopted from https://stackoverflow.com/questions/25179693/how-to-check-whether-the-point-is-in-the-tetrahedron-or-not/60745339
- init_Flux(Flux)[source]
load the interior flux calculation and save it to the
AcousticsSimulationclass.- Parameters:
Flux (edg_acoustics.Flux) – the flux object.
- init_TimeIntegrator(time_integrator)[source]
load the time integrator to be used to and save it to the
AcousticsSimulationclass.- Parameters:
time_integrator (edg_acoustics.TimeIntegrator)
- RHS_operator(P, Vx, Vy, Vz, BCvar)[source]
Compute the right-hand side of the linear acoustic equations with the DG discretization.
- Parameters:
P (numpy.ndarray) – Pressure field.
Vx (numpy.ndarray) – Velocity field in the x-direction.
Vy (numpy.ndarray) – Velocity field in the y-direction.
Vz (numpy.ndarray) – Velocity field in the z-direction.
BCvar (list[dict]) – see
edg_acoustics.AbsorbBC.BCvar.
- Returns:
RHS_P (numpy.ndarray) – Right-hand side values for pressure.
RHS_Vx (numpy.ndarray) – Right-hand side values for velocity in the x-direction.
RHS_Vy (numpy.ndarray) – Right-hand side values for velocity in the y-direction.
RHS_Vz (numpy.ndarray) – Right-hand side values for velocity in the z-direction.
BCvar (list[dict]) – updated boundary condition variables.
- time_integration(**kwargs)[source]
Perform time integration for the acoustics simulation. Additional keyword arguments are optional and can vary:
- Parameters:
n_time_steps (int) – number of time steps to be performed.
total_time (float) – total simulation time to be performed, determines the number of time steps given the current time step.
delta_step (int) – print solution every delta_step time steps.
save_step (int) – save solution every save_step time steps.
format (str) – the format of the file to save the results. Can be either ‘mat’ or ‘npy’. The default format is ‘mat’.
- Returns:
prec (numpy.ndarray) – Pressure field at the microphone locations.
- class edg_acoustics.Mesh(filename, BC_labels, Nx=Nx_default, **kwargs)[source]
Mesh data structure generated from common mesh file formats or geometry files of GMSH.
If the mesh file is in .msh format, the mesh data is read and stored in the object. If the mesh file is in .geo format, the mesh is generated iteratively with a characteristic length determined by the bisection method such that the PPW is between 8 and 10.
Data structure containing mesh definition. Mesh data is obtained from data stored in common mesh file format by Gmsh format (versions 2.2, 4.0, and 4.1
.msh). Since mesh reading relies on meshio, other mesh generators can be made available in the future (e.g., DOLFIN XML (.xml), Netgen (.vol,.vol.gz).Meshdefines the domain discretisation and is used to define the finite element approximation of the solution to the acoustic wave propagation.- Parameters:
- Raises:
ValueError – If BC_labels[‘my_label’] is not present in the mesh, an error is raised. If a label is present in the mesh but not in BC_labels, an error is raised.
- BC_labels
a dictionary containing the human readable label of each boundary condition and the associated lable number used in the mesh generator. BC_labels[‘my_label’] returns the label number of the label named ‘my_label’. If BC_labels[‘my_label’] is not present in the mesh, an error is raised. If a label is present in the mesh but not in BC_labels, an error is raised.
- BC_triangles
a dictionary containing the list of triangles that have a certain boundary condition. self.BC_triangles[‘BC_label’] is a numpy.ndarray with the nodes of each triangle where boundary condition of type ‘BC_label’ is to be implemented. The nodes defining each triangle in the numpy.ndarray are stored per row.
- Type:
- EToE
an
[4, N_tets]array containing the information of which elements are neighbors of an element, i.e., EToE[j, i] returns the index of the jth neighbor of element i. The definition of j-th neighbor follows the mesh generator’s convention.- Type:
- EToF
an
[4 x N_tets]array containing the information of which face is shared between the element and its neighbor, i.e., EToF[j, i] returns the face index of the j-th neighbor of element i. Face indices follow the same convention as neighbor indices.- Type:
- EToV
An
[4 x self.N_tets]array containing the 4 indices of the vertices of theN_tetstetrahedra that make up the mesh.- Type:
- filename
the file (pathlike) to read geometry or mesh data from. Current supported format includes Gmsh (format versions 2.2, 4.0, and 4.1, .msh) or .geo files.
Since mesh reading relies on meshio, the following formats can be made available in the future: Abaqus (.inp), ANSYS msh (.msh), AVS-UCD (.avs), CGNS (.cgns), DOLFIN XML (.xml), Exodus (.e, .exo), FLAC3D (.f3grid), H5M (.h5m), Kratos/MDPA (.mdpa), Medit (.mesh, .meshb), MED/Salome (.med), Nastran (bulk data, .bdf, .fem, .nas), Netgen (.vol, .vol.gz), Neuroglancer precomputed format, Gmsh (format versions 2.2, 4.0, and 4.1, .msh), OBJ (.obj), OFF (.off), PERMAS (.post, .post.gz, .dato, .dato.gz), PLY (.ply), STL (.stl), Tecplot .dat, TetGen .node/.ele, SVG (2D output only) (.svg), SU2 (.su2), UGRID (.ugrid), VTK (.vtk), VTU (.vtu), WKT (TIN) (.wkt), XDMF (.xdmf, .xmf).
- Type:
- N_BC_triangles
the number of triangles on the boundary of the domain associated to each boundary label in self.BC_labels. For example self.N_BC_triangles[‘my_label’] returns the number of boundary triangles associated to lable ‘my_label’. ‘my_label’ must be a key of self.BC_labels.
- vertices
An
[3 x self.N_vertices]array containing the 3 coordinates of theN_verticesvertices that make up the mesh.- Type:
Example
An element of this class can be initialized in the following way
>>> import edg_acoustics >>> BC_labels = {'slip': 11, 'impedance1': 13, 'impedance2': 14, 'impedance3': 15} >>> filename = "../data/tests/mesh/CoarseMesh.msh" >>> mesh = edg_acoustics.Mesh(filename, BC_labels) >>> mesh.N_BC_triangles {'slip': 5347, 'impedance1': 400, 'impedance2': 3576, 'impedance3': 3294}
- static create_mesh_from_geo_file(geo_file, freq_max, Nx)[source]
Create a mesh from a geo file. The mesh is generated iteratively with a characteristic length by the bisection method such that the PPW is between 8 and 10.
- static compute_PPW(msh_file, wavelength, Np, volume)[source]
Compute the points per wavelength (PPW) of the mesh.
- Parameters:
msh_file (str) – The path to the .msh file.
wavelength (float) – The wavelength of the maximum frequency.
Np (int) – see
edg_acoustics.AcousticsSimulation.Np.volume (float) – The volume of the space.
- Returns:
PPW (float) – The points per wavelength of the mesh.
N_tetra (int) – The number of tetrahedra in the mesh.
- static mesh_geo_file(geo_file, length_of_mesh)[source]
Generate a mesh with a characteristic length of length_of_mesh. Save the mesh to a .msh file, keeping the same name as the .geo file.
- init_from_mesh_file(filename, BC_labels)[source]
Initializes the
Meshobject from a mesh file. The mesh file is read and the relevant data is stored
- __eq__(other)[source]
Equality operator for the
Meshclass.- Parameters:
other (edg_acoustics.Mesh) – The mesh object to compare with.
- Returns:
bool – True if the two mesh objects are equal, False otherwise.
- static compute_element_connectivity(EToV)[source]
Computes element connectivity.
Given a mesh made up of
N_tetstetrahedra, compute the element connectivity. Element connectivity contains the information of the index of the neighbor over each of the four faces of the element and the face index of the neighbor shared. This information is returned in two arrays: EToE and EToF.- Parameters:
EToV (numpy.ndarray) – see
EToV.- Returns:
- class edg_acoustics.BoundaryCondition[source]
Bases:
abc.ABCAbstract base class for boundary conditions.
- static init_ADEvariables(BCpara, BCnode)[source]
Initiate ADE variables, normal velocity, characteristic waves (outgoing and incoming).
- Parameters:
BCnode (list[dict]) – see
edg_acoustics.AcousticsSimulation.BCnode.BCpara (list[dict]) – see
edg_acoustics.AbsorbBC.BCpara.
- Returns:
BCvar (list [dict]) – see
edg_acoustics.AbsorbBC.BCvar.
- static check_BCpara(BCnode, BCpara, freq_max)[source]
Check if BCpara is compatible with AcousticsSimulation.BCnode and satisfies physical admissibility condition.
Given an acoustics simulation data structure with a set of boundary conditions specified in acoustics_simulation.BCnode, check if the list of boundary conditions specification and parameters are compatible. By compatible we mean that all boundary conditions (keys) in BCpara exist in acoustics_simulation.BCnode, and vice-versa. Also, to satisfy the causality and reality conditions, multi-pole model parameters \(\zeta_i\) (stored in first row of numpy.array BCpara[BC_label] need to be positive. To satisfy the passivity condition, the magnitude of the reflection coefficient from the multi-pole model need to be smaller than 1, that is, \(|R(\omega)|\leq 1\), where
\[R(\omega)\approx{R}_\infty+\sum_{k=1}^{S}\frac{A_k}{\zeta_k+\mathrm{i}\omega}+ \sum_{l=1}^{T} \frac{1}{2}\Big( \frac{B_l-\mathrm{i}C_l}{\alpha_l-\mathrm{i}\beta_l+\mathrm{i}\omega}+\frac{B_l+\mathrm{i}C_l}{\alpha_l+\mathrm{i}\beta_l+\mathrm{i}\omega} \Big)\]- Parameters:
BCnode (list[dict]) – see
edg_acoustics.AcousticsSimulation.BCnode.BCpara (list[dict]) – see
edg_acoustics.AbsorbBC.BCpara.freq_max (float) – maximum resolvable frequency of the simulation. <default>:
edg_acoustics.boundary_condition.FREQ_MAX
- Raises:
AssertionError – If BCpara.[index][‘label’] is not present in the acoustics_simulation.BCnode.[index][‘label’], an error is raised. If a label is present in the acoustics_simulation.BCnode.[index][‘label’] but not in BCpara, an error is raised. If the labels in BCpara and BCnode are not the same, an error is raised.
AssertionError – If the reflection coefficient is not smaller than 1, an error is raised.
AssertionError – If the number of BC types is not the same in the BC_labels and BC_para, an error is raised.
AssertionError – If the causality and reality conditions are not met, an error is raised.
- static compute_Re(omega, paras)[source]
Computes the reflection coefficient given the passed parameter of the multi-pole model at the frequencies of omega.
- Parameters:
omega (numpy.ndarray) – angular frequency.
paras (dict) – see
edg_acoustics.boundary_condition.AbsorbBC.BCpara.
- Returns:
Re (numpy.ndarray) – reflection coefficient at the frequencies of omega.
- class edg_acoustics.AbsorbBC(BCnode, BCpara, freq_max=FREQ_MAX)[source]
Bases:
BoundaryConditionSetup absorptive boundary condition of a DG acoustics simulation for a specific scenario.
AbsorbBCis used to load the boundary condition parameters, and to initiate the ADE variables.- Parameters:
BCnode (list[dict]) – see
edg_acoustics.AcousticsSimulation.BCnode.freq_max (float) – maximum resolvable frequency of the simulation. <default>:
edg_acoustics.boundary_condition.FREQ_MAX
- BCpara
a list of boundary conditon parameters from the multi-pole model. Each element is a dictionary with keys (values) [‘label’(int),’RI’(float),’RP’(numpy.ndarray),’CP’(numpy.ndarray)]. ‘RI’ refers to the limit value of the reflection coefficient as the frequency approaches infinity, i.e., \(R_\infty\). ‘RP’ refers to real pole pairs, i.e., \(A\) (stored in 1st row), \(\zeta\) (stored in 2nd row). ‘CP’ refers to complex pole pairs, i.e., \(B\) (stored in 1st row), \(C\) (stored in 2nd row), \(\alpha\) (stored in 3rd row), \(\beta\) (stored in 4th row).
More details about the multi-pole model parameters and boundary condition can be found in reference https://doi.org/10.1121/10.0001128.
BCpara[:][‘label’] must contain the same integer elements as acoustics_simulation.BCnode[:][‘label’], i.e., all boundary conditions in the simulation must have an associated boundary condition parameters.
- class edg_acoustics.InitialCondition[source]
Bases:
abc.ABCSetup initial condition of a DG acoustics simulation for a specific scenario.
InitialConditionis used to setup initial condition.- abstract Pinit(xyz)[source]
Setup initial condition for pressure.
- Parameters:
xyz (numpy.ndarray)
- abstract VXinit(xyz)[source]
Setup initial condition for velocity in x-direction.
- Parameters:
xyz (numpy.ndarray)
- abstract VYinit(xyz)[source]
Setup initial condition for velocity in y-direction.
- Parameters:
xyz (numpy.ndarray)
- abstract VZinit(xyz)[source]
Setup initial condition for velocity in z-direction.
- Parameters:
xyz (numpy.ndarray)
- class edg_acoustics.Monopole_IC(source_xyz, frequency)[source]
Bases:
InitialConditionSetup a monopole source for a specific scenario.
Monopole_ICis used to setup monopple source initial condition.- Parameters:
xyz (numpy.ndarray) – see
edg_acoustics.AcousticsSimulation.xyz.source_xyz (numpy.ndarray) – an (3,) array containing the physical coordinates of the monopole source.
halfwidth (float) – half-bandwidth of the initial Gaussian pulse.
frequency (float)
- source_xyz
an (3,) array containing the physical coordinates of the monopole source.
- Type:
- static solve_halfwidth(frequency)[source]
Solve halfwidth of the initial Gaussian pulse, given a frequency, using linear interpolation. Avoids root choosing issue with analytical spectra of Gaussian pulse.
- Parameters:
frequency (float) – frequency of the monopole source.
- Returns:
halfwidth (float) – halfwidth of the initial Gaussian pulse.
- Pinit(xyz)[source]
Setup initial condition for pressure.
- Parameters:
xyz (numpy.ndarray)
- VXinit(xyz)[source]
Setup initial condition for velocity in x-direction.
- Parameters:
xyz (numpy.ndarray)
- VYinit(xyz)[source]
Setup initial condition for velocity in y-direction.
- Parameters:
xyz (numpy.ndarray)
- VZinit(xyz)[source]
Setup initial condition for velocity in z-direction.
- Parameters:
xyz (numpy.ndarray)
- class edg_acoustics.UpwindFlux(rho0, c0, n_xyz)[source]
Bases:
FluxCalculation of upwind fluxes. :param rho0: The reference density. :type rho0: float :param c0: The reference speed of sound. :type c0: float :param n_xyz: The array representing the normal vector of the face. :type n_xyz: numpy.ndarray
- Parameters:
rho0 (float)
c0 (float)
n_xyz (numpy.ndarray)
- rho0
see
edg_acoustics.AcousticsSimulation.rho0.- Type:
- c0
see
edg_acoustics.AcousticsSimulation.c0.- Type:
- n_xyz
see
edg_acoustics.AcousticsSimulation.n_xyz.- Type:
- cn1s
precalculated constant for the calculation of the flux of velocity in x-direction.
- Type:
- cn2s
precalculated constant for the calculation of the flux of velocity in y-direction.
- Type:
- cn3s
precalculated constant for the calculation of the flux of velocity in z-direction.
- Type:
- cn1n2
precalculated constant for the calculation of the flux of velocity.
- Type:
- cn1n3
precalculated constant for the calculation of the flux of velocity.
- Type:
- cn2n3
precalculated constant for the calculation of the flux of velocity.
- Type:
- n1rho
precalculated constant for the calculation of the flux of velocity in x-direction.
- Type:
- n2rho
precalculated constant for the calculation of the flux of velocity in y-direction.
- Type:
- n3rho
precalculated constant for the calculation of the flux of velocity in z-direction.
- Type:
- csn1rho
precalculated constant for the calculation of the pressure flux.
- Type:
- csn2rho
precalculated constant for the calculation of the pressure flux.
- Type:
- csn3rho
precalculated constant for the calculation of the pressure flux.
- Type:
- FluxP(dvx, dvy, dvz, dp)[source]
This method calculates the pressure flux using the given input arrays.
- Parameters:
dvx (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the x-direction.
dvy (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the y-direction.
dvz (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the z-direction.
dp (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in pressure.
- Returns:
numpy.ndarray – The calculated pressure flux.
- FluxVx(dvx, dvy, dvz, dp)[source]
This method calculates the flux of velocity in x-direction using the given input arrays.
- Parameters:
dvx (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the x-direction.
dvy (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the y-direction.
dvz (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the z-direction.
dp (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in pressure.
- Returns:
numpy.ndarray – The calculated flux of velocity in x-direction.
- FluxVy(dvx, dvy, dvz, dp)[source]
This method calculates the flux of velocity in y-direction using the given input arrays.
- Parameters:
dvx (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the x-direction.
dvy (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the y-direction.
dvz (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the z-direction.
dp (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in pressure.
- Returns:
numpy.ndarray – The calculated flux of velocity in y-direction.
- FluxVz(dvx, dvy, dvz, dp)[source]
This method calculates the flux of velocity in z-direction using the given input arrays.
- Parameters:
dvx (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the x-direction.
dvy (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the y-direction.
dvz (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in velocity in the z-direction.
dp (numpy.ndarray) – The array representing jump values across the faces of neighboring elements in pressure.
- Returns:
numpy.ndarray – The calculated flux of velocity in z-direction.
- class edg_acoustics.TimeIntegrator(L_operator, dtscale, CFL=CFL_Default)[source]
Bases:
abc.ABCBase class for time integrators.
- Parameters:
L_operator (Callable[[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, list]]) – the function in AcousticSimulation that enables the computation of Lq, given q = [P, Vx, Vy, Vz]
dtscale (float) – the time step scale based on the mesh size measure
CFL (float) – the CFL number
- L_operator
the function in AcousticSimulation that enables the computation of Lq, given q = [P, Vx, Vy, Vz]
- Type:
Callable[[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, list]]
- class edg_acoustics.TSI_TI(L_operator, dtscale, CFL=CFL_Default, **kwargs)[source]
Bases:
TimeIntegratorClass for time integrator of Taylor-series time integration scheme.
TSI_TIis used to evolve the pressure and velocity at the time T to time T + dt, based on the Taylor-series time integration scheme.- Parameters:
L_operator (Callable[[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, list]]) – the function in AcousticSimulation that enables the computation of Lq, given q = [P, Vx, Vy, Vz]
dtscale (float) – the time step scale based on the mesh size measure
CFL (float) – the CFL number
Nt (int) – the order of the time integration scheme.
- L_operator
the function in AcousticSimulation that enables the computation of Lq, given q = [P, Vx, Vy, Vz]
- Type:
Callable[[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, list]]
- step_dt(P, Vx, Vy, Vz, BC)[source]
Takes the pressure, velocity at the time T and evolves/updates them to time T + dt.
- Parameters:
P (numpy.ndarray) – the pressure at the time T, will be updated to the pressure at the time T + dt
Vx (numpy.ndarray) – the x-component of the velocity at the time T, will be updated to the x-component of the velocity at the time T + dt
Vy (numpy.ndarray) – the y-component of the velocity at the time T, will be updated to the y-component of the velocity at the time T + dt
Vz (numpy.ndarray) – the z-component of the velocity at the time T, will be updated to the z-component of the velocity at the time T + dt
BC (edg_acoustics.AbsorbBC) – the boundary condition object
- class edg_acoustics.Monopole_postprocessor(sim, delta_step, sampling_freq=Sampling_Freq)[source]
Postprocessor for monopole source simulation results.
Monopole_postprocessoris used to postprocess the simulation results of a monopole source, mainly to correct the source spectrum. Furthermore, it collects simulation results to be saved.- Parameters:
sim (edg_acoustics.AcousticsSimulation) – The acoustic simulation object instance of
edg_acoustics.AcousticsSimulation.delta_step (float) – Factor by which the simulation results are saved.
sampling_freq (float) – The desired sampling frequency. Default is set to Sampling_Freq = 44100 Hz.
- sim
The acoustic simulation object instance of
edg_acoustics.AcousticsSimulation.
- sampling_freq
The desired sampling frequency. Default is set to Sampling_Freq = 44100 Hz.
- Type:
- IRold
The impulse response at the receiver locations.
- Type:
- IRnew
The resampled impulse response at the receiver locations.
- Type:
- TR_original
The original transfer function at the receiver locations.
- Type:
- TR_free
The free transfer function at the receiver locations.
- Type:
- TR
The corrected transfer function at the receiver locations.
- Type:
- freqs
The frequency vector.
- Type:
- apply_resample()[source]
Resamples the impulse response to the desired sampling frequency.
- Returns:
IRnew (numpy.ndarray) – see
IRnew.