Page MenuHomec4science

helmholtz_box_scattering_problem_engine.py
No OneTemporary

File Metadata

Created
Fri, May 24, 14:18

helmholtz_box_scattering_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 rrompy.hfengines.base.scattering_problem_engine import \
ScatteringProblemEngine
__all__ = ['HelmholtzBoxScatteringProblemEngine']
class HelmholtzBoxScatteringProblemEngine(ScatteringProblemEngine):
"""
Solver for scattering problem outside a box with parametric wavenumber.
- \Delta u - omega^2 * n^2 * u = 0 in \Omega
u = 0 on \Gamma_D
\partial_nu - i k u = 0 on \Gamma_R
with exact solution a transmitted plane wave.
"""
def __init__(self, R:float, kappa:float, theta:float, n:int,
degree_threshold : int = np.inf, verbosity : int = 10):
super().__init__(degree_threshold = degree_threshold,
verbosity = verbosity)
self.omega = kappa
import mshr
scatterer = mshr.Polygon([fen.Point(-1, -.5), fen.Point(1, -.5),
fen.Point(1, .5), fen.Point(.8, .5),
fen.Point(.8, -.3), fen.Point(-.8, -.3),
fen.Point(-.8, .5), fen.Point(-1, .5),])
mesh = mshr.generate_mesh(mshr.Circle(fen.Point(0, 0), R)-scatterer, n)
self.V = fen.FunctionSpace(mesh, "P", 3)
self.DirichletBoundary = (lambda x, on_boundary:
on_boundary and (x[0]**2+x[1]**2)**.5 < .95 * R)
self.RobinBoundary = "REST"
c, s = np.cos(theta), np.sin(theta)
x, y = fen.SpatialCoordinate(mesh)[:]
u0R = - fen.cos(kappa * (c * x + s * y))
u0I = - fen.sin(kappa * (c * x + s * y))
self.DirichletDatum = [u0R, u0I]

Event Timeline