Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F100228160
mixed_problem_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, Jan 29, 04:44
Size
6 KB
Mime Type
text/x-python
Expires
Fri, Jan 31, 04:44 (2 d)
Engine
blob
Format
Raw Data
Handle
23926739
Attached To
R6746 RationalROMPy
mixed_problem_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
numpy
as
np
import
fenics
as
fen
from
scipy.sparse
import
csr_matrix
from
matplotlib
import
pyplot
as
plt
from
rrompy.utilities.base.types
import
Np1D
,
strLst
from
.problem_engine_base
import
ProblemEngineBase
from
rrompy.utilities.base
import
purgeList
,
getNewFilename
,
verbosityDepth
__all__
=
[
'MixedProblemEngineBase'
]
class
MixedProblemEngineBase
(
ProblemEngineBase
):
"""
Generic solver for parametric mixed problems.
Attributes:
verbosity: Verbosity level.
BCManager: Boundary condition manager.
V: Real mixed FE space.
u: Generic mixed trial functions for variational form evaluation.
v: Generic mixed test functions for variational form evaluation.
As: Scipy sparse array representation (in CSC format) of As.
bs: Numpy array representation of bs.
energyNormMatrix: Scipy sparse matrix representing inner product over
V.
bsmu: Mu value of last bs evaluation.
liftDirichletDatamu: Mu value of last Dirichlet datum evaluation.
liftedDirichletDatum: Dofs of Dirichlet datum lifting.
mu0BC: Mu value of last Dirichlet datum lifting.
degree_threshold: Threshold for ufl expression interpolation degree.
"""
nAs
,
nbs
=
1
,
1
functional
=
lambda
self
,
u
:
0.
def
__init__
(
self
,
degree_threshold
:
int
=
np
.
inf
,
verbosity
:
int
=
10
):
super
()
.
__init__
(
degree_threshold
=
degree_threshold
,
verbosity
=
verbosity
)
V1
=
fen
.
FiniteElement
(
"P"
,
fen
.
triangle
,
1
)
V2
=
fen
.
FiniteElement
(
"R"
,
fen
.
triangle
,
0
)
element
=
fen
.
MixedElement
([
V1
,
V2
])
self
.
V
=
fen
.
FunctionSpace
(
fen
.
UnitSquareMesh
(
10
,
10
),
element
)
self
.
primalIndices
=
[
0
]
def
buildEnergyNormForm
(
self
):
# L2 primal
"""
Build sparse matrix (in CSR format) representative of scalar product.
"""
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"INIT"
,
"Assembling energy matrix."
,
end
=
""
)
a
=
0.
for
j
in
self
.
primalIndices
:
a
=
a
+
fen
.
dot
(
self
.
u
[
j
],
self
.
v
[
j
])
normMatFen
=
fen
.
assemble
(
a
*
fen
.
dx
)
normMat
=
fen
.
as_backend_type
(
normMatFen
)
.
mat
()
self
.
energyNormMatrix
=
csr_matrix
(
normMat
.
getValuesCSR
()[::
-
1
],
shape
=
normMat
.
size
)
if
self
.
verbosity
>=
20
:
verbosityDepth
(
"DEL"
,
" Done."
,
inline
=
True
)
def
plot
(
self
,
u
:
Np1D
,
name
:
strLst
=
"u"
,
save
:
str
=
None
,
what
:
strLst
=
'all'
,
saveFormat
:
str
=
"eps"
,
saveDPI
:
int
=
100
,
**
figspecs
):
"""
Do some nice plots of the complex-valued function with given dofs.
Args:
u: numpy complex array with function dofs.
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'.
saveFormat(optional): Format for saved plot(s). Defaults to "eps".
saveDPI(optional): DPI for saved plot(s). Defaults to 100.
figspecs(optional key args): Optional arguments for matplotlib
figure creation.
"""
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'
]
=
(
13.
*
len
(
what
)
/
4
,
3
)
if
isinstance
(
name
,
(
str
,)):
name
=
[
"{}_sub{}"
.
format
(
name
,
j
+
1
)
\
for
j
in
range
(
self
.
V
.
num_sub_spaces
())]
for
j
in
range
(
self
.
V
.
num_sub_spaces
()):
subplotcode
=
100
+
len
(
what
)
*
10
Vj
=
self
.
V
.
sub
(
j
)
dofs
=
Vj
.
dofmap
()
.
dofs
()
uj
=
u
[
dofs
]
plt
.
figure
(
**
figspecs
)
Vj
=
Vj
.
collapse
()
plt
.
jet
()
if
'ABS'
in
what
:
uAb
=
fen
.
Function
(
Vj
)
uAb
.
vector
()[:]
=
np
.
array
(
np
.
abs
(
uj
),
dtype
=
float
)
subplotcode
=
subplotcode
+
1
plt
.
subplot
(
subplotcode
)
p
=
fen
.
plot
(
uAb
,
title
=
"|{}|"
.
format
(
name
[
j
]))
plt
.
colorbar
(
p
)
if
'PHASE'
in
what
:
uPh
=
fen
.
Function
(
Vj
)
uPh
.
vector
()[:]
=
np
.
array
(
np
.
angle
(
uj
),
dtype
=
float
)
subplotcode
=
subplotcode
+
1
plt
.
subplot
(
subplotcode
)
p
=
fen
.
plot
(
uPh
,
title
=
"phase({})"
.
format
(
name
[
j
]))
plt
.
colorbar
(
p
)
if
'REAL'
in
what
:
uRe
=
fen
.
Function
(
Vj
)
uRe
.
vector
()[:]
=
np
.
array
(
np
.
real
(
uj
),
dtype
=
float
)
subplotcode
=
subplotcode
+
1
plt
.
subplot
(
subplotcode
)
p
=
fen
.
plot
(
uRe
,
title
=
"Re({})"
.
format
(
name
[
j
]))
plt
.
colorbar
(
p
)
if
'IMAG'
in
what
:
uIm
=
fen
.
Function
(
Vj
)
uIm
.
vector
()[:]
=
np
.
array
(
np
.
imag
(
uj
),
dtype
=
float
)
subplotcode
=
subplotcode
+
1
plt
.
subplot
(
subplotcode
)
p
=
fen
.
plot
(
uIm
,
title
=
"Im({})"
.
format
(
name
[
j
]))
plt
.
colorbar
(
p
)
if
save
is
not
None
:
save
=
save
.
strip
()
plt
.
savefig
(
getNewFilename
(
"{}_sub{}_fig_"
.
format
(
save
,
j
+
1
),
saveFormat
),
format
=
saveFormat
,
dpi
=
saveDPI
)
plt
.
show
()
plt
.
close
()
Event Timeline
Log In to Comment