Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F62822610
vector_fenics_engine_base.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
Wed, May 15, 21:53
Size
6 KB
Mime Type
text/x-python
Expires
Fri, May 17, 21:53 (2 d)
Engine
blob
Format
Raw Data
Handle
17694114
Attached To
R6746 RationalROMPy
vector_fenics_engine_base.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
fenics
as
fen
import
numpy
as
np
from
matplotlib
import
pyplot
as
plt
from
.fenics_engine_base
import
FenicsEngineBase
from
rrompy.utilities.base.types
import
Np1D
,
List
,
strLst
,
Tuple
,
FigHandle
from
rrompy.utilities.base.data_structures
import
purgeList
,
getNewFilename
from
rrompy.solver.fenics
import
fenplot
__all__
=
[
'VectorFenicsEngineBase'
]
class
VectorFenicsEngineBase
(
FenicsEngineBase
):
"""Generic solver for parametric vector fenics problems."""
def
__init__
(
self
,
degree_threshold
:
int
=
np
.
inf
,
verbosity
:
int
=
10
,
timestamp
:
bool
=
True
):
super
()
.
__init__
(
degree_threshold
=
degree_threshold
,
verbosity
=
verbosity
,
timestamp
=
timestamp
)
self
.
V
=
fen
.
VectorFunctionSpace
(
fen
.
UnitSquareMesh
(
10
,
10
),
"P"
,
1
)
def
plot
(
self
,
u
:
Np1D
,
warping
:
List
[
callable
]
=
None
,
is_state
:
bool
=
False
,
name
:
str
=
"u"
,
save
:
str
=
None
,
what
:
strLst
=
'all'
,
forceNewFile
:
bool
=
True
,
saveFormat
:
str
=
"eps"
,
saveDPI
:
int
=
100
,
show
:
bool
=
True
,
colorMap
:
str
=
"jet"
,
fenplotArgs
:
dict
=
{},
**
figspecs
)
->
Tuple
[
List
[
FigHandle
],
List
[
str
]]:
"""
Do some nice plots of the complex-valued function with given dofs.
Args:
u: numpy complex array with function dofs.
warping(optional): Domain warping functions.
is_state(optional): whether given u is value before multiplication
by c. Defaults to False.
name(optional): Name to be shown as title of the plots. Defaults to
'u'.
save(optional): Where to save plot(s). Defaults to None, i.e. no
saving.
what(optional): Which plots to do. If list, can contain 'ABS',
'PHASE', 'REAL', 'IMAG'. If str, same plus wildcard 'ALL'.
Defaults to 'ALL'.
forceNewFile(optional): Whether to create new output file.
saveFormat(optional): Format for saved plot(s). Defaults to "eps".
saveDPI(optional): DPI for saved plot(s). Defaults to.
show(optional): Whether to show figure. Defaults to True.
colorMap(optional): Pyplot colormap. Defaults to 'jet'.
fenplotArgs(optional): Optional arguments for fenplot.
figspecs(optional key args): Optional arguments for matplotlib
figure creation.
Returns:
List of output filenames and list of figure handles.
"""
if
not
is_state
and
not
self
.
isCEye
:
return
super
()
.
plot
(
u
,
warping
,
False
,
name
,
save
,
what
,
forceNewFile
,
saveFormat
,
saveDPI
,
show
,
fenplotArgs
,
**
figspecs
)
if
isinstance
(
what
,
(
str
,)):
if
what
.
upper
()
==
'ALL'
:
what
=
[
'ABS'
,
'PHASE'
,
'REAL'
,
'IMAG'
]
else
:
what
=
[
what
]
what
=
purgeList
(
what
,
[
'ABS'
,
'PHASE'
,
'REAL'
,
'IMAG'
],
listname
=
self
.
name
()
+
".what"
,
baselevel
=
1
)
if
len
(
what
)
==
0
:
return
if
'figsize'
not
in
figspecs
.
keys
():
figspecs
[
'figsize'
]
=
plt
.
figaspect
(
1.
/
len
(
what
))
figs
=
[
None
]
*
self
.
V
.
num_sub_spaces
()
fileOut
=
[
None
]
*
self
.
V
.
num_sub_spaces
()
for
j
in
range
(
self
.
V
.
num_sub_spaces
()):
II
=
self
.
V
.
sub
(
j
)
.
dofmap
()
.
dofs
()
Vj
=
self
.
V
.
sub
(
j
)
.
collapse
()
subplotidx
=
0
figs
[
j
]
=
plt
.
figure
(
**
figspecs
)
plt
.
set_cmap
(
colorMap
)
if
'ABS'
in
what
:
uAb
=
fen
.
Function
(
Vj
)
uAb
.
vector
()
.
set_local
(
np
.
abs
(
u
[
II
]))
subplotidx
=
subplotidx
+
1
ax
=
figs
[
j
]
.
add_subplot
(
1
,
len
(
what
),
subplotidx
)
p
=
fenplot
(
uAb
,
warping
=
warping
,
title
=
"|{}_comp{}|"
.
format
(
name
,
j
),
**
fenplotArgs
)
if
self
.
V
.
mesh
()
.
geometric_dimension
()
>
1
:
figs
[
j
]
.
colorbar
(
p
,
ax
=
ax
)
if
'PHASE'
in
what
:
uPh
=
fen
.
Function
(
Vj
)
uPh
.
vector
()
.
set_local
(
np
.
angle
(
u
[
II
]))
subplotidx
=
subplotidx
+
1
ax
=
figs
[
j
]
.
add_subplot
(
1
,
len
(
what
),
subplotidx
)
p
=
fenplot
(
uPh
,
warping
=
warping
,
title
=
"phase({}_comp{})"
.
format
(
name
,
j
),
**
fenplotArgs
)
if
self
.
V
.
mesh
()
.
geometric_dimension
()
>
1
:
figs
[
j
]
.
colorbar
(
p
,
ax
=
ax
)
if
'REAL'
in
what
:
uRe
=
fen
.
Function
(
Vj
)
uRe
.
vector
()
.
set_local
(
np
.
real
(
u
[
II
]))
subplotidx
=
subplotidx
+
1
ax
=
figs
[
j
]
.
add_subplot
(
1
,
len
(
what
),
subplotidx
)
p
=
fenplot
(
uRe
,
warping
=
warping
,
title
=
"Re({}_comp{})"
.
format
(
name
,
j
),
**
fenplotArgs
)
if
self
.
V
.
mesh
()
.
geometric_dimension
()
>
1
:
figs
[
j
]
.
colorbar
(
p
,
ax
=
ax
)
if
'IMAG'
in
what
:
uIm
=
fen
.
Function
(
Vj
)
uIm
.
vector
()
.
set_local
(
np
.
imag
(
u
[
II
]))
subplotidx
=
subplotidx
+
1
ax
=
figs
[
j
]
.
add_subplot
(
1
,
len
(
what
),
subplotidx
)
p
=
fenplot
(
uIm
,
warping
=
warping
,
title
=
"Im({}_comp{})"
.
format
(
name
,
j
),
**
fenplotArgs
)
if
self
.
V
.
mesh
()
.
geometric_dimension
()
>
1
:
figs
[
j
]
.
colorbar
(
p
,
ax
=
ax
)
plt
.
tight_layout
()
if
save
is
not
None
:
save
=
save
.
strip
()
if
forceNewFile
:
fileOut
=
getNewFilename
(
"{}_comp{}_fig_"
.
format
(
save
,
j
),
saveFormat
)
else
:
fileOut
=
"{}_comp{}_fig.{}"
.
format
(
save
,
j
,
saveFormat
)
figs
[
j
]
.
savefig
(
fileOut
,
format
=
saveFormat
,
dpi
=
saveDPI
)
if
show
:
plt
.
show
()
if
fileOut
[
0
]
is
None
:
return
figs
return
figs
,
fileOut
Event Timeline
Log In to Comment