Page MenuHomec4science

reduced_basis_pivoted.py
No OneTemporary

File Metadata

Created
Sun, Apr 28, 13:32

reduced_basis_pivoted.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_approximant import GenericPivotedApproximant
from rrompy.hfengines.base import MarginalProxyEngine
from rrompy.reduction_methods.trained_model import (TrainedModelPivotedData,
TrainedModelPivotedReducedBasis as tModel)
from rrompy.reduction_methods.base.reduced_basis_utils import \
projectAffineDecomposition
from rrompy.utilities.base.types import (Np1D, Np2D, List, ListAny, Tuple,
DictAny, HFEng, paramVal, sampList)
from rrompy.utilities.base import verbosityManager as vbMng, freepar as fp
from rrompy.utilities.numerical import customPInv
from rrompy.utilities.exception_manager import (RROMPyWarning, RROMPyException,
RROMPyAssert)
__all__ = ['ReducedBasisPivoted']
class ReducedBasisPivoted(GenericPivotedApproximant):
"""
ROM pivoted RB approximant 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;
- '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;
- 'samplerMarginal': marginal sample point generator;
- 'R': rank for pivot Galerkin projection; defaults to prod(S);
- 'PODTolerance': tolerance for pivot snapshots POD; defaults to
-1;
- 'polybasisMarginal': type of polynomial basis for marginal
interpolation; allowed values include 'MONOMIAL', 'CHEBYSHEV'
and 'LEGENDRE'; defaults to 'MONOMIAL';
- 'MMarginal': degree of marginal interpolant; defaults to 0;
- 'radialBasisMarginal': radial basis family for marginal
interpolant; defaults to 0, i.e. no radial basis;
- 'radialBasisWeightsMarginal': radial basis weights for marginal
interpolant; defaults to 0, i.e. identity;
- 'interpRcondMarginal': tolerance for marginal interpolation;
defaults to None.
Defaults to empty dict.
homogeneized(optional): Whether to homogeneize Dirichlet BCs. 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.
homogeneized: Whether to homogeneize Dirichlet BCs.
approxRadius: Dummy radius of approximant (i.e. distance from mu0 to
farthest sample point).
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;
- 'R': rank for Galerkin projection;
- 'PODTolerance': tolerance for snapshots POD;
- 'polybasisMarginal': type of polynomial basis for marginal
interpolation;
- 'MMarginal': degree of marginal interpolant;
- 'radialBasisMarginal': radial basis family for marginal
interpolant;
- 'radialBasisWeightsMarginal': radial basis weights for marginal
interpolant;
- 'interpRcondMarginal': tolerance for marginal interpolation.
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;
- 'samplerMarginal': marginal sample point generator.
POD: Whether to compute POD of snapshots.
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.
samplerMarginal: Marginal sample point generator.
R: Rank for Galerkin projection.
PODTolerance: Tolerance for pivot snapshots POD.
polybasisMarginal: Type of polynomial basis for marginal interpolation.
MMarginal: Degree of marginal interpolant.
radialBasisMarginal: Radial basis family for marginal interpolant.
radialBasisWeightsMarginal: Radial basis weights for marginal
interpolant.
interpRcondMarginal: Tolerance for marginal interpolation.
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.
As: List of sparse matrices (in CSC format) representing coefficients
of linear system matrix.
bs: List of numpy vectors representing coefficients of linear system
RHS.
ARBs: List of sparse matrices (in CSC format) representing coefficients
of compressed linear system matrix.
bRBs: List of numpy vectors representing coefficients of compressed
linear system RHS.
"""
def __init__(self, HFEngine:HFEng, mu0 : paramVal = None,
directionPivot : ListAny = [0],
approxParameters : DictAny = {}, homogeneized : bool = False,
verbosity : int = 10, timestamp : bool = True):
self._preInit()
self._addParametersToList(["R", "PODTolerance"], [1, -1])
super().__init__(HFEngine = HFEngine, mu0 = mu0,
directionPivot = directionPivot,
approxParameters = approxParameters,
homogeneized = homogeneized,
verbosity = verbosity, timestamp = timestamp)
self._postInit()
@property
def S(self):
"""Value of S."""
return self._S
@S.setter
def S(self, S):
GenericPivotedApproximant.S.fset(self, S)
if not hasattr(self, "_R"):
self._R = np.prod(self.S) * np.prod(self.SMarginal)
@property
def R(self):
"""Value of R. Its assignment may change S."""
return self._R
@R.setter
def R(self, R):
if R < 0: raise RROMPyException("R must be non-negative.")
self._R = R
self._approxParameters["R"] = self.R
if (hasattr(self, "_S") and hasattr(self, "_SMarginal")
and np.prod(self.S) * np.prod(self.SMarginal) < self.R):
RROMPyWarning(("Prescribed S and SMarginal are too small. "
"Decreasing R."))
self.R = np.prod(self.S) * np.prod(self.SMarginal)
@property
def PODTolerance(self):
"""Value of PODTolerance."""
return self._PODTolerance
@PODTolerance.setter
def PODTolerance(self, PODTolerance):
self._PODTolerance = PODTolerance
self._approxParameters["PODTolerance"] = self.PODTolerance
def _setupProjectionMatrix(self):
"""Compute projection matrix."""
RROMPyAssert(self._mode, message = "Cannot setup numerator.")
vbMng(self, "INIT", "Starting computation of projection matrix.", 7)
if self.POD:
U, s, _ = np.linalg.svd(self.samplingEngine.RPODCoalesced)
s = s ** 2.
else:
Gramian = self.HFEngine.innerProduct(
self.samplingEngine.samplesCoalesced,
self.samplingEngine.samplesCoalesced)
U, s, _ = np.linalg.svd(Gramian)
nsamples = self.samplingEngine.samplesCoalesced.shape[1]
snorm = np.cumsum(s[::-1]) / np.sum(s)
nPODTrunc = min(nsamples - np.argmax(snorm > self.PODTolerance),
self.R)
pMat = self.samplingEngine.samplesCoalesced.dot(U[:, : nPODTrunc])
vbMng(self, "MAIN",
"Assembling {}x{} projection matrix from {} samples.".format(
*(pMat.shape), nsamples), 5)
vbMng(self, "DEL", "Done computing projection matrix.", 7)
return pMat
def _setupAffineBlocks(self):
"""Compute list of marginalized affine blocks of system."""
hasAs = hasattr(self, "AsList") and self.AsList is not None
hasbs = hasattr(self, "bsList") and self.bsList is not None
if hasAs and hasbs: return
vbMng(self, "INIT", "Computing affine blocks of system.", 10)
mu0Eff = self.mu0.data[0, self.directionPivot]
if not hasAs: self.AsList = [None] * len(self.musMarginal)
if not hasbs: self.bsList = [None] * len(self.musMarginal)
for k, muM in enumerate(self.musMarginal):
muEff = [fp] * self.npar
for jj, kk in enumerate(self.directionMarginal):
muEff[kk] = muM(0, jj)
MEnginek = MarginalProxyEngine(self.HFEngine, muEff)
if not hasAs:
self.AsList[k] = MEnginek.affineLinearSystemA(mu0Eff,
self.scaleFactorPivot)
if not hasbs:
self.bsList[k] = MEnginek.affineLinearSystemb(mu0Eff,
self.scaleFactorPivot,
self.homogeneized)
vbMng(self, "DEL", "Done computing affine blocks.", 10)
def setupApprox(self):
"""Compute RB projection matrix."""
if self.checkComputedApprox():
return
RROMPyAssert(self._mode, message = "Cannot setup approximant.")
vbMng(self, "INIT", "Setting up {}.".format(self.name()), 5)
self.computeScaleFactor()
self.computeSnapshots()
self._setupAffineBlocks()
pMat = self._setupProjectionMatrix()
ARBsList, bRBsList, pList = self.assembleReducedSystemMarginal(pMat)
if self.trainedModel is None:
self.trainedModel = tModel()
self.trainedModel.verbosity = self.verbosity
self.trainedModel.timestamp = self.timestamp
data = TrainedModelPivotedData(self.trainedModel.name(), self.mu0,
pList, self.scaleFactor,
self.HFEngine.rescalingExp,
self.directionPivot)
data.musPivot = copy(self.musPivot)
data.musMarginal = copy(self.musMarginal)
self.trainedModel.data = data
else:
self.trainedModel = self.trainedModel
self.trainedModel.data.projMat = copy(pList)
self.trainedModel.data.marginalInterp = self._setupMarginalInterp()
self.trainedModel.data.ARBsList = ARBsList
self.trainedModel.data.bRBsList = bRBsList
self.trainedModel.data.approxParameters = copy(self.approxParameters)
vbMng(self, "DEL", "Done setting up approximant.", 5)
def assembleReducedSystemMarginal(self, pMat : sampList = None)\
-> Tuple[List[List[Np2D]],
List[List[Np1D]], List[Np2D]]:
"""Build affine blocks of RB linear system through projections."""
if pMat is None:
self.setupApprox()
ARBsList = self.trainedModel.data.ARBsList
bRBsList = self.trainedModel.data.bRBsList
else:
vbMng(self, "INIT", "Projecting affine terms of HF model.", 10)
ARBsList = [None] * len(self.musMarginal)
bRBsList = [None] * len(self.musMarginal)
projList = [None] * len(self.musMarginal)
for k, (As, bs, sample) in enumerate(zip(self.AsList, self.bsList,
self.samplingEngine.samples)):
compLocal = self.HFEngine.innerProduct(sample, pMat)
projList[k] = sample.dot(customPInv(compLocal))
ARBsList[k], bRBsList[k] = projectAffineDecomposition(As, bs,
projList[k])
vbMng(self, "DEL", "Done projecting affine terms.", 10)
return ARBsList, bRBsList, projList

Event Timeline