Page MenuHomec4science

trained_model_rb.py
No OneTemporary

File Metadata

Created
Fri, May 3, 14:57

trained_model_rb.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/>.
#
import numpy as np
from scipy.linalg import eigvals
from .trained_model import TrainedModel
from rrompy.utilities.base.types import Np1D, paramList, sampList
from rrompy.utilities.base import verbosityDepth
from rrompy.utilities.exception_manager import RROMPyWarning, RROMPyAssert
from rrompy.parameter import checkParameterList
from rrompy.sampling import emptySampleList
__all__ = ['TrainedModelRB']
class TrainedModelRB(TrainedModel):
"""
ROM approximant evaluation for RB approximant.
Attributes:
Data: dictionary with all that can be pickled.
"""
def getApproxReduced(self, mu:paramList) -> sampList:
"""
Evaluate reduced representation of approximant at arbitrary parameter.
Args:
mu: Target parameter.
"""
mus, wasPar = checkParameterList(mu, self.data.npar)
if (not hasattr(self, "lastSolvedAppReduced")
or self.lastSolvedAppReduced != mus):
if self.verbosity >= 5:
verbosityDepth("INIT", ("Computing RB solution at mu = "
"{}.").format(mus),
timestamp = self.timestamp)
thetaAs, thetabs = self.data.thetaAs, self.data.thetabs
ARBs, bRBs = self.data.ARBs, self.data.bRBs
self.uAppReduced = emptySampleList()
self.uAppReduced.reset((ARBs[0].shape[0], len(mu)),
self.data.projMat.dtype)
for i, mu in enumerate(mus):
if self.verbosity >= 10:
verbosityDepth("INIT", ("Assembling reduced model for mu "
"= {}.").format(mu),
timestamp = self.timestamp)
ARBmu = eval(thetaAs[0]) * ARBs[0]
bRBmu = eval(thetabs[0]) * bRBs[0]
for j in range(1, len(ARBs)):
ARBmu += eval(thetaAs[j]) * ARBs[j]
for j in range(1, len(bRBs)):
bRBmu += eval(thetabs[j]) * bRBs[j]
if self.verbosity >= 10:
verbosityDepth("DEL", "Done assembling reduced model.",
timestamp = self.timestamp)
if self.verbosity >= 5:
verbosityDepth("INIT", ("Solving reduced model for mu = "
"{}.").format(mu),
timestamp = self.timestamp)
self.uAppReduced[i] = np.linalg.solve(ARBmu, bRBmu)
if self.verbosity >= 5:
verbosityDepth("DEL", "Done solving reduced model.",
timestamp = self.timestamp)
if self.verbosity >= 5:
verbosityDepth("DEL", "Done computing RB solution.",
timestamp = self.timestamp)
self.lastSolvedAppReduced = mus
if wasPar: return self.uAppReduced[0]
return self.uAppReduced
def getPoles(self) -> Np1D:
"""
Obtain approximant poles.
Returns:
Numpy complex vector of poles.
"""
RROMPyAssert(self.data.npar, 1, "Number of parameters")
RROMPyWarning(("Impossible to compute poles in general affine "
"parameter dependence. Results subject to "
"interpretation/rescaling, or possibly completely "
"wrong."))
ARBs = self.data.ARBs
R = ARBs[0].shape[0]
if len(ARBs) < 2:
return
A = np.eye(R * (len(ARBs) - 1), dtype = np.complex)
B = np.zeros_like(A)
A[: R, : R] = - ARBs[0]
for j in range(len(ARBs) - 1):
Aj = ARBs[j + 1]
B[: R, j * R : (j + 1) * R] = Aj
II = np.arange(R, R * (len(ARBs) - 1))
B[II, II - R] = 1.
return np.power(eigvals(A, B)
+ self.data.mu0(0) ** self.data.rescalingExp,
1. / self.data.rescalingExp)

Event Timeline