Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F60950308
helmholtz_square_simplified_domain_problem_engine.py
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Fri, May 3, 13:51
Size
4 KB
Mime Type
text/x-python
Expires
Sun, May 5, 13:51 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
17444319
Attached To
R6746 RationalROMPy
helmholtz_square_simplified_domain_problem_engine.py
View Options
# 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.utilities.base.types
import
paramVal
from
rrompy.solver.fenics
import
fenZERO
from
rrompy.hfengines.linear_problem.helmholtz_problem_engine
import
(
HelmholtzProblemEngine
)
from
rrompy.utilities.base
import
verbosityManager
as
vbMng
from
rrompy.solver.fenics
import
fenics2Sparse
class
HelmholtzSquareSimplifiedDomainProblemEngine
(
HelmholtzProblemEngine
):
"""
Solver for square Helmholtz problems with parametric laplacian.
- \dxx u - mu_2**2 \dyy u - mu_1**2 * u = f in \Omega_mu = [0,\pi]^2
u = 0 on \partial\Omega
"""
def
__init__
(
self
,
kappa
:
float
,
theta
:
float
,
n
:
int
,
mu0
:
paramVal
=
[
12.
**
.
5
,
1.
],
degree_threshold
:
int
=
np
.
inf
,
verbosity
:
int
=
10
,
timestamp
:
bool
=
True
):
super
()
.
__init__
(
mu0
=
mu0
,
degree_threshold
=
degree_threshold
,
verbosity
=
verbosity
,
timestamp
=
timestamp
)
self
.
_affinePoly
=
True
self
.
nAs
=
3
self
.
npar
=
2
self
.
rescalingExp
=
[
2.
,
2.
]
pi
=
np
.
pi
mesh
=
fen
.
RectangleMesh
(
fen
.
Point
(
0
,
0
),
fen
.
Point
(
pi
,
pi
),
3
*
n
,
3
*
n
)
self
.
V
=
fen
.
FunctionSpace
(
mesh
,
"P"
,
1
)
c
,
s
=
np
.
cos
(
theta
),
np
.
sin
(
theta
)
x
,
y
=
fen
.
SpatialCoordinate
(
mesh
)[:]
self
.
forcingTerm
=
[
fen
.
cos
(
kappa
*
(
c
*
x
+
s
*
y
)),
fen
.
sin
(
kappa
*
(
c
*
x
+
s
*
y
))]
def
buildA
(
self
):
"""Build terms of operator of linear system."""
if
self
.
thAs
[
0
]
is
None
:
self
.
thAs
=
self
.
getMonomialWeights
(
self
.
nAs
)
if
self
.
As
[
0
]
is
None
:
vbMng
(
self
,
"INIT"
,
"Assembling operator term A0."
,
20
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
a0Re
=
fen
.
dot
(
self
.
u
.
dx
(
0
),
self
.
v
.
dx
(
0
))
*
fen
.
dx
self
.
As
[
0
]
=
fenics2Sparse
(
a0Re
,
{},
DirichletBC0
,
1
)
vbMng
(
self
,
"DEL"
,
"Done assembling operator term."
,
20
)
if
self
.
As
[
1
]
is
None
:
vbMng
(
self
,
"INIT"
,
"Assembling operator term A1."
,
20
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
nRe
,
nIm
=
self
.
refractionIndex
n2Re
,
n2Im
=
nRe
*
nRe
-
nIm
*
nIm
,
2
*
nRe
*
nIm
parsRe
=
self
.
iterReduceQuadratureDegree
(
zip
([
n2Re
],
[
"refractionIndexSquaredReal"
]))
parsIm
=
self
.
iterReduceQuadratureDegree
(
zip
([
n2Im
],
[
"refractionIndexSquaredImag"
]))
a1Re
=
-
n2Re
*
fen
.
dot
(
self
.
u
,
self
.
v
)
*
fen
.
dx
a1Im
=
-
n2Im
*
fen
.
dot
(
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
)
if
self
.
As
[
2
]
is
None
:
vbMng
(
self
,
"INIT"
,
"Assembling operator term A2."
,
20
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
a2Re
=
fen
.
dot
(
self
.
u
.
dx
(
1
),
self
.
v
.
dx
(
1
))
*
fen
.
dx
self
.
As
[
2
]
=
fenics2Sparse
(
a2Re
,
{},
DirichletBC0
,
0
)
vbMng
(
self
,
"DEL"
,
"Done assembling operator term."
,
20
)
Event Timeline
Log In to Comment