Skip to content

High-Fidelity Solver

SchroedingerEquation is a high-fidelity (HF), Schrödinger-equation solver for local, complex interactions.

By default, rose will provide HF solution using scipy.integrate.solve_ivp. For details about providing your own solutions, see Basis documentation.

SchroedingerEquation(interaction, rk_tols=[1e-09, 1e-09], s_0=None, domain=None)

Solver for the single-channel, reduced, radial Schrödinger equation using scipy.integrate https://docs.scipy.org/doc/scipy/reference/generated/scipy.integrate.solve_ivp.html

All other solvers for the single-channel, reduced, radial Schrödinger equation inherit from this class

Solves the Shrödinger equation for local, complex potentials.

Parameters:

Name Type Description Default
interaction Interaction required
rk_tols list

2-element list of numbers specifying tolerances for the Runge-Kutta solver: the relative tolerance and the absolute tolerance

[1e-09, 1e-09]
s_0 float)
None
domain ndarray)
None

Returns:

Name Type Description
solver SchroedingerEquation

instance of SchroedingerEquation

Source code in src/rose/schroedinger.py
 47
 48
 49
 50
 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
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
def __init__(
    self,
    interaction: Interaction,
    rk_tols: list = [1e-9, 1e-9],
    s_0=None,
    domain=None,
):
    r"""Solves the Shrödinger equation for local, complex potentials.

    Parameters:
        interaction (Interaction): See [Interaction documentation](interaction.md).
        rk_tols (list): 2-element list of numbers specifying tolerances for the
            Runge-Kutta solver: the relative tolerance and the  absolute tolerance
        s_0 (float) :
        domain (ndarray) :

    Returns:
        solver (SchroedingerEquation): instance of `SchroedingerEquation`

    """
    if domain is None:
        domain = np.array([self.DEFAULT_S_MIN, self.DEFAULT_S_MAX])

    if s_0 is None:
        s_0 = domain[1] - np.pi
        assert s_0 > 0

    assert domain[1] > s_0

    self.s_0 = s_0
    self.domain = domain.copy()
    self.interaction = interaction
    self.rk_tols = rk_tols

    if self.interaction is not None:
        if self.interaction.k_c == 0:
            self.eta = 0
        else:
            # There is Coulomb, but we haven't (yet) worked out how to emulate
            # across energies, so we can precompute H+ and H- stuff.
            self.eta = self.interaction.eta(None)

        self.Hm = H_minus(self.s_0, self.interaction.ell, self.eta)
        self.Hp = H_plus(self.s_0, self.interaction.ell, self.eta)
        self.Hmp = H_minus_prime(self.s_0, self.interaction.ell, self.eta)
        self.Hpp = H_plus_prime(self.s_0, self.interaction.ell, self.eta)

        self.domain[0], self.init_cond = self.initial_conditions(
            self.eta,
            self.PHI_THRESHOLD,
            self.interaction.ell,
            self.domain[0],
        )

        assert self.domain[0] < self.s_0

delta(alpha, **kwargs)

Calculates the $\ell$-th partial wave phase shift

Parameters:

Name Type Description Default
alpha ndarray

parameter vector

required

Returns:

Name Type Description
delta float

phase shift extracted from the reduced, radial wave function

Source code in src/rose/schroedinger.py
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
def delta(
    self,
    alpha: np.array,
    **kwargs,
):
    r"""Calculates the $\ell$-th partial wave phase shift

    Parameters:
        alpha (ndarray): parameter vector

    Returns:
        delta (float): phase shift extracted from the reduced, radial
            wave function

    """
    sl = self.smatrix(alpha, **kwargs)

    return np.log(sl) / 2j

initial_conditions(eta, phi_threshold, l, rho_0=None)

Returns:

Type Description

initial_conditions (tuple) : initial conditions [phi, phi'] at rho_0

Parameters:

Name Type Description Default
eta float

sommerfield param

required
phi_threshold float

minimum $\phi$ value; The wave function is considered zero below this value.

required
rho_0 float

starting point for the solver

None
Source code in src/rose/schroedinger.py
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
def initial_conditions(self, eta: float, phi_threshold: float, l: int, rho_0=None):
    r"""
    Returns:
        initial_conditions (tuple) : initial conditions [phi, phi'] at rho_0

    Parameters:
        eta (float): sommerfield param
        phi_threshold (float): minimum $\phi$ value; The wave function is
            considered zero below this value.
        rho_0 (float): starting point for the solver
    """

    C_l = Gamow_factor(l, eta)
    min_rho_0 = (phi_threshold / C_l) ** (1 / (l + 1))

    if min_rho_0 > rho_0:
        rho_0 = min_rho_0

    phi_0 = C_l * rho_0 ** (l + 1)
    phi_prime_0 = C_l * (l + 1) * rho_0**l

    if self.interaction.is_complex:
        initial_conditions = np.array([phi_0 * (1 + 0j), phi_prime_0 * (1 + 0j)])
    else:
        initial_conditions = np.array([phi_0, phi_prime_0])

    return rho_0, initial_conditions

phi(alpha, s_mesh, **kwargs)

Computes the reduced, radial wave function $\phi$ (or $u$) on s_mesh using the Runge-Kutta method

Parameters:

Name Type Description Default
alpha ndarray

parameter vector

required
s_mesh ndarray

values of $s$ at which $\phi$ is calculated

required
kwargs dict)

