simulations.phonons#

Classes:

Phonon(S, force_recalc, **kwargs)

Base class for phonon simulations in a linear chain of masses and springs.

PhononNum(S, force_recalc, **kwargs)

Numerical model to simulate coherent acoustic phonons.

PhononAna(S, force_recalc, **kwargs)

Analytical model to simulate coherent acoustic phonons.

class udkm1Dsim.simulations.phonons.Phonon(S, force_recalc, **kwargs)[source]#

Bases: Simulation

Base class for phonon simulations in a linear chain of masses and springs.

Parameters:
  • S (Structure) – sample to do simulations with.

  • force_recalc (boolean) – force recalculation of results.

Keyword Arguments:
  • save_data (boolean) – true to save simulation results.

  • cache_dir (str) – path to cached data.

  • disp_messages (boolean) – true to display messages from within the simulations.

  • progress_bar (boolean) – enable tqdm progress bar.

  • only_heat (boolean) – true when including only thermal expansion without coherent phonon dynamics.

Attributes:
  • S (Structure) – sample structure to calculate simulations on.

  • force_recalc (boolean) – force recalculation of results.

  • save_data (boolean) – true to save simulation results.

  • cache_dir (str) – path to cached data.

  • disp_messages (boolean) – true to display messages from within the simulations.

  • progress_bar (boolean) – enable tqdm progress bar.

  • only_heat (boolean) – true when including only thermal expansion without coherent phonon dynamics.

Methods:

get_hash(delays, temp_map, delta_temp_map, ...)

Calculates an unique hash given by the delays, and temp_map and delta_temp_map as well as the sample structure hash for relevant lattice parameters.

get_all_strains_per_unique_layer(strain_map)

Determines all values of the strain per unique layer.

get_reduced_strains_per_unique_layer(strain_map)

Calculates all strains per unique layer that are given by the input strain_map, but with a reduced number.

check_temp_maps(temp_map, delta_temp_map, delays)

Check temperature profiles for correct dimensions.

calc_sticks_from_temp_map(temp_map, ...)

Calculates the sticks to insert into the layer springs which model the external force (thermal stress).

conv_with_function(y, x, handle)

Convolutes the array \(y(x)\) with a function given by the handle on the argument array \(x\).

disp_message(message)

Wrapper to display messages for that class.

save(full_filename, data, *args)

Save data to file.

get_hash(delays, temp_map, delta_temp_map, **kwargs)[source]#

Calculates an unique hash given by the delays, and temp_map and delta_temp_map as well as the sample structure hash for relevant lattice parameters.

Parameters:
  • delays (ndarray[float]) – delay grid for the simulation.

  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

  • **kwargs (float, optional) – optional parameters.

Returns:

hash (str) – unique hash.

get_all_strains_per_unique_layer(strain_map)[source]#

Determines all values of the strain per unique layer.

Parameters:

strain_map (ndarray[float]) – spatio-temporal strain profile.

Returns:

