Page MenuHomec4science

linear_elasticity_helmholtz_problem_engine.py
No OneTemporary

File Metadata

Created
Thu, May 16, 03:02

linear_elasticity_helmholtz_problem_engine.py

# Copyright (C) 2018 by the RROMPy authors
#
# This file is part of RROMPy.
#
# RROMPy is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# RROMPy is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with RROMPy. If not, see <http://www.gnu.org/licenses/>.
#
import numpy as np
import fenics as fen
from .linear_elasticity_problem_engine import LinearElasticityProblemEngine
from rrompy.utilities.base.types import List, ScOp, paramVal
from rrompy.solver.fenics import (fenZERO, fenZEROS, fenONE, L2NormMatrix,
elasticNormMatrix, elasticDualNormMatrix)
from rrompy.utilities.base import verbosityManager as vbMng
from rrompy.utilities.poly_fitting.polynomial import (
hashDerivativeToIdx as hashD)
from rrompy.solver.fenics import fenics2Sparse
__all__ = ['LinearElasticityHelmholtzProblemEngine']
class LinearElasticityHelmholtzProblemEngine(LinearElasticityProblemEngine):
"""
Solver for generic linear elasticity Helmholtz problems with parametric
wavenumber.
- div(lambda_ * div(u) * I + 2 * mu_ * epsilon(u))
- rho_ * mu^2 * u = f in \Omega
u = u0 on \Gamma_D
\partial_nu = g1 on \Gamma_N
\partial_nu + h u = g2 on \Gamma_R
Attributes:
verbosity: Verbosity level.
BCManager: Boundary condition manager.
V: Real vector FE space.
u: Generic vector trial functions for variational form evaluation.
v: Generic vector test functions for variational form evaluation.
As: Scipy sparse array representation (in CSC format) of As.
bs: Numpy array representation of bs.
energyNormMatrix: Scipy sparse matrix representing inner product over
V.
energyNormDualMatrix: Scipy sparse matrix representing inner product
over V'.
dualityMatrix: Scipy sparse matrix representing duality V-V'.
energyNormPartialDualMatrix: Scipy sparse matrix representing dual
inner product between Riesz representers V-V.
liftedDirichletDatum: Dofs of Dirichlet datum lifting.
mu0BC: Mu value of last Dirichlet datum lifting.
degree_threshold: Threshold for ufl expression interpolation degree.
omega: Value of omega.
lambda_: Value of lambda_.
mu_: Value of mu_.
forcingTerm: Value of f.
DirichletDatum: Value of u0.
NeumannDatum: Value of g1.
RobinDatumG: Value of g2.
RobinDatumH: Value of h.
DirichletBoundary: Function handle to \Gamma_D.
NeumannBoundary: Function handle to \Gamma_N.
RobinBoundary: Function handle to \Gamma_R.
ds: Boundary measure 2-tuple (resp. for Neumann and Robin boundaries).
A0: Scipy sparse array representation (in CSC format) of A0.
A1: Scipy sparse array representation (in CSC format) of A1.
b0: Numpy array representation of b0.
dsToBeSet: Whether ds needs to be set.
"""
def __init__(self, mu0 : paramVal = [0.], degree_threshold : int = np.inf,
verbosity : int = 10, timestamp : bool = True):
super().__init__(mu0 = [mu0], degree_threshold = degree_threshold,
verbosity = verbosity, timestamp = timestamp)
self.nAs = 2
self.omega = np.abs(self.mu0(0, 0))
self.rho_ = fenONE
self.rescalingExp = [2.]
@property
def rho_(self):
"""Value of rho_."""
return self._rho_
@rho_.setter
def rho_(self, rho_):
self.resetAs()
if not isinstance(rho_, (list, tuple,)):
rho_ = [rho_, fenZERO]
self._rho_ = rho_
def buildEnergyNormForm(self): # energy + omega norm
"""
Build sparse matrix (in CSR format) representative of scalar product.
"""
vbMng(self, "INIT", "Assembling energy matrix.", 20)
self.energyNormMatrix = elasticNormMatrix(
self.V, self.lambda_[0], self.mu_[0],
np.abs(self.omega)**2 * self.rho_[0])
vbMng(self, "DEL", "Done assembling energy matrix.", 20)
def buildEnergyNormDualForm(self):
"""
Build sparse matrix (in CSR format) representative of dual scalar
product.
"""
vbMng(self, "INIT", "Assembling energy dual matrix.", 20)
self.energyNormDualMatrix = elasticDualNormMatrix(
self.V, self.lambda_[0], self.mu_[0],
np.abs(self.omega)**2 * self.rho_[0],
compressRank = self._energyDualNormCompress)
vbMng(self, "DEL", "Done assembling energy dual matrix.", 20)
def buildEnergyNormPartialDualForm(self):
"""
Build sparse matrix (in CSR format) representative of dual scalar
product without duality.
"""
vbMng(self, "INIT", "Assembling energy partial dual matrix.", 20)
self.energyNormPartialDualMatrix = elasticDualNormMatrix(
self.V, self.lambda_[0], self.mu_[0],
np.abs(self.omega)**2 * self.rho_[0],
compressRank = self._energyDualNormCompress,
duality = False)
vbMng(self, "DEL", "Done assembling energy partial dual matrix.", 20)
def A(self, mu : paramVal = [], der : List[int] = 0) -> ScOp:
"""Assemble (derivative of) operator of linear system."""
mu = self.checkParameter(mu)
if not hasattr(der, "__len__"): der = [der] * self.npar
derI = hashD(der)
self.autoSetDS()
if derI <= 0 and self.As[0] is None:
vbMng(self, "INIT", "Assembling operator term A0.", 20)
DirichletBC0 = fen.DirichletBC(self.V,
fenZEROS(self.V.mesh().topology().dim()),
self.DirichletBoundary)
lambda_Re, lambda_Im = self.lambda_
mu_Re, mu_Im = self.mu_
hRe, hIm = self.RobinDatumH
termNames = ["lambda_", "mu_", "RobinDatumH"]
parsRe = self.iterReduceQuadratureDegree(zip(
[lambda_Re, mu_Re, hRe],
[x + "Real" for x in termNames]))
parsIm = self.iterReduceQuadratureDegree(zip(
[lambda_Im, mu_Re, hIm],
[x + "Imag" for x in termNames]))
epsilon = lambda u: 0.5 * (fen.grad(u) + fen.nabla_grad(u))
sigma = lambda u, l_, m_: (
l_ * fen.div(u) * fen.Identity(u.geometric_dimension())
+ 2. * m_ * epsilon(u))
a0Re = (fen.inner(sigma(self.u, lambda_Re, mu_Re),
epsilon(self.v)) * fen.dx
+ hRe * fen.inner(self.u, self.v) * self.ds(1))
a0Im = (fen.inner(sigma(self.u, lambda_Im, mu_Im),
epsilon(self.v)) * fen.dx
+ hIm * fen.inner(self.u, self.v) * self.ds(1))
self.As[0] = (fenics2Sparse(a0Re, parsRe, DirichletBC0, 1)
+ 1.j * fenics2Sparse(a0Im, parsIm, DirichletBC0, 0))
vbMng(self, "DEL", "Done assembling operator term.", 20)
if derI <= 1 and self.As[1] is None:
vbMng(self, "INIT", "Assembling operator term A1.", 20)
DirichletBC0 = fen.DirichletBC(self.V,
fenZEROS(self.V.mesh().topology().dim()),
self.DirichletBoundary)
rho_Re, rho_Im = self.rho_
parsRe = self.iterReduceQuadratureDegree(zip([rho_Re],
["rho_Real"]))
parsIm = self.iterReduceQuadratureDegree(zip([rho_Im],
["rho_Imag"]))
a1Re = - rho_Re * fen.inner(self.u, self.v) * fen.dx
a1Im = - rho_Im * fen.inner(self.u, self.v) * fen.dx
self.As[1] = (fenics2Sparse(a1Re, parsRe, DirichletBC0, 0)
+ 1.j * fenics2Sparse(a1Im, parsIm, DirichletBC0, 0))
vbMng(self, "DEL", "Done assembling operator term.", 20)
return self._assembleA(mu, der, derI)

Event Timeline