Quickstart/Cheat-Sheet¶
Since this documentation is written in a jupyter-notebook we will import a little ipython helper function to display file with syntax highlighting.
[1]:
from glotaran.utils.ipython import display_file
To start using pyglotaran
in your project, you have to import it first. In addition we need to import some extra components for later use.
[2]:
from glotaran.analysis.optimize import optimize
from glotaran.io import load_model
from glotaran.io import load_parameters
from glotaran.io import save_dataset
from glotaran.io.prepare_dataset import prepare_time_trace_dataset
from glotaran.project.scheme import Scheme
Let us get some example data to analyze:
[3]:
from glotaran.examples.sequential import dataset
dataset
[3]:
<xarray.Dataset> Dimensions: (time: 2100, spectral: 72) Coordinates: * time (time) float64 -1.0 -0.99 -0.98 -0.97 ... 19.96 19.97 19.98 19.99 * spectral (spectral) float64 600.0 601.4 602.8 604.2 ... 696.6 698.0 699.4 Data variables: data (time, spectral) float64 0.004406 0.002114 ... 1.692 1.547
- time: 2100
- spectral: 72
- time(time)float64-1.0 -0.99 -0.98 ... 19.98 19.99
array([-1. , -0.99, -0.98, ..., 19.97, 19.98, 19.99])
- spectral(spectral)float64600.0 601.4 602.8 ... 698.0 699.4
array([600. , 601.4, 602.8, 604.2, 605.6, 607. , 608.4, 609.8, 611.2, 612.6, 614. , 615.4, 616.8, 618.2, 619.6, 621. , 622.4, 623.8, 625.2, 626.6, 628. , 629.4, 630.8, 632.2, 633.6, 635. , 636.4, 637.8, 639.2, 640.6, 642. , 643.4, 644.8, 646.2, 647.6, 649. , 650.4, 651.8, 653.2, 654.6, 656. , 657.4, 658.8, 660.2, 661.6, 663. , 664.4, 665.8, 667.2, 668.6, 670. , 671.4, 672.8, 674.2, 675.6, 677. , 678.4, 679.8, 681.2, 682.6, 684. , 685.4, 686.8, 688.2, 689.6, 691. , 692.4, 693.8, 695.2, 696.6, 698. , 699.4])
- data(time, spectral)float640.004406 0.002114 ... 1.692 1.547
array([[ 4.40640583e-03, 2.11430457e-03, 6.16223200e-03, ..., 2.46116898e-03, 7.52202480e-03, 1.68531044e-03], [ 9.76765902e-03, 9.24351743e-03, 6.57267519e-03, ..., 8.57267681e-04, 8.62908518e-03, 1.33994183e-03], [ 6.19691597e-04, 3.36874081e-03, -7.72069381e-03, ..., 2.90945698e-03, 2.68457302e-03, 6.29536990e-03], ..., [ 1.46403952e+00, 1.63522605e+00, 1.80513953e+00, ..., 1.89378532e+00, 1.69785761e+00, 1.53055474e+00], [ 1.47573943e+00, 1.64152650e+00, 1.81016020e+00, ..., 1.88131920e+00, 1.69261786e+00, 1.52690121e+00], [ 1.46267451e+00, 1.64286647e+00, 1.80474098e+00, ..., 1.89704089e+00, 1.69249504e+00, 1.54737863e+00]])
Like all data in pyglotaran
, the dataset is a xarray.Dataset. You can find more information about the xarray
library the xarray hompage.
The loaded dataset is a simulated sequential model.
Plotting raw data¶
Now we lets plot some time traces.
[4]:
plot_data = dataset.data.sel(spectral=[620, 630, 650], method="nearest")
plot_data.plot.line(x="time", aspect=2, size=5);
We can also plot spectra at different times.
[5]:
plot_data = dataset.data.sel(time=[1, 10, 20], method="nearest")
plot_data.plot.line(x="spectral", aspect=2, size=5);
Preparing data¶
To get an idea about how to model your data, you should inspect the singular value decomposition. Pyglotaran has a function to calculate it (among other things).
[6]:
dataset = prepare_time_trace_dataset(dataset)
dataset
[6]:
<xarray.Dataset> Dimensions: (time: 2100, spectral: 72, left_singular_value_index: 72, singular_value_index: 72, right_singular_value_index: 72) Coordinates: * time (time) float64 -1.0 -0.99 -0.98 ... 19.98 19.99 * spectral (spectral) float64 600.0 601.4 ... 698.0 699.4 Dimensions without coordinates: left_singular_value_index, singular_value_index, right_singular_value_index Data variables: data (time, spectral) float64 0.004406 ... 1.547 data_left_singular_vectors (time, left_singular_value_index) float64 -1... data_singular_values (singular_value_index) float64 4.62e+03 ... ... data_right_singular_vectors (right_singular_value_index, spectral) float64 ...
- time: 2100
- spectral: 72
- left_singular_value_index: 72
- singular_value_index: 72
- right_singular_value_index: 72
- time(time)float64-1.0 -0.99 -0.98 ... 19.98 19.99
array([-1. , -0.99, -0.98, ..., 19.97, 19.98, 19.99])
- spectral(spectral)float64600.0 601.4 602.8 ... 698.0 699.4
array([600. , 601.4, 602.8, 604.2, 605.6, 607. , 608.4, 609.8, 611.2, 612.6, 614. , 615.4, 616.8, 618.2, 619.6, 621. , 622.4, 623.8, 625.2, 626.6, 628. , 629.4, 630.8, 632.2, 633.6, 635. , 636.4, 637.8, 639.2, 640.6, 642. , 643.4, 644.8, 646.2, 647.6, 649. , 650.4, 651.8, 653.2, 654.6, 656. , 657.4, 658.8, 660.2, 661.6, 663. , 664.4, 665.8, 667.2, 668.6, 670. , 671.4, 672.8, 674.2, 675.6, 677. , 678.4, 679.8, 681.2, 682.6, 684. , 685.4, 686.8, 688.2, 689.6, 691. , 692.4, 693.8, 695.2, 696.6, 698. , 699.4])
- data(time, spectral)float640.004406 0.002114 ... 1.692 1.547
array([[ 4.40640583e-03, 2.11430457e-03, 6.16223200e-03, ..., 2.46116898e-03, 7.52202480e-03, 1.68531044e-03], [ 9.76765902e-03, 9.24351743e-03, 6.57267519e-03, ..., 8.57267681e-04, 8.62908518e-03, 1.33994183e-03], [ 6.19691597e-04, 3.36874081e-03, -7.72069381e-03, ..., 2.90945698e-03, 2.68457302e-03, 6.29536990e-03], ..., [ 1.46403952e+00, 1.63522605e+00, 1.80513953e+00, ..., 1.89378532e+00, 1.69785761e+00, 1.53055474e+00], [ 1.47573943e+00, 1.64152650e+00, 1.81016020e+00, ..., 1.88131920e+00, 1.69261786e+00, 1.52690121e+00], [ 1.46267451e+00, 1.64286647e+00, 1.80474098e+00, ..., 1.89704089e+00, 1.69249504e+00, 1.54737863e+00]])
- data_left_singular_vectors(time, left_singular_value_index)float64-1.634e-06 3.552e-06 ... 0.03732
array([[-1.63397956e-06, 3.55229012e-06, 4.41236917e-05, ..., 1.23961893e-03, 2.88783261e-03, 1.29354289e-03], [-3.30423988e-07, -7.57717220e-06, 8.89692346e-05, ..., 1.36331283e-02, 9.30819939e-04, -1.19431864e-02], [ 2.98036845e-06, 1.35547896e-05, -7.14712926e-05, ..., -3.70481659e-02, 1.72472720e-02, 7.19463249e-03], ..., [-1.22425909e-02, -8.34920848e-03, 1.19990316e-02, ..., -3.57640918e-02, -2.27705903e-02, -2.73824090e-02], [-1.22391467e-02, -8.36765748e-03, 1.19006059e-02, ..., 6.77543051e-03, 6.93989656e-03, -1.93649218e-02], [-1.22194140e-02, -8.34615736e-03, 1.19060123e-02, ..., 3.12538398e-02, -5.16025808e-03, 3.73217258e-02]])
- data_singular_values(singular_value_index)float644.62e+03 1.126e+03 ... 0.3779
array([4.62000286e+03, 1.12633314e+03, 1.79751200e+02, 5.38553033e-01, 5.33920885e-01, 5.31937048e-01, 5.29579871e-01, 5.19738874e-01, 5.17405429e-01, 5.15384273e-01, 5.11816784e-01, 5.10256642e-01, 5.06771426e-01, 5.05634516e-01, 5.03659448e-01, 5.00646955e-01, 4.98573398e-01, 4.95635803e-01, 4.95387767e-01, 4.94395652e-01, 4.92134982e-01, 4.88835430e-01, 4.86206387e-01, 4.83984293e-01, 4.83464313e-01, 4.80980077e-01, 4.78212858e-01, 4.76454159e-01, 4.74176335e-01, 4.72614360e-01, 4.70756311e-01, 4.69095471e-01, 4.67107317e-01, 4.63198901e-01, 4.60971354e-01, 4.59380889e-01, 4.58614036e-01, 4.56862846e-01, 4.55518792e-01, 4.51288190e-01, 4.47766935e-01, 4.46043108e-01, 4.44993781e-01, 4.44097253e-01, 4.42271835e-01, 4.39281362e-01, 4.38732177e-01, 4.37956030e-01, 4.35605797e-01, 4.34005547e-01, 4.29969086e-01, 4.28473762e-01, 4.27661255e-01, 4.23085704e-01, 4.22913643e-01, 4.21339517e-01, 4.19968576e-01, 4.17595229e-01, 4.15093640e-01, 4.12109834e-01, 4.11674520e-01, 4.06690100e-01, 4.05821646e-01, 4.03755195e-01, 4.02826971e-01, 3.97008816e-01, 3.96541028e-01, 3.91919042e-01, 3.89826145e-01, 3.86443967e-01, 3.81254593e-01, 3.77851662e-01])
- data_right_singular_vectors(right_singular_value_index, spectral)float64-0.03526 -0.03908 ... 0.05706
array([[-0.03526232, -0.03907575, -0.04313042, ..., -0.02800184, -0.02528602, -0.02276965], [ 0.09065977, 0.09931153, 0.10801094, ..., -0.03120991, -0.02818644, -0.02537568], [ 0.15354044, 0.16821143, 0.18260036, ..., 0.0071333 , 0.00645198, 0.00583452], ..., [-0.03987942, -0.05484248, -0.07426609, ..., 0.12492518, -0.12920086, 0.00441674], [ 0.22772252, 0.00561649, -0.16994312, ..., 0.03688477, 0.09298589, 0.12265546], [-0.02928189, 0.12690273, -0.05653684, ..., -0.08193683, -0.04207367, 0.0570617 ]])
First, take a look at the first 10 singular values:
[7]:
plot_data = dataset.data_singular_values.sel(singular_value_index=range(0, 10))
plot_data.plot(yscale="log", marker="o", linewidth=0, aspect=2, size=5);
Working with models¶
To analyze our data, we need to create a model.
Create a file called model.yaml
in your working directory and fill it with the following:
[8]:
display_file("model.yaml", syntax="yaml")
[8]:
type: kinetic-spectrum
initial_concentration:
input:
compartments: [s1, s2, s3]
parameters: [input.1, input.0, input.0]
k_matrix:
k1:
matrix:
(s2, s1): kinetic.1
(s3, s2): kinetic.2
(s3, s3): kinetic.3
megacomplex:
m1:
k_matrix: [k1]
irf:
irf1:
type: gaussian
center: irf.center
width: irf.width
dataset:
dataset1:
initial_concentration: input
megacomplex: [m1]
irf: irf1
Now you can load the model file.
[9]:
model = load_model("model.yaml")
You can check your model for problems with model.validate
.
[10]:
model.validate()
[10]:
'Your model is valid.'
Working with parameters¶
Now define some starting parameters. Create a file called parameters.yaml
with the following content.
[11]:
display_file("parameters.yaml", syntax="yaml")
[11]:
input:
- ['1', 1, {'vary': False, 'non-negative': False}]
- ['0', 0, {'vary': False, 'non-negative': False}]
kinetic: [
0.5,
0.3,
0.1,
]
irf:
- ['center', 0.3]
- ['width', 0.1]
[12]:
parameters = load_parameters("parameters.yaml")
You can model.validate
also to check for missing parameters.
[13]:
model.validate(parameters=parameters)
[13]:
'Your model is valid.'
Since not all problems in the model can be detected automatically it is wise to visually inspect the model. For this purpose, you can just print the model.
[14]:
model
[14]:
Model¶
Type: kinetic-spectrum
Initial Concentration¶
input:
Label: input
Compartments: [‘s1’, ‘s2’, ‘s3’]
Parameters: [input.1, input.0, input.0]
Exclude From Normalize: []
K Matrix¶
k1:
Label: k1
Matrix:
(‘s2’, ‘s1’): kinetic.1
(‘s3’, ‘s2’): kinetic.2
(‘s3’, ‘s3’): kinetic.3
Irf¶
irf1 (gaussian):
Label: irf1
Type: gaussian
Center: irf.center
Width: irf.width
Normalize: True
Backsweep: False
Dataset¶
dataset1:
Label: dataset1
Megacomplex: [‘m1’]
Initial Concentration: input
Irf: irf1
Megacomplex¶
m1 (None):
Label: m1
K Matrix: [‘k1’]
The same way you should inspect your parameters.
[15]:
parameters
[15]:
input:
Label
Value
StdErr
Min
Max
Vary
Non-Negative
Expr
1
1
0
-inf
inf
False
False
None
0
0
0
-inf
inf
False
False
None
irf:
Label
Value
StdErr
Min
Max
Vary
Non-Negative
Expr
center
0.3
0
-inf
inf
True
False
None
width
0.1
0
-inf
inf
True
False
None
kinetic:
Label
Value
StdErr
Min
Max
Vary
Non-Negative
Expr
1
0.5
0
-inf
inf
True
False
None
2
0.3
0
-inf
inf
True
False
None
3
0.1
0
-inf
inf
True
False
None
Optimizing data¶
Now we have everything together to optimize our parameters. First we import optimize.
[16]:
scheme = Scheme(model, parameters, {"dataset1": dataset})
result = optimize(scheme)
result
Iteration Total nfev Cost Cost reduction Step norm Optimality
0 1 7.5779e+00 1.01e+02
1 2 7.5777e+00 1.69e-04 5.40e-05 5.85e-02
2 3 7.5777e+00 3.06e-11 1.16e-08 4.27e-06
`ftol` termination condition is satisfied.
Function evaluations 3, initial cost 7.5779e+00, final cost 7.5777e+00, first-order optimality 4.27e-06.
[16]:
Optimization Result |
|
---|---|
Number of residual evaluation |
3 |
Number of variables |
5 |
Number of datapoints |
151200 |
Degrees of freedom |
151195 |
Chi Square |
1.52e+01 |
Reduced Chi Square |
1.00e-04 |
Root Mean Square Error (RMSE) |
1.00e-02 |
Model¶
Type: kinetic-spectrum
Initial Concentration¶
input:
Label: input
Compartments: [‘s1’, ‘s2’, ‘s3’]
Parameters: [input.1: 1.00000e+00 (fixed), input.0: 0.00000e+00 (fixed), input.0: 0.00000e+00 (fixed)]
Exclude From Normalize: []
K Matrix¶
k1:
Label: k1
Matrix:
(‘s2’, ‘s1’): kinetic.1: 5.00042e-01 (StdErr: 7e-05 ,initial: 5.00000e-01)
(‘s3’, ‘s2’): kinetic.2: 2.99967e-01 (StdErr: 4e-05 ,initial: 3.00000e-01)
(‘s3’, ‘s3’): kinetic.3: 1.00002e-01 (StdErr: 5e-06 ,initial: 1.00000e-01)
Irf¶
irf1 (gaussian):
Label: irf1
Type: gaussian
Center: irf.center: 2.99997e-01 (StdErr: 5e-06 ,initial: 3.00000e-01)
Width: irf.width: 1.00004e-01 (StdErr: 7e-06 ,initial: 1.00000e-01)
Normalize: True
Backsweep: False
Dataset¶
dataset1:
Label: dataset1
Megacomplex: [‘m1’]
Initial Concentration: input
Irf: irf1
Megacomplex¶
m1 (None):
Label: m1
K Matrix: [‘k1’]
[17]:
result.optimized_parameters
[17]:
input:
Label
Value
StdErr
Min
Max
Vary
Non-Negative
Expr
1
1
0
-inf
inf
False
False
None
0
0
0
-inf
inf
False
False
None
irf:
Label
Value
StdErr
Min
Max
Vary
Non-Negative
Expr
center
0.299997
5.01701e-06
-inf
inf
True
False
None
width
0.100004
6.71203e-06
-inf
inf
True
False
None
kinetic:
Label
Value
StdErr
Min
Max
Vary
Non-Negative
Expr
1
0.500042
7.26493e-05
-inf
inf
True
False
None
2
0.299967
4.19612e-05
-inf
inf
True
False
None
3
0.100002
4.78672e-06
-inf
inf
True
False
None
You can get the resulting data for your dataset with result.get_dataset
.
[18]:
result_dataset = result.data["dataset1"]
result_dataset
[18]:
<xarray.Dataset> Dimensions: (time: 2100, spectral: 72, left_singular_value_index: 72, singular_value_index: 72, right_singular_value_index: 72, clp_label: 3, species: 3, component: 3, to_species: 3, from_species: 3) Coordinates: * time (time) float64 -1.0 ... 19.99 * spectral (spectral) float64 600.0 ... 699.4 * clp_label (clp_label) <U2 's1' 's2' 's3' * species (species) <U2 's1' 's2' 's3' rate (component) float64 -0.5 -0.3 -0.1 lifetime (component) float64 -2.0 ... -10.0 * to_species (to_species) <U2 's1' 's2' 's3' * from_species (from_species) <U2 's1' 's2' 's3' Dimensions without coordinates: left_singular_value_index, singular_value_index, right_singular_value_index, component Data variables: (12/23) data (time, spectral) float64 0.0044... data_left_singular_vectors (time, left_singular_value_index) float64 ... data_singular_values (singular_value_index) float64 ... data_right_singular_vectors (spectral, right_singular_value_index) float64 ... matrix (time, clp_label) float64 6.135... clp (spectral, clp_label) float64 1... ... ... decay_associated_spectra (spectral, component) float64 2... a_matrix (component, species) float64 1.... k_matrix (to_species, from_species) float64 ... k_matrix_reduced (to_species, from_species) float64 ... irf_center float64 0.3 irf_width float64 0.1 Attributes: root_mean_square_error: 0.010011690198544799 weighted_root_mean_square_error: 0.010011690198544799
- time: 2100
- spectral: 72
- left_singular_value_index: 72
- singular_value_index: 72
- right_singular_value_index: 72
- clp_label: 3
- species: 3
- component: 3
- to_species: 3
- from_species: 3
- time(time)float64-1.0 -0.99 -0.98 ... 19.98 19.99
array([-1. , -0.99, -0.98, ..., 19.97, 19.98, 19.99])
- spectral(spectral)float64600.0 601.4 602.8 ... 698.0 699.4
array([600. , 601.4, 602.8, 604.2, 605.6, 607. , 608.4, 609.8, 611.2, 612.6, 614. , 615.4, 616.8, 618.2, 619.6, 621. , 622.4, 623.8, 625.2, 626.6, 628. , 629.4, 630.8, 632.2, 633.6, 635. , 636.4, 637.8, 639.2, 640.6, 642. , 643.4, 644.8, 646.2, 647.6, 649. , 650.4, 651.8, 653.2, 654.6, 656. , 657.4, 658.8, 660.2, 661.6, 663. , 664.4, 665.8, 667.2, 668.6, 670. , 671.4, 672.8, 674.2, 675.6, 677. , 678.4, 679.8, 681.2, 682.6, 684. , 685.4, 686.8, 688.2, 689.6, 691. , 692.4, 693.8, 695.2, 696.6, 698. , 699.4])
- clp_label(clp_label)<U2's1' 's2' 's3'
array(['s1', 's2', 's3'], dtype='<U2')
- species(species)<U2's1' 's2' 's3'
array(['s1', 's2', 's3'], dtype='<U2')
- rate(component)float64-0.5 -0.3 -0.1
array([-0.50004222, -0.29996672, -0.10000202])
- lifetime(component)float64-2.0 -3.334 -10.0
array([-1.99983115, -3.33370315, -9.99979789])
- to_species(to_species)<U2's1' 's2' 's3'
array(['s1', 's2', 's3'], dtype='<U2')
- from_species(from_species)<U2's1' 's2' 's3'
array(['s1', 's2', 's3'], dtype='<U2')
- data(time, spectral)float640.004406 0.002114 ... 1.692 1.547
array([[ 4.40640583e-03, 2.11430457e-03, 6.16223200e-03, ..., 2.46116898e-03, 7.52202480e-03, 1.68531044e-03], [ 9.76765902e-03, 9.24351743e-03, 6.57267519e-03, ..., 8.57267681e-04, 8.62908518e-03, 1.33994183e-03], [ 6.19691597e-04, 3.36874081e-03, -7.72069381e-03, ..., 2.90945698e-03, 2.68457302e-03, 6.29536990e-03], ..., [ 1.46403952e+00, 1.63522605e+00, 1.80513953e+00, ..., 1.89378532e+00, 1.69785761e+00, 1.53055474e+00], [ 1.47573943e+00, 1.64152650e+00, 1.81016020e+00, ..., 1.88131920e+00, 1.69261786e+00, 1.52690121e+00], [ 1.46267451e+00, 1.64286647e+00, 1.80474098e+00, ..., 1.89704089e+00, 1.69249504e+00, 1.54737863e+00]])
- data_left_singular_vectors(time, left_singular_value_index)float64-1.634e-06 3.552e-06 ... 0.03732
array([[-1.63397956e-06, 3.55229012e-06, 4.41236917e-05, ..., 1.23961893e-03, 2.88783261e-03, 1.29354289e-03], [-3.30423988e-07, -7.57717220e-06, 8.89692346e-05, ..., 1.36331283e-02, 9.30819939e-04, -1.19431864e-02], [ 2.98036845e-06, 1.35547896e-05, -7.14712926e-05, ..., -3.70481659e-02, 1.72472720e-02, 7.19463249e-03], ..., [-1.22425909e-02, -8.34920848e-03, 1.19990316e-02, ..., -3.57640918e-02, -2.27705903e-02, -2.73824090e-02], [-1.22391467e-02, -8.36765748e-03, 1.19006059e-02, ..., 6.77543051e-03, 6.93989656e-03, -1.93649218e-02], [-1.22194140e-02, -8.34615736e-03, 1.19060123e-02, ..., 3.12538398e-02, -5.16025808e-03, 3.73217258e-02]])
- data_singular_values(singular_value_index)float644.62e+03 1.126e+03 ... 0.3779
array([4.62000286e+03, 1.12633314e+03, 1.79751200e+02, 5.38553033e-01, 5.33920885e-01, 5.31937048e-01, 5.29579871e-01, 5.19738874e-01, 5.17405429e-01, 5.15384273e-01, 5.11816784e-01, 5.10256642e-01, 5.06771426e-01, 5.05634516e-01, 5.03659448e-01, 5.00646955e-01, 4.98573398e-01, 4.95635803e-01, 4.95387767e-01, 4.94395652e-01, 4.92134982e-01, 4.88835430e-01, 4.86206387e-01, 4.83984293e-01, 4.83464313e-01, 4.80980077e-01, 4.78212858e-01, 4.76454159e-01, 4.74176335e-01, 4.72614360e-01, 4.70756311e-01, 4.69095471e-01, 4.67107317e-01, 4.63198901e-01, 4.60971354e-01, 4.59380889e-01, 4.58614036e-01, 4.56862846e-01, 4.55518792e-01, 4.51288190e-01, 4.47766935e-01, 4.46043108e-01, 4.44993781e-01, 4.44097253e-01, 4.42271835e-01, 4.39281362e-01, 4.38732177e-01, 4.37956030e-01, 4.35605797e-01, 4.34005547e-01, 4.29969086e-01, 4.28473762e-01, 4.27661255e-01, 4.23085704e-01, 4.22913643e-01, 4.21339517e-01, 4.19968576e-01, 4.17595229e-01, 4.15093640e-01, 4.12109834e-01, 4.11674520e-01, 4.06690100e-01, 4.05821646e-01, 4.03755195e-01, 4.02826971e-01, 3.97008816e-01, 3.96541028e-01, 3.91919042e-01, 3.89826145e-01, 3.86443967e-01, 3.81254593e-01, 3.77851662e-01])
- data_right_singular_vectors(spectral, right_singular_value_index)float64-0.03526 0.09066 ... 0.1227 0.05706
array([[-0.03526232, 0.09065977, 0.15354044, ..., -0.03987942, 0.22772252, -0.02928189], [-0.03907575, 0.09931153, 0.16821143, ..., -0.05484248, 0.00561649, 0.12690273], [-0.04313042, 0.10801094, 0.18260036, ..., -0.07426609, -0.16994312, -0.05653684], ..., [-0.02800184, -0.03120991, 0.0071333 , ..., 0.12492518, 0.03688477, -0.08193683], [-0.02528602, -0.02818644, 0.00645198, ..., -0.12920086, 0.09298589, -0.04207367], [-0.02276965, -0.02537568, 0.00583452, ..., 0.00441674, 0.12265546, 0.0570617 ]])
- matrix(time, clp_label)float646.135e-39 2.328e-41 ... 0.2516
array([[6.13514971e-39, 2.32768320e-41, 5.27597823e-44], [2.25681381e-38, 8.62709393e-41, 1.97007373e-43], [8.21956225e-38, 3.16600667e-40, 7.28437234e-43], ..., [5.35667661e-05, 6.71343722e-03, 2.52111561e-01], [5.32995783e-05, 6.69359610e-03, 2.51879669e-01], [5.30337231e-05, 6.67381308e-03, 2.51647949e-01]])
- clp(spectral, clp_label)float6415.0 0.04572 ... -0.0004073 6.108
array([[ 1.49960096e+01, 4.57160153e-02, 5.83132627e+00], [ 1.64713223e+01, 7.03188190e-02, 6.48802777e+00], [ 1.79688090e+01, 1.16465927e-01, 7.19409224e+00], [ 1.94633659e+01, 1.99530065e-01, 7.95095980e+00], [ 2.09455086e+01, 3.21273823e-01, 8.76451789e+00], [ 2.23863588e+01, 5.11762284e-01, 9.62994279e+00], [ 2.37600274e+01, 7.90207471e-01, 1.05498251e+01], [ 2.50528063e+01, 1.18244101e+00, 1.15222796e+01], [ 2.62331418e+01, 1.72746991e+00, 1.25466574e+01], [ 2.72823535e+01, 2.45531646e+00, 1.36219055e+01], [ 2.81842207e+01, 3.39429109e+00, 1.47430216e+01], [ 2.89216181e+01, 4.56429885e+00, 1.59094107e+01], [ 2.94699962e+01, 5.98118837e+00, 1.71149777e+01], [ 2.98326345e+01, 7.62327066e+00, 1.83565855e+01], [ 2.99927718e+01, 9.45147870e+00, 1.96310784e+01], [ 2.99480980e+01, 1.14099554e+01, 2.09291222e+01], [ 2.97025693e+01, 1.34017794e+01, 2.22472462e+01], [ 2.92599502e+01, 1.53235772e+01, 2.35749868e+01], [ 2.86273395e+01, 1.70460852e+01, 2.49091663e+01], [ 2.78168673e+01, 1.84631534e+01, 2.62360624e+01], ... [ 2.36709456e-01, 2.07649279e-03, 2.68030913e+01], [ 1.85896999e-01, -2.03327125e-03, 2.54799369e+01], [ 1.41486682e-01, -2.21098318e-04, 2.41471086e+01], [ 1.08603055e-01, 7.14355787e-04, 2.28148344e+01], [ 7.97808202e-02, 6.69212614e-04, 2.14928918e+01], [ 6.43138191e-02, -2.41084834e-03, 2.01861752e+01], [ 4.72773084e-02, -2.23838162e-03, 1.89013993e+01], [ 3.61819950e-02, -1.43225142e-03, 1.76445992e+01], [ 2.34120589e-02, 4.98988560e-04, 1.64220417e+01], [ 1.44045062e-02, 1.97349219e-03, 1.52376376e+01], [ 1.32559656e-02, 1.56797378e-03, 1.40954062e+01], [ 7.77258158e-03, 2.03043426e-03, 1.30008235e+01], [ 7.35747123e-03, -1.46528736e-03, 1.19565818e+01], [ 5.64128771e-03, -2.65545284e-04, 1.09595307e+01], [ 1.88974995e-03, 1.93339968e-03, 1.00165852e+01], [ 2.89379540e-03, -1.86583680e-03, 9.12891060e+00], [ 1.10595484e-03, 1.14779160e-03, 8.29273849e+00], [ 6.02513002e-06, 6.64000324e-04, 7.51101054e+00], [-1.06834010e-04, 1.66389161e-04, 6.78276776e+00], [ 9.41671831e-04, -4.07315711e-04, 6.10778224e+00]])
- weighted_residual(time, spectral)float640.004406 0.002114 ... 0.01037
array([[ 0.00440641, 0.0021143 , 0.00616223, ..., 0.00246117, 0.00752202, 0.00168531], [ 0.00976766, 0.00924352, 0.00657268, ..., 0.00085727, 0.00862909, 0.00133994], [ 0.00061969, 0.00336874, -0.00772069, ..., 0.00290946, 0.00268457, 0.00629537], ..., [-0.00721545, -0.00183516, -0.01031871, ..., 0.00016827, -0.01215767, -0.00928509], [ 0.00584162, 0.00597562, -0.00362267, ..., -0.01055609, -0.01582454, -0.01152228], [-0.00586718, 0.00882476, -0.0073678 , ..., 0.00690606, -0.01437566, 0.01037043]])
- residual(time, spectral)float640.004406 0.002114 ... 0.01037
array([[ 0.00440641, 0.0021143 , 0.00616223, ..., 0.00246117, 0.00752202, 0.00168531], [ 0.00976766, 0.00924352, 0.00657268, ..., 0.00085727, 0.00862909, 0.00133994], [ 0.00061969, 0.00336874, -0.00772069, ..., 0.00290946, 0.00268457, 0.00629537], ..., [-0.00721545, -0.00183516, -0.01031871, ..., 0.00016827, -0.01215767, -0.00928509], [ 0.00584162, 0.00597562, -0.00362267, ..., -0.01055609, -0.01582454, -0.01152228], [-0.00586718, 0.00882476, -0.0073678 , ..., 0.00690606, -0.01437566, 0.01037043]])
- weighted_residual_left_singular_vectors(time, left_singular_value_index)float64-0.005691 0.001499 ... 0.03928
array([[-0.00569091, 0.00149892, -0.01204452, ..., 0.02164489, 0.00866273, -0.00613372], [-0.01267347, -0.03282439, 0.00308261, ..., 0.00266544, 0.00612207, -0.01044496], [ 0.0090834 , 0.03708844, 0.01981182, ..., 0.00080255, 0.00748311, 0.00485063], ..., [-0.01575108, 0.02379694, 0.00015748, ..., 0.01991437, -0.02535916, -0.02495013], [ 0.03181592, -0.00922295, 0.02997191, ..., 0.02202831, 0.01296788, -0.02364961], [ 0.01063036, 0.03058216, -0.00600193, ..., -0.02131906, -0.01005235, 0.03928321]])
- weighted_residual_singular_values(singular_value_index)float640.5388 0.5349 ... 0.3805 0.3771
array([0.53875738, 0.53492355, 0.53234792, 0.53141546, 0.52350347, 0.52014971, 0.51763156, 0.515837 , 0.51175972, 0.50770305, 0.50653803, 0.50445833, 0.50335449, 0.50117685, 0.49901635, 0.49574294, 0.49512834, 0.49332981, 0.49130765, 0.4902113 , 0.48767831, 0.48464475, 0.48250013, 0.4818346 , 0.47785268, 0.47639144, 0.47378499, 0.47323422, 0.47059672, 0.46866317, 0.46627588, 0.46372974, 0.46117045, 0.4600454 , 0.45882264, 0.45770538, 0.45284726, 0.45145812, 0.44986437, 0.44857042, 0.44628129, 0.44476373, 0.44335807, 0.44182479, 0.43920228, 0.43825261, 0.43795698, 0.4355028 , 0.43383009, 0.43183149, 0.42896918, 0.42754143, 0.42721465, 0.42298028, 0.42238096, 0.42107737, 0.41773087, 0.41623072, 0.41101751, 0.40972659, 0.40814534, 0.40562001, 0.40310644, 0.40218224, 0.40070383, 0.39649309, 0.39289089, 0.38978767, 0.38597177, 0.38568146, 0.38053178, 0.37705321])
- weighted_residual_right_singular_vectors(right_singular_value_index, spectral)float64-0.1852 0.03689 ... 0.03781
array([[-0.18518233, 0.03689312, 0.14841771, ..., -0.17127057, -0.16493907, 0.11690938], [ 0.08846118, 0.06741675, -0.15494277, ..., 0.26800598, 0.00627419, 0.11649434], [-0.04971632, 0.01471108, 0.07888437, ..., -0.02673401, 0.09243867, 0.06280086], ..., [-0.01033244, -0.07789872, -0.01835027, ..., -0.0295117 , -0.00388295, -0.1014005 ], [ 0.22232868, -0.01267628, -0.17212122, ..., 0.07502039, 0.07038176, 0.08953982], [-0.04669009, 0.09773565, -0.04477255, ..., -0.06418638, -0.03701197, 0.03781061]])
- residual_left_singular_vectors(time, left_singular_value_index)float64-0.005691 0.001499 ... 0.03928
array([[-0.00569091, 0.00149892, -0.01204452, ..., 0.02164489, 0.00866273, -0.00613372], [-0.01267347, -0.03282439, 0.00308261, ..., 0.00266544, 0.00612207, -0.01044496], [ 0.0090834 , 0.03708844, 0.01981182, ..., 0.00080255, 0.00748311, 0.00485063], ..., [-0.01575108, 0.02379694, 0.00015748, ..., 0.01991437, -0.02535916, -0.02495013], [ 0.03181592, -0.00922295, 0.02997191, ..., 0.02202831, 0.01296788, -0.02364961], [ 0.01063036, 0.03058216, -0.00600193, ..., -0.02131906, -0.01005235, 0.03928321]])
- residual_singular_values(singular_value_index)float640.5388 0.5349 ... 0.3805 0.3771
array([0.53875738, 0.53492355, 0.53234792, 0.53141546, 0.52350347, 0.52014971, 0.51763156, 0.515837 , 0.51175972, 0.50770305, 0.50653803, 0.50445833, 0.50335449, 0.50117685, 0.49901635, 0.49574294, 0.49512834, 0.49332981, 0.49130765, 0.4902113 , 0.48767831, 0.48464475, 0.48250013, 0.4818346 , 0.47785268, 0.47639144, 0.47378499, 0.47323422, 0.47059672, 0.46866317, 0.46627588, 0.46372974, 0.46117045, 0.4600454 , 0.45882264, 0.45770538, 0.45284726, 0.45145812, 0.44986437, 0.44857042, 0.44628129, 0.44476373, 0.44335807, 0.44182479, 0.43920228, 0.43825261, 0.43795698, 0.4355028 , 0.43383009, 0.43183149, 0.42896918, 0.42754143, 0.42721465, 0.42298028, 0.42238096, 0.42107737, 0.41773087, 0.41623072, 0.41101751, 0.40972659, 0.40814534, 0.40562001, 0.40310644, 0.40218224, 0.40070383, 0.39649309, 0.39289089, 0.38978767, 0.38597177, 0.38568146, 0.38053178, 0.37705321])
- residual_right_singular_vectors(right_singular_value_index, spectral)float64-0.1852 0.03689 ... 0.03781
array([[-0.18518233, 0.03689312, 0.14841771, ..., -0.17127057, -0.16493907, 0.11690938], [ 0.08846118, 0.06741675, -0.15494277, ..., 0.26800598, 0.00627419, 0.11649434], [-0.04971632, 0.01471108, 0.07888437, ..., -0.02673401, 0.09243867, 0.06280086], ..., [-0.01033244, -0.07789872, -0.01835027, ..., -0.0295117 , -0.00388295, -0.1014005 ], [ 0.22232868, -0.01267628, -0.17212122, ..., 0.07502039, 0.07038176, 0.08953982], [-0.04669009, 0.09773565, -0.04477255, ..., -0.06418638, -0.03701197, 0.03781061]])
- fitted_data(time, spectral)float64-1.909e-14 -2.657e-14 ... 1.537
array([[-1.90897645e-14, -2.65681574e-14, 1.32464352e-13, ..., -2.74259782e-15, 3.27515792e-15, -3.31310500e-15], [-2.56912547e-14, -6.54996890e-14, -1.90342533e-14, ..., -1.99276359e-16, -5.50479801e-14, -3.61514204e-14], [-4.94179350e-15, -1.14864715e-14, -2.19364457e-14, ..., -4.79317107e-14, -5.52596163e-15, -2.29694735e-14], ..., [ 1.47125497e+00, 1.63706121e+00, 1.81545824e+00, ..., 1.89361705e+00, 1.71001528e+00, 1.53983983e+00], [ 1.46989781e+00, 1.63555089e+00, 1.81378287e+00, ..., 1.89187529e+00, 1.70844240e+00, 1.53842349e+00], [ 1.46854169e+00, 1.63404171e+00, 1.81210878e+00, ..., 1.89013483e+00, 1.70687070e+00, 1.53700820e+00]])
- species_associated_spectra(spectral, species)float6415.0 0.04572 ... -0.0004073 6.108
array([[ 1.49960096e+01, 4.57160153e-02, 5.83132627e+00], [ 1.64713223e+01, 7.03188190e-02, 6.48802777e+00], [ 1.79688090e+01, 1.16465927e-01, 7.19409224e+00], [ 1.94633659e+01, 1.99530065e-01, 7.95095980e+00], [ 2.09455086e+01, 3.21273823e-01, 8.76451789e+00], [ 2.23863588e+01, 5.11762284e-01, 9.62994279e+00], [ 2.37600274e+01, 7.90207471e-01, 1.05498251e+01], [ 2.50528063e+01, 1.18244101e+00, 1.15222796e+01], [ 2.62331418e+01, 1.72746991e+00, 1.25466574e+01], [ 2.72823535e+01, 2.45531646e+00, 1.36219055e+01], [ 2.81842207e+01, 3.39429109e+00, 1.47430216e+01], [ 2.89216181e+01, 4.56429885e+00, 1.59094107e+01], [ 2.94699962e+01, 5.98118837e+00, 1.71149777e+01], [ 2.98326345e+01, 7.62327066e+00, 1.83565855e+01], [ 2.99927718e+01, 9.45147870e+00, 1.96310784e+01], [ 2.99480980e+01, 1.14099554e+01, 2.09291222e+01], [ 2.97025693e+01, 1.34017794e+01, 2.22472462e+01], [ 2.92599502e+01, 1.53235772e+01, 2.35749868e+01], [ 2.86273395e+01, 1.70460852e+01, 2.49091663e+01], [ 2.78168673e+01, 1.84631534e+01, 2.62360624e+01], ... [ 2.36709456e-01, 2.07649279e-03, 2.68030913e+01], [ 1.85896999e-01, -2.03327125e-03, 2.54799369e+01], [ 1.41486682e-01, -2.21098318e-04, 2.41471086e+01], [ 1.08603055e-01, 7.14355787e-04, 2.28148344e+01], [ 7.97808202e-02, 6.69212614e-04, 2.14928918e+01], [ 6.43138191e-02, -2.41084834e-03, 2.01861752e+01], [ 4.72773084e-02, -2.23838162e-03, 1.89013993e+01], [ 3.61819950e-02, -1.43225142e-03, 1.76445992e+01], [ 2.34120589e-02, 4.98988560e-04, 1.64220417e+01], [ 1.44045062e-02, 1.97349219e-03, 1.52376376e+01], [ 1.32559656e-02, 1.56797378e-03, 1.40954062e+01], [ 7.77258158e-03, 2.03043426e-03, 1.30008235e+01], [ 7.35747123e-03, -1.46528736e-03, 1.19565818e+01], [ 5.64128771e-03, -2.65545284e-04, 1.09595307e+01], [ 1.88974995e-03, 1.93339968e-03, 1.00165852e+01], [ 2.89379540e-03, -1.86583680e-03, 9.12891060e+00], [ 1.10595484e-03, 1.14779160e-03, 8.29273849e+00], [ 6.02513002e-06, 6.64000324e-04, 7.51101054e+00], [-1.06834010e-04, 1.66389161e-04, 6.78276776e+00], [ 9.41671831e-04, -4.07315711e-04, 6.10778224e+00]])
- species_concentration(time, species)float646.135e-39 2.328e-41 ... 0.2516
array([[6.13514971e-39, 2.32768320e-41, 5.27597823e-44], [2.25681381e-38, 8.62709393e-41, 1.97007373e-43], [8.21956225e-38, 3.16600667e-40, 7.28437234e-43], ..., [5.35667661e-05, 6.71343722e-03, 2.52111561e-01], [5.32995783e-05, 6.69359610e-03, 2.51879669e-01], [5.30337231e-05, 6.67381308e-03, 2.51647949e-01]])
- decay_associated_spectra(spectral, component)float6425.81 -21.75 10.93 ... -22.9 11.45
array([[ 25.80997627, -21.74824494, 10.93427828], [ 28.45449439, -24.14882668, 12.1656546 ], [ 31.15985044, -26.68063243, 13.48959104], [ 33.86521883, -29.31064096, 14.90878801], [ 36.56774354, -32.05651989, 16.43428498], [ 39.15436552, -34.82504372, 18.05703704], [ 41.55603454, -37.57790888, 19.78190173], [ 43.69094962, -40.24348756, 21.60534425], [ 45.42885181, -42.72185779, 23.52614777], [ 46.6740538 , -44.93403812, 25.54233786], [ 47.33020465, -46.79051863, 27.64453465], [ 47.32931636, -48.23932085, 29.8316226 ], [ 46.59580641, -49.21798286, 32.0921726 ], [ 45.18128231, -49.76895034, 34.42030256], [ 43.16070742, -49.97803066, 36.81009504], [ 40.65388112, -49.94983095, 39.24404781], [ 37.90048714, -49.91357036, 41.71565256], [ 35.14303906, -50.08837836, 44.2052895 ], [ 32.70574519, -50.78540584, 46.70700013], [ 30.84031584, -52.21850202, 49.19505351], ... [ 50.46197358, -100.48354945, 50.25828532], [ 47.94176898, -95.53311983, 47.77724785], [ 45.39503684, -90.53162133, 45.27807117], [ 42.86306079, -85.53439104, 42.7799333 ], [ 40.3569588 , -80.57834604, 40.30116806], [ 37.9003317 , -75.68697109, 37.85095321], [ 35.47512403, -70.86972575, 35.44187903], [ 33.10670203, -66.15578209, 33.08526205], [ 30.79796607, -61.56740725, 30.79285325], [ 28.56563549, -57.12321654, 28.57198555], [ 26.42489658, -52.84183625, 26.43019563], [ 24.36694981, -48.73692869, 24.37775147], [ 22.41830561, -44.83064948, 22.41970134], [ 20.54506276, -41.08955931, 20.55013784], [ 18.76868422, -37.54882193, 18.78202746], [ 17.11563299, -34.23029434, 17.11755515], [ 15.53928119, -31.08783002, 15.54965479], [ 14.07438959, -28.15822585, 14.08384228], [ 12.71075375, -25.42918004, 12.71831945], [ 11.44827644, -22.89999368, 11.45265891]])
- a_matrix(component, species)float641.0 -2.499 1.874 ... 0.0 0.0 1.875
array([[ 1. , -2.49926767, 1.8740545 ], [ 0. , 2.49926767, -3.74914737], [ 0. , 0. , 1.87509287]])
- k_matrix(to_species, from_species)float64-0.5 0.0 0.0 0.5 ... 0.0 0.3 -0.1
array([[-0.50004222, 0. , 0. ], [ 0.50004222, -0.29996672, 0. ], [ 0. , 0.29996672, -0.10000202]])
- k_matrix_reduced(to_species, from_species)float640.0 0.0 0.0 0.5 0.0 0.0 0.0 0.3 0.1
array([[0. , 0. , 0. ], [0.50004222, 0. , 0. ], [0. , 0.29996672, 0.10000202]])
- irf_center()float640.3
array(0.2999967)
- irf_width()float640.1
array(0.10000371)
- root_mean_square_error :
- 0.010011690198544799
- weighted_root_mean_square_error :
- 0.010011690198544799
Visualize the Result¶
The resulting data can be visualized the same way as the dataset. To judge the quality of the fit, you should look at first left and right singular vectors of the residual.
[19]:
residual_left = result_dataset.residual_left_singular_vectors.sel(left_singular_value_index=0)
residual_right = result_dataset.residual_right_singular_vectors.sel(right_singular_value_index=0)
residual_left.plot.line(x="time", aspect=2, size=5)
residual_right.plot.line(x="spectral", aspect=2, size=5);
Finally, you can save your result.
[20]:
save_dataset(result_dataset, "dataset1.nc")