Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F61952027
trained_model_reduced_basis.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 10, 00:42
Size
5 KB
Mime Type
text/x-python
Expires
Sun, May 12, 00:42 (2 d)
Engine
blob
Format
Raw Data
Handle
17583125
Attached To
R6746 RationalROMPy
trained_model_reduced_basis.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
from
rrompy.reduction_methods.base.trained_model
import
TrainedModel
from
rrompy.utilities.base.types
import
(
Np1D
,
ListAny
,
paramVal
,
paramList
,
sampList
)
from
rrompy.utilities.base
import
verbosityManager
as
vbMng
,
freepar
as
fp
from
rrompy.utilities.numerical
import
(
eigvalsNonlinearDense
,
marginalizePolyList
)
from
rrompy.utilities.expression
import
expressionEvaluator
from
rrompy.utilities.exception_manager
import
RROMPyException
,
RROMPyWarning
from
rrompy.parameter
import
checkParameter
,
checkParameterList
from
rrompy.sampling
import
emptySampleList
__all__
=
[
'TrainedModelReducedBasis'
]
class
TrainedModelReducedBasis
(
TrainedModel
):
"""
ROM approximant evaluation for RB approximant.
Attributes:
Data: dictionary with all that can be pickled.
"""
def
reset
(
self
):
super
()
.
reset
()
if
hasattr
(
self
,
"data"
)
and
hasattr
(
self
.
data
,
"lastSetupMu"
):
self
.
data
.
lastSetupMu
=
None
def
assembleReducedModel
(
self
,
mu
:
paramVal
):
mu
=
checkParameter
(
mu
,
self
.
data
.
npar
)
if
not
(
hasattr
(
self
.
data
,
"lastSetupMu"
)
and
self
.
data
.
lastSetupMu
==
mu
):
vbMng
(
self
,
"INIT"
,
"Assembling reduced model for mu = {}."
\
.
format
(
mu
),
17
)
muEff
=
mu
**
self
.
data
.
rescalingExp
self
.
data
.
ARBmu
,
self
.
data
.
bRBmu
=
0.
,
0.
for
thA
,
ARB
in
zip
(
self
.
data
.
thAs
,
self
.
data
.
ARBs
):
self
.
data
.
ARBmu
=
(
expressionEvaluator
(
thA
[
0
],
muEff
)
*
ARB
+
self
.
data
.
ARBmu
)
for
thb
,
bRB
in
zip
(
self
.
data
.
thbs
,
self
.
data
.
bRBs
):
self
.
data
.
bRBmu
=
(
expressionEvaluator
(
thb
[
0
],
muEff
)
*
bRB
+
self
.
data
.
bRBmu
)
vbMng
(
self
,
"DEL"
,
"Done assembling reduced model."
,
17
)
self
.
data
.
lastSetupMu
=
mu
def
getApproxReduced
(
self
,
mu
:
paramList
=
[])
->
sampList
:
"""
Evaluate reduced representation of approximant at arbitrary parameter.
Args:
mu: Target parameter.
"""
mu
=
checkParameterList
(
mu
,
self
.
data
.
npar
)[
0
]
if
(
not
hasattr
(
self
,
"lastSolvedApproxReduced"
)
or
self
.
lastSolvedApproxReduced
!=
mu
):
vbMng
(
self
,
"INIT"
,
"Computing RB solution at mu = {}."
.
format
(
mu
),
12
)
self
.
uApproxReduced
=
emptySampleList
()
for
i
in
range
(
len
(
mu
)):
self
.
assembleReducedModel
(
mu
[
i
])
vbMng
(
self
,
"INIT"
,
"Solving reduced model for mu = {}."
.
format
(
mu
[
i
]),
15
)
uAppR
=
np
.
linalg
.
solve
(
self
.
data
.
ARBmu
,
self
.
data
.
bRBmu
)
if
i
==
0
:
#self.data.ARBs[0].shape[-1], len(mu)
self
.
uApproxReduced
.
reset
((
len
(
uAppR
),
len
(
mu
)),
dtype
=
uAppR
.
dtype
)
self
.
uApproxReduced
[
i
]
=
uAppR
vbMng
(
self
,
"DEL"
,
"Done solving reduced model."
,
15
)
vbMng
(
self
,
"DEL"
,
"Done computing RB solution."
,
12
)
self
.
lastSolvedApproxReduced
=
mu
return
self
.
uApproxReduced
def
getPoles
(
self
,
marginalVals
:
ListAny
=
[
fp
],
jSupp
:
int
=
1
,
**
kwargs
)
->
Np1D
:
"""
Obtain approximant poles.
Returns:
Numpy complex vector of poles.
"""
if
not
self
.
data
.
affinePoly
:
RROMPyWarning
((
"Unable to compute approximate poles due "
"to parametric dependence (detected non-"
"polynomial). Change HFEngine.affinePoly to True "
"if necessary."
))
return
if
not
hasattr
(
marginalVals
,
"__len__"
):
marginalVals
=
[
marginalVals
]
mVals
=
list
(
marginalVals
)
try
:
rDim
=
mVals
.
index
(
fp
)
if
rDim
<
len
(
mVals
)
-
1
and
fp
in
mVals
[
rDim
+
1
:]:
raise
except
:
raise
RROMPyException
((
"Exactly 1 'freepar' entry in "
"marginalVals must be provided."
))
ARBs
=
self
.
data
.
ARBs
if
self
.
data
.
npar
>
1
:
mVals
[
rDim
]
=
self
.
data
.
mu0
(
rDim
)
mVals
=
checkParameter
(
mVals
)
.
data
.
flatten
()
mVals
[
rDim
]
=
fp
ARBs
=
marginalizePolyList
(
ARBs
,
mVals
,
"auto"
)
ev
=
eigvalsNonlinearDense
(
ARBs
,
jSupp
=
jSupp
,
**
kwargs
)
return
np
.
power
(
ev
,
1.
/
self
.
data
.
rescalingExp
[
rDim
])
Event Timeline
Log In to Comment