Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F62179106
helmholtz_square_bubble_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
Sat, May 11, 10:48
Size
11 KB
Mime Type
text/x-python
Expires
Mon, May 13, 10:48 (2 d)
Engine
blob
Format
Raw Data
Handle
17615305
Attached To
R6746 RationalROMPy
helmholtz_square_bubble_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
scipy.sparse
as
scsp
import
fenics
as
fen
from
rrompy.utilities.base.types
import
Np1D
,
ScOp
,
Tuple
,
FenExpr
from
rrompy.utilities.base.fenics
import
fenZERO
from
.helmholtz_problem_engine
import
HelmholtzProblemEngine
from
rrompy.utilities.base
import
verbosityDepth
__all__
=
[
'HelmholtzSquareBubbleDomainProblemEngine'
]
class
HelmholtzSquareBubbleDomainProblemEngine
(
HelmholtzProblemEngine
):
"""
Solver for square bubble Helmholtz problems with parametric domain heigth.
- \Delta u - kappa^2 * u = f in \Omega_mu = [0,\pi] x [0,\mu\pi]
u = 0 on \Gamma_mu = \partial\Omega_mu
with exact solution square bubble times plane wave.
"""
nAs
,
nbs
=
3
,
20
rescalingExp
=
1.
def
__init__
(
self
,
kappa
:
float
,
theta
:
float
,
n
:
int
,
mu0
:
np
.
complex
=
1.
,
degree_threshold
:
int
=
np
.
inf
,
verbosity
:
int
=
10
,
timestamp
:
bool
=
True
):
super
()
.
__init__
(
degree_threshold
=
degree_threshold
,
verbosity
=
verbosity
,
timestamp
=
timestamp
)
self
.
omega
=
kappa
self
.
kappa
=
kappa
self
.
theta
=
theta
self
.
mu0
=
mu0
self
.
forcingTermMu
=
np
.
nan
mesh
=
fen
.
RectangleMesh
(
fen
.
Point
(
0
,
0
),
fen
.
Point
(
np
.
pi
,
np
.
pi
),
n
,
n
)
self
.
V
=
fen
.
FunctionSpace
(
mesh
,
"P"
,
3
)
def
buildEnergyNormForm
(
self
):
# H1
"""
Build sparse matrix (in CSR format) representative of scalar product.
"""
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"INIT"
,
"Assembling energy matrix."
,
timestamp
=
self
.
timestamp
)
mudx
=
np
.
abs
(
self
.
mu0
)
*
fen
.
dot
(
self
.
u
.
dx
(
0
),
self
.
v
.
dx
(
0
))
*
fen
.
dx
muM
=
np
.
abs
(
self
.
mu0
)
*
fen
.
dot
(
self
.
u
,
self
.
v
)
*
fen
.
dx
imudy
=
1.
/
np
.
abs
(
self
.
mu0
)
*
(
fen
.
dot
(
self
.
u
.
dx
(
1
),
self
.
v
.
dx
(
1
))
*
fen
.
dx
)
normMatFen
=
fen
.
assemble
(
mudx
+
imudy
+
muM
)
normMat
=
fen
.
as_backend_type
(
normMatFen
)
.
mat
()
self
.
energyNormMatrix
=
scsp
.
csr_matrix
(
normMat
.
getValuesCSR
()[::
-
1
],
shape
=
normMat
.
size
)
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"DEL"
,
"Done assembling energy matrix."
,
timestamp
=
self
.
timestamp
)
def
getForcingTerm
(
self
,
mu
:
complex
)
->
Tuple
[
FenExpr
,
FenExpr
]:
"""Compute forcing term."""
if
not
np
.
isclose
(
mu
,
self
.
forcingTermMu
):
if
self
.
verbosity
>=
25
:
verbosityDepth
(
"INIT"
,
(
"Assembling base expression for "
"forcing term."
),
timestamp
=
self
.
timestamp
)
pi
=
np
.
pi
c
,
s
=
np
.
cos
(
self
.
theta
),
np
.
sin
(
self
.
theta
)
x
,
y
=
fen
.
SpatialCoordinate
(
self
.
V
.
mesh
())[:]
muR
,
muI
=
np
.
real
(
mu
),
np
.
imag
(
mu
)
mu2R
,
mu2I
=
np
.
real
(
mu
**
2.
),
np
.
imag
(
mu
**
2.
)
C
=
16.
/
pi
**
4.
bR
=
C
*
(
2
*
(
x
*
(
pi
-
x
)
+
y
*
(
pi
-
y
))
+
(
self
.
kappa
*
s
)
**
2.
*
(
mu2R
-
1.
)
*
x
*
(
pi
-
x
)
*
y
*
(
pi
-
y
))
bI
=
C
*
(
2
*
self
.
kappa
*
(
c
*
(
pi
-
2
*
x
)
*
y
*
(
pi
-
y
)
+
s
*
x
*
(
pi
-
x
)
*
(
pi
-
2
*
y
))
+
(
self
.
kappa
*
s
)
**
2.
*
mu2I
*
x
*
(
pi
-
x
)
*
y
*
(
pi
-
y
))
wR
=
(
fen
.
cos
(
self
.
kappa
*
(
c
*
x
+
s
*
muR
*
y
))
*
fen
.
exp
(
self
.
kappa
*
s
*
muI
*
y
))
wI
=
(
fen
.
sin
(
self
.
kappa
*
(
c
*
x
+
s
*
muR
*
y
))
*
fen
.
exp
(
self
.
kappa
*
s
*
muI
*
y
))
self
.
forcingTerm
=
[
bR
*
wR
+
bI
*
wI
,
bI
*
wR
-
bR
*
wI
]
self
.
forcingTermMu
=
mu
if
self
.
verbosity
>=
25
:
verbosityDepth
(
"DEL"
,
"Done assembling base expression."
,
timestamp
=
self
.
timestamp
)
return
self
.
forcingTerm
def
getExtraFactorB
(
self
,
mu
:
complex
,
der
:
int
)
->
Tuple
[
FenExpr
,
FenExpr
]:
"""Compute extra expression in RHS."""
def
getPowMinusj
(
x
,
power
):
powR
=
x
**
power
powI
=
fenZERO
if
power
%
2
==
1
:
powR
,
powI
=
powI
,
powR
if
(
power
+
3
)
%
4
<
2
:
powR
,
powI
=
-
powR
,
-
powI
return
powR
,
powI
if
self
.
verbosity
>=
25
:
verbosityDepth
(
"INIT"
,
(
"Assembling auxiliary expression for "
"forcing term derivative."
),
timestamp
=
self
.
timestamp
)
from
math
import
factorial
as
fact
y
=
fen
.
SpatialCoordinate
(
self
.
V
.
mesh
())[
1
]
powR
,
powI
=
[(
self
.
kappa
*
np
.
sin
(
self
.
theta
))
**
der
*
k
\
for
k
in
getPowMinusj
(
y
,
der
)]
mu2R
,
mu2I
=
np
.
real
(
mu
**
2.
),
np
.
imag
(
mu
**
2.
)
exprR
=
mu2R
*
powR
-
mu2I
*
powI
exprI
=
mu2I
*
powR
+
mu2R
*
powI
if
der
>=
1
:
muR
,
muI
=
np
.
real
(
2.
*
mu
),
np
.
imag
(
2.
*
mu
)
powR
,
powI
=
[(
self
.
kappa
*
np
.
sin
(
self
.
theta
))
**
(
der
-
1
)
*
k
\
*
der
for
k
in
getPowMinusj
(
y
,
der
-
1
)]
exprR
+=
muR
*
powR
-
muI
*
powI
exprI
+=
muI
*
powR
+
muR
*
powI
if
der
>=
2
:
powR
,
powI
=
[(
self
.
kappa
*
np
.
sin
(
self
.
theta
))
**
(
der
-
2
)
*
k
\
*
der
*
(
der
-
1
)
for
k
in
getPowMinusj
(
y
,
der
-
2
)]
exprR
+=
powR
exprI
+=
powI
fac
=
fact
(
der
)
if
self
.
verbosity
>=
25
:
verbosityDepth
(
"DEL"
,
"Done assembling auxiliary expression."
,
timestamp
=
self
.
timestamp
)
return
[
exprR
/
fac
,
exprI
/
fac
]
def
A
(
self
,
mu
:
complex
,
der
:
int
=
0
)
->
ScOp
:
"""Assemble (derivative of) operator of linear system."""
Anull
=
self
.
checkAInBounds
(
der
)
if
Anull
is
not
None
:
return
Anull
self
.
autoSetDS
()
if
der
<=
0
and
self
.
As
[
0
]
is
None
:
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"INIT"
,
"Assembling operator term A0."
,
timestamp
=
self
.
timestamp
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
a0Re
=
fen
.
dot
(
self
.
u
.
dx
(
1
),
self
.
v
.
dx
(
1
))
*
fen
.
dx
A0Re
=
fen
.
assemble
(
a0Re
)
DirichletBC0
.
apply
(
A0Re
)
A0ReMat
=
fen
.
as_backend_type
(
A0Re
)
.
mat
()
A0Rer
,
A0Rec
,
A0Rev
=
A0ReMat
.
getValuesCSR
()
self
.
As
[
0
]
=
scsp
.
csr_matrix
((
A0Rev
,
A0Rec
,
A0Rer
),
shape
=
A0ReMat
.
size
,
dtype
=
np
.
complex
)
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"DEL"
,
"Done assembling operator term."
,
timestamp
=
self
.
timestamp
)
if
der
<=
2
and
self
.
As
[
2
]
is
None
:
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"INIT"
,
"Assembling operator term A2."
,
timestamp
=
self
.
timestamp
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
nRe
,
nIm
=
self
.
refractionIndex
n2Re
,
n2Im
=
nRe
*
nRe
-
nIm
*
nIm
,
2
*
nRe
*
nIm
k2Re
,
k2Im
=
np
.
real
(
self
.
omega
**
2
),
np
.
imag
(
self
.
omega
**
2
)
k2n2Re
=
k2Re
*
n2Re
-
k2Im
*
n2Im
k2n2Im
=
k2Re
*
n2Im
+
k2Im
*
n2Re
parsRe
=
self
.
iterReduceQuadratureDegree
(
zip
([
k2n2Re
],
[
"kappaSquaredRefractionIndexSquaredReal"
]))
parsIm
=
self
.
iterReduceQuadratureDegree
(
zip
([
k2n2Im
],
[
"kappaSquaredRefractionIndexSquaredImag"
]))
a2Re
=
(
fen
.
dot
(
self
.
u
.
dx
(
0
),
self
.
v
.
dx
(
0
))
-
k2n2Re
*
fen
.
dot
(
self
.
u
,
self
.
v
))
*
fen
.
dx
a2Im
=
-
k2n2Im
*
fen
.
dot
(
self
.
u
,
self
.
v
)
*
fen
.
dx
A2Re
=
fen
.
assemble
(
a2Re
,
form_compiler_parameters
=
parsRe
)
A2Im
=
fen
.
assemble
(
a2Im
,
form_compiler_parameters
=
parsIm
)
DirichletBC0
.
zero
(
A2Re
)
DirichletBC0
.
zero
(
A2Im
)
A2ReMat
=
fen
.
as_backend_type
(
A2Re
)
.
mat
()
A2ImMat
=
fen
.
as_backend_type
(
A2Im
)
.
mat
()
A2Rer
,
A2Rec
,
A2Rev
=
A2ReMat
.
getValuesCSR
()
A2Imr
,
A2Imc
,
A2Imv
=
A2ImMat
.
getValuesCSR
()
self
.
As
[
2
]
=
(
scsp
.
csr_matrix
((
A2Rev
,
A2Rec
,
A2Rer
),
shape
=
A2ReMat
.
size
)
+
1.j
*
scsp
.
csr_matrix
((
A2Imv
,
A2Imc
,
A2Imr
),
shape
=
A2ImMat
.
size
))
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"DEL"
,
"Done assembling operator term."
,
timestamp
=
self
.
timestamp
)
if
der
==
0
:
return
self
.
As
[
0
]
+
mu
**
2
*
self
.
As
[
2
]
if
der
==
1
:
return
2.
*
mu
*
self
.
As
[
2
]
return
self
.
As
[
2
]
def
b
(
self
,
mu
:
complex
,
der
:
int
=
0
,
homogeneized
:
bool
=
False
)
->
Np1D
:
"""Assemble (derivative of) RHS of linear system."""
bnull
=
self
.
checkbInBounds
(
der
,
homogeneized
)
if
bnull
is
not
None
:
return
bnull
if
homogeneized
and
not
np
.
isclose
(
self
.
mu0BC
,
mu
):
self
.
u0BC
=
self
.
liftDirichletData
(
mu
)
if
not
np
.
isclose
(
self
.
bsmu
,
mu
):
self
.
bsmu
=
mu
self
.
resetbs
()
b
=
self
.
bsH
[
der
]
if
homogeneized
else
self
.
bs
[
der
]
if
b
is
None
:
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"INIT"
,
(
"Assembling forcing term "
"b{}."
)
.
format
(
der
),
timestamp
=
self
.
timestamp
)
if
der
<
self
.
nbs
:
fRe
,
fIm
=
self
.
getForcingTerm
(
mu
)
cRe
,
cIm
=
self
.
getExtraFactorB
(
mu
,
der
)
cfRe
=
cRe
*
fRe
-
cIm
*
fIm
cfIm
=
cRe
*
fIm
+
cIm
*
fRe
else
:
cfRe
,
cfIm
=
fenZERO
,
fenZERO
parsRe
=
self
.
iterReduceQuadratureDegree
(
zip
([
cfRe
],
[
"forcingTermDer{}Real"
.
format
(
der
)]))
parsIm
=
self
.
iterReduceQuadratureDegree
(
zip
([
cfIm
],
[
"forcingTermDer{}Imag"
.
format
(
der
)]))
L0Re
=
fen
.
dot
(
cfRe
,
self
.
v
)
*
fen
.
dx
L0Im
=
fen
.
dot
(
cfIm
,
self
.
v
)
*
fen
.
dx
b0Re
=
fen
.
assemble
(
L0Re
,
form_compiler_parameters
=
parsRe
)
b0Im
=
fen
.
assemble
(
L0Im
,
form_compiler_parameters
=
parsIm
)
if
homogeneized
:
Ader
=
self
.
A
(
mu
,
der
)
b0Re
[:]
-=
np
.
real
(
Ader
.
dot
(
self
.
u0BC
))
b0Im
[:]
-=
np
.
imag
(
Ader
.
dot
(
self
.
u0BC
))
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
DirichletBC0
.
apply
(
b0Re
)
DirichletBC0
.
apply
(
b0Im
)
b
=
np
.
array
(
b0Re
[:]
+
1.j
*
b0Im
[:],
dtype
=
np
.
complex
)
if
homogeneized
:
self
.
bsH
[
der
]
=
b
else
:
self
.
bs
[
der
]
=
b
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"DEL"
,
"Done assembling forcing term."
,
timestamp
=
self
.
timestamp
)
return
b
Event Timeline
Log In to Comment