Skip to content

Reduced-Basis Emulator

A ReducedBasisEmulator builds the reduced basis for a given interaction and partial wave.

It instantiates a Basis and precomputes the integrals it can with the basis states and operators. When the user asks for an emulated phase shift, the emulator computes the $\hat{\phi}$ coefficients.

ReducedBasisEmulator(interaction, basis, s_0=6 * np.pi, initialize_emulator=True)

A ReducedBasisEmulator (RBE) uses the specified interaction and theta_train to generate solutions to the Schrödinger equation at a specific energy (energy) and partial wave (l).

Using the Galerkin projection method, a linear combination of those solutions (or a PCA basis of them) is found at some arbitrary point in parameter space, theta.

Trains a reduced-basis emulator based on the provided interaction and basis.

Parameters:

Name Type Description Default
interaction Interaction

local interaction

required
basis Basis required
s_0 float

$s$ point where the phase shift is extracted

6 * pi
initialize_emulator

whether to construct matrix elements for emulation now, or leave them for later. Initialization is required for emulation of phase shifts or wavefunctions, but not for the exact phase shifts. This argument is typically always true, unless no emulation is desired (e.g. in ScatteringAmplitudeEmulator.HIFI_solver)

True

Returns:

Name Type Description
rbe ReducedBasisEmulator

$\ell$-specific, reduced basis emulator

Source code in src/rose/reduced_basis_emulator.py
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
def __init__(
    self,
    interaction: Interaction,
    basis: Basis,
    s_0: float = 6 * np.pi,
    initialize_emulator=True,
):
    r"""Trains a reduced-basis emulator based on the provided interaction and basis.

    Parameters:
        interaction (Interaction): local interaction
        basis (Basis): see [Basis documentation](basis.md)
        s_0 (float): $s$ point where the phase shift is extracted
        initialize_emulator : whether to construct matrix elements for emulation now, or
            leave them for later. Initialization is required for emulation of phase shifts or
            wavefunctions, but not for the exact phase shifts. This argument is typically
            always true, unless no emulation is desired (e.g. in
            ScatteringAmplitudeEmulator.HIFI_solver)

    Returns:
        rbe (ReducedBasisEmulator): $\ell$-specific, reduced basis emulator

    """
    self.interaction = interaction
    self.basis = basis
    self.l = self.basis.l

    if isinstance(self.basis, CustomBasis):
        self.basis.solver = basis.solver

    self.s_mesh = np.copy(basis.rho_mesh)
    self.s_0 = s_0

    if initialize_emulator:
        self.initialize_emulator()

coefficients(theta)

Calculated the coefficients in the $\hat{\phi}$ expansion.

Parameters:

Name Type Description Default
theta ndarray

point in parameters space

required

Returns:

Name Type Description
coefficients ndarray

coefficients in $\hat{\phi}$ expansion

Source code in src/rose/reduced_basis_emulator.py
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
def coefficients(self, theta: np.array):
    r"""Calculated the coefficients in the $\hat{\phi}$ expansion.

    Parameters:
        theta (ndarray): point in parameters space

    Returns:
        coefficients (ndarray): coefficients in $\hat{\phi}$ expansion

    """
    invk, beta = self.interaction.coefficients(theta)

    A_utilde = np.einsum("i,ijk", beta, self.A_2)
    A = self.A_13 + A_utilde + invk * self.A_3_coulomb

    b_utilde = beta @ self.b_2
    b = self.b_13 + b_utilde + invk * self.b_3_coulomb

    return np.linalg.solve(A, b)

emulate_phase_shift(theta)

Calculate the emulated phase shift.

Parameters:

Name Type Description Default
theta ndarray

point in parameters space

required

Returns:

Name Type Description
delta float | complex

emulated phase shift (extracted at $s_0$)

Source code in src/rose/reduced_basis_emulator.py
247
248
249
250
251
252
253
254
255
256
257
def emulate_phase_shift(self, theta: np.array):
    r"""Calculate the emulated phase shift.

    Parameters:
        theta (ndarray): point in parameters space

    Returns:
        delta (float | complex): emulated phase shift (extracted at $s_0$)

    """
    return np.log(self.S_matrix_element(theta)) / 2j

