Page MenuHomec4science

rational_interpolant_pivoted_greedy.py
No OneTemporary

File Metadata

Created
Wed, May 1, 10:30

rational_interpolant_pivoted_greedy.py

# 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 .generic_pivoted_greedy_approximant import GenericPivotedGreedyApproximant
from rrompy.utilities.numerical import dot
from rrompy.reduction_methods.standard import RationalInterpolant
from rrompy.reduction_methods.pivoted import (RationalInterpolantPivoted,
PODGlobal)
from rrompy.utilities.base.types import paramList
from rrompy.utilities.base import verbosityManager as vbMng
from rrompy.utilities.exception_manager import RROMPyAssert
from rrompy.parameter import checkParameterList, emptyParameterList
__all__ = ['RationalInterpolantPivotedGreedy']
class RationalInterpolantPivotedGreedy(GenericPivotedGreedyApproximant,
RationalInterpolantPivoted):
"""
ROM pivoted greedy rational interpolant computation for parametric
problems.
Args:
HFEngine: HF problem solver.
mu0(optional): Default parameter. Defaults to 0.
directionPivot(optional): Pivot components. 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;
- 'matchingWeight': weight for pole matching optimization; defaults
to 1;
- 'cutOffTolerance': tolerance for ignoring parasitic poles;
defaults to np.inf;
- 'matchingWeightError': weight for pole matching optimization in
error estimation; defaults to 0;
- 'cutOffToleranceError': tolerance for ignoring parasitic poles
in error estimation; defaults to np.inf;
- 'S': total number of pivot samples current approximant relies
upon;
- 'samplerPivot': pivot sample point generator;
- 'SMarginal': number of starting marginal samples;
- 'samplerMarginalGrid': marginal sample point generator via sparse
grid;
- 'polybasis': type of polynomial basis for pivot interpolation;
defaults to 'MONOMIAL';
- 'polybasisMarginal': type of polynomial basis for marginal
interpolation; allowed values include 'MONOMIAL', 'CHEBYSHEV'
and 'LEGENDRE'; defaults to 'MONOMIAL';
- 'M': degree of rational interpolant numerator; defaults to 0;
- 'N': degree of rational interpolant denominator; defaults to 0;
- 'MMarginal': degree of marginal interpolant; defaults to 0;
- 'greedyTolMarginal': uniform error tolerance for marginal greedy
algorithm; defaults to 1e-1;
- 'maxIterMarginal': maximum number of marginal greedy steps;
defaults to 1e2;
- 'polydegreetypeMarginal': type of polynomial degree for marginal;
defaults to 'TOTAL';
- 'radialDirectionalWeights': radial basis weights for pivot
numerator; defaults to 1;
- 'radialDirectionalWeightsMarginal': radial basis weights for
marginal interpolant; defaults to 1;
- 'nNearestNeighbor': number of pivot nearest neighbors considered
if polybasis allows; defaults to -1;
- 'nNearestNeighborMarginal': number of marginal nearest neighbors
considered if polybasisMarginal allows; defaults to -1;
- 'interpRcond': tolerance for pivot interpolation; defaults to
None;
- 'interpRcondMarginal': tolerance for marginal interpolation;
defaults to None;
- 'robustTol': tolerance for robust rational denominator
management; defaults to 0;
- 'correctorForce': whether corrector should forcefully delete bad
poles; defaults to False;
- 'correctorTol': tolerance for corrector step; defaults to 0.,
i.e. no bad poles;
- 'correctorMaxIter': maximum number of corrector iterations;
defaults to 1.
Defaults to empty dict.
approx_state(optional): Whether to approximate state. Defaults to
False.
verbosity(optional): Verbosity level. Defaults to 10.
Attributes:
HFEngine: HF problem solver.
mu0: Default parameter.
directionPivot: Pivot components.
mus: Array of snapshot parameters.
musMarginal: Array of marginal 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;
- 'matchingWeight': weight for pole matching optimization;
- 'cutOffTolerance': tolerance for ignoring parasitic poles;
- 'matchingWeightError': weight for pole matching optimization in
error estimation;
- 'cutOffToleranceError': tolerance for ignoring parasitic poles
in error estimation;
- 'polybasis': type of polynomial basis for pivot interpolation;
- 'polybasisMarginal': type of polynomial basis for marginal
interpolation;
- 'M': degree of rational interpolant numerator;
- 'N': degree of rational interpolant denominator;
- 'MMarginal': degree of marginal interpolant;
- 'greedyTolMarginal': uniform error tolerance for marginal greedy
algorithm;
- 'maxIterMarginal': maximum number of marginal greedy steps;
- 'polydegreetypeMarginal': type of polynomial degree for marginal;
- 'radialDirectionalWeights': radial basis weights for pivot
numerator;
- 'radialDirectionalWeightsMarginal': radial basis weights for
marginal interpolant;
- 'nNearestNeighbor': number of pivot nearest neighbors considered
if polybasis allows;
- 'nNearestNeighborMarginal': number of marginal nearest neighbors
considered if polybasisMarginal allows;
- 'interpRcond': tolerance for pivot interpolation;
- 'interpRcondMarginal': tolerance for marginal interpolation;
- 'robustTol': tolerance for robust rational denominator
management;
- 'correctorForce': whether corrector should forcefully delete bad
poles;
- 'correctorTol': tolerance for corrector step;
- 'correctorMaxIter': maximum number of corrector iterations.
parameterListCritical: Recognized keys of critical approximant
parameters:
- 'S': total number of pivot samples current approximant relies
upon;
- 'samplerPivot': pivot sample point generator;
- 'SMarginal': total number of marginal samples current approximant
relies upon;
- 'samplerMarginalGrid': marginal sample point generator via sparse
grid.
approx_state: Whether to approximate state.
verbosity: Verbosity level.
POD: Whether to compute POD of snapshots.
matchingWeight: Weight for pole matching optimization.
cutOffTolerance: Tolerance for ignoring parasitic poles.
matchingWeightError: Weight for pole matching optimization in error
estimation.
cutOffToleranceError: Tolerance for ignoring parasitic poles in error
estimation.
S: Total number of pivot samples current approximant relies upon.
samplerPivot: Pivot sample point generator.
SMarginal: Total number of marginal samples current approximant relies
upon.
samplerMarginalGrid: Marginal sample point generator via sparse grid.
polybasis: Type of polynomial basis for pivot interpolation.
polybasisMarginal: Type of polynomial basis for marginal interpolation.
M: Degree of rational interpolant numerator.
N: Degree of rational interpolant denominator.
MMarginal: Degree of marginal interpolant.
greedyTolMarginal: Uniform error tolerance for marginal greedy
algorithm.
maxIterMarginal: Maximum number of marginal greedy steps.
polydegreetypeMarginal: Type of polynomial degree for marginal.
radialDirectionalWeights: Radial basis weights for pivot numerator.
radialDirectionalWeightsMarginal: Radial basis weights for marginal
interpolant.
nNearestNeighbor: Number of pivot nearest neighbors considered if
polybasis allows.
nNearestNeighborMarginal: Number of marginal nearest neighbors
considered if polybasisMarginal allows.
interpRcond: Tolerance for pivot interpolation.
interpRcondMarginal: Tolerance for marginal interpolation.
robustTol: Tolerance for robust rational denominator management.
correctorForce: Whether corrector should forcefully delete bad poles.
correctorTol: Tolerance for corrector step.
correctorMaxIter: Maximum number of corrector iterations.
muBoundsPivot: list of bounds for pivot parameter values.
muBoundsMarginal: list of bounds for marginal 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.
"""
def _finalizeSnapshots(self):
self.samplingEngine = self._samplingEngineOld
for muM, sEN in zip(self.musMargLoc, self.samplingEngs):
self.samplingEngine.samples += [sEN.samples]
self.samplingEngine.nsamples += [sEN.nsamples]
self.samplingEngine.mus += [sEN.mus]
self.samplingEngine.musMarginal.append(muM)
self.samplingEngine._derIdxs += [[(0,) * self.npar]
for _ in range(sEN.nsamples)]
if self.POD:
self.samplingEngine.RPOD += [sEN.RPOD]
self.samplingEngine.samples_full += [copy(sEN.samples_full)]
if self.POD == PODGlobal:
self.samplingEngine.coalesceSamples(self.interpRcondMarginal)
else:
self.samplingEngine.coalesceSamples()
def computeSnapshots(self):
"""Compute snapshots of solution map."""
RROMPyAssert(self._mode,
message = "Cannot start snapshot computation.")
vbMng(self, "INIT", "Starting computation of snapshots.", 5)
self.musPivot = self.samplerPivot.generatePoints(self.S)
self.mus = emptyParameterList()
self.mus.reset((self.S, self.HFEngine.npar))
self.samplingEngine.resetHistory()
for k in range(self.S):
self.mus.data[k, self.directionPivot] = self.musPivot[k].data
self.mus.data[k, self.directionMarginal] = self.musMargLoc[-1].data
self.samplingEngine.iterSample(self.mus)
vbMng(self, "DEL", "Done computing snapshots.", 5)
self._m_selfmus = copy(self.mus)
self._mus = self.musPivot
self._m_mu0 = copy(self.mu0)
self._m_HFErescalingExp = copy(self.HFEngine.rescalingExp)
self._mu0 = checkParameterList(self.mu0(self.directionPivot), 1)[0]
self.HFEngine.rescalingExp = [self.HFEngine.rescalingExp[
self.directionPivot[0]]]
def setupApproxPivoted(self, mus:paramList) -> int:
if self.checkComputedApprox(): return -1
RROMPyAssert(self._mode, message = "Cannot setup approximant.")
vbMng(self, "INIT", "Setting up {}.". format(self.name()), 10)
self.computeScaleFactor()
if self.trainedModel is None:
self.trainedModel = self.tModelType()
self.trainedModel.verbosity = self.verbosity
self.trainedModel.timestamp = self.timestamp
datadict = {"mu0": self.mu0, "projMat": np.zeros((0, 0)),
"scaleFactor": self.scaleFactor,
"rescalingExp": self.HFEngine.rescalingExp,
"directionPivot": self.directionPivot}
self.trainedModel.data = self.initializeModelData(datadict)[0]
self.trainedModel.data.Qs, self.trainedModel.data.Ps = [], []
_trainedModelOld = copy(self.trainedModel)
self._scaleFactorOldPivot = copy(self.scaleFactor)
self.scaleFactor = self.scaleFactorPivot
self._temporaryPivot = 1
self._samplingEngineOld = copy(self.samplingEngine)
self.musMargLoc, self.samplingEngs = [], [None] * len(mus)
Qs, Ps = [None] * len(mus), [None] * len(mus)
self.verbosity -= 15
for j, mu in enumerate(mus):
RationalInterpolant.setupSampling(self)
self.trainedModel = None
self.musMargLoc += [mu]
RationalInterpolant.setupApprox(self)
self._mu0 = self._m_mu0
self._mus = self._m_selfmus
self.HFEngine.rescalingExp = self._m_HFErescalingExp
del self._m_mu0, self._m_selfmus, self._m_HFErescalingExp
self.samplingEngs[j] = copy(self.samplingEngine)
Qs[j] = copy(self.trainedModel.data.Q)
Ps[j] = copy(self.trainedModel.data.P)
self.scaleFactor = self._scaleFactorOldPivot
del self._scaleFactorOldPivot, self._temporaryPivot
self._finalizeSnapshots()
del self._samplingEngineOld, self.musMargLoc, self.samplingEngs
self._mus = self.samplingEngine.musCoalesced
self.trainedModel = _trainedModelOld
self.trainedModel.data.mus = copy(self.mus)
self.trainedModel.data.musMarginal = copy(self.musMarginal)
padRight = (self.samplingEngine.nsamplesTot
- self.trainedModel.data.projMat.shape[1])
for j in range(len(self.trainedModel.data.Ps)):
nsj = self.samplingEngine.nsamples[j]
self.trainedModel.data.Ps[j].pad(0, padRight)
self.trainedModel.data.HIs[j].pad(0, padRight)
padLeft = self.trainedModel.data.projMat.shape[1]
for j in range(len(mus)):
nsj = self.samplingEngine.nsamples[j]
if self.POD == PODGlobal:
rRightj = self.samplingEngine.RPODCPart[:,
padLeft : padLeft + nsj]
Ps[j].postmultiplyTensorize(rRightj.T)
else:
padRight -= nsj
Ps[j].pad(padLeft, padRight)
padLeft += nsj
pMat = self.samplingEngine.samplesCoalesced.data
pMatEff = dot(self.HFEngine.C, pMat) if self.approx_state else pMat
self.trainedModel.data.projMat = pMatEff
self.trainedModel.data.Qs += Qs
self.trainedModel.data.Ps += Ps
self.verbosity += 15
vbMng(self, "DEL", "Done setting up approximant.", 10)
return 0

Event Timeline