Skip to content

Commit

Permalink
Extend formal integral to line_interaction_type macroatom (#856)
Browse files Browse the repository at this point in the history
* Extend formal integral to line_interaction_type macroatom

* Simplify indexing
  • Loading branch information
chvogl authored and wkerzendorf committed Jul 31, 2018
1 parent 4622143 commit 29f0193
Show file tree
Hide file tree
Showing 2 changed files with 44 additions and 1 deletion.
10 changes: 10 additions & 0 deletions tardis/io/atomic.py
Original file line number Diff line number Diff line change
Expand Up @@ -367,10 +367,20 @@ def prepare_atom_data(
self.macro_atom_data['destination_level_number']
])

tmp_macro_source_level_idx = pd.MultiIndex.from_arrays([
self.macro_atom_data['atomic_number'],
self.macro_atom_data['ion_number'],
self.macro_atom_data['source_level_number']
])

self.macro_atom_data.loc[:, 'destination_level_idx'] = self.macro_atom_references.loc[
tmp_macro_destination_level_idx, "references_idx"
].astype(np.int64).values # why it is named `destination_level_idx` ?! It is reference index

self.macro_atom_data.loc[:, 'source_level_idx'] = self.macro_atom_references.loc[
tmp_macro_source_level_idx, "references_idx"
].astype(np.int64).values

elif line_interaction_type == 'downbranch':
# Sets all the destination levels to -1 to indicate that they
# are not used in downbranch calculations
Expand Down
35 changes: 34 additions & 1 deletion tardis/montecarlo/formal_integral.py
Original file line number Diff line number Diff line change
@@ -1,6 +1,7 @@
import warnings
import numpy as np
import pandas as pd
import scipy.sparse as sp
from astropy import units as u
from astropy import constants as const

Expand Down Expand Up @@ -44,10 +45,11 @@ def raise_or_return(message):
'FormalIntegrator.'
)

if not self.runner.line_interaction_type == 'downbranch':
if not self.runner.line_interaction_type in ['downbranch', 'macroatom']:
return raise_or_return(
'The FormalIntegrator currently only works for '
'line_interaction_type == "downbranch"'
'and line_interaction_type == "macroatom"'
)

return True
Expand Down Expand Up @@ -101,6 +103,19 @@ def make_source_function(self):
plasma = self.plasma
runner = self.runner
atomic_data = self.plasma.atomic_data
macro_ref = atomic_data.macro_atom_references
macro_data = atomic_data.macro_atom_data

no_lvls = len(atomic_data.levels)
no_shells = len(model.w)

if runner.line_interaction_type == 'macroatom':
internal_jump_mask = (macro_data.transition_type >= 0).values
ma_int_data = macro_data[internal_jump_mask]
internal = plasma.transition_probabilities[internal_jump_mask]

source_level_idx = ma_int_data.source_level_idx.values
destination_level_idx = ma_int_data.destination_level_idx.values

Edotlu_norm_factor = (1 / (runner.time_of_simulation * model.volume))
exptau = 1 - np.exp(- plasma.tau_sobolevs)
Expand All @@ -118,6 +133,22 @@ def make_source_function(self):
upper_level_index = atomic_data.lines.index.droplevel('level_number_lower')
e_dot_lu = pd.DataFrame(Edotlu, index=upper_level_index)
e_dot_u = e_dot_lu.groupby(level=[0, 1, 2]).sum()
e_dot_u_src_idx = macro_ref.loc[e_dot_u.index].references_idx.values

if runner.line_interaction_type == 'macroatom':
C_frame = pd.DataFrame(
columns=np.arange(no_shells), index=macro_ref.index
)
q_indices = (source_level_idx, destination_level_idx)
for shell in xrange(no_shells):
Q = sp.coo_matrix(
(internal[shell], q_indices), shape=(no_lvls, no_lvls)
)
inv_N = sp.identity(no_lvls) - Q
e_dot_u_vec = np.zeros(no_lvls)
e_dot_u_vec[e_dot_u_src_idx] = e_dot_u[shell].values
C_frame[shell] = sp.linalg.spsolve(inv_N.T, e_dot_u_vec)

e_dot_u.index.names = ['atomic_number', 'ion_number', 'source_level_number'] # To make the q_ul e_dot_u product work, could be cleaner
transitions = atomic_data.macro_atom_data[atomic_data.macro_atom_data.transition_type == -1].copy()
transitions_index = transitions.set_index(['atomic_number', 'ion_number', 'source_level_number']).index.copy()
Expand All @@ -126,6 +157,8 @@ def make_source_function(self):
t = model.time_explosion.value
lines = atomic_data.lines.set_index('line_id')
wave = lines.wavelength_cm.loc[transitions.transition_line_id].values.reshape(-1,1)
if runner.line_interaction_type == 'macroatom':
e_dot_u = C_frame.loc[e_dot_u.index]
att_S_ul = (wave * (q_ul * e_dot_u) * t / (4 * np.pi))

result = pd.DataFrame(att_S_ul.as_matrix(), index=transitions.transition_line_id.values)
Expand Down

0 comments on commit 29f0193

Please sign in to comment.