Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F61082039
rational_interpolant.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 4, 10:24
Size
26 KB
Mime Type
text/x-python
Expires
Mon, May 6, 10:24 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
17462407
Attached To
R6746 RationalROMPy
rational_interpolant.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/>.
#
from
copy
import
deepcopy
as
copy
import
numpy
as
np
from
rrompy.reduction_methods.base
import
checkRobustTolerance
from
.generic_standard_approximant
import
GenericStandardApproximant
from
rrompy.utilities.poly_fitting.polynomial
import
(
polybases
as
ppb
,
polyfitname
,
polyvander
as
pvP
,
polyvanderTotal
as
pvTP
,
PolynomialInterpolator
as
PI
)
from
rrompy.utilities.poly_fitting.radial_basis
import
(
polybases
as
rbpb
,
RadialBasisInterpolator
as
RBI
)
from
rrompy.utilities.poly_fitting.moving_least_squares
import
(
polybases
as
mlspb
,
MovingLeastSquaresInterpolator
as
MLSI
)
from
rrompy.reduction_methods.trained_model
import
(
TrainedModelRational
as
tModel
)
from
rrompy.reduction_methods.trained_model
import
TrainedModelData
from
rrompy.utilities.base.types
import
(
Np1D
,
Np2D
,
HFEng
,
DictAny
,
Tuple
,
List
,
paramVal
,
sampList
)
from
rrompy.utilities.base
import
verbosityManager
as
vbMng
from
rrompy.utilities.numerical
import
(
multifactorial
,
customPInv
,
dot
,
fullDegreeN
,
totalDegreeN
,
degreeTotalToFull
,
fullDegreeMaxMask
,
totalDegreeMaxMask
,
nextDerivativeIndices
,
hashDerivativeToIdx
as
hashD
,
hashIdxToDerivative
as
hashI
)
from
rrompy.utilities.exception_manager
import
(
RROMPyException
,
RROMPyAssert
,
RROMPyWarning
)
__all__
=
[
'RationalInterpolant'
]
class
RationalInterpolant
(
GenericStandardApproximant
):
"""
ROM rational interpolant computation for parametric problems.
Args:
HFEngine: HF problem solver.
mu0(optional): Default parameter. Defaults to 0.
approxParameters(optional): Dictionary containing values for main
parameters of approximant. Recognized keys are:
- 'POD': whether to compute POD of snapshots; defaults to True;
- 'S': total number of samples current approximant relies upon;
- 'sampler': sample point generator;
- 'polybasis': type of polynomial basis for interpolation; defaults
to 'MONOMIAL';
- 'M': degree of rational interpolant numerator; defaults to 0;
- 'N': degree of rational interpolant denominator; defaults to 0;
- 'polydegreetype': type of polynomial degree; defaults to 'TOTAL';
- 'radialDirectionalWeights': radial basis weights for interpolant
numerator; defaults to 0, i.e. identity;
- 'nNearestNeighbor': mumber of nearest neighbors considered in
numerator if polybasis allows; defaults to -1;
- 'interpRcond': tolerance for interpolation; defaults to None;
- 'robustTol': tolerance for robust rational denominator
management; defaults to 0.
Defaults to empty dict.
verbosity(optional): Verbosity level. Defaults to 10.
Attributes:
HFEngine: HF problem solver.
mu0: Default parameter.
mus: Array of snapshot parameters.
approxParameters: Dictionary containing values for main parameters of
approximant. Recognized keys are in parameterList.
parameterListSoft: Recognized keys of soft approximant parameters:
- 'POD': whether to compute POD of snapshots;
- 'polybasis': type of polynomial basis for interpolation;
- 'M': degree of rational interpolant numerator;
- 'N': degree of rational interpolant denominator;
- 'polydegreetype': type of polynomial degree;
- 'radialDirectionalWeights': radial basis weights for interpolant
numerator;
- 'nNearestNeighbor': mumber of nearest neighbors considered in
numerator if polybasis allows;
- 'interpRcond': tolerance for interpolation via numpy.polyfit;
- 'robustTol': tolerance for robust rational denominator
management.
parameterListCritical: Recognized keys of critical approximant
parameters:
- 'S': total number of samples current approximant relies upon;
- 'sampler': sample point generator.
POD: Whether to compute POD of snapshots.
S: Number of solution snapshots over which current approximant is
based upon.
sampler: Sample point generator.
polybasis: type of polynomial basis for interpolation.
M: Numerator degree of approximant.
N: Denominator degree of approximant.
polydegreetype: Type of polynomial degree.
radialDirectionalWeights: Radial basis weights for interpolant
numerator.
nNearestNeighbor: Number of nearest neighbors considered in numerator
if polybasis allows.
interpRcond: Tolerance for interpolation via numpy.polyfit.
robustTol: Tolerance for robust rational denominator management.
muBounds: list of bounds for parameter values.
samplingEngine: Sampling engine.
uHF: High fidelity solution(s) with parameter(s) lastSolvedHF as
sampleList.
lastSolvedHF: Parameter(s) corresponding to last computed high fidelity
solution(s) as parameterList.
uApproxReduced: Reduced approximate solution(s) with parameter(s)
lastSolvedApprox as sampleList.
lastSolvedApproxReduced: Parameter(s) corresponding to last computed
reduced approximate solution(s) as parameterList.
uApprox: Approximate solution(s) with parameter(s) lastSolvedApprox as
sampleList.
lastSolvedApprox: Parameter(s) corresponding to last computed
approximate solution(s) as parameterList.
Q: Numpy 1D vector containing complex coefficients of approximant
denominator.
P: Numpy 2D vector whose columns are FE dofs of coefficients of
approximant numerator.
"""
def
__init__
(
self
,
HFEngine
:
HFEng
,
mu0
:
paramVal
=
None
,
approxParameters
:
DictAny
=
{},
verbosity
:
int
=
10
,
timestamp
:
bool
=
True
):
self
.
_preInit
()
self
.
_addParametersToList
([
"polybasis"
,
"M"
,
"N"
,
"polydegreetype"
,
"radialDirectionalWeights"
,
"nNearestNeighbor"
,
"interpRcond"
,
"robustTol"
],
[
"MONOMIAL"
,
0
,
0
,
"TOTAL"
,
1
,
-
1
,
-
1
,
0
])
super
()
.
__init__
(
HFEngine
=
HFEngine
,
mu0
=
mu0
,
approxParameters
=
approxParameters
,
verbosity
=
verbosity
,
timestamp
=
timestamp
)
self
.
catchInstability
=
False
self
.
_postInit
()
@property
def
polybasis
(
self
):
"""Value of polybasis."""
return
self
.
_polybasis
@polybasis.setter
def
polybasis
(
self
,
polybasis
):
try
:
polybasis
=
polybasis
.
upper
()
.
strip
()
.
replace
(
" "
,
""
)
if
polybasis
not
in
ppb
+
rbpb
+
mlspb
:
raise
RROMPyException
(
"Prescribed polybasis not recognized."
)
self
.
_polybasis
=
polybasis
except
:
RROMPyWarning
((
"Prescribed polybasis not recognized. Overriding "
"to 'MONOMIAL'."
))
self
.
_polybasis
=
"MONOMIAL"
self
.
_approxParameters
[
"polybasis"
]
=
self
.
polybasis
@property
def
polybasis0
(
self
):
if
"_"
in
self
.
polybasis
:
return
self
.
polybasis
.
split
(
"_"
)[
0
]
return
self
.
polybasis
@property
def
interpRcond
(
self
):
"""Value of interpRcond."""
return
self
.
_interpRcond
@interpRcond.setter
def
interpRcond
(
self
,
interpRcond
):
self
.
_interpRcond
=
interpRcond
self
.
_approxParameters
[
"interpRcond"
]
=
self
.
interpRcond
@property
def
radialDirectionalWeights
(
self
):
"""Value of radialDirectionalWeights."""
return
self
.
_radialDirectionalWeights
@radialDirectionalWeights.setter
def
radialDirectionalWeights
(
self
,
radialDirectionalWeights
):
self
.
_radialDirectionalWeights
=
radialDirectionalWeights
self
.
_approxParameters
[
"radialDirectionalWeights"
]
=
(
self
.
radialDirectionalWeights
)
@property
def
nNearestNeighbor
(
self
):
"""Value of nNearestNeighbor."""
return
self
.
_nNearestNeighbor
@nNearestNeighbor.setter
def
nNearestNeighbor
(
self
,
nNearestNeighbor
):
self
.
_nNearestNeighbor
=
nNearestNeighbor
self
.
_approxParameters
[
"nNearestNeighbor"
]
=
self
.
nNearestNeighbor
@property
def
M
(
self
):
"""Value of M."""
return
self
.
_M
@M.setter
def
M
(
self
,
M
):
if
M
<
0
:
raise
RROMPyException
(
"M must be non-negative."
)
self
.
_M
=
M
self
.
_approxParameters
[
"M"
]
=
self
.
M
@property
def
N
(
self
):
"""Value of N."""
return
self
.
_N
@N.setter
def
N
(
self
,
N
):
if
N
<
0
:
raise
RROMPyException
(
"N must be non-negative."
)
self
.
_N
=
N
self
.
_approxParameters
[
"N"
]
=
self
.
N
@property
def
polydegreetype
(
self
):
"""Value of polydegreetype."""
return
self
.
_polydegreetype
@polydegreetype.setter
def
polydegreetype
(
self
,
polydegreetype
):
try
:
polydegreetype
=
polydegreetype
.
upper
()
.
strip
()
.
replace
(
" "
,
""
)
if
polydegreetype
not
in
[
"TOTAL"
,
"FULL"
]:
raise
RROMPyException
((
"Prescribed polydegreetype not "
"recognized."
))
self
.
_polydegreetype
=
polydegreetype
except
:
RROMPyWarning
((
"Prescribed polydegreetype not recognized. "
"Overriding to 'TOTAL'."
))
self
.
_polydegreetype
=
"TOTAL"
self
.
_approxParameters
[
"polydegreetype"
]
=
self
.
polydegreetype
@property
def
robustTol
(
self
):
"""Value of tolerance for robust rational denominator management."""
return
self
.
_robustTol
@robustTol.setter
def
robustTol
(
self
,
robustTol
):
if
robustTol
<
0.
:
RROMPyWarning
((
"Overriding prescribed negative robustness "
"tolerance to 0."
))
robustTol
=
0.
self
.
_robustTol
=
robustTol
self
.
_approxParameters
[
"robustTol"
]
=
self
.
robustTol
def
resetSamples
(
self
):
"""Reset samples."""
super
()
.
resetSamples
()
self
.
_musUniqueCN
=
None
self
.
_derIdxs
=
None
self
.
_reorder
=
None
def
_setupInterpolationIndices
(
self
):
"""Setup parameters for polyvander."""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot setup interpolation indices."
)
if
self
.
_musUniqueCN
is
None
or
len
(
self
.
_reorder
)
!=
len
(
self
.
mus
):
self
.
_musUniqueCN
,
musIdxsTo
,
musIdxs
,
musCount
=
(
self
.
trainedModel
.
centerNormalize
(
self
.
mus
)
.
unique
(
return_index
=
True
,
return_inverse
=
True
,
return_counts
=
True
))
self
.
_musUnique
=
self
.
mus
[
musIdxsTo
]
self
.
_derIdxs
=
[
None
]
*
len
(
self
.
_musUniqueCN
)
self
.
_reorder
=
np
.
empty
(
len
(
musIdxs
),
dtype
=
int
)
filled
=
0
for
j
,
cnt
in
enumerate
(
musCount
):
self
.
_derIdxs
[
j
]
=
nextDerivativeIndices
([],
self
.
mus
.
shape
[
1
],
cnt
)
jIdx
=
np
.
nonzero
(
musIdxs
==
j
)[
0
]
self
.
_reorder
[
jIdx
]
=
np
.
arange
(
filled
,
filled
+
cnt
)
filled
+=
cnt
def
_setupDenominator
(
self
):
"""Compute rational denominator."""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot setup denominator."
)
vbMng
(
self
,
"INIT"
,
"Starting computation of denominator."
,
7
)
while
self
.
N
>
0
:
invD
,
fitinv
=
self
.
_computeInterpolantInverseBlocks
()
if
self
.
POD
:
ev
,
eV
=
self
.
findeveVGQR
(
self
.
samplingEngine
.
RPOD
,
invD
)
else
:
ev
,
eV
=
self
.
findeveVGExplicit
(
self
.
samplingEngine
.
samples
,
invD
)
nevBad
=
checkRobustTolerance
(
ev
,
self
.
robustTol
)
if
nevBad
<=
1
:
break
if
self
.
catchInstability
:
raise
RROMPyException
((
"Instability in denominator "
"computation: eigenproblem is poorly "
"conditioned."
))
RROMPyWarning
((
"Smallest {} eigenvalues below tolerance. Reducing "
"N by 1."
)
.
format
(
nevBad
))
self
.
N
=
self
.
N
-
1
if
self
.
N
<=
0
:
self
.
_N
=
0
eV
=
np
.
ones
((
1
,
1
))
q
=
PI
()
q
.
npar
=
self
.
npar
q
.
polybasis
=
self
.
polybasis0
if
self
.
polydegreetype
==
"TOTAL"
:
q
.
coeffs
=
degreeTotalToFull
(
tuple
([
self
.
N
+
1
]
*
self
.
npar
),
self
.
npar
,
eV
[:,
0
])
else
:
q
.
coeffs
=
eV
[:,
0
]
.
reshape
([
self
.
N
+
1
]
*
self
.
npar
)
vbMng
(
self
,
"DEL"
,
"Done computing denominator."
,
7
)
return
q
,
fitinv
def
_setupNumerator
(
self
):
"""Compute rational numerator."""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot setup numerator."
)
vbMng
(
self
,
"INIT"
,
"Starting computation of numerator."
,
7
)
Qevaldiag
=
np
.
zeros
((
len
(
self
.
mus
),
len
(
self
.
mus
)),
dtype
=
np
.
complex
)
verb
=
self
.
trainedModel
.
verbosity
self
.
trainedModel
.
verbosity
=
0
self
.
_setupInterpolationIndices
()
idxGlob
=
0
for
j
,
derIdxs
in
enumerate
(
self
.
_derIdxs
):
nder
=
len
(
derIdxs
)
idxLoc
=
np
.
arange
(
len
(
self
.
mus
))[(
self
.
_reorder
>=
idxGlob
)
*
(
self
.
_reorder
<
idxGlob
+
nder
)]
idxGlob
+=
nder
Qval
=
[
0
]
*
nder
for
der
in
range
(
nder
):
derIdx
=
hashI
(
der
,
self
.
npar
)
Qval
[
der
]
=
(
self
.
trainedModel
.
getQVal
(
self
.
_musUnique
[
j
],
derIdx
,
scl
=
np
.
power
(
self
.
scaleFactor
,
-
1.
))
/
multifactorial
(
derIdx
))
for
derU
,
derUIdx
in
enumerate
(
derIdxs
):
for
derQ
,
derQIdx
in
enumerate
(
derIdxs
):
diffIdx
=
[
x
-
y
for
(
x
,
y
)
in
zip
(
derUIdx
,
derQIdx
)]
if
all
([
x
>=
0
for
x
in
diffIdx
]):
diffj
=
hashD
(
diffIdx
)
Qevaldiag
[
idxLoc
[
derU
],
idxLoc
[
derQ
]]
=
Qval
[
diffj
]
if
self
.
POD
:
Qevaldiag
=
Qevaldiag
.
dot
(
self
.
samplingEngine
.
RPOD
.
T
)
self
.
trainedModel
.
verbosity
=
verb
cfun
=
totalDegreeN
if
self
.
polydegreetype
==
"TOTAL"
else
fullDegreeN
M
=
copy
(
self
.
M
)
while
len
(
self
.
mus
)
<
cfun
(
M
,
self
.
npar
):
M
-=
1
if
M
<
self
.
M
:
RROMPyWarning
((
"M too large compared to S. Reducing M by "
"{}"
)
.
format
(
self
.
M
-
M
))
self
.
M
=
M
while
self
.
M
>=
0
:
if
self
.
polybasis
in
ppb
:
p
=
PI
()
wellCond
,
msg
=
p
.
setupByInterpolation
(
self
.
_musUniqueCN
,
Qevaldiag
,
self
.
M
,
self
.
polybasis
,
self
.
verbosity
>=
5
,
self
.
polydegreetype
==
"TOTAL"
,
{
"derIdxs"
:
self
.
_derIdxs
,
"reorder"
:
self
.
_reorder
,
"scl"
:
np
.
power
(
self
.
scaleFactor
,
-
1.
)},
{
"rcond"
:
self
.
interpRcond
})
elif
self
.
polybasis
in
rbpb
:
p
=
RBI
()
wellCond
,
msg
=
p
.
setupByInterpolation
(
self
.
_musUniqueCN
,
Qevaldiag
,
self
.
M
,
self
.
polybasis
,
self
.
radialDirectionalWeights
,
self
.
verbosity
>=
5
,
self
.
polydegreetype
==
"TOTAL"
,
{
"derIdxs"
:
self
.
_derIdxs
,
"reorder"
:
self
.
_reorder
,
"scl"
:
np
.
power
(
self
.
scaleFactor
,
-
1.
),
"nNearestNeighbor"
:
self
.
nNearestNeighbor
},
{
"rcond"
:
self
.
interpRcond
})
else
:
# if self.polybasis in mlspb:
p
=
MLSI
()
wellCond
,
msg
=
p
.
setupByInterpolation
(
self
.
_musUniqueCN
,
Qevaldiag
,
self
.
M
,
self
.
polybasis
,
self
.
radialDirectionalWeights
,
self
.
verbosity
>=
5
,
self
.
polydegreetype
==
"TOTAL"
,
{
"derIdxs"
:
self
.
_derIdxs
,
"reorder"
:
self
.
_reorder
,
"scl"
:
np
.
power
(
self
.
scaleFactor
,
-
1.
),
"nNearestNeighbor"
:
self
.
nNearestNeighbor
})
vbMng
(
self
,
"MAIN"
,
msg
,
5
)
if
wellCond
:
break
if
self
.
catchInstability
:
raise
RROMPyException
((
"Instability in numerator computation: "
"polyfit is poorly conditioned."
))
RROMPyWarning
(
"Polyfit is poorly conditioned. Reducing M by 1."
)
self
.
M
=
self
.
M
-
1
if
self
.
M
<
0
:
raise
RROMPyException
((
"Instability in computation of numerator. "
"Aborting."
))
vbMng
(
self
,
"DEL"
,
"Done computing numerator."
,
7
)
return
p
def
setupApprox
(
self
):
"""
Compute rational interpolant.
SVD-based robust eigenvalue management.
"""
if
self
.
checkComputedApprox
():
return
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot setup approximant."
)
vbMng
(
self
,
"INIT"
,
"Setting up {}."
.
format
(
self
.
name
()),
5
)
self
.
computeSnapshots
()
if
self
.
trainedModel
is
None
:
self
.
trainedModel
=
tModel
()
self
.
trainedModel
.
verbosity
=
self
.
verbosity
self
.
trainedModel
.
timestamp
=
self
.
timestamp
data
=
TrainedModelData
(
self
.
trainedModel
.
name
(),
self
.
mu0
,
self
.
samplingEngine
.
samples
,
self
.
scaleFactor
,
self
.
HFEngine
.
rescalingExp
)
data
.
mus
=
copy
(
self
.
mus
)
self
.
trainedModel
.
data
=
data
else
:
self
.
trainedModel
=
self
.
trainedModel
self
.
trainedModel
.
data
.
projMat
=
copy
(
self
.
samplingEngine
.
samples
)
if
self
.
N
>
0
:
Q
=
self
.
_setupDenominator
()[
0
]
else
:
Q
=
PI
()
Q
.
coeffs
=
np
.
ones
(
tuple
([
1
]
*
self
.
npar
),
dtype
=
np
.
complex
)
Q
.
npar
=
self
.
npar
Q
.
polybasis
=
self
.
polybasis
self
.
trainedModel
.
data
.
Q
=
Q
self
.
trainedModel
.
data
.
P
=
self
.
_setupNumerator
()
self
.
trainedModel
.
data
.
approxParameters
=
copy
(
self
.
approxParameters
)
vbMng
(
self
,
"DEL"
,
"Done setting up approximant."
,
5
)
def
_computeInterpolantInverseBlocks
(
self
)
->
Tuple
[
List
[
Np2D
],
Np2D
]:
"""
Compute inverse factors for minimal interpolant target functional.
"""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot solve eigenvalue problem."
)
self
.
_setupInterpolationIndices
()
cfun
=
totalDegreeN
if
self
.
polydegreetype
==
"TOTAL"
else
fullDegreeN
N
=
copy
(
self
.
N
)
while
len
(
self
.
mus
)
<
cfun
(
N
,
self
.
npar
):
N
-=
1
if
N
<
self
.
N
:
RROMPyWarning
((
"N too large compared to S. Reducing N by "
"{}"
)
.
format
(
self
.
N
-
N
))
self
.
N
=
N
while
self
.
N
>=
0
:
if
self
.
polydegreetype
==
"TOTAL"
:
TE
,
_
,
argIdxs
=
pvTP
(
self
.
_musUniqueCN
,
self
.
N
,
self
.
polybasis0
,
self
.
_derIdxs
,
self
.
_reorder
,
scl
=
np
.
power
(
self
.
scaleFactor
,
-
1.
))
TE
=
TE
[:,
argIdxs
]
idxsB
=
totalDegreeMaxMask
(
self
.
N
,
self
.
npar
)
else
:
#if self.polydegreetype == "FULL":
TE
=
pvP
(
self
.
_musUniqueCN
,
[
self
.
N
]
*
self
.
npar
,
self
.
polybasis0
,
self
.
_derIdxs
,
self
.
_reorder
,
scl
=
np
.
power
(
self
.
scaleFactor
,
-
1.
))
idxsB
=
fullDegreeMaxMask
(
self
.
N
,
self
.
npar
)
fitOut
=
customPInv
(
TE
,
rcond
=
self
.
interpRcond
,
full
=
True
)
vbMng
(
self
,
"MAIN"
,
(
"Fitting {} samples with degree {} through {}... "
"Conditioning of pseudoinverse system: {:.4e}."
)
.
format
(
TE
.
shape
[
0
],
self
.
N
,
polyfitname
(
self
.
polybasis0
),
fitOut
[
1
][
1
][
0
]
/
fitOut
[
1
][
1
][
-
1
]),
5
)
if
fitOut
[
1
][
0
]
==
TE
.
shape
[
1
]:
fitinv
=
fitOut
[
0
][
idxsB
,
:]
break
if
self
.
catchInstability
:
raise
RROMPyException
((
"Instability in denominator "
"computation: polyfit is poorly "
"conditioned."
))
RROMPyWarning
(
"Polyfit is poorly conditioned. Reducing N by 1."
)
self
.
N
=
self
.
N
-
1
if
self
.
polydegreetype
==
"TOTAL"
:
TN
,
_
,
argIdxs
=
pvTP
(
self
.
_musUniqueCN
,
self
.
N
,
self
.
polybasis0
,
self
.
_derIdxs
,
self
.
_reorder
,
scl
=
np
.
power
(
self
.
scaleFactor
,
-
1.
))
TN
=
TN
[:,
argIdxs
]
else
:
#if self.polydegreetype == "FULL":
TN
=
pvP
(
self
.
_musUniqueCN
,
[
self
.
N
]
*
self
.
npar
,
self
.
polybasis0
,
self
.
_derIdxs
,
self
.
_reorder
,
scl
=
np
.
power
(
self
.
scaleFactor
,
-
1.
))
invD
=
[
None
]
*
(
len
(
idxsB
))
for
k
in
range
(
len
(
idxsB
)):
pseudoInv
=
np
.
diag
(
fitinv
[
k
,
:])
idxGlob
=
0
for
j
,
derIdxs
in
enumerate
(
self
.
_derIdxs
):
nder
=
len
(
derIdxs
)
idxGlob
+=
nder
if
nder
>
1
:
idxLoc
=
np
.
arange
(
len
(
self
.
mus
))[
(
self
.
_reorder
>=
idxGlob
-
nder
)
*
(
self
.
_reorder
<
idxGlob
)]
invLoc
=
fitinv
[
k
,
idxLoc
]
pseudoInv
[
np
.
ix_
(
idxLoc
,
idxLoc
)]
=
0.
for
diffj
,
diffjIdx
in
enumerate
(
derIdxs
):
for
derQ
,
derQIdx
in
enumerate
(
derIdxs
):
derUIdx
=
[
x
-
y
for
(
x
,
y
)
in
zip
(
diffjIdx
,
derQIdx
)]
if
all
([
x
>=
0
for
x
in
derUIdx
]):
derU
=
hashD
(
derUIdx
)
pseudoInv
[
idxLoc
[
derU
],
idxLoc
[
derQ
]]
=
(
invLoc
[
diffj
])
invD
[
k
]
=
dot
(
pseudoInv
,
TN
)
return
invD
,
fitinv
def
findeveVGExplicit
(
self
,
sampleE
:
sampList
,
invD
:
List
[
Np2D
])
->
Tuple
[
Np1D
,
Np2D
]:
"""
Compute explicitly eigenvalues and eigenvectors of rational denominator
matrix.
"""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot solve eigenvalue problem."
)
nEnd
=
invD
[
0
]
.
shape
[
1
]
eWidth
=
len
(
invD
)
vbMng
(
self
,
"INIT"
,
"Building gramian matrix."
,
10
)
gramian
=
self
.
HFEngine
.
innerProduct
(
sampleE
,
sampleE
)
G
=
np
.
zeros
((
nEnd
,
nEnd
),
dtype
=
np
.
complex
)
for
k
in
range
(
eWidth
):
G
+=
dot
(
dot
(
gramian
,
invD
[
k
])
.
T
,
invD
[
k
]
.
conj
())
.
T
vbMng
(
self
,
"DEL"
,
"Done building gramian."
,
10
)
vbMng
(
self
,
"INIT"
,
"Solving eigenvalue problem for gramian matrix."
,
7
)
ev
,
eV
=
np
.
linalg
.
eigh
(
G
)
vbMng
(
self
,
"MAIN"
,
(
"Solved eigenvalue problem of size {} with condition number "
"{:.4e}."
)
.
format
(
nEnd
,
ev
[
-
1
]
/
ev
[
0
]),
5
)
vbMng
(
self
,
"DEL"
,
"Done solving eigenvalue problem."
,
7
)
return
ev
,
eV
def
findeveVGQR
(
self
,
RPODE
:
Np2D
,
invD
:
List
[
Np2D
])
->
Tuple
[
Np1D
,
Np2D
]:
"""
Compute eigenvalues and eigenvectors of rational denominator matrix
through SVD of R factor.
"""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot solve eigenvalue problem."
)
nEnd
=
invD
[
0
]
.
shape
[
1
]
S
=
RPODE
.
shape
[
0
]
eWidth
=
len
(
invD
)
vbMng
(
self
,
"INIT"
,
"Building half-gramian matrix stack."
,
10
)
Rstack
=
np
.
zeros
((
S
*
eWidth
,
nEnd
),
dtype
=
np
.
complex
)
for
k
in
range
(
eWidth
):
Rstack
[
k
*
S
:
(
k
+
1
)
*
S
,
:]
=
dot
(
RPODE
,
invD
[
k
])
vbMng
(
self
,
"DEL"
,
"Done building half-gramian."
,
10
)
vbMng
(
self
,
"INIT"
,
"Solving svd for square root of gramian matrix."
,
7
)
_
,
s
,
eV
=
np
.
linalg
.
svd
(
Rstack
,
full_matrices
=
False
)
ev
=
s
[::
-
1
]
eV
=
eV
[::
-
1
,
:]
.
T
.
conj
()
vbMng
(
self
,
"MAIN"
,
(
"Solved svd problem of size {} x {} with condition number "
"{:.4e}."
)
.
format
(
*
Rstack
.
shape
,
s
[
0
]
/
s
[
-
1
]),
5
)
vbMng
(
self
,
"DEL"
,
"Done solving svd."
,
7
)
return
ev
,
eV
def
getResidues
(
self
,
*
args
,
**
kwargs
)
->
Np1D
:
"""
Obtain approximant residues.
Returns:
Matrix with residues as columns.
"""
return
self
.
trainedModel
.
getResidues
(
*
args
,
**
kwargs
)
Event Timeline
Log In to Comment