strains (list[ndarray[float]] – all strains per unique layer.

get_reduced_strains_per_unique_layer(strain_map, N=100)[source]#

Calculates all strains per unique layer that are given by the input strain_map, but with a reduced number. The reduction is done by equally spacing the strains between the min and max strain with a given number \(N\), which can be also an array of the \(len(N) = L\), where \(L\) is the number of unique layers.

Parameters:
  • strain_map (ndarray[float]) – spatio-temporal strain profile.

  • N (int, optional) – number of reduced strains. Defaults to 100.

Returns:

strains (list[ndarray[float]] – reduced strains per unique layer.

check_temp_maps(temp_map, delta_temp_map, delays)[source]#

Check temperature profiles for correct dimensions.

Parameters:
  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

  • delays (ndarray[float]) – delay grid for the simulation.

Returns:

(tuple)

  • temp_map (ndarray[float]) - checked spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) - checked spatio-temporal differential temperature profile.

calc_sticks_from_temp_map(temp_map, delta_temp_map)[source]#

Calculates the sticks to insert into the layer springs which model the external force (thermal stress). The length of \(l_i\) of the \(i\)-th spacer stick is calculated from the temperature-dependent linear thermal expansion \(\alpha(T)\) of the layer:

\[\alpha(T) = \frac{1}{L} \frac{d L}{d T}\]

which results after integration in

\[l = \Delta L = L_1 \exp(A(T_2) - A(T_1)) - L_1\]

where \(A(T)\) is the integrated lin. therm. expansion coefficient in respect to the temperature \(T\). The indices 1 and 2 indicate the initial and final state.

Parameters:
  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

Returns:

(tuple)

  • sticks (ndarray[float]) - summed spacer sticks.

  • sticks_sub_systems (ndarray[float]) - spacer sticks per sub-system.

static conv_with_function(y, x, handle)#

Convolutes the array \(y(x)\) with a function given by the handle on the argument array \(x\).

Parameters:
  • y (ndarray[float]) – y data.

  • x (ndarray[float]) – x data.

  • handle (@lamdba) – convolution function.

Returns:

y_conv (ndarray[float]) – convoluted data.

disp_message(message)#

Wrapper to display messages for that class.

Parameters:

message (str) – message to display.

save(full_filename, data, *args)#

Save data to file. The variable name can be handed as variable argument.

Parameters:
  • full_filename (str) – full file name to data file.

  • data (ndarray) – actual data to save.

  • *args (str, optional) – variable name within the data file.

class udkm1Dsim.simulations.phonons.PhononNum(S, force_recalc, **kwargs)[source]#

Bases: Phonon

Numerical model to simulate coherent acoustic phonons.

Parameters:
  • S (Structure) – sample to do simulations with.

  • force_recalc (boolean) – force recalculation of results.

Keyword Arguments:
  • save_data (boolean) – true to save simulation results.

  • cache_dir (str) – path to cached data.

  • disp_messages (boolean) – true to display messages from within the simulations.

  • progress_bar (boolean) – enable tqdm progress bar.

  • only_heat (boolean) – true when including only thermal expansion without coherent phonon dynamics.

Attributes:
  • S (Structure) – sample structure to calculate simulations on.

  • force_recalc (boolean) – force recalculation of results.

  • save_data (boolean) – true to save simulation results.

  • cache_dir (str) – path to cached data.

  • disp_messages (boolean) – true to display messages from within the simulations.

  • progress_bar (boolean) – enable tqdm progress bar.

  • only_heat (boolean) – true when including only thermal expansion without coherent phonon dynamics.

  • ode_options (dict) – options for scipy solve_ivp ode solver.

References

Methods:

get_strain_map(delays, temp_map, delta_temp_map)

Returns a strain profile for the sample structure for given temperature profile.

calc_strain_map(delays, temp_map, delta_temp_map)

Calculates the strain_map of the sample structure for a given temp_map and delta_temp_map and delay array.

ode_func(t, X, delays, force_from_heat, ...)

Provides the according ode function for the ode solver which has to be solved.

calc_force_from_spring(d_X1, d_X2, spring_consts)

Calculates the force \(F_i^{spring}\) acting on each mass due to the displacement between the left and right site of that mass.

calc_force_from_heat(sticks, spring_consts)

Calculates the force acting on each mass due to the heat expansion, which is modeled by spacer sticks.

calc_force_from_damping(v, damping, masses)

Calculates the force acting on each mass in a linear spring due to damping (\(\gamma_i\)) according to the shift velocity difference \(v_{i}-v_{i-1}\) with \(v_i(t) = \dot{x}_i(t)\):

calc_sticks_from_temp_map(temp_map, ...)

Calculates the sticks to insert into the layer springs which model the external force (thermal stress).

check_temp_maps(temp_map, delta_temp_map, delays)

Check temperature profiles for correct dimensions.

conv_with_function(y, x, handle)

Convolutes the array \(y(x)\) with a function given by the handle on the argument array \(x\).

disp_message(message)

Wrapper to display messages for that class.

get_all_strains_per_unique_layer(strain_map)

Determines all values of the strain per unique layer.

get_hash(delays, temp_map, delta_temp_map, ...)

Calculates an unique hash given by the delays, and temp_map and delta_temp_map as well as the sample structure hash for relevant lattice parameters.

get_reduced_strains_per_unique_layer(strain_map)

Calculates all strains per unique layer that are given by the input strain_map, but with a reduced number.

save(full_filename, data, *args)

Save data to file.

get_strain_map(delays, temp_map, delta_temp_map)[source]#

Returns a strain profile for the sample structure for given temperature profile. The result can be saved using an unique hash of the sample and the simulation parameters in order to reuse it.

Parameters:
  • delays (ndarray[Quantity]) – delays range of simulation [s].

  • temp_map (ndarray[float]) – spatio-temporal temperature map.

  • delta_temp_map (ndarray[float]) – spatio-temporal differential

  • map. (temperature) –

Returns:

strain_map (ndarray[float]) – spatio-temporal strain profile.

calc_strain_map(delays, temp_map, delta_temp_map)[source]#

Calculates the strain_map of the sample structure for a given temp_map and delta_temp_map and delay array. Further details are given in Ref. [7]. The coupled differential equations are solved for each oscillator in a linear chain of masses and springs:

\[\begin{split}m_i\ddot{x}_i = & -k_i(x_i-x_{i-1})-k_{i+1}(x_i-x_{i+1}) \\ & + m_i\gamma_i(\dot{x}_i-\dot{x}_{i-1}) + F_i^{heat}(t)\end{split}\]

where \(x_i(t) = z_{i}(t) - z_i^0\) is the shift of each layer. \(m_i\) is the mass and \(k_i = m_i\, v_i^2/c_i^2\) is the spring constant of each layer. Furthermore an empirical damping term \(F_i^{damp} = \gamma_i(\dot{x}_i-\dot{x}_{i-1})\) is introduced and the external force (thermal stress) \(F_i^{heat}(t)\). The thermal stresses are modelled as spacer sticks which are calculated from the linear thermal expansion coefficients. The equation of motion can be reformulated as:

\[m_i\ddot{x}_i = F_i^{spring} + F_i^{damp} + F_i^{heat}(t)\]

The numerical solution also allows for non-harmonic inter-atomic potentials of up to the order \(M\). Accordingly \(k_i = (k_i^1 \ldots k_i^{M-1})\) can be an array accounting for higher orders of the potential which is in the harmonic case purely quadratic (\(k_i = k_i^1\)). The resulting force from the displacement of the springs

\[F_i^{spring} = -k_i(x_i-x_{i-1})-k_{i+1}(x_i-x_{i+1})\]

includes:

\[k_i(x_i-x_{i-1}) = \sum_{j=1}^{M-1} k_i^j (x_i-x_{i-1})^j\]
Parameters:
  • delays (ndarray[Quantity]) – delays range of simulation [s].

  • temp_map (ndarray[float]) – spatio-temporal temperature map.

  • delta_temp_map (ndarray[float]) – spatio-temporal differential temperature map.

Returns:

(tuple)

  • strain_map (ndarray[float]) - spatio-temporal strain profile.

  • sticks_sub_systems (ndarray[float]) - spacer sticks per sub-system.

  • velocities (ndarray[float]) - spatio-temporal velocity profile.

static ode_func(t, X, delays, force_from_heat, damping, spring_consts, masses, L, pbar=None, state=None)[source]#

Provides the according ode function for the ode solver which has to be solved. The ode function has the input \(t\) and \(X(t)\) and calculates the temporal derivative \(\dot{X}(t)\) where the vector

\[X(t) = [x(t) \; \dot{x}(t)] \quad \mbox{and } \quad \dot{X}(t) = [\dot{x}(t) \; \ddot{x}(t)] .\]

\(x(t)\) is the actual shift of each layer.

Parameters:
  • t (ndarray[float]) – internal time steps of the ode solver.

  • X (ndarray[float]) – internal variable of the ode solver.

  • delays (ndarray[float]) – delays range of simulation [s].

  • force_from_heat (ndarray[float]) – force due to thermal expansion.

  • damping (ndarray[float]) – phonon damping.

  • spring_consts (ndarray[float]) – spring constants of masses.

  • masses (ndarray[float]) – masses of layers.

  • L (int) – number of layers.

  • pbar (tqdm) – tqdm progressbar.

  • state (list[float]) – state variables for progress bar.

Returns:

X_prime (ndarray[float]) – velocities and accelerations of masses.

static calc_force_from_spring(d_X1, d_X2, spring_consts)[source]#

Calculates the force \(F_i^{spring}\) acting on each mass due to the displacement between the left and right site of that mass.

\[F_i^{spring} = -k_i(x_i-x_{i-1})-k_{i+1}(x_i-x_{i+1})\]

We introduce-higher order inter-atomic potentials by

\[k_i(x_i-x_{i-1}) = \sum_{j=1}^{M} k_i^j (x_i-x_{i-1})^j\]

where \(M\) is the order of the spring constants.

Parameters:
  • d_X1 (ndarray[float]) – left displacements.

  • d_X2 (ndarray[float]) – right displacements.

  • spring_consts (ndarray[float]) – spring constants of masses.

Returns:

F (ndarray[float]) – force from springs.

static calc_force_from_heat(sticks, spring_consts)[source]#

Calculates the force acting on each mass due to the heat expansion, which is modeled by spacer sticks.

Parameters:
  • sticks (ndarray[float]) – spacer sticks.

  • spring_consts (ndarray[float]) – spring constants of masses.

Returns:

F (ndarray[float]) – force from thermal expansion.

static calc_force_from_damping(v, damping, masses)[source]#

Calculates the force acting on each mass in a linear spring due to damping (\(\gamma_i\)) according to the shift velocity difference \(v_{i}-v_{i-1}\) with \(v_i(t) = \dot{x}_i(t)\):

\[F_i^{damp} = \gamma_i(\dot{x}_i-\dot{x}_{i-1})\]
Parameters:
  • v (ndarray[float]) – velocity of masses.

  • damping (ndarray[float]) – phonon damping.

  • masses (ndarray[float]) – masses.

Returns:

F (ndarray[float]) – force from damping.

calc_sticks_from_temp_map(temp_map, delta_temp_map)#

Calculates the sticks to insert into the layer springs which model the external force (thermal stress). The length of \(l_i\) of the \(i\)-th spacer stick is calculated from the temperature-dependent linear thermal expansion \(\alpha(T)\) of the layer:

\[\alpha(T) = \frac{1}{L} \frac{d L}{d T}\]

which results after integration in

\[l = \Delta L = L_1 \exp(A(T_2) - A(T_1)) - L_1\]

where \(A(T)\) is the integrated lin. therm. expansion coefficient in respect to the temperature \(T\). The indices 1 and 2 indicate the initial and final state.

Parameters:
  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

Returns:

(tuple)

  • sticks (ndarray[float]) - summed spacer sticks.

  • sticks_sub_systems (ndarray[float]) - spacer sticks per sub-system.

check_temp_maps(temp_map, delta_temp_map, delays)#

Check temperature profiles for correct dimensions.

Parameters:
  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

  • delays (ndarray[float]) – delay grid for the simulation.

Returns:

(tuple)

  • temp_map (ndarray[float]) - checked spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) - checked spatio-temporal differential temperature profile.

