elastodynamicsx.solvers.eigensolver
Module contents
- class elastodynamicsx.solvers.eigensolver.EigenmodesSolver(comm: Comm, M: Mat, C: Mat, K: Mat, **kwargs)[source]
Bases:
EPS
Convenience class inhereted from SLEPc.EPS, with methods and default parameters that are relevant for computing the resonances of an elastic component.
- Parameters:
comm – The MPI communicator
M – The mass matrix
C – The damping matrix. C=None means no dissipation. C!=None is not supported yet.
K – The stiffness matrix
- Keyword Arguments:
nev – The number of eigenvalues to be computed
Example
# ### # Free resonances of an elastic cube # ### from mpi4py import MPI from dolfinx import mesh, fem from elastodynamicsx.solvers import EigenmodesSolver from elastodynamicsx.pde import material, PDE domain = mesh.create_box(MPI.COMM_WORLD, [[0,0,0], [1,1,1]], [10,10,10]) V = dolfinx.fem.functionspace(domain, ("CG", 1, (3,))) rho, lambda_, mu = 1, 2, 1 mat = material(V, rho, lambda_, mu) pde = PDE(V, materials=[mat]) nev = 6 + 6 # the first 6 resonances are rigid body motion eps = EigenmodesSolver(V.mesh.comm, pde.M(), None, pde.K(), nev=nev) eps.solve() eps.plot(V) freqs = eps.getEigenfrequencies() print('First resonance frequencies:', freqs)
- getEigenmodes(which='all') List[Vec] [source]
Returns the desired modeshapes
- Parameters:
which – ‘all’, or an integer, or a list of integers, or a slice object
Example
getEigenmodes() # returns all computed eigenmodes getEigenmodes(3) # returns mode number 4 getEigenmodes([3,5]) # returns modes number 4 and 6 getEigenmodes(slice(0,None,2)) # returns even modes
- getModalBasis() ModalBasis [source]
elastodynamicsx.solvers.eigensolver.eigenmodessolver
- class elastodynamicsx.solvers.eigensolver.eigenmodessolver.EigenmodesSolver(comm: Comm, M: Mat, C: Mat, K: Mat, **kwargs)[source]
Bases:
EPS
Convenience class inhereted from SLEPc.EPS, with methods and default parameters that are relevant for computing the resonances of an elastic component.
- Parameters:
comm – The MPI communicator
M – The mass matrix
C – The damping matrix. C=None means no dissipation. C!=None is not supported yet.
K – The stiffness matrix
- Keyword Arguments:
nev – The number of eigenvalues to be computed
Example
# ### # Free resonances of an elastic cube # ### from mpi4py import MPI from dolfinx import mesh, fem from elastodynamicsx.solvers import EigenmodesSolver from elastodynamicsx.pde import material, PDE domain = mesh.create_box(MPI.COMM_WORLD, [[0,0,0], [1,1,1]], [10,10,10]) V = dolfinx.fem.functionspace(domain, ("CG", 1, (3,))) rho, lambda_, mu = 1, 2, 1 mat = material(V, rho, lambda_, mu) pde = PDE(V, materials=[mat]) nev = 6 + 6 # the first 6 resonances are rigid body motion eps = EigenmodesSolver(V.mesh.comm, pde.M(), None, pde.K(), nev=nev) eps.solve() eps.plot(V) freqs = eps.getEigenfrequencies() print('First resonance frequencies:', freqs)
- getEigenmodes(which='all') List[Vec] [source]
Returns the desired modeshapes
- Parameters:
which – ‘all’, or an integer, or a list of integers, or a slice object
Example
getEigenmodes() # returns all computed eigenmodes getEigenmodes(3) # returns mode number 4 getEigenmodes([3,5]) # returns modes number 4 and 6 getEigenmodes(slice(0,None,2)) # returns even modes
- getModalBasis() ModalBasis [source]