passed to scipy.integrate.solve_ivp

{}

Returns:

Name Type Description
phi ndarray

reduced, radial wave function

Source code in src/rose/schroedinger.py
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
def phi(
    self,
    alpha: np.array,
    s_mesh: np.array,
    **kwargs,
):
    r"""Computes the reduced, radial wave function $\phi$ (or $u$) on `s_mesh` using the
    Runge-Kutta method

    Parameters:
        alpha (ndarray): parameter vector
        s_mesh (ndarray): values of $s$ at which $\phi$ is calculated
        kwargs (dict) : passed to scipy.integrate.solve_ivp

    Returns:
        phi (ndarray): reduced, radial wave function

    """
    solution = np.zeros_like(s_mesh, dtype=np.complex128)
    mask = np.logical_and(s_mesh >= self.domain[0], s_mesh <= self.domain[1])
    phi = self.solve_se(alpha, **kwargs)
    solution[mask] = phi(s_mesh)[0][mask]
    return solution

rmatrix(alpha, **kwargs)

Calculates the $\ell$-th partial wave R-matrix element at the specified energy. using the Runge-Kutta method for integrating the Radial SE. kwargs are passed to solve_se.

Parameters:

Name Type Description Default
alpha ndarray

parameter vector

required
kwargs dict)

passed to scipy.integrate.solve_ivp

{}

Returns:

Type Description

rl (float) : r-matrix element, or logarithmic derivative of wavefunction at the channel radius; s_0

Source code in src/rose/schroedinger.py
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
def rmatrix(
    self,
    alpha: np.array,
    **kwargs,
):
    r"""Calculates the $\ell$-th partial wave R-matrix element at the specified energy.
        using the Runge-Kutta method for integrating the Radial SE. kwargs are passed to
        solve_se.

    Parameters:
        alpha (ndarray): parameter vector
        kwargs (dict) : passed to scipy.integrate.solve_ivp

    Returns:
        rl (float)  : r-matrix element, or logarithmic derivative of wavefunction at the channel
            radius; s_0

    """
    solution = self.solve_se(alpha, **kwargs)
    u = solution(self.s_0)
    rl = 1 / self.s_0 * (u[0] / u[1])
    return rl

solve_se(alpha, **kwargs)

Solves the reduced, radial Schrödinger equation using the builtin in Runge-Kutta solver in scipy.integrate.solve_ivp

Parameters:

Name Type Description Default
alpha ndarray

parameter vector

required
kwargs dict)

passed to scipy.integrate.solve_ivp

{}

Returns:

Type Description

sol (scipy.integrate.OdeSolution) : the radial wavefunction

and its first derivative; see https://docs.scipy.org/doc/scipy/reference/generated/scipy.integrate.OdeSolution.html#scipy.integrate.OdeSolution

Source code in src/rose/schroedinger.py
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
def solve_se(
    self,
    alpha: np.array,
    **kwargs,
):
    r"""Solves the reduced, radial Schrödinger equation using the builtin in Runge-Kutta
        solver in scipy.integrate.solve_ivp

    Parameters:
        alpha (ndarray): parameter vector
        kwargs (dict) : passed to scipy.integrate.solve_ivp

    Returns:
        sol  (scipy.integrate.OdeSolution) : the radial wavefunction
        and its first derivative; see https://docs.scipy.org/doc/scipy/reference/generated/scipy.integrate.OdeSolution.html#scipy.integrate.OdeSolution
    """

    args = self.interaction.bundle_gcoeff_args(alpha)
    sol = solve_ivp(
        lambda s, phi: np.array(
            [
                phi[1],
                -1 * g_coeff(s, *args) * phi[0],
            ]
        ),
        self.domain,
        self.init_cond,
        rtol=self.rk_tols[0],
        atol=self.rk_tols[1],
        dense_output=True,
        **kwargs,
    )

    return sol.sol