Manual
FEM machines
Acoustics: volume
FinEtoolsAcoustics.FEMMAcoustModule.FEMMAcoust
— TypeFEMMAcoust{ID<:IntegDomain, M} <: AbstractFEMM
Type for linear acoustics finite element modeling machine.
FinEtoolsAcoustics.FEMMAcoustModule.acousticstiffness
— Functionacousticstiffness(self::FEMMAcoust, assembler::A, geom::NodalField, P::NodalField{T}) where {T<:Number, A<:AbstractSysmatAssembler}
Compute the acoustic "stiffness" matrix.
Arguments
self
= acoustics modelassembler
= matrix assemblergeom
= geometry fieldP
= acoustic (perturbation) pressure field
Return a matrix.
The acoustic "stiffness" matrix is by convention called stiffness, however its mechanical meaning is quite different. (It has to do with kinetic energy.) It is the matrix $\mathbf{K}_a$ in this matrix ODE system for the acoustic pressure:
\[\mathbf{M}_a \mathbf{\ddot{p}} + \mathbf{K}_a \mathbf{{p}} = \mathbf{{f}}\]
The bilinear-form function bilform_diffusion
is used to compute the matrix.
acousticstiffness(self::FEMMAcoustNICE, assembler::A, geom::NodalField, P::NodalField{T}) where {T<:Number, A<:AbstractSysmatAssembler}
Compute the acoustic mass matrix.
Arguments
self
= acoustics modelassembler
= matrix assemblergeom
= geometry fieldP
= acoustic (perturbation) pressure field
Return a matrix.
FinEtoolsAcoustics.FEMMAcoustModule.acousticmass
— Functionacousticmass(self::FEMMAcoust, assembler::A,
geom::NodalField,
Pddot::NodalField{T}) where {T<:Number,
A<:AbstractSysmatAssembler}
Compute the acoustic mass matrix.
Arguments
self
= acoustics modelassembler
= matrix assemblergeom
= geometry fieldPddot
= second order rate of the acoustic (perturbation) pressure field
The acoustic "mass" matrix is by convention called mass, however its mechanical meaning is quite different. (It has to do with potential energy.) It is the matrix $\mathbf{M}_a$ in this matrix ODE system for the acoustic pressure:
\[\mathbf{M}_a \mathbf{\ddot{p}} + \mathbf{K}_a \mathbf{{p}} = \mathbf{{f}}\]
The bilinear-form function bilform_dot
is used to compute the matrix.
FinEtoolsAcoustics.FEMMAcoustModule.inspectintegpoints
— Functioninspectintegpoints(self::FEMMAcoust,
geom::NodalField{GFT},
P::NodalField{T},
temp::NodalField{FT},
felist::VecOrMat{IntT},
inspector::F,
idat,
quantity = :gradient;
context...) where {T <: Number, GFT, FT, IntT, F <: Function}
Inspect integration point quantities.
Arguments
geom
- reference geometry fieldP
- pressure fieldtemp
- temperature field (ignored)felist
- indexes of the finite elements that are to be inspected: The fes to be included are:fes[felist]
.context
- struct: see the update!() method of the material.inspector
- function with the signatureidat = inspector(idat, j, conn, x, out, loc);
whereidat
- a structure or an array that the inspector may use to maintain some state, for instance gradient,j
is the element number,conn
is the element connectivity,out
is the output of theupdate!()
method,loc
is the location of the integration point in the reference configuration.
Output
The updated inspector data is returned.
Acoustics: surface
FinEtoolsAcoustics.FEMMAcoustSurfModule.FEMMAcoustSurf
— TypeFEMMAcoustSurf{ID<:IntegDomain, M} <: AbstractFEMM
Class for linear acoustics finite element modeling machine.
FinEtoolsAcoustics.FEMMAcoustSurfModule.acousticABC
— FunctionacousticABC(self::FEMMAcoustSurf, assembler::A,
geom::NodalField,
Pdot::NodalField{T}) where {T<:Number, A<:AbstractSysmatAssembler}
Compute the acoustic ABC (Absorbing Boundary Condition) matrix.
Arguments
self
= acoustics modelassembler
= matrix assembler; must be able to assemble unsymmetric matrixgeom
= geometry fieldPdot
= rate of the acoustic (perturbation) pressure field
We assume here that the impedance of this boundary is $ho c$.
FinEtoolsAcoustics.FEMMAcoustSurfModule.acousticrobin
— Functionacousticrobin(
self::FEMMAcoustSurf,
assembler::A,
geom::NodalField,
Pdot::NodalField{T},
impedance
) where {T<:Number,A<:AbstractSysmatAssembler}
Compute the acoustic "Robin boundary condition" (damping) matrix.
Arguments
self
= acoustics modelassembler
= matrix assembler; must be able to assemble unsymmetric matrixgeom
= geometry fieldPdot
= rate of the acoustic (perturbation) pressure fieldimpedance
= acoustic impedance of the boundary
We assume here that the impedance of this boundary is $ho c$.
FinEtoolsAcoustics.FEMMAcoustSurfModule.pressure2resultantforce
— Functionpressure2resultantforce(self::FEMMAcoustSurf, assembler::A,
geom::NodalField,
P::NodalField{T},
Force::Field) where {T<:Number, A<:AbstractSysmatAssembler}
Compute the rectangular coupling matrix that transcribes given pressure on the surface into the resultant force acting on the surface.
Arguments
self
= acoustics modelassembler
= matrix assembler; must be able to assemble unsymmetric matrixgeom
= geometry fieldP
= acoustic (perturbation) pressure fieldForce
= field for the force resultant
FinEtoolsAcoustics.FEMMAcoustSurfModule.pressure2resultanttorque
— Functionpressure2resultanttorque(self::FEMMAcoustSurf, assembler::A,
geom::NodalField,
P::NodalField{T},
Torque::GeneralField, CG::FFltVec) where {T<:Number, A<:AbstractSysmatAssembler}
Compute the rectangular coupling matrix that transcribes given pressure on the surface into the resultant torque acting on the surface with respect to the CG.
Arguments
self
= acoustics modelassembler
= matrix assembler; must be able to assemble unsymmetric matrixgeom
= geometry fieldP
= acoustic (perturbation) pressure fieldTorque
= field for the torque resultant
FinEtoolsAcoustics.FEMMAcoustSurfModule.acousticcouplingpanels
— Functionacousticcouplingpanels(self::FEMMAcoustSurf, geom::NodalField, u::NodalField{T}) where {T}
Compute the acoustic pressure-structure coupling matrix.
The acoustic pressure-nodal force matrix transforms the pressure distributed along the surface to forces acting on the nodes of the finite element model. Its transpose transforms displacements (or velocities, or accelerations) into the normal component of the displacement (or velocity, or acceleration) along the surface.
Arguments
geom
=geometry fieldu
= displacement field
n
= outer normal (pointing into the acoustic medium).- The pressures along the surface are assumed constant (uniform) along each finite element –- panel. The panel pressures are assumed to be given the same numbers as the serial numbers of the finite elements in the set.
Algorithms
Acoustics
FinEtoolsAcoustics.AlgoAcoustModule.steadystate
— Functionsteadystate(modeldata::FDataDict)
Steady-state acoustics solver.
modeldata
= dictionary with string keys
"fens"
= finite element node set"regions"
= array of region dictionaries"essential_bcs"
= array of essential boundary condition dictionaries"ABCs"
= array of absorbing boundary condition dictionaries"flux_bcs"
= array of flux boundary condition dictionaries
For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains items:
"femm"
= finite element mmodel machine (mandatory);
For essential boundary conditions (optional) each dictionary would hold
"pressure"
= fixed (prescribed) pressure (scalar), or a function with signature function T = f(x) If not given, zero pressure assumed."node_list"
= list of nodes on the boundary to which the condition applies (mandatory)
For absorbing boundary conditions (optional) each dictionary may hold
"femm"
= finite element mmodel machine (mandatory).
For flux boundary conditions (optional) each dictionary would hold
"femm"
= finite element mmodel machine (mandatory);"normal_flux"
= normal component of the flux through the boundary (scalar), which is the normal derivative of the pressure.
Output
modeldata
= the dictionary is augmented with
"geom"
= the nodal field that is the geometry"P"
= the nodal field that is the computed pressure (in the general a complex-number field)
Material models for acoustics
FinEtoolsAcoustics.MatAcoustFluidModule.MatAcoustFluid
— TypeMatAcoustFluid <: AbstractMat
Type for acoustic fluid material.
FinEtoolsAcoustics.MatAcoustFluidModule.bulkmodulus
— Functionbulkmodulus(self::MatAcoustFluid)
Return the bulk modulus.
Modules
FinEtoolsAcoustics.FinEtoolsAcoustics
— ModuleFinEtools (C) 2017-2024, Petr Krysl
Finite Element tools. Julia implementation of the finite element method for continuum mechanics. Package for acoustic problems.
FinEtoolsAcoustics.FEMMAcoustModule
— ModuleFEMMAcoustModule
Module for operations on interiors of domains to construct system matrices and system vectors for linear acoustics.
FinEtoolsAcoustics.FEMMAcoustSurfModule
— ModuleFEMMAcoustSurfModule
Module for operations on boundaries of domains to construct system matrices and system vectors for linear acoustics.
FinEtoolsAcoustics.FEMMAcoustNICEModule.FEMMAcoustNICE
— TypeFEMMAcoustNICE{S<:AbstractFESet, F<:Function, M} <: AbstractFEMM
Type for linear acoustics finite element modeling machine.
FinEtoolsAcoustics.AlgoAcoustModule
— ModuleAlgoAcoustModule
Module for linear acoustics algorithms.
FinEtoolsAcoustics.MatAcoustFluidModule
— ModuleMatAcoustFluidModule
Module for acoustic-fluid material.