Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F60776180
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
Thu, May 2, 13:06
Size
33 KB
Mime Type
text/x-python
Expires
Sat, May 4, 13:06 (2 d)
Engine
blob
Format
Raw Data
Handle
17411129
Attached To
R6746 RationalROMPy
rational_interpolant.py
View Options
# Copyright (C) 2018-2020 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
scipy.linalg
import
eig
from
collections.abc
import
Iterable
from
.generic_standard_approximant
import
GenericStandardApproximant
from
rrompy.utilities.poly_fitting.polynomial
import
(
polybases
as
ppb
,
polyfitname
,
polyvander
as
pvP
,
polyTimes
,
PolynomialInterpolator
as
PI
,
PolynomialInterpolatorNodal
as
PIN
)
from
rrompy.utilities.poly_fitting.heaviside
import
rational2heaviside
from
rrompy.utilities.poly_fitting.radial_basis
import
(
polybases
as
rbpb
,
RadialBasisInterpolator
as
RBI
)
from
rrompy.utilities.base.types
import
(
Np1D
,
Np2D
,
Tuple
,
List
,
paramList
,
interpEng
)
from
rrompy.utilities.base
import
verbosityManager
as
vbMng
from
rrompy.utilities.numerical
import
pseudoInverse
,
dot
,
baseDistanceMatrix
from
rrompy.utilities.numerical.factorials
import
multifactorial
from
rrompy.utilities.numerical.hash_derivative
import
(
nextDerivativeIndices
,
hashDerivativeToIdx
as
hashD
,
hashIdxToDerivative
as
hashI
)
from
rrompy.utilities.numerical.degree
import
(
reduceDegreeN
,
degreeTotalToFull
,
fullDegreeMaxMask
,
totalDegreeMaxMask
)
from
rrompy.utilities.exception_manager
import
(
RROMPyException
,
RROMPyAssert
,
RROMPyWarning
)
__all__
=
[
'RationalInterpolant'
]
def
polyTimesTable
(
P
:
interpEng
,
mus
:
Np1D
,
reorder
:
List
[
int
],
derIdxs
:
List
[
List
[
List
[
int
]]],
scl
:
Np1D
=
None
)
->
Np2D
:
"""Table of polynomial products."""
if
not
isinstance
(
P
,
PI
):
raise
RROMPyException
((
"Polynomial to evaluate must be a polynomial "
"interpolator."
))
Pvals
=
[[
0.
]
*
len
(
derIdx
)
for
derIdx
in
derIdxs
]
for
j
,
derIdx
in
enumerate
(
derIdxs
):
nder
=
len
(
derIdx
)
for
der
in
range
(
nder
):
derI
=
hashI
(
der
,
P
.
npar
)
Pvals
[
j
][
der
]
=
P
([
mus
[
j
]],
derI
,
scl
)
/
multifactorial
(
derI
)
return
blockDiagDer
(
Pvals
,
reorder
,
derIdxs
)
def
vanderInvTable
(
vanInv
:
Np2D
,
idxs
:
List
[
int
],
reorder
:
List
[
int
],
derIdxs
:
List
[
List
[
List
[
int
]]])
->
Np2D
:
"""Table of Vandermonde pseudo-inverse."""
S
=
len
(
reorder
)
Ts
=
[
None
]
*
len
(
idxs
)
for
k
in
range
(
len
(
idxs
)):
invLocs
=
[
None
]
*
len
(
derIdxs
)
idxGlob
=
0
for
j
,
derIdx
in
enumerate
(
derIdxs
):
nder
=
len
(
derIdx
)
idxGlob
+=
nder
idxLoc
=
np
.
arange
(
S
)[(
reorder
>=
idxGlob
-
nder
)
*
(
reorder
<
idxGlob
)]
invLocs
[
j
]
=
vanInv
[
k
,
idxLoc
]
Ts
[
k
]
=
blockDiagDer
(
invLocs
,
reorder
,
derIdxs
,
[
2
,
1
,
0
])
return
Ts
def
blockDiagDer
(
vals
:
List
[
Np1D
],
reorder
:
List
[
int
],
derIdxs
:
List
[
List
[
List
[
int
]]],
permute
:
List
[
int
]
=
None
)
->
Np2D
:
"""Table of derivative values for point confluence."""
S
=
len
(
reorder
)
T
=
np
.
zeros
((
S
,
S
),
dtype
=
np
.
complex
)
if
permute
is
None
:
permute
=
[
0
,
1
,
2
]
idxGlob
=
0
for
j
,
derIdx
in
enumerate
(
derIdxs
):
nder
=
len
(
derIdx
)
idxGlob
+=
nder
idxLoc
=
np
.
arange
(
S
)[(
reorder
>=
idxGlob
-
nder
)
*
(
reorder
<
idxGlob
)]
val
=
vals
[
j
]
for
derI
,
derIdxI
in
enumerate
(
derIdx
):
for
derJ
,
derIdxJ
in
enumerate
(
derIdx
):
diffIdx
=
[
x
-
y
for
(
x
,
y
)
in
zip
(
derIdxI
,
derIdxJ
)]
if
all
([
x
>=
0
for
x
in
diffIdx
]):
diffj
=
hashD
(
diffIdx
)
i1
,
i2
,
i3
=
np
.
array
([
derI
,
derJ
,
diffj
])[
permute
]
T
[
idxLoc
[
i1
],
idxLoc
[
i2
]]
=
val
[
i3
]
return
T
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': kind of snapshots orthogonalization; allowed values
include 0, 1/2, and 1; defaults to 1, i.e. POD;
- 'scaleFactorDer': scaling factors for derivative computation;
defaults to 'AUTO';
- '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
'AUTO', i.e. maximum allowed;
- 'N': degree of rational interpolant denominator; defaults to
'AUTO', i.e. maximum allowed;
- 'polydegreetype': type of polynomial degree; defaults to 'TOTAL';
- 'radialDirectionalWeights': radial basis weights for interpolant
numerator; defaults to 1;
- 'radialDirectionalWeightsAdapt': bounds for adaptive rescaling of
radial basis weights; defaults to [-1, -1];
- 'functionalSolve': strategy for minimization of denominator
functional; allowed values include 'NORM', 'DOMINANT',
'BARYCENTRIC[_NORM]', and 'BARYCENTRIC_AVERAGE' (check pdf in
main folder for explanation); defaults to 'NORM';
- 'interpTol': tolerance for interpolation; defaults to None;
- 'QTol': 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': kind of snapshots orthogonalization;
- 'scaleFactorDer': scaling factors for derivative computation;
- '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;
- 'radialDirectionalWeightsAdapt': bounds for adaptive rescaling of
radial basis weights;
- 'functionalSolve': strategy for minimization of denominator
functional;
- 'interpTol': tolerance for interpolation via numpy.polyfit;
- 'QTol': 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.
verbosity: Verbosity level.
POD: Kind of snapshots orthogonalization.
scaleFactorDer: Scaling factors for derivative computation.
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.
radialDirectionalWeightsAdapt: Bounds for adaptive rescaling of radial
basis weights.
functionalSolve: Strategy for minimization of denominator functional.
interpTol: Tolerance for interpolation via numpy.polyfit.
QTol: 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.
"""
_allowedFunctionalSolveKinds
=
[
"NORM"
,
"DOMINANT"
,
"BARYCENTRIC_NORM"
,
"BARYCENTRIC_AVERAGE"
]
def
__init__
(
self
,
*
args
,
**
kwargs
):
self
.
_preInit
()
self
.
_addParametersToList
([
"polybasis"
,
"M"
,
"N"
,
"polydegreetype"
,
"radialDirectionalWeights"
,
"radialDirectionalWeightsAdapt"
,
"functionalSolve"
,
"interpTol"
,
"QTol"
],
[
"MONOMIAL"
,
"AUTO"
,
"AUTO"
,
"TOTAL"
,
1.
,
[
-
1.
,
-
1.
],
"NORM"
,
-
1
,
0.
])
super
()
.
__init__
(
*
args
,
**
kwargs
)
self
.
_postInit
()
@property
def
tModelType
(
self
):
from
.trained_model.trained_model_rational
import
TrainedModelRational
return
TrainedModelRational
@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
:
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
functionalSolve
(
self
):
"""Value of functionalSolve."""
return
self
.
_functionalSolve
@functionalSolve.setter
def
functionalSolve
(
self
,
functionalSolve
):
try
:
functionalSolve
=
functionalSolve
.
upper
()
.
strip
()
.
replace
(
" "
,
""
)
if
functionalSolve
==
"BARYCENTRIC"
:
functionalSolve
+=
"_NORM"
if
functionalSolve
not
in
self
.
_allowedFunctionalSolveKinds
:
raise
RROMPyException
((
"Prescribed functionalSolve not "
"recognized."
))
self
.
_functionalSolve
=
functionalSolve
except
:
RROMPyWarning
((
"Prescribed functionalSolve not recognized. "
"Overriding to 'NORM'."
))
self
.
_functionalSolve
=
"NORM"
self
.
_approxParameters
[
"functionalSolve"
]
=
self
.
functionalSolve
@property
def
interpTol
(
self
):
"""Value of interpTol."""
return
self
.
_interpTol
@interpTol.setter
def
interpTol
(
self
,
interpTol
):
self
.
_interpTol
=
interpTol
self
.
_approxParameters
[
"interpTol"
]
=
self
.
interpTol
@property
def
radialDirectionalWeights
(
self
):
"""Value of radialDirectionalWeights."""
return
self
.
_radialDirectionalWeights
@radialDirectionalWeights.setter
def
radialDirectionalWeights
(
self
,
radialDirectionalWeights
):
if
isinstance
(
radialDirectionalWeights
,
Iterable
):
radialDirectionalWeights
=
list
(
radialDirectionalWeights
)
else
:
radialDirectionalWeights
=
[
radialDirectionalWeights
]
self
.
_radialDirectionalWeights
=
radialDirectionalWeights
self
.
_approxParameters
[
"radialDirectionalWeights"
]
=
(
self
.
radialDirectionalWeights
)
@property
def
radialDirectionalWeightsAdapt
(
self
):
"""Value of radialDirectionalWeightsAdapt."""
return
self
.
_radialDirectionalWeightsAdapt
@radialDirectionalWeightsAdapt.setter
def
radialDirectionalWeightsAdapt
(
self
,
radialDirectionalWeightsAdapt
):
self
.
_radialDirectionalWeightsAdapt
=
radialDirectionalWeightsAdapt
self
.
_approxParameters
[
"radialDirectionalWeightsAdapt"
]
=
(
self
.
radialDirectionalWeightsAdapt
)
@property
def
M
(
self
):
"""Value of M."""
return
self
.
_M
@M.setter
def
M
(
self
,
M
):
if
isinstance
(
M
,
str
):
M
=
M
.
strip
()
.
replace
(
" "
,
""
)
if
"-"
not
in
M
:
M
=
M
+
"-0"
self
.
_M_isauto
,
self
.
_M_shift
=
True
,
int
(
M
.
split
(
"-"
)[
-
1
])
M
=
0
if
M
<
0
:
raise
RROMPyException
(
"M must be non-negative."
)
self
.
_M
=
M
self
.
_approxParameters
[
"M"
]
=
self
.
M
def
_setMAuto
(
self
):
self
.
M
=
max
(
0
,
reduceDegreeN
(
self
.
S
,
self
.
S
,
self
.
npar
,
self
.
polydegreetype
)
-
self
.
_M_shift
)
vbMng
(
self
,
"MAIN"
,
"Automatically setting M to {}."
.
format
(
self
.
M
),
25
)
@property
def
N
(
self
):
"""Value of N."""
return
self
.
_N
@N.setter
def
N
(
self
,
N
):
if
isinstance
(
N
,
str
):
N
=
N
.
strip
()
.
replace
(
" "
,
""
)
if
"-"
not
in
N
:
N
=
N
+
"-0"
self
.
_N_isauto
,
self
.
_N_shift
=
True
,
int
(
N
.
split
(
"-"
)[
-
1
])
N
=
0
if
N
<
0
:
raise
RROMPyException
(
"N must be non-negative."
)
self
.
_N
=
N
self
.
_approxParameters
[
"N"
]
=
self
.
N
def
_setNAuto
(
self
):
self
.
N
=
max
(
0
,
reduceDegreeN
(
self
.
S
,
self
.
S
,
self
.
npar
,
self
.
polydegreetype
)
-
self
.
_N_shift
)
vbMng
(
self
,
"MAIN"
,
"Automatically setting N to {}."
.
format
(
self
.
N
),
25
)
@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
QTol
(
self
):
"""Value of tolerance for robust rational denominator management."""
return
self
.
_QTol
@QTol.setter
def
QTol
(
self
,
QTol
):
if
QTol
<
0.
:
RROMPyWarning
((
"Overriding prescribed negative robustness "
"tolerance to 0."
))
QTol
=
0.
self
.
_QTol
=
QTol
self
.
_approxParameters
[
"QTol"
]
=
self
.
QTol
def
resetSamples
(
self
):
"""Reset samples."""
super
()
.
resetSamples
()
self
.
_musUniqueCN
=
None
self
.
_derIdxs
=
None
self
.
_reorder
=
None
def
_setupInterpolationIndices
(
self
):
"""Setup parameters for polyvander."""
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
)
if
hasattr
(
self
,
"_N_isauto"
):
self
.
_setNAuto
()
else
:
N
=
reduceDegreeN
(
self
.
N
,
self
.
S
,
self
.
npar
,
self
.
polydegreetype
)
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
.
functionalSolve
!=
"NORM"
and
self
.
npar
>
1
:
RROMPyWarning
((
"Strategy for functional optimization must be "
"'NORM' for more than one parameter. "
"Overriding to 'NORM'."
))
self
.
functionalSolve
=
"NORM"
if
(
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
and
self
.
N
+
1
<
self
.
S
):
RROMPyWarning
((
"Barycentric strategy cannot be applied with "
"Least Squares. Overriding to 'NORM'."
))
self
.
functionalSolve
=
"NORM"
if
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
:
invD
,
TN
=
None
,
None
self
.
_setupInterpolationIndices
()
if
len
(
self
.
_musUnique
)
!=
self
.
S
:
RROMPyWarning
((
"Barycentric functional optimization "
"cannot be applied to repeated samples. "
"Overriding to 'NORM'."
))
self
.
functionalSolve
=
"NORM"
if
self
.
functionalSolve
[:
11
]
!=
"BARYCENTRIC"
:
invD
,
TN
=
self
.
_computeInterpolantInverseBlocks
()
if
self
.
POD
==
1
:
sampleE
=
self
.
samplingEngine
.
Rscale
Rscaling
=
None
elif
self
.
POD
==
1
/
2
:
sampleE
=
self
.
samplingEngine
.
samples_normal
Rscaling
=
self
.
samplingEngine
.
Rscale
else
:
sampleE
=
self
.
samplingEngine
.
samples
Rscaling
=
None
ev
,
eV
=
self
.
findeveVG
(
sampleE
,
invD
,
TN
,
Rscaling
)
if
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
:
break
nevBad
=
np
.
sum
(
np
.
abs
(
ev
/
ev
[
-
1
])
<
self
.
QTol
)
if
not
nevBad
:
break
if
self
.
npar
==
1
:
dN
=
nevBad
else
:
#if self.npar > 1 and self.functionalSolve == "NORM":
dN
=
self
.
N
-
reduceDegreeN
(
self
.
N
,
len
(
eV
)
-
nevBad
,
self
.
npar
,
self
.
polydegreetype
)
vbMng
(
self
,
"MAIN"
,
(
"Smallest {} eigenvalue{} below tolerance. Reducing N by "
"{}."
)
.
format
(
nevBad
,
"s"
*
(
nevBad
>
1
),
dN
),
10
)
self
.
N
=
self
.
N
-
dN
if
hasattr
(
self
,
"_gram"
):
del
self
.
_gram
if
self
.
N
<=
0
:
self
.
N
,
eV
=
0
,
np
.
ones
((
1
,)
*
self
.
npar
,
dtype
=
np
.
complex
)
if
self
.
N
>
0
and
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
:
q
=
PIN
()
q
.
polybasis
,
q
.
nodes
=
self
.
polybasis0
,
eV
else
:
q
=
PI
()
q
.
npar
,
q
.
polybasis
=
self
.
npar
,
self
.
polybasis0
if
self
.
polydegreetype
==
"TOTAL"
:
q
.
coeffs
=
degreeTotalToFull
(
tuple
([
self
.
N
+
1
]
*
self
.
npar
),
self
.
npar
,
eV
)
else
:
q
.
coeffs
=
eV
.
reshape
([
self
.
N
+
1
]
*
self
.
npar
)
vbMng
(
self
,
"DEL"
,
"Done computing denominator."
,
7
)
return
q
def
_setupNumerator
(
self
):
"""Compute rational numerator."""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot setup numerator."
)
vbMng
(
self
,
"INIT"
,
"Starting computation of numerator."
,
7
)
self
.
_setupInterpolationIndices
()
Qevaldiag
=
polyTimesTable
(
self
.
trainedModel
.
data
.
Q
,
self
.
_musUniqueCN
,
self
.
_reorder
,
self
.
_derIdxs
,
self
.
scaleFactorRel
)
if
self
.
POD
==
1
:
Qevaldiag
=
Qevaldiag
.
dot
(
self
.
samplingEngine
.
Rscale
.
T
)
elif
self
.
POD
==
1
/
2
:
Qevaldiag
=
Qevaldiag
*
self
.
samplingEngine
.
Rscale
if
hasattr
(
self
,
"_M_isauto"
):
self
.
_setMAuto
()
M
=
self
.
M
else
:
M
=
reduceDegreeN
(
self
.
M
,
self
.
S
,
self
.
npar
,
self
.
polydegreetype
)
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
:
pParRest
=
[
self
.
M
,
self
.
polybasis
,
self
.
verbosity
>=
5
,
self
.
polydegreetype
==
"TOTAL"
,
{
"derIdxs"
:
self
.
_derIdxs
,
"reorder"
:
self
.
_reorder
,
"scl"
:
self
.
scaleFactorRel
}]
if
self
.
polybasis
in
ppb
:
p
=
PI
()
else
:
self
.
computeScaleFactor
()
rDWEff
=
np
.
array
([
w
*
f
for
w
,
f
in
zip
(
self
.
radialDirectionalWeights
,
self
.
scaleFactor
)])
pParRest
=
pParRest
[:
2
]
+
[
rDWEff
]
+
pParRest
[
2
:]
pParRest
[
-
1
][
"optimizeScalingBounds"
]
=
(
self
.
radialDirectionalWeightsAdapt
)
p
=
RBI
()
if
self
.
polybasis
in
ppb
+
rbpb
:
pParRest
+=
[{
"rcond"
:
self
.
interpTol
}]
wellCond
,
msg
=
p
.
setupByInterpolation
(
self
.
_musUniqueCN
,
Qevaldiag
,
*
pParRest
)
vbMng
(
self
,
"MAIN"
,
msg
,
5
)
if
wellCond
:
break
vbMng
(
self
,
"MAIN"
,
(
"Polyfit is poorly conditioned. Reducing M "
"by 1."
),
10
)
self
.
M
=
self
.
M
-
1
if
self
.
M
<
0
:
raise
RROMPyException
((
"Instability in computation of numerator. "
"Aborting."
))
self
.
M
=
M
vbMng
(
self
,
"DEL"
,
"Done computing numerator."
,
7
)
return
p
def
setupApprox
(
self
)
->
int
:
"""Compute rational interpolant."""
if
self
.
checkComputedApprox
():
return
-
1
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot setup approximant."
)
vbMng
(
self
,
"INIT"
,
"Setting up {}."
.
format
(
self
.
name
()),
5
)
self
.
computeSnapshots
()
self
.
_setupTrainedModel
(
self
.
samplingEngine
.
projectionMatrix
)
self
.
_setupRational
(
self
.
_setupDenominator
())
self
.
trainedModel
.
data
.
approxParameters
=
copy
(
self
.
approxParameters
)
vbMng
(
self
,
"DEL"
,
"Done setting up approximant."
,
5
)
return
0
def
_setupRational
(
self
,
Q
:
interpEng
,
P
:
interpEng
=
None
):
vbMng
(
self
,
"INIT"
,
"Starting approximant finalization."
,
5
)
self
.
trainedModel
.
data
.
Q
=
Q
if
P
is
None
:
P
=
self
.
_setupNumerator
()
while
self
.
N
>
0
and
self
.
npar
==
1
:
if
self
.
HFEngine
.
_ignoreResidues
:
pls
=
Q
.
roots
()
cfs
,
projMat
=
None
,
None
else
:
cfs
,
pls
,
_
=
rational2heaviside
(
P
,
Q
)
cfs
=
cfs
[:
self
.
N
]
.
T
if
self
.
POD
!=
1
:
projMat
=
self
.
samplingEngine
.
projectionMatrix
else
:
projMat
=
None
foci
=
self
.
sampler
.
normalFoci
()
plsA
=
self
.
mapParameterList
(
self
.
mapParameterList
(
self
.
mu0
)(
0
,
0
)
+
self
.
scaleFactor
*
pls
,
"B"
)(
0
)
idxBad
=
self
.
HFEngine
.
flagBadPolesResiduesAbsolute
(
plsA
,
cfs
,
projMat
)
if
not
self
.
HFEngine
.
_ignoreResidues
:
cfs
[:,
idxBad
]
=
0.
idxBad
+=
self
.
HFEngine
.
flagBadPolesResiduesRelative
(
pls
,
cfs
,
projMat
,
foci
)
idxBad
=
idxBad
>
0
if
not
np
.
any
(
idxBad
):
break
vbMng
(
self
,
"MAIN"
,
"Removing {} spurious pole{} out of {}."
.
format
(
np
.
sum
(
idxBad
),
"s"
*
(
np
.
sum
(
idxBad
)
>
1
),
self
.
N
),
10
)
if
isinstance
(
Q
,
PIN
):
Q
.
nodes
=
Q
.
nodes
[
idxBad
==
False
]
else
:
Q
=
PI
()
Q
.
npar
=
self
.
npar
Q
.
polybasis
=
self
.
polybasis0
Q
.
coeffs
=
np
.
ones
(
1
,
dtype
=
np
.
complex
)
for
pl
in
pls
[
idxBad
==
False
]:
Q
.
coeffs
=
polyTimes
(
Q
.
coeffs
,
[
-
pl
,
1.
],
Pbasis
=
Q
.
polybasis
,
Rbasis
=
Q
.
polybasis
)
Q
.
coeffs
/=
np
.
linalg
.
norm
(
Q
.
coeffs
)
self
.
trainedModel
.
data
.
Q
=
Q
self
.
N
=
Q
.
deg
[
0
]
P
=
self
.
_setupNumerator
()
self
.
trainedModel
.
data
.
P
=
P
vbMng
(
self
,
"DEL"
,
"Terminated approximant finalization."
,
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
()
pvPPar
=
[
self
.
polybasis0
,
self
.
_derIdxs
,
self
.
_reorder
,
self
.
scaleFactorRel
]
full
=
self
.
N
+
1
==
self
.
S
==
len
(
self
.
_musUniqueCN
)
if
full
:
mus
=
self
.
_musUniqueCN
[
self
.
_reorder
]
dist
=
baseDistanceMatrix
(
mus
,
magnitude
=
False
)[
...
,
0
]
dist
[
np
.
arange
(
self
.
N
+
1
),
np
.
arange
(
self
.
N
+
1
)]
=
multifactorial
([
self
.
N
])
fitinvE
=
np
.
prod
(
dist
,
axis
=
1
)
**
-
1
vbMng
(
self
,
"MAIN"
,
(
"Evaluating quasi-Lagrangian basis of degree {} at {} "
"sample points."
)
.
format
(
self
.
N
,
self
.
N
+
1
),
5
)
invD
=
[
np
.
diag
(
fitinvE
)]
TN
=
pvP
(
self
.
_musUniqueCN
,
self
.
N
,
*
pvPPar
)
else
:
while
self
.
N
>=
0
:
if
self
.
polydegreetype
==
"TOTAL"
:
Neff
=
self
.
N
idxsB
=
totalDegreeMaxMask
(
self
.
N
,
self
.
npar
)
else
:
#if self.polydegreetype == "FULL":
Neff
=
[
self
.
N
]
*
self
.
npar
idxsB
=
fullDegreeMaxMask
(
self
.
N
,
self
.
npar
)
TN
=
pvP
(
self
.
_musUniqueCN
,
Neff
,
*
pvPPar
)
fitOut
=
pseudoInverse
(
TN
,
rcond
=
self
.
interpTol
,
full
=
True
)
vbMng
(
self
,
"MAIN"
,
(
"Fitting {} samples with degree {} through {}... "
"Conditioning of pseudoinverse system: {:.4e}."
)
.
format
(
TN
.
shape
[
0
],
self
.
N
,
polyfitname
(
self
.
polybasis0
),
fitOut
[
1
][
1
][
0
]
/
fitOut
[
1
][
1
][
-
1
]),
5
)
if
fitOut
[
1
][
0
]
==
TN
.
shape
[
1
]:
fitinv
=
fitOut
[
0
][
idxsB
,
:]
break
vbMng
(
self
,
"MAIN"
,
"Polyfit is poorly conditioned. Reducing N by 1."
,
10
)
self
.
N
=
self
.
N
-
1
if
self
.
N
<
0
:
raise
RROMPyException
((
"Instability in computation of "
"denominator. Aborting."
))
invD
=
vanderInvTable
(
fitinv
,
idxsB
,
self
.
_reorder
,
self
.
_derIdxs
)
return
invD
,
TN
def
findeveVG
(
self
,
sampleE
:
Np2D
,
invD
:
List
[
Np2D
],
TN
:
Np2D
,
Rscaling
:
Np1D
=
None
)
->
Tuple
[
Np1D
,
Np2D
]:
"""
Compute eigenvalues and eigenvectors of rational denominator matrix, or
of its right chol factor if POD.
"""
RROMPyAssert
(
self
.
_mode
,
message
=
"Cannot solve spectral problem."
)
if
self
.
POD
==
1
:
if
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
:
Rstack
=
sampleE
else
:
vbMng
(
self
,
"INIT"
,
"Building generalized half-gramian."
,
10
)
S
,
eWidth
=
sampleE
.
shape
[
0
],
len
(
invD
)
Rstack
=
np
.
zeros
((
S
*
eWidth
,
TN
.
shape
[
1
]),
dtype
=
np
.
complex
)
for
k
in
range
(
eWidth
):
Rstack
[
k
*
S
:
(
k
+
1
)
*
S
,
:]
=
dot
(
sampleE
,
dot
(
invD
[
k
],
TN
))
vbMng
(
self
,
"DEL"
,
"Done building half-gramian."
,
10
)
_
,
s
,
Vh
=
np
.
linalg
.
svd
(
Rstack
,
full_matrices
=
False
)
evG
,
eVG
=
s
[::
-
1
],
Vh
[::
-
1
]
.
T
.
conj
()
evExp
,
probKind
=
-
2.
,
"svd "
else
:
if
not
hasattr
(
self
,
"_gram"
):
vbMng
(
self
,
"INIT"
,
"Building gramian matrix."
,
10
)
self
.
_gram
=
self
.
HFEngine
.
innerProduct
(
sampleE
,
sampleE
,
is_state
=
True
)
if
Rscaling
is
not
None
:
self
.
_gram
=
(
self
.
_gram
.
T
*
Rscaling
.
conj
())
.
T
*
Rscaling
vbMng
(
self
,
"DEL"
,
"Done building gramian."
,
10
)
if
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
:
G
=
self
.
_gram
else
:
vbMng
(
self
,
"INIT"
,
"Building generalized gramian."
,
10
)
G
=
np
.
zeros
((
TN
.
shape
[
1
],)
*
2
,
dtype
=
np
.
complex
)
for
k
in
range
(
len
(
invD
)):
iDkN
=
dot
(
invD
[
k
],
TN
)
G
+=
dot
(
dot
(
self
.
_gram
,
iDkN
)
.
T
,
iDkN
.
conj
())
.
T
vbMng
(
self
,
"DEL"
,
"Done building gramian."
,
10
)
evG
,
eVG
=
np
.
linalg
.
eigh
(
G
)
evExp
,
probKind
=
-
1.
,
"eigen"
if
(
self
.
functionalSolve
in
[
"NORM"
,
"BARYCENTRIC_NORM"
]
or
np
.
sum
(
np
.
abs
(
evG
)
<
np
.
finfo
(
float
)
.
eps
*
np
.
abs
(
evG
[
-
1
])
*
len
(
evG
))
==
1
):
eV
=
eVG
[:,
0
]
elif
self
.
functionalSolve
==
"BARYCENTRIC_AVERAGE"
:
eV
=
eVG
.
dot
(
evG
**
evExp
*
np
.
sum
(
eVG
,
axis
=
0
)
.
conj
())
else
:
eV
=
eVG
.
dot
(
evG
**
evExp
*
eVG
[
0
]
.
conj
())
vbMng
(
self
,
"MAIN"
,
(
"Solved {}problem of size {} with condition number "
"{:.4e}."
)
.
format
(
probKind
,
len
(
evG
)
-
1
,
evG
[
-
1
]
/
evG
[
1
]),
5
)
if
self
.
functionalSolve
[:
11
]
==
"BARYCENTRIC"
:
S
,
mus
=
len
(
eV
),
self
.
_musUniqueCN
[
self
.
_reorder
]
.
flatten
()
arrow
=
np
.
zeros
((
S
+
1
,)
*
2
,
dtype
=
np
.
complex
)
arrow
[
1
:,
0
]
=
1.
arrow
[
0
,
1
:]
=
eV
arrow
[
np
.
arange
(
1
,
S
+
1
),
np
.
arange
(
1
,
S
+
1
)]
=
mus
active
=
np
.
eye
(
S
+
1
)
active
[
0
,
0
]
=
0.
poles
,
qTm1
=
eig
(
arrow
,
active
)
eVgood
=
np
.
isinf
(
poles
)
+
np
.
isnan
(
poles
)
==
False
poles
=
poles
[
eVgood
]
self
.
N
=
len
(
poles
)
if
self
.
QTol
>
0
:
# compare optimal score with self.N poles to those obtained
# by removing one of the poles
qTm1
=
qTm1
[
1
:,
eVgood
]
.
conj
()
**
-
1.
dists
=
mus
.
reshape
(
-
1
,
1
)
-
mus
dists
[
np
.
arange
(
S
),
np
.
arange
(
S
)]
=
multifactorial
([
self
.
N
])
dists
=
np
.
prod
(
dists
,
axis
=
1
)
.
conj
()
**
-
1.
qComp
=
np
.
empty
((
self
.
N
+
1
,
S
),
dtype
=
np
.
complex
)
qComp
[
0
]
=
dists
*
np
.
prod
(
qTm1
,
axis
=
1
)
for
j
in
range
(
self
.
N
):
qTmj
=
np
.
prod
(
qTm1
[:,
np
.
arange
(
self
.
N
)
!=
j
],
axis
=
1
)
qComp
[
j
+
1
]
=
dists
*
qTmj
Lqs
=
qComp
.
dot
(
eVG
)
scores
=
np
.
real
(
np
.
sum
(
Lqs
*
evG
**
-
evExp
*
Lqs
.
conj
(),
axis
=
1
))
evBad
=
scores
[
1
:]
<
self
.
QTol
*
scores
[
0
]
nevBad
=
np
.
sum
(
evBad
)
if
nevBad
:
vbMng
(
self
,
"MAIN"
,
(
"Suboptimal pole{} detected. Reducing N by "
"{}."
)
.
format
(
"s"
*
(
nevBad
>
1
),
nevBad
),
10
)
self
.
N
=
self
.
N
-
nevBad
poles
=
poles
[
evBad
==
False
]
eV
=
poles
return
evG
[
1
:],
eV
def
getResidues
(
self
,
*
args
,
**
kwargs
)
->
Tuple
[
paramList
,
Np2D
]:
"""
Obtain approximant residues.
Returns:
Matrix with residues as columns.
"""
return
self
.
trainedModel
.
getResidues
(
*
args
,
**
kwargs
)
Event Timeline
Log In to Comment