static conv_with_function(y, x, handle)#

Convolutes the array \(y(x)\) with a function given by the handle on the argument array \(x\).

Parameters:
  • y (ndarray[float]) – y data.

  • x (ndarray[float]) – x data.

  • handle (@lamdba) – convolution function.

Returns:

y_conv (ndarray[float]) – convoluted data.

disp_message(message)#

Wrapper to display messages for that class.

Parameters:

message (str) – message to display.

get_all_strains_per_unique_layer(strain_map)#

Determines all values of the strain per unique layer.

Parameters:

strain_map (ndarray[float]) – spatio-temporal strain profile.

Returns:

strains (list[ndarray[float]] – all strains per unique layer.

get_hash(delays, temp_map, delta_temp_map, **kwargs)#

Calculates an unique hash given by the delays, and temp_map and delta_temp_map as well as the sample structure hash for relevant lattice parameters.

Parameters:
  • delays (ndarray[float]) – delay grid for the simulation.

  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

  • **kwargs (float, optional) – optional parameters.

Returns:

hash (str) – unique hash.

get_reduced_strains_per_unique_layer(strain_map, N=100)#

Calculates all strains per unique layer that are given by the input strain_map, but with a reduced number. The reduction is done by equally spacing the strains between the min and max strain with a given number \(N\), which can be also an array of the \(len(N) = L\), where \(L\) is the number of unique layers.

Parameters:
  • strain_map (ndarray[float]) – spatio-temporal strain profile.

  • N (int, optional) – number of reduced strains. Defaults to 100.

Returns:

strains (list[ndarray[float]] – reduced strains per unique layer.

save(full_filename, data, *args)#

Save data to file. The variable name can be handed as variable argument.

Parameters:
  • full_filename (str) – full file name to data file.

  • data (ndarray) – actual data to save.

  • *args (str, optional) – variable name within the data file.

class udkm1Dsim.simulations.phonons.PhononAna(S, force_recalc, **kwargs)[source]#

Bases: Phonon

Analytical model to simulate coherent acoustic phonons.

Parameters:
  • S (Structure) – sample to do simulations with.

  • force_recalc (boolean) – force recalculation of results.

Keyword Arguments:
  • save_data (boolean) – true to save simulation results.

  • cache_dir (str) – path to cached data.

  • disp_messages (boolean) – true to display messages from within the simulations.

  • progress_bar (boolean) – enable tqdm progress bar.

  • only_heat (boolean) – true when including only thermal expansion without coherent phonon dynamics.

Attributes:
  • S (Structure) – sample structure to calculate simulations on.

  • force_recalc (boolean) – force recalculation of results.

  • save_data (boolean) – true to save simulation results.

  • cache_dir (str) – path to cached data.

  • disp_messages (boolean) – true to display messages from within the simulations.

  • progress_bar (boolean) – enable tqdm progress bar.

  • only_heat (boolean) – true when including only thermal expansion without coherent phonon dynamics.

References

Methods:

get_strain_map(delays, temp_map, delta_temp_map)

Returns a strain profile for the sample structure for given temperature profile.

calc_strain_map(delays, temp_map, delta_temp_map)

Calculates the strain_map of the sample structure for a given temp_map and delta_temp_map and delay array.

solve_eigenproblem()

Creates the real and symmetric \(K\) matrix (\(L \times L\)) of spring constants \(k_i\) and masses \(m_i\) and calculates the eigenvectors \(\Xi_j\) and eigenfrequencies \(\omega_j\) for the matrix which are used to calculate the strain_map of the structure.

get_energy_per_eigenmode(A, B)

Returns the sorted energy per Eigenmode of the coherent phonons of the 1D sample.

calc_sticks_from_temp_map(temp_map, ...)

Calculates the sticks to insert into the layer springs which model the external force (thermal stress).

check_temp_maps(temp_map, delta_temp_map, delays)

Check temperature profiles for correct dimensions.

conv_with_function(y, x, handle)

Convolutes the array \(y(x)\) with a function given by the handle on the argument array \(x\).

disp_message(message)

Wrapper to display messages for that class.

get_all_strains_per_unique_layer(strain_map)

Determines all values of the strain per unique layer.

get_hash(delays, temp_map, delta_temp_map, ...)

Calculates an unique hash given by the delays, and temp_map and delta_temp_map as well as the sample structure hash for relevant lattice parameters.

get_reduced_strains_per_unique_layer(strain_map)

Calculates all strains per unique layer that are given by the input strain_map, but with a reduced number.

save(full_filename, data, *args)

Save data to file.

get_strain_map(delays, temp_map, delta_temp_map)[source]#

Returns a strain profile for the sample structure for given temperature profile. The result can be saved using an unique hash of the sample and the simulation parameters in order to reuse it.

Parameters:
  • delays (ndarray[Quantity]) – delays range of simulation [s].

  • temp_map (ndarray[float]) – spatio-temporal temperature map.

  • delta_temp_map (ndarray[float]) – spatio-temporal differential

  • map. (temperature) –

Returns:

(tuple)

  • strain_map (ndarray[float]) - spatio-temporal strain profile.

  • A (ndarray[float]) - coefficient vector A of general solution.

  • B (ndarray[float]) - coefficient vector B of general solution.

calc_strain_map(delays, temp_map, delta_temp_map)[source]#

Calculates the strain_map of the sample structure for a given temp_map and delta_temp_map and delay array. Further details are given in Ref. [8]. Within the linear chain of \(N\) masses (\(m_i\)) at position \(z_i\) coupled with spring constants \(k_i\) one can formulate the differential equation of motion as follow:

\[m_i\ddot{x}_i = -k_i(x_i-x_{i-1})-k_{i+1}(x_i-x_{i+1}) + F_i^{heat}(t)\]

Since we only consider nearest-neighbor interaction one can write:

\[\ddot{x}_i = \sum_{n=1}^N \kappa_{i,n} x_n = \Delta_i(t)\]

Here \(x_i(t) = z_i(t)-z_i^0\) is the shift of each layer, \(F_i^{heat}(t)\) is the external force (thermal stress) of each layer and \(\kappa_{i,i} = -(k_i + k_{i+1})/m_i\), and \(\kappa_{i,i+1} = \kappa_{i+1,i} = k_{i+1}/m_i\).

\(k_i = m_i\, v_i^2/c_i^2\) is the spring constant and \(c_i\) and \(v_i\) are the thickness and longitudinal sound velocity of each layer respectively. One can rewrite the homogeneous differential equation in matrix form to obtain the general solution

\[\frac{d^2}{dt^2} X = K \, X\]

Here \(X = (x_1 \ldots x_N)\) and \(K\) is the tri-diagonal matrix of \(\kappa\) which is real and symmetric. The differential equation can be solved with the ansatz:

\[X(t) = \sum_j \Xi_j \, (A_j \cos(\omega_j \, t) + B_j \sin(\omega_j \, t))\]

where \(\Xi_j = (\xi_1^j \ldots \xi_N^j)\) are the eigenvectors of the matrix \(K\). Thus by solving the Eigenproblem for \(K\) one gets the eigenvecotrs \(\Xi_j\) and the eigenfrequencies \(\omega_j\). From the initial conditions

\[X(0) = \sum_j \Xi_j \, A_j = \Xi \, A \qquad V(0) = \dot{X}(0) = \sum_j \Xi_j \, \omega_j\, B_j = \Xi \, \omega \, B\]

one can determine the real coefficient vectors \(A\) and \(B\) in order to calculate \(X(t)\) and \(V(t)\) using the ansatz:

\[A = \Xi \setminus X(0) \qquad B = (\Xi \setminus V(0)) / \omega\]

The external force is implemented as spacer sticks which are inserted into the springs and hence the layers have a new equilibrium positions \(z_i(\infty) = z_i^\infty\). Thus we can do a coordination transformation:

\[z_i(t) = z_i^0 + x_i(t) = z_i^\infty + x_i^\infty(t)\]

and

\[x_i^\infty(t) = z_i^0 - z_i^\infty + x_i(t)\]

with the initial condition \(x_i(0) = 0\) it becomes

\[x_i^\infty(0) = z_i^0 - z_i^\infty = \sum_{j = i+1}^N l_j\]

\(x_i^\infty(0)\) is the new initial condition after the excitation where \(l_i\) is the length of the \(i\)-th spacer stick. The spacer sticks are calculated from the temperature change and the linear thermal expansion coefficients. The actual strain \(\epsilon_i(t)\) of each layer is calculates as follows:

\[\epsilon_i(t) = [ \Delta x_i(t) + l_i) ] / c_i\]

with \(\Delta x_i = x_i - x_{i-1}\). The stick \(l_i\) have to be added here, because \(x_i\) has been transformed into the new coordinate system \(x_i^\infty\).

Parameters:
  • delays (ndarray[Quantity]) – delays range of simulation [s].

  • temp_map (ndarray[float]) – spatio-temporal temperature map.

  • delta_temp_map (ndarray[float]) – spatio-temporal differential temperature map.

Returns:

(tuple)

  • strain_map (ndarray[float]) - spatio-temporal strain profile.

  • A (ndarray[float]) - coefficient vector A of general solution.

  • B (ndarray[float]) - coefficient vector B of general solution.

solve_eigenproblem()[source]#

Creates the real and symmetric \(K\) matrix (\(L \times L\)) of spring constants \(k_i\) and masses \(m_i\) and calculates the eigenvectors \(\Xi_j\) and eigenfrequencies \(\omega_j\) for the matrix which are used to calculate the strain_map of the structure. If the result has been save to file, load it from there.

Returns:

(tuple)

  • Xi (ndarray[float]) - eigenvectors.

  • omega (ndarray[float]) - eigenfrequencies.

get_energy_per_eigenmode(A, B)[source]#

Returns the sorted energy per Eigenmode of the coherent phonons of the 1D sample.

\[E_j = \frac{1}{2} (A^2_j + B^2_j)\, \omega_j^2\, m_j \, \| \Xi_j\|^2\]

Frequencies are in [Hz] and energy per mode in [J].

Parameters:
  • A (ndarray[float]) – coefficient vector A of general solution.

  • B (ndarray[float]) – coefficient vector B of general solution.

Returns:

(tuple)

  • omega (ndarray[float]) - eigenfrequencies.

  • E (ndarray[float]) - energy per eigenmode.

calc_sticks_from_temp_map(temp_map, delta_temp_map)#

Calculates the sticks to insert into the layer springs which model the external force (thermal stress). The length of \(l_i\) of the \(i\)-th spacer stick is calculated from the temperature-dependent linear thermal expansion \(\alpha(T)\) of the layer:

\[\alpha(T) = \frac{1}{L} \frac{d L}{d T}\]

which results after integration in

\[l = \Delta L = L_1 \exp(A(T_2) - A(T_1)) - L_1\]

where \(A(T)\) is the integrated lin. therm. expansion coefficient in respect to the temperature \(T\). The indices 1 and 2 indicate the initial and final state.

Parameters:
  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

Returns:

(tuple)

  • sticks (ndarray[float]) - summed spacer sticks.

  • sticks_sub_systems (ndarray[float]) - spacer sticks per sub-system.

check_temp_maps(temp_map, delta_temp_map, delays)#

Check temperature profiles for correct dimensions.

Parameters:
  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

  • delays (ndarray[float]) – delay grid for the simulation.

Returns:

(tuple)

  • temp_map (ndarray[float]) - checked spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) - checked spatio-temporal differential temperature profile.

