Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F61006792
membrane_fracture_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, 22:54
Size
7 KB
Mime Type
text/x-python
Expires
Sun, May 5, 22:54 (2 d)
Engine
blob
Format
Raw Data
Handle
17452336
Attached To
R6746 RationalROMPy
membrane_fracture_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
import
mshr
,
ufl
from
rrompy.utilities.base.types
import
paramVal
from
rrompy.solver.fenics
import
fenZERO
,
fenONE
from
rrompy.hfengines.linear_problem.helmholtz_problem_engine
import
(
HelmholtzProblemEngine
)
from
rrompy.utilities.base
import
verbosityManager
as
vbMng
from
rrompy.solver.fenics
import
fenics2Sparse
__all__
=
[
'MembraneFractureEngine'
]
class
MembraneFractureEngine
(
HelmholtzProblemEngine
):
def
__init__
(
self
,
mu0
:
paramVal
=
[
20.
**
.
5
,
.
6
],
H
:
float
=
1.
,
L
:
float
=
.
75
,
delta
:
float
=
.
05
,
n
:
int
=
50
,
degree_threshold
:
int
=
np
.
inf
,
homogeneized
:
bool
=
False
,
verbosity
:
int
=
10
,
timestamp
:
bool
=
True
):
super
()
.
__init__
(
mu0
=
mu0
,
degree_threshold
=
degree_threshold
,
homogeneized
=
homogeneized
,
verbosity
=
verbosity
,
timestamp
=
timestamp
)
self
.
_affinePoly
=
False
self
.
nAs
=
5
self
.
npar
=
2
self
.
H
=
H
self
.
rescalingExp
=
[
2.
,
1.
]
domain
=
(
mshr
.
Rectangle
(
fen
.
Point
(
0.
,
-
H
/
2.
),
fen
.
Point
(
2.
*
L
+
delta
,
H
/
2.
))
-
mshr
.
Rectangle
(
fen
.
Point
(
L
,
0.
),
fen
.
Point
(
L
+
delta
,
H
/
2.
)))
mesh
=
mshr
.
generate_mesh
(
domain
,
n
)
self
.
V
=
fen
.
FunctionSpace
(
mesh
,
"P"
,
1
)
self
.
NeumannBoundary
=
lambda
x
,
on_b
:
(
on_b
and
x
[
1
]
>=
-
H
/
4.
and
x
[
0
]
>=
L
and
x
[
0
]
<=
L
+
delta
)
self
.
DirichletBoundary
=
"REST"
x
,
y
=
fen
.
SpatialCoordinate
(
mesh
)[:]
self
.
_belowIndicator
=
ufl
.
conditional
(
ufl
.
le
(
y
,
0.
),
fenONE
,
fenZERO
)
self
.
_aboveIndicator
=
fenONE
-
self
.
_belowIndicator
self
.
DirichletDatum
=
[
fen
.
exp
(
-
10.
*
(
H
/
2.
+
y
)
/
H
-
.
5
*
((
x
-
.
6
*
L
)
/
(
.
1
*
L
))
**
2.
)
*
self
.
_belowIndicator
,
fenZERO
]
def
buildA
(
self
):
"""Build terms of operator of linear system."""
if
self
.
As
[
0
]
is
None
:
vbMng
(
self
,
"INIT"
,
"Assembling operator term A0."
,
20
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
a0Re
=
fenZERO
*
fen
.
dot
(
self
.
u
,
self
.
v
)
*
fen
.
dx
self
.
As
[
0
]
=
fenics2Sparse
(
a0Re
,
{},
DirichletBC0
,
1
)
self
.
thAs
[
0
]
=
self
.
getMonomialSingleWeight
([
0
,
0
])
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
)
a1Re
=
(
self
.
H
**
3
/
4.
*
self
.
_aboveIndicator
*
fen
.
dot
(
self
.
u
.
dx
(
1
),
self
.
v
.
dx
(
1
))
*
fen
.
dx
)
self
.
As
[
1
]
=
fenics2Sparse
(
a1Re
,
{},
DirichletBC0
,
0
)
self
.
thAs
[
1
]
=
[(
'x'
,
'()'
,
1
,
'*'
,
-
2.
,
'+'
,
self
.
H
),
(
0.
,),
(
-
2.
,),
None
]
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
=
self
.
H
**
2
*
fen
.
dot
(
self
.
u
.
dx
(
0
),
self
.
v
.
dx
(
0
))
*
fen
.
dx
self
.
As
[
2
]
=
fenics2Sparse
(
a2Re
,
{},
DirichletBC0
,
0
)
self
.
thAs
[
2
]
=
[(
'x'
,
'()'
,
1
,
'*'
,
-
1.
,
'+'
,
self
.
H
,
'*'
,
(
'x'
,
'()'
,
1
),
'**'
,
2.
),
(
0.
,),
(
'x'
,
'()'
,
1
,
'**'
,
2.
,
'*'
,
4.
,
'-'
,
(
'x'
,
'()'
,
1
,
'*'
,
6.
*
self
.
H
),
'+'
,
2.
*
self
.
H
**
2.
,
'*'
,
(
'x'
,
'()'
,
1
)),
(
0.
,),
(
0.
,),
(
'x'
,
'()'
,
1
,
'**'
,
2.
,
'*'
,
6.
,
'-'
,
(
'x'
,
'()'
,
1
,
'*'
,
6.
*
self
.
H
),
'+'
,
self
.
H
**
2.
),
(
0.
,),
(
0.
,),
(
0.
,),
(
'x'
,
'()'
,
1
,
'*'
,
4.
,
'-'
,
2.
*
self
.
H
),
(
0.
,),
(
0.
,),
(
0.
,),
(
0.
,),
(
1.
,),
None
]
vbMng
(
self
,
"DEL"
,
"Done assembling operator term."
,
20
)
if
self
.
As
[
3
]
is
None
:
vbMng
(
self
,
"INIT"
,
"Assembling operator term A3."
,
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"
]))
a3Re
=
-
n2Re
*
fen
.
dot
(
self
.
u
,
self
.
v
)
*
fen
.
dx
a3Im
=
-
n2Im
*
fen
.
dot
(
self
.
u
,
self
.
v
)
*
fen
.
dx
self
.
As
[
3
]
=
(
fenics2Sparse
(
a3Re
,
parsRe
,
DirichletBC0
,
0
)
+
1.j
*
fenics2Sparse
(
a3Im
,
parsIm
,
DirichletBC0
,
0
))
self
.
thAs
[
3
]
=
[(
'x'
,
'()'
,
1
,
'*'
,
-
1.
,
'+'
,
self
.
H
,
'*'
,
(
'x'
,
'()'
,
1
),
'**'
,
2.
,
'*'
,
(
'x'
,
'()'
,
0
)),
(
'x'
,
'()'
,
1
,
'*'
,
-
1.
,
'+'
,
self
.
H
,
'*'
,
(
'x'
,
'()'
,
1
),
'**'
,
2.
),
(
2.
*
self
.
H
**
2.
,
'-'
,
(
'x'
,
'()'
,
1
,
'*'
,
6.
*
self
.
H
),
'+'
,
(
'x'
,
'()'
,
1
,
'**'
,
2.
,
'*'
,
4.
),
'*'
,
(
'x'
,
'()'
,
1
),
'*'
,
(
'x'
,
'()'
,
0
)),
(
0.
,),
(
'x'
,
'()'
,
1
,
'**'
,
2.
,
'*'
,
4.
,
'-'
,
(
'x'
,
'()'
,
1
,
'*'
,
6.
*
self
.
H
),
'+'
,
2.
*
self
.
H
**
2.
,
'*'
,
(
'x'
,
'()'
,
1
)),
(
'x'
,
'()'
,
1
,
'**'
,
2.
,
'*'
,
6.
,
'-'
,
(
'x'
,
'()'
,
1
,
'*'
,
6.
*
self
.
H
),
'+'
,
self
.
H
**
2.
,
'*'
,
(
'x'
,
'()'
,
0
)),
(
0.
,),
(
0.
,),
(
'x'
,
'()'
,
1
,
'**'
,
2.
,
'*'
,
6.
,
'-'
,
(
'x'
,
'()'
,
1
,
'*'
,
6.
*
self
.
H
),
'+'
,
self
.
H
**
2.
),
(
'x'
,
'()'
,
1
,
'*'
,
4.
,
'-'
,
2.
*
self
.
H
,
'*'
,
(
'x'
,
'()'
,
0
)),
(
0.
,),
(
0.
,),
(
0.
,),
(
'x'
,
'()'
,
1
,
'*'
,
4.
,
'-'
,
2.
*
self
.
H
),
(
'x'
,
'()'
,
0
),
(
0.
,),
(
0.
,),
(
0.
,),
(
0.
,),
(
1.
,),
None
]
vbMng
(
self
,
"DEL"
,
"Done assembling operator term."
,
20
)
if
self
.
As
[
4
]
is
None
:
vbMng
(
self
,
"INIT"
,
"Assembling operator term A4."
,
20
)
DirichletBC0
=
fen
.
DirichletBC
(
self
.
V
,
fenZERO
,
self
.
DirichletBoundary
)
a4Re
=
.
25
*
self
.
H
**
2
*
fen
.
dot
(
self
.
u
.
dx
(
1
),
self
.
v
.
dx
(
1
))
*
fen
.
dx
self
.
As
[
4
]
=
fenics2Sparse
(
a4Re
,
{},
DirichletBC0
,
0
)
self
.
thAs
[
4
]
=
self
.
getMonomialSingleWeight
([
0
,
2
])
vbMng
(
self
,
"DEL"
,
"Done assembling operator term."
,
20
)
Event Timeline
Log In to Comment