emulate_wave_function(theta)

Calculate the emulated wave function.

Parameters:

Name Type Description Default
theta ndarray

point in parameters space

required

Returns:

Name Type Description
wave_function ndarray

$\hat{\phi}$ on the default $s$ mesh

Source code in src/rose/reduced_basis_emulator.py
234
235
236
237
238
239
240
241
242
243
244
245
def emulate_wave_function(self, theta: np.array):
    r"""Calculate the emulated wave function.

    Parameters:
        theta (ndarray): point in parameters space

    Returns:
        wave_function (ndarray): $\hat{\phi}$ on the default $s$ mesh

    """
    x = self.coefficients(theta)
    return self.basis.phi_hat(x)

from_train(interaction, theta_train, n_basis=4, use_svd=True, s_mesh=DEFAULT_RHO_MESH, s_0=6 * np.pi, base_solver=SchroedingerEquation(None)) classmethod

Trains a reduced-basis emulator based on the provided interaction and training space.

Parameters:

Name Type Description Default
interaction Interaction

local interaction

required
theta_train ndarray

training points in parameter space; shape = (n_points, n_parameters)

required
n_basis int

number of basis vectors for $\hat{\phi}$ expansion

4
use_svd bool

Use principal components of training wave functions?

True
s_mesh ndarray

$s$ (or $\rho$) grid on which wave functions are evaluated

DEFAULT_RHO_MESH
s_0 float

$s$ point where the phase shift is extracted

6 * pi
base_solver

the solver used in the basis.

SchroedingerEquation(None)

Returns:

Name Type Description
rbe ReducedBasisEmulator

$\ell$-specific, reduced basis emulator

Source code in src/rose/reduced_basis_emulator.py
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
@classmethod
def from_train(
    cls,
    interaction: Interaction,
    theta_train: np.array,
    n_basis: int = 4,
    use_svd: bool = True,
    s_mesh: np.array = DEFAULT_RHO_MESH,
    s_0: float = 6 * np.pi,
    base_solver: SchroedingerEquation = SchroedingerEquation(None),
):
    r"""Trains a reduced-basis emulator based on the provided interaction and training space.

    Parameters:
        interaction (Interaction): local interaction
        theta_train (ndarray): training points in parameter space;
            shape = (n_points, n_parameters)
        n_basis (int): number of basis vectors for $\hat{\phi}$ expansion
        use_svd (bool): Use principal components of training wave functions?
        s_mesh (ndarray): $s$ (or $\rho$) grid on which wave functions are evaluated
        s_0 (float): $s$ point where the phase shift is extracted
        base_solver : the solver used in the basis.

    Returns:
        rbe (ReducedBasisEmulator): $\ell$-specific, reduced basis emulator

    """

    basis = RelativeBasis(
        base_solver.clone_for_new_interaction(interaction),
        theta_train,
        s_mesh,
        n_basis,
        interaction.ell,
        use_svd,
    )
    return cls(interaction, basis, s_0=s_0)

load(obj, filename) classmethod

Loads a previously trained emulator.

Parameters:

Name Type Description Default
filename string

name of file

required

Returns:

Name Type Description
emulator ReducedBasisEmulator

previously trainined ReducedBasisEmulator

Source code in src/rose/reduced_basis_emulator.py
36
37
38
39
40
41
42
43
44
45
46
47
48
49
@classmethod
def load(obj, filename):
    r"""Loads a previously trained emulator.

    Parameters:
        filename (string): name of file

    Returns:
        emulator (ReducedBasisEmulator): previously trainined `ReducedBasisEmulator`

    """
    with open(filename, "rb") as f:
        rbe = pickle.load(f)
    return rbe

save(filename)

Write the current emulator to file.

Parameters:

Name Type Description Default
filename string

name of the file where the emulator is saved

required

Returns:

Name Type Description
_ None

nothing

Source code in src/rose/reduced_basis_emulator.py
271
272
273
274
275
276
277
278
279
280
281
def save(self, filename):
    r"""Write the current emulator to file.

    Parameters:
        filename (string): name of the file where the emulator is saved

    Returns:
        _ (None): nothing
    """
    with open(filename, "wb") as f:
        pickle.dump(self, f)