elastodynamicsx.solvers.frequencydomain

Module contents

class elastodynamicsx.solvers.frequencydomain.FrequencyDomainSolver(comm: Comm, M: Mat, C: Mat, K: Mat, b: Vec, b_update_function: Callable | None = None, **kwargs)[source]

Bases: object

Class for solving frequency domain problems.

Parameters:
  • comm – The MPI communicator

  • M – The mass matrix

  • C – The damping matrix

  • K – The stiffness matrix

  • b – The load vector

  • b_update_function – A function that updates the load vector (in-place) The function must take b,omega as parameters. e.g.: b_update_function = lambda b,omega: b[:]=omega If set to None, the call is ignored.

Keyword Arguments:

petsc_options – Options that are passed to the linear algebra backend PETSc. For available choices for the ‘petsc_options’ kwarg, see the PETSc documentation.

Example

from mpi4py import MPI

from dolfinx import mesh, fem
import ufl

from elastodynamicsx.solvers import FrequencyDomainSolver
from elastodynamicsx.pde import material, BoundaryCondition, PDE

# domain
length, height = 10, 10
Nx, Ny = 10, 10
domain = mesh.create_rectangle(MPI.COMM_WORLD, [[0,0], [length,height]], [Nx,Ny])
V      = dolfinx.fem.functionspace(domain, ("Lagrange", 1, (2,)))

# material
rho, lambda_, mu = 1, 2, 1
mat = material(V, rho, lambda_, mu)

# absorbing boundary condition
Z_N, Z_T = mat.Z_N, mat.Z_T  # P and S mechanical impedances
bcs = [ BoundaryCondition(V, 'Dashpot', (Z_N, Z_T)) ]

# gaussian source term
F0     = fem.Constant(domain, PETSc.ScalarType([1, 0]))  # polarization
R0     = 0.1  # radius
x0, y0 = length/2, height/2  # center
x      = ufl.SpatialCoordinate(domain)
gaussianBF = F0 * ufl.exp(-((x[0]-x0)**2 + (x[1]-y0)**2) /2 /R0**2) / (2 * 3.141596*R0**2)
bf         = BodyForce(V, gaussianBF)

# PDE
pde = PDE(V, materials=[mat], bodyforces=[bf], bcs=bcs)

# solve
fdsolver = FrequencyDomainSolver(V.mesh.comm, pde.M(), pde.C(), pde.K(), pde.b())
omega    = 1.0
u        = fem.Function(V, name='solution')
fdsolver.solve(omega=omega, out=u.x.petsc_vec)
default_petsc_options = {'ksp_type': 'preonly', 'pc_type': 'lu'}
solve(omega: float | ndarray, out: Vec | None = None, callbacks: List[Callable] = [], **kwargs) Vec[source]

Solve the linear problem

Parameters:
  • omega – The angular frequency (scalar or array)

  • out – The solution (displacement field) to the last solve. If None a new PETSc.Vec is created

  • callbacks – If omega is an array, list of callback functions to be called after each solve (e.g. plot, store solution, …). Ignored if omega is a scalar.

Keyword Arguments:

live_plotter – a plotter object that can refresh through a live_plotter.live_plotter_update_function(i, out) function

Returns:

out

property M: Mat

The mass matrix

property C: Mat

The damping matrix

property K: Mat

The stiffness matrix

property b: Vec

The load vector

elastodynamicsx.solvers.frequencydomain.frequencydomainsolver

class elastodynamicsx.solvers.frequencydomain.frequencydomainsolver.FrequencyDomainSolver(comm: Comm, M: Mat, C: Mat, K: Mat, b: Vec, b_update_function: Callable | None = None, **kwargs)[source]

Bases: object

Class for solving frequency domain problems.

Parameters:
  • comm – The MPI communicator

  • M – The mass matrix

  • C – The damping matrix

  • K – The stiffness matrix

  • b – The load vector

  • b_update_function – A function that updates the load vector (in-place) The function must take b,omega as parameters. e.g.: b_update_function = lambda b,omega: b[:]=omega If set to None, the call is ignored.

Keyword Arguments:

petsc_options

Options that are passed to the linear algebra backend PETSc. For available choices for the ‘petsc_options’ kwarg, see the PETSc documentation.

Example

from mpi4py import MPI

from dolfinx import mesh, fem
import ufl

from elastodynamicsx.solvers import FrequencyDomainSolver
from elastodynamicsx.pde import material, BoundaryCondition, PDE

# domain
length, height = 10, 10
Nx, Ny = 10, 10
domain = mesh.create_rectangle(MPI.COMM_WORLD, [[0,0], [length,height]], [Nx,Ny])
V      = dolfinx.fem.functionspace(domain, ("Lagrange", 1, (2,)))

# material
rho, lambda_, mu = 1, 2, 1
mat = material(V, rho, lambda_, mu)

# absorbing boundary condition
Z_N, Z_T = mat.Z_N, mat.Z_T  # P and S mechanical impedances
bcs = [ BoundaryCondition(V, 'Dashpot', (Z_N, Z_T)) ]

# gaussian source term
F0     = fem.Constant(domain, PETSc.ScalarType([1, 0]))  # polarization
R0     = 0.1  # radius
x0, y0 = length/2, height/2  # center
x      = ufl.SpatialCoordinate(domain)
gaussianBF = F0 * ufl.exp(-((x[0]-x0)**2 + (x[1]-y0)**2) /2 /R0**2) / (2 * 3.141596*R0**2)
bf         = BodyForce(V, gaussianBF)

# PDE
pde = PDE(V, materials=[mat], bodyforces=[bf], bcs=bcs)

# solve
fdsolver = FrequencyDomainSolver(V.mesh.comm, pde.M(), pde.C(), pde.K(), pde.b())
omega    = 1.0
u        = fem.Function(V, name='solution')
fdsolver.solve(omega=omega, out=u.x.petsc_vec)
default_petsc_options = {'ksp_type': 'preonly', 'pc_type': 'lu'}
solve(omega: float | ndarray, out: Vec | None = None, callbacks: List[Callable] = [], **kwargs) Vec[source]

Solve the linear problem

Parameters:
  • omega – The angular frequency (scalar or array)

  • out – The solution (displacement field) to the last solve. If None a new PETSc.Vec is created

  • callbacks – If omega is an array, list of callback functions to be called after each solve (e.g. plot, store solution, …). Ignored if omega is a scalar.

Keyword Arguments:

live_plotter – a plotter object that can refresh through a live_plotter.live_plotter_update_function(i, out) function

Returns:

out

property M: Mat

The mass matrix

property C: Mat

The damping matrix

property K: Mat

The stiffness matrix

property b: Vec

The load vector