static conv_with_function(y, x, handle)#

Convolutes the array \(y(x)\) with a function given by the handle on the argument array \(x\).

Parameters:
  • y (ndarray[float]) – y data.

  • x (ndarray[float]) – x data.

  • handle (@lamdba) – convolution function.

Returns:

y_conv (ndarray[float]) – convoluted data.

disp_message(message)#

Wrapper to display messages for that class.

Parameters:

message (str) – message to display.

get_all_strains_per_unique_layer(strain_map)#

Determines all values of the strain per unique layer.

Parameters:

strain_map (ndarray[float]) – spatio-temporal strain profile.

Returns:

strains (list[ndarray[float]] – all strains per unique layer.

get_hash(delays, temp_map, delta_temp_map, **kwargs)#

Calculates an unique hash given by the delays, and temp_map and delta_temp_map as well as the sample structure hash for relevant lattice parameters.

Parameters:
  • delays (ndarray[float]) – delay grid for the simulation.

  • temp_map (ndarray[float]) – spatio-temporal temperature profile.

  • delta_temp_map (ndarray[float]) – spatio-temporal temperature difference profile.

  • **kwargs (float, optional) – optional parameters.

Returns:

hash (str) – unique hash.

get_reduced_strains_per_unique_layer(strain_map, N=100)#

Calculates all strains per unique layer that are given by the input strain_map, but with a reduced number. The reduction is done by equally spacing the strains between the min and max strain with a given number \(N\), which can be also an array of the \(len(N) = L\), where \(L\) is the number of unique layers.

Parameters:
  • strain_map (ndarray[float]) – spatio-temporal strain profile.

  • N (int, optional) – number of reduced strains. Defaults to 100.

Returns:

strains (list[ndarray[float]] – reduced strains per unique layer.

save(full_filename, data, *args)#

Save data to file. The variable name can be handed as variable argument.

Parameters:
  • full_filename (str) – full file name to data file.

  • data (ndarray) – actual data to save.

  • *args (str, optional) – variable name within the data file.