/
lightcurve.html
3623 lines (3164 loc) · 448 KB
/
lightcurve.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
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
102
103
104
105
106
107
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
135
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
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
<!DOCTYPE html>
<html>
<head>
<meta charset="utf-8" />
<meta name="viewport" content="width=device-width, initial-scale=1.0" />
<title>lightkurve.lightcurve — Lightkurve </title>
<!-- Loaded before other Sphinx assets -->
<link href="../../_static/styles/theme.css?digest=1999514e3f237ded88cf" rel="stylesheet">
<link href="../../_static/styles/pydata-sphinx-theme.css?digest=1999514e3f237ded88cf" rel="stylesheet">
<link rel="stylesheet"
href="../../_static/vendor/fontawesome/5.13.0/css/all.min.css">
<link rel="preload" as="font" type="font/woff2" crossorigin
href="../../_static/vendor/fontawesome/5.13.0/webfonts/fa-solid-900.woff2">
<link rel="preload" as="font" type="font/woff2" crossorigin
href="../../_static/vendor/fontawesome/5.13.0/webfonts/fa-brands-400.woff2">
<link rel="stylesheet" type="text/css" href="../../_static/pygments.css" />
<link rel="stylesheet" type="text/css" href="../../_static/css/custom.css" />
<!-- Pre-loaded scripts that we'll load fully later -->
<link rel="preload" as="script" href="../../_static/scripts/pydata-sphinx-theme.js?digest=1999514e3f237ded88cf">
<script data-url_root="../../" id="documentation_options" src="../../_static/documentation_options.js"></script>
<script src="../../_static/jquery.js"></script>
<script src="../../_static/underscore.js"></script>
<script src="../../_static/doctools.js"></script>
<link rel="index" title="Index" href="../../genindex.html" />
<link rel="search" title="Search" href="../../search.html" />
<meta name="viewport" content="width=device-width, initial-scale=1" />
<meta name="docsearch:language" content="None">
<!-- Google Analytics -->
<script async="" src="https://www.google-analytics.com/analytics.js"></script>
<script>
window.ga = window.ga || function () {
(ga.q = ga.q || []).push(arguments) };
ga.l = +new Date;
ga('create', 'UA-69171-9', 'auto');
ga('set', 'anonymizeIp', true);
ga('send', 'pageview');
</script>
</head>
<body data-spy="scroll" data-target="#bd-toc-nav" data-offset="60">
<div class="container-fluid" id="banner"></div>
<nav class="navbar navbar-light navbar-expand-lg bg-light fixed-top bd-navbar" id="navbar-main"><div class="container-xl">
<div id="navbar-start">
<a class="navbar-brand" href="../../index.html">
<p class="title">Lightkurve v2.4</p>
</a>
</div>
<button class="navbar-toggler" type="button" data-toggle="collapse" data-target="#navbar-collapsible" aria-controls="navbar-collapsible" aria-expanded="false" aria-label="Toggle navigation">
<span class="navbar-toggler-icon"></span>
</button>
<div id="navbar-collapsible" class="col-lg-9 collapse navbar-collapse">
<div id="navbar-center" class="mr-auto">
<div class="navbar-center-item">
<ul id="navbar-main-elements" class="navbar-nav">
<li class="toctree-l1 nav-item">
<a class="reference internal nav-link" href="../../whats-new-v2.html">
What's new?
</a>
</li>
<li class="toctree-l1 nav-item">
<a class="reference internal nav-link" href="../../quickstart.html">
Quickstart
</a>
</li>
<li class="toctree-l1 nav-item">
<a class="reference internal nav-link" href="../../tutorials/index.html">
Tutorials
</a>
</li>
<li class="toctree-l1 nav-item">
<a class="reference internal nav-link" href="../../reference/index.html">
API
</a>
</li>
<li class="toctree-l1 nav-item">
<a class="reference internal nav-link" href="../../about/index.html">
About Lightkurve
</a>
</li>
<li class="toctree-l1 nav-item">
<a class="reference internal nav-link" href="../../development/index.html">
Developing for Lightkurve
</a>
</li>
</ul>
</div>
</div>
<div id="navbar-end">
<div class="navbar-end-item">
<ul id="navbar-icon-links" class="navbar-nav" aria-label="Icon Links">
<li class="nav-item">
<a class="nav-link" href="https://github.com/lightkurve/lightkurve" rel="noopener" target="_blank" title="GitHub"><span><i class="fab fa-github-square"></i></span>
<label class="sr-only">GitHub</label></a>
</li>
</ul>
</div>
</div>
</div>
</div>
</nav>
<div class="container-xl">
<div class="row">
<!-- Only show if we have sidebars configured, else just a small margin -->
<div class="col-12 col-md-3 bd-sidebar">
<div class="sidebar-start-items"><form class="bd-search d-flex align-items-center" action="../../search.html" method="get">
<i class="icon fas fa-search"></i>
<input type="search" class="form-control" name="q" id="search-input" placeholder="Search the docs ..." aria-label="Search the docs ..." autocomplete="off" >
</form><nav class="bd-links" id="bd-docs-nav" aria-label="Main navigation">
<div class="bd-toc-item active">
</div>
</nav>
</div>
<div class="sidebar-end-items">
</div>
</div>
<div class="d-none d-xl-block col-xl-2 bd-toc">
</div>
<main class="col-12 col-md-9 col-xl-7 py-md-5 pl-md-5 pr-md-4 bd-content" role="main">
<div>
<h1>Source code for lightkurve.lightcurve</h1><div class="highlight"><pre>
<span></span><span class="sd">"""Defines LightCurve, KeplerLightCurve, and TessLightCurve."""</span>
<span class="kn">import</span> <span class="nn">os</span>
<span class="kn">import</span> <span class="nn">datetime</span>
<span class="kn">import</span> <span class="nn">logging</span>
<span class="kn">import</span> <span class="nn">warnings</span>
<span class="kn">import</span> <span class="nn">collections</span>
<span class="kn">from</span> <span class="nn">collections.abc</span> <span class="kn">import</span> <span class="n">Sequence</span>
<span class="kn">import</span> <span class="nn">numpy</span> <span class="k">as</span> <span class="nn">np</span>
<span class="kn">from</span> <span class="nn">scipy.signal</span> <span class="kn">import</span> <span class="n">savgol_filter</span>
<span class="kn">from</span> <span class="nn">scipy.interpolate</span> <span class="kn">import</span> <span class="n">interp1d</span>
<span class="kn">import</span> <span class="nn">matplotlib</span>
<span class="kn">from</span> <span class="nn">matplotlib</span> <span class="kn">import</span> <span class="n">pyplot</span> <span class="k">as</span> <span class="n">plt</span>
<span class="kn">from</span> <span class="nn">copy</span> <span class="kn">import</span> <span class="n">deepcopy</span>
<span class="kn">from</span> <span class="nn">astropy.table</span> <span class="kn">import</span> <span class="n">Table</span><span class="p">,</span> <span class="n">Column</span><span class="p">,</span> <span class="n">MaskedColumn</span>
<span class="kn">from</span> <span class="nn">astropy.io</span> <span class="kn">import</span> <span class="n">fits</span>
<span class="kn">from</span> <span class="nn">astropy.time</span> <span class="kn">import</span> <span class="n">TimeBase</span><span class="p">,</span> <span class="n">Time</span><span class="p">,</span> <span class="n">TimeDelta</span>
<span class="kn">from</span> <span class="nn">astropy</span> <span class="kn">import</span> <span class="n">units</span> <span class="k">as</span> <span class="n">u</span>
<span class="kn">from</span> <span class="nn">astropy.units</span> <span class="kn">import</span> <span class="n">Quantity</span>
<span class="kn">from</span> <span class="nn">astropy.timeseries</span> <span class="kn">import</span> <span class="n">TimeSeries</span><span class="p">,</span> <span class="n">aggregate_downsample</span>
<span class="kn">from</span> <span class="nn">astropy.table</span> <span class="kn">import</span> <span class="n">vstack</span>
<span class="kn">from</span> <span class="nn">astropy.stats</span> <span class="kn">import</span> <span class="n">calculate_bin_edges</span>
<span class="kn">from</span> <span class="nn">astropy.utils.decorators</span> <span class="kn">import</span> <span class="n">deprecated</span><span class="p">,</span> <span class="n">deprecated_renamed_argument</span>
<span class="kn">from</span> <span class="nn">astropy.utils.exceptions</span> <span class="kn">import</span> <span class="n">AstropyUserWarning</span>
<span class="kn">from</span> <span class="nn">astropy.utils.masked</span> <span class="kn">import</span> <span class="n">Masked</span>
<span class="kn">from</span> <span class="nn">.</span> <span class="kn">import</span> <span class="n">PACKAGEDIR</span><span class="p">,</span> <span class="n">MPLSTYLE</span>
<span class="kn">from</span> <span class="nn">.utils</span> <span class="kn">import</span> <span class="p">(</span>
<span class="n">running_mean</span><span class="p">,</span>
<span class="n">bkjd_to_astropy_time</span><span class="p">,</span>
<span class="n">btjd_to_astropy_time</span><span class="p">,</span>
<span class="n">validate_method</span><span class="p">,</span>
<span class="n">_query_solar_system_objects</span><span class="p">,</span>
<span class="n">finalize_notebook_url</span>
<span class="p">)</span>
<span class="kn">from</span> <span class="nn">.utils</span> <span class="kn">import</span> <span class="n">LightkurveWarning</span><span class="p">,</span> <span class="n">LightkurveDeprecationWarning</span>
<span class="n">__all__</span> <span class="o">=</span> <span class="p">[</span><span class="s2">"LightCurve"</span><span class="p">,</span> <span class="s2">"KeplerLightCurve"</span><span class="p">,</span> <span class="s2">"TessLightCurve"</span><span class="p">,</span> <span class="s2">"FoldedLightCurve"</span><span class="p">]</span>
<span class="n">log</span> <span class="o">=</span> <span class="n">logging</span><span class="o">.</span><span class="n">getLogger</span><span class="p">(</span><span class="vm">__name__</span><span class="p">)</span>
<span class="n">_HAS_VAR_BINS</span> <span class="o">=</span> <span class="s1">'time_bin_end'</span> <span class="ow">in</span> <span class="n">aggregate_downsample</span><span class="o">.</span><span class="vm">__kwdefaults__</span>
<span class="k">def</span> <span class="nf">_to_unitless_day</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">data</span><span class="p">,</span> <span class="n">Quantity</span><span class="p">):</span>
<span class="k">return</span> <span class="n">data</span><span class="o">.</span><span class="n">to</span><span class="p">(</span><span class="n">u</span><span class="o">.</span><span class="n">day</span><span class="p">)</span><span class="o">.</span><span class="n">value</span>
<span class="k">elif</span> <span class="ow">not</span> <span class="n">np</span><span class="o">.</span><span class="n">isscalar</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="k">return</span> <span class="n">np</span><span class="o">.</span><span class="n">asarray</span><span class="p">([</span><span class="n">_to_unitless_day</span><span class="p">(</span><span class="n">item</span><span class="p">)</span> <span class="k">for</span> <span class="n">item</span> <span class="ow">in</span> <span class="n">data</span><span class="p">])</span><span class="o">.</span><span class="n">flatten</span><span class="p">()</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">return</span> <span class="n">data</span>
<span class="k">def</span> <span class="nf">_is_dict_like</span><span class="p">(</span><span class="n">data1</span><span class="p">):</span>
<span class="k">return</span> <span class="nb">hasattr</span><span class="p">(</span><span class="n">data1</span><span class="p">,</span> <span class="s2">"keys"</span><span class="p">)</span> <span class="ow">and</span> <span class="nb">callable</span><span class="p">(</span><span class="nb">getattr</span><span class="p">(</span><span class="n">data1</span><span class="p">,</span> <span class="s2">"keys"</span><span class="p">))</span>
<span class="k">def</span> <span class="nf">_is_list_like</span><span class="p">(</span><span class="n">data1</span><span class="p">):</span>
<span class="c1"># https://stackoverflow.com/a/37842328</span>
<span class="k">return</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">data1</span><span class="p">,</span> <span class="n">Sequence</span><span class="p">)</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">data1</span><span class="p">,</span> <span class="nb">str</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">_is_np_structured_array</span><span class="p">(</span><span class="n">data1</span><span class="p">):</span>
<span class="k">return</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">data1</span><span class="p">,</span> <span class="n">np</span><span class="o">.</span><span class="n">ndarray</span><span class="p">)</span> <span class="ow">and</span> <span class="n">data1</span><span class="o">.</span><span class="n">dtype</span><span class="o">.</span><span class="n">names</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span>
<div class="viewcode-block" id="LightCurve"><a class="viewcode-back" href="../../reference/api/lightkurve.LightCurve.html#lightkurve.LightCurve">[docs]</a><span class="k">class</span> <span class="nc">LightCurve</span><span class="p">(</span><span class="n">TimeSeries</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""</span>
<span class="sd"> Subclass of AstroPy `~astropy.table.Table` guaranteed to have *time*, *flux*, and *flux_err* columns.</span>
<span class="sd"> Compared to the generic `~astropy.timeseries.TimeSeries` class, `LightCurve`</span>
<span class="sd"> ensures that each object has `time`, `flux`, and `flux_err` columns.</span>
<span class="sd"> These three columns are special for two reasons:</span>
<span class="sd"> 1. they are the key columns upon which all light curve operations operate;</span>
<span class="sd"> 2. they are always present (though they may be populated with ``NaN`` values).</span>
<span class="sd"> `LightCurve` objects also provide user-friendly attribute access to</span>
<span class="sd"> columns and meta data.</span>
<span class="sd"> Parameters</span>
<span class="sd"> ----------</span>
<span class="sd"> data : numpy ndarray, dict, list, `~astropy.table.Table`, or table-like object, optional</span>
<span class="sd"> Data to initialize time series. This does not need to contain the times</span>
<span class="sd"> or fluxes, which can be provided separately, but if it does contain the</span>
<span class="sd"> times and fluxes they should be in columns called ``'time'``,</span>
<span class="sd"> ``'flux'``, and ``'flux_err'`` to be automatically recognized.</span>
<span class="sd"> time : `~astropy.time.Time` or iterable</span>
<span class="sd"> Time values. They can either be given directly as a</span>
<span class="sd"> `~astropy.time.Time` array or as any iterable that initializes the</span>
<span class="sd"> `~astropy.time.Time` class.</span>
<span class="sd"> flux : `~astropy.units.Quantity` or iterable</span>
<span class="sd"> Flux values for every time point.</span>
<span class="sd"> flux_err : `~astropy.units.Quantity` or iterable</span>
<span class="sd"> Uncertainty on each flux data point.</span>
<span class="sd"> **kwargs : dict</span>
<span class="sd"> Additional keyword arguments are passed to `~astropy.table.QTable`.</span>
<span class="sd"> Attributes</span>
<span class="sd"> ----------</span>
<span class="sd"> meta : `dict`</span>
<span class="sd"> meta data associated with the lightcurve. The header of the underlying FITS file (if applicable)</span>
<span class="sd"> is store in this dictionary. By convention, keys in this dictionary are usually in uppercase.</span>
<span class="sd"> Notes</span>
<span class="sd"> -----</span>
<span class="sd"> *Attribute access*: You can access a column or a ``meta`` value directly as an attribute.</span>
<span class="sd"> >>> lc.flux # shortcut for lc['flux'] # doctest: +SKIP</span>
<span class="sd"> >>> lc.sector # shortcut for lc.meta['SECTOR'] # doctest: +SKIP</span>
<span class="sd"> >>> lc.flux = lc.flux * 1.05 # update the values of a column. # doctest: +SKIP</span>
<span class="sd"> In case the given name is both a column name and a key in ``meta``, the column will be returned.</span>
<span class="sd"> Note that you *cannot* create a new column using the attribute interface. If you do so,</span>
<span class="sd"> a new attribute is created instead, and a warning is raised.</span>
<span class="sd"> If you do create such attributes on purpose, please note that the attributes are not carried</span>
<span class="sd"> over when the lightcurve object is copied, or a new lightcurve object is derived</span>
<span class="sd"> based on a copy, e.g., ``normalize()``.</span>
<span class="sd"> Examples</span>
<span class="sd"> --------</span>
<span class="sd"> >>> import lightkurve as lk</span>
<span class="sd"> >>> lc = lk.LightCurve(time=[1, 2, 3, 4], flux=[0.98, 1.02, 1.03, 0.97])</span>
<span class="sd"> >>> lc.time</span>
<span class="sd"> <Time object: scale='tdb' format='jd' value=[1. 2. 3. 4.]></span>
<span class="sd"> >>> lc.flux</span>
<span class="sd"> <Quantity [0.98, 1.02, 1.03, 0.97]></span>
<span class="sd"> >>> lc.bin(time_bin_size=2, time_bin_start=0.5).flux</span>
<span class="sd"> <Quantity [1., 1.]></span>
<span class="sd"> """</span>
<span class="c1"># The constructor of the `TimeSeries` base class will enforce the presence</span>
<span class="c1"># of these columns:</span>
<span class="n">_required_columns</span> <span class="o">=</span> <span class="p">[</span><span class="s2">"time"</span><span class="p">,</span> <span class="s2">"flux"</span><span class="p">,</span> <span class="s2">"flux_err"</span><span class="p">]</span>
<span class="c1"># The following keywords were removed in Lightkurve v2.0.</span>
<span class="c1"># Their use will trigger a warning.</span>
<span class="n">_deprecated_keywords</span> <span class="o">=</span> <span class="p">(</span>
<span class="s2">"targetid"</span><span class="p">,</span>
<span class="s2">"label"</span><span class="p">,</span>
<span class="s2">"time_format"</span><span class="p">,</span>
<span class="s2">"time_scale"</span><span class="p">,</span>
<span class="s2">"flux_unit"</span><span class="p">,</span>
<span class="p">)</span>
<span class="n">_deprecated_column_keywords</span> <span class="o">=</span> <span class="p">[</span>
<span class="s2">"centroid_col"</span><span class="p">,</span>
<span class="s2">"centroid_row"</span><span class="p">,</span>
<span class="s2">"cadenceno"</span><span class="p">,</span>
<span class="s2">"quality"</span><span class="p">,</span>
<span class="p">]</span>
<span class="c1"># If an iterable is passed for ``time``, we will initialize an AstroPy</span>
<span class="c1"># ``Time`` object using the following format and scale:</span>
<span class="n">_default_time_format</span> <span class="o">=</span> <span class="s2">"jd"</span>
<span class="n">_default_time_scale</span> <span class="o">=</span> <span class="s2">"tdb"</span>
<span class="c1"># To emulate pandas, we do not support creating new columns or meta data</span>
<span class="c1"># fields via attribute assignment, and raise a warning in __setattr__ when</span>
<span class="c1"># a new attribute is created. We need to relax this warning during the</span>
<span class="c1"># initial construction of the object using `_new_attributes_relax`.</span>
<span class="n">_new_attributes_relax</span> <span class="o">=</span> <span class="kc">True</span>
<span class="c1"># cf. issue #925</span>
<span class="n">__array_priority__</span> <span class="o">=</span> <span class="mi">100_000</span>
<div class="viewcode-block" id="LightCurve.__init__"><a class="viewcode-back" href="../../reference/api/lightkurve.LightCurve.html#lightkurve.LightCurve.__init__">[docs]</a> <span class="k">def</span> <span class="fm">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">data</span><span class="o">=</span><span class="kc">None</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="n">time</span><span class="o">=</span><span class="kc">None</span><span class="p">,</span> <span class="n">flux</span><span class="o">=</span><span class="kc">None</span><span class="p">,</span> <span class="n">flux_err</span><span class="o">=</span><span class="kc">None</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span>
<span class="c1"># the ` {has,get,set}_time_in_data()`: helpers to handle `data` of different types</span>
<span class="c1"># in some cases, they also need to access kwargs["names"] as well</span>
<span class="k">def</span> <span class="nf">get_time_idx_in</span><span class="p">(</span><span class="n">names</span><span class="p">):</span>
<span class="n">time_indices</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">argwhere</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">asarray</span><span class="p">(</span><span class="n">names</span><span class="p">)</span> <span class="o">==</span> <span class="s2">"time"</span><span class="p">)</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">time_indices</span><span class="p">)</span> <span class="o">></span> <span class="mi">0</span><span class="p">:</span>
<span class="k">return</span> <span class="n">time_indices</span><span class="p">[</span><span class="mi">0</span><span class="p">][</span><span class="mi">0</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">return</span> <span class="kc">None</span>
<span class="k">def</span> <span class="nf">get_time_in_data_list</span><span class="p">():</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">data</span><span class="p">)</span> <span class="o"><</span> <span class="mi">1</span><span class="p">:</span>
<span class="k">return</span> <span class="kc">None</span>
<span class="n">names</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"names"</span><span class="p">)</span>
<span class="k">if</span> <span class="n">names</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="c1"># the first item MUST be time if no names specified</span>
<span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">data</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">TimeBase</span><span class="p">):</span> <span class="c1"># Time or TimeDelta</span>
<span class="k">return</span> <span class="n">data</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">return</span> <span class="kc">None</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">time_idx</span> <span class="o">=</span> <span class="n">get_time_idx_in</span><span class="p">(</span><span class="n">names</span><span class="p">)</span>
<span class="k">if</span> <span class="n">time_idx</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span><span class="p">:</span>
<span class="k">return</span> <span class="n">data</span><span class="p">[</span><span class="n">time_idx</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">return</span> <span class="kc">None</span>
<span class="k">def</span> <span class="nf">set_time_in_data_list</span><span class="p">(</span><span class="n">value</span><span class="p">):</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">data</span><span class="p">)</span> <span class="o"><</span> <span class="mi">1</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"data should be non-empty"</span><span class="p">)</span>
<span class="n">names</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"names"</span><span class="p">)</span>
<span class="k">if</span> <span class="n">names</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="c1"># the first item MUST be time if no names specified</span>
<span class="c1"># this is to support base Table's select columns</span>
<span class="c1"># in __getitem__()</span>
<span class="c1"># https://github.com/astropy/astropy/blob/326435449ad8d859f1abf36800c3fb88d49c27ea/astropy/table/table.py#L1888</span>
<span class="n">data</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span> <span class="o">=</span> <span class="n">value</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">time_idx</span> <span class="o">=</span> <span class="n">get_time_idx_in</span><span class="p">(</span><span class="n">names</span><span class="p">)</span>
<span class="k">if</span> <span class="n">time_idx</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">data</span><span class="p">[</span><span class="n">time_idx</span><span class="p">]</span> <span class="o">=</span> <span class="n">value</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"data should have time column"</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">get_time_in_data_np_structured_array</span><span class="p">():</span>
<span class="k">if</span> <span class="n">data</span><span class="o">.</span><span class="n">dtype</span><span class="o">.</span><span class="n">names</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span> <span class="c1"># no labeled filed, not a structured array</span>
<span class="k">return</span> <span class="kc">None</span>
<span class="k">if</span> <span class="s2">"time"</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">data</span><span class="o">.</span><span class="n">dtype</span><span class="o">.</span><span class="n">names</span><span class="p">:</span>
<span class="k">return</span> <span class="kc">None</span>
<span class="k">return</span> <span class="n">data</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span>
<span class="k">def</span> <span class="nf">remove_time_from_data_np_structured_array</span><span class="p">():</span>
<span class="k">if</span> <span class="n">data</span><span class="o">.</span><span class="n">dtype</span><span class="o">.</span><span class="n">names</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"data should be a numpy structured array"</span><span class="p">)</span>
<span class="k">if</span> <span class="s2">"time"</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">data</span><span class="o">.</span><span class="n">dtype</span><span class="o">.</span><span class="n">names</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"data should have a time field"</span><span class="p">)</span>
<span class="n">filtered_names</span> <span class="o">=</span> <span class="p">[</span><span class="n">n</span> <span class="k">for</span> <span class="n">n</span> <span class="ow">in</span> <span class="n">data</span><span class="o">.</span><span class="n">dtype</span><span class="o">.</span><span class="n">names</span> <span class="k">if</span> <span class="n">n</span> <span class="o">!=</span> <span class="s2">"time"</span><span class="p">]</span>
<span class="k">return</span> <span class="n">data</span><span class="p">[</span><span class="n">filtered_names</span><span class="p">]</span>
<span class="k">def</span> <span class="nf">has_time_in_data</span><span class="p">():</span>
<span class="w"> </span><span class="sd">"""Check if the data has a column with the name"""</span>
<span class="k">if</span> <span class="n">data</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="k">return</span> <span class="kc">False</span>
<span class="k">elif</span> <span class="n">_is_dict_like</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># data is a dict-like object with keys</span>
<span class="k">return</span> <span class="s2">"time"</span> <span class="ow">in</span> <span class="n">data</span><span class="o">.</span><span class="n">keys</span><span class="p">()</span>
<span class="k">elif</span> <span class="n">_is_list_like</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># case data is a list-like object (a list of columns, etc.)</span>
<span class="k">return</span> <span class="n">get_time_in_data_list</span><span class="p">()</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span>
<span class="k">elif</span> <span class="n">_is_np_structured_array</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># case numpy structured array (supported by base TimeSeries)</span>
<span class="c1"># https://numpy.org/doc/stable/user/basics.rec.html</span>
<span class="k">return</span> <span class="n">get_time_in_data_np_structured_array</span><span class="p">()</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span><span class="sa">f</span><span class="s2">"Unsupported type for time in data: </span><span class="si">{</span><span class="nb">type</span><span class="p">(</span><span class="n">data</span><span class="p">)</span><span class="si">}</span><span class="s2">"</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">get_time_in_data</span><span class="p">():</span>
<span class="k">if</span> <span class="n">_is_dict_like</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># data is a dict-like object with keys</span>
<span class="k">return</span> <span class="n">data</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span>
<span class="k">elif</span> <span class="n">_is_list_like</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="k">return</span> <span class="n">get_time_in_data_list</span><span class="p">()</span>
<span class="k">elif</span> <span class="n">_is_np_structured_array</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="k">return</span> <span class="n">get_time_in_data_np_structured_array</span><span class="p">()</span>
<span class="k">else</span><span class="p">:</span>
<span class="c1"># should never reach here. It'd have been caught by `has_time_in()``</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"Unsupported type for time in data"</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">set_time_in_data</span><span class="p">(</span><span class="n">value</span><span class="p">):</span>
<span class="k">if</span> <span class="n">_is_dict_like</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># data is a dict-like object with keys</span>
<span class="n">data</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span> <span class="o">=</span> <span class="n">value</span>
<span class="k">elif</span> <span class="n">_is_list_like</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="n">set_time_in_data_list</span><span class="p">(</span><span class="n">value</span><span class="p">)</span>
<span class="k">elif</span> <span class="n">_is_np_structured_array</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># astropy Time cannot be assigned to a column in np structured array</span>
<span class="c1"># we have special codepath handling it outside this function</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"Setting Time instances to np structured array is not supported"</span><span class="p">)</span>
<span class="k">else</span><span class="p">:</span>
<span class="c1"># should never reach here. It'd have been caught by `has_time_in()``</span>
<span class="k">raise</span> <span class="ne">AssertionError</span><span class="p">(</span><span class="s2">"Unsupported type for time in data"</span><span class="p">)</span>
<span class="c1"># Delay checking for required columns until the end</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_required_columns_relax</span> <span class="o">=</span> <span class="kc">True</span>
<span class="c1"># Lightkurve v1.x supported passing time, flux, and flux_err as</span>
<span class="c1"># positional arguments. We support it here for backwards compatibility.</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="ow">in</span> <span class="p">[</span><span class="mi">1</span><span class="p">,</span> <span class="mi">2</span><span class="p">]:</span>
<span class="n">warnings</span><span class="o">.</span><span class="n">warn</span><span class="p">(</span>
<span class="s2">"passing flux as a positional argument is deprecated"</span>
<span class="s2">", please use ``flux=...`` instead."</span><span class="p">,</span>
<span class="n">LightkurveDeprecationWarning</span><span class="p">,</span>
<span class="p">)</span>
<span class="n">time</span> <span class="o">=</span> <span class="n">data</span>
<span class="n">flux</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span>
<span class="n">data</span> <span class="o">=</span> <span class="kc">None</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">2</span><span class="p">:</span>
<span class="n">flux_err</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span>
<span class="c1"># For backwards compatibility with Lightkurve v1.x,</span>
<span class="c1"># we support passing deprecated keywords via **kwargs.</span>
<span class="n">deprecated_kws</span> <span class="o">=</span> <span class="p">{}</span>
<span class="k">for</span> <span class="n">kw</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">_deprecated_keywords</span><span class="p">:</span>
<span class="k">if</span> <span class="n">kw</span> <span class="ow">in</span> <span class="n">kwargs</span><span class="p">:</span>
<span class="n">deprecated_kws</span><span class="p">[</span><span class="n">kw</span><span class="p">]</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="n">kw</span><span class="p">)</span>
<span class="n">deprecated_column_kws</span> <span class="o">=</span> <span class="p">{}</span>
<span class="k">for</span> <span class="n">kw</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">_deprecated_column_keywords</span><span class="p">:</span>
<span class="k">if</span> <span class="n">kw</span> <span class="ow">in</span> <span class="n">kwargs</span><span class="p">:</span>
<span class="n">deprecated_column_kws</span><span class="p">[</span><span class="n">kw</span><span class="p">]</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="n">kw</span><span class="p">)</span>
<span class="c1"># If `time` is passed as keyword argument, we populate it with integer numbers</span>
<span class="k">if</span> <span class="n">data</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">or</span> <span class="ow">not</span> <span class="n">has_time_in_data</span><span class="p">():</span>
<span class="k">if</span> <span class="n">time</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">and</span> <span class="n">flux</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">time</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">arange</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="n">flux</span><span class="p">))</span>
<span class="c1"># We are tolerant of missing time format</span>
<span class="k">if</span> <span class="n">time</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">time</span><span class="p">,</span> <span class="p">(</span><span class="n">Time</span><span class="p">,</span> <span class="n">TimeDelta</span><span class="p">)):</span>
<span class="c1"># Lightkurve v1.x supported specifying the time_format</span>
<span class="c1"># as a constructor kwarg</span>
<span class="n">time</span> <span class="o">=</span> <span class="n">Time</span><span class="p">(</span>
<span class="n">time</span><span class="p">,</span>
<span class="nb">format</span><span class="o">=</span><span class="n">deprecated_kws</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"time_format"</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">_default_time_format</span><span class="p">),</span>
<span class="n">scale</span><span class="o">=</span><span class="n">deprecated_kws</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"time_scale"</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">_default_time_scale</span><span class="p">),</span>
<span class="p">)</span>
<span class="c1"># Also be tolerant of missing time format if time is passed via `data`</span>
<span class="k">if</span> <span class="n">data</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span> <span class="ow">and</span> <span class="n">has_time_in_data</span><span class="p">():</span>
<span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">get_time_in_data</span><span class="p">(),</span> <span class="p">(</span><span class="n">Time</span><span class="p">,</span> <span class="n">TimeDelta</span><span class="p">)):</span>
<span class="n">tmp_time</span> <span class="o">=</span> <span class="n">Time</span><span class="p">(</span>
<span class="n">get_time_in_data</span><span class="p">(),</span>
<span class="nb">format</span><span class="o">=</span><span class="n">deprecated_kws</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"time_format"</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">_default_time_format</span><span class="p">),</span>
<span class="n">scale</span><span class="o">=</span><span class="n">deprecated_kws</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"time_scale"</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">_default_time_scale</span><span class="p">),</span>
<span class="p">)</span>
<span class="k">if</span> <span class="n">_is_np_structured_array</span><span class="p">(</span><span class="n">data</span><span class="p">):</span>
<span class="c1"># special case for np structured array</span>
<span class="c1"># one cannot set a `Time` instance to it</span>
<span class="c1"># so we set the time to the `time` param, and take it out of data</span>
<span class="n">time</span> <span class="o">=</span> <span class="n">tmp_time</span>
<span class="n">data</span> <span class="o">=</span> <span class="n">remove_time_from_data_np_structured_array</span><span class="p">()</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">set_time_in_data</span><span class="p">(</span><span class="n">tmp_time</span><span class="p">)</span>
<span class="c1"># Allow overriding the required columns</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_required_columns</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s2">"_required_columns"</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">_required_columns</span><span class="p">)</span>
<span class="c1"># Call the SampledTimeSeries constructor.</span>
<span class="c1"># Disable required columns for now; we'll check those later.</span>
<span class="n">tmp</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">_required_columns</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_required_columns</span> <span class="o">=</span> <span class="p">[]</span>
<span class="nb">super</span><span class="p">()</span><span class="o">.</span><span class="fm">__init__</span><span class="p">(</span><span class="n">data</span><span class="o">=</span><span class="n">data</span><span class="p">,</span> <span class="n">time</span><span class="o">=</span><span class="n">time</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">)</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_required_columns</span> <span class="o">=</span> <span class="n">tmp</span>
<span class="c1"># For some operations, an empty time series needs to be created, then</span>
<span class="c1"># columns added one by one. We should check that when columns are added</span>
<span class="c1"># manually, time is added first and is of the right type.</span>
<span class="k">if</span> <span class="n">data</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">and</span> <span class="n">time</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">and</span> <span class="n">flux</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">and</span> <span class="n">flux_err</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_required_columns_relax</span> <span class="o">=</span> <span class="kc">True</span>
<span class="k">return</span>
<span class="c1"># Load `time`, `flux`, and `flux_err` from the table as local variable names</span>
<span class="n">time</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span> <span class="c1"># super().__init__() guarantees this is a column</span>
<span class="k">if</span> <span class="s2">"flux"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">colnames</span><span class="p">:</span>
<span class="k">if</span> <span class="n">flux</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">flux</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span>
<span class="sa">f</span><span class="s2">"'flux' has been given both in the `data` table and as a keyword argument"</span>
<span class="p">)</span>
<span class="k">if</span> <span class="s2">"flux_err"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">colnames</span><span class="p">:</span>
<span class="k">if</span> <span class="n">flux_err</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">flux_err</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span>
<span class="sa">f</span><span class="s2">"'flux_err' has been given both in the `data` table and as a keyword argument"</span>
<span class="p">)</span>
<span class="c1"># Ensure `flux` and `flux_err` are populated with NaNs if missing</span>
<span class="k">if</span> <span class="n">flux</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">and</span> <span class="n">time</span> <span class="ow">is</span> <span class="ow">not</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">flux</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">empty</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="n">time</span><span class="p">))</span>
<span class="n">flux</span><span class="p">[:]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">nan</span>
<span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">flux</span><span class="p">,</span> <span class="n">Quantity</span><span class="p">):</span>
<span class="n">flux</span> <span class="o">=</span> <span class="n">Quantity</span><span class="p">(</span><span class="n">flux</span><span class="p">,</span> <span class="n">deprecated_kws</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s2">"flux_unit"</span><span class="p">))</span>
<span class="k">if</span> <span class="n">flux_err</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">flux_err</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">empty</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="n">flux</span><span class="p">))</span>
<span class="n">flux_err</span><span class="p">[:]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">nan</span>
<span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">flux_err</span><span class="p">,</span> <span class="n">Quantity</span><span class="p">):</span>
<span class="n">flux_err</span> <span class="o">=</span> <span class="n">Quantity</span><span class="p">(</span><span class="n">flux_err</span><span class="p">,</span> <span class="n">flux</span><span class="o">.</span><span class="n">unit</span><span class="p">)</span>
<span class="c1"># Backwards compatibility with Lightkurve v1.x</span>
<span class="c1"># Ensure attributes are set if passed via deprecated kwargs</span>
<span class="k">for</span> <span class="n">kw</span> <span class="ow">in</span> <span class="n">deprecated_kws</span><span class="p">:</span>
<span class="k">if</span> <span class="n">kw</span> <span class="ow">not</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="p">:</span>
<span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="p">[</span><span class="n">kw</span><span class="o">.</span><span class="n">upper</span><span class="p">()]</span> <span class="o">=</span> <span class="n">deprecated_kws</span><span class="p">[</span><span class="n">kw</span><span class="p">]</span>
<span class="c1"># Ensure all required columns are in the right order</span>
<span class="k">with</span> <span class="bp">self</span><span class="o">.</span><span class="n">_delay_required_column_checks</span><span class="p">():</span>
<span class="k">for</span> <span class="n">idx</span><span class="p">,</span> <span class="n">col</span> <span class="ow">in</span> <span class="nb">enumerate</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_required_columns</span><span class="p">):</span>
<span class="k">if</span> <span class="n">col</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">colnames</span><span class="p">:</span>
<span class="bp">self</span><span class="o">.</span><span class="n">remove_column</span><span class="p">(</span><span class="n">col</span><span class="p">)</span>
<span class="bp">self</span><span class="o">.</span><span class="n">add_column</span><span class="p">(</span><span class="nb">locals</span><span class="p">()[</span><span class="n">col</span><span class="p">],</span> <span class="n">index</span><span class="o">=</span><span class="n">idx</span><span class="p">,</span> <span class="n">name</span><span class="o">=</span><span class="n">col</span><span class="p">)</span>
<span class="c1"># Ensure columns are set if passed via deprecated kwargs</span>
<span class="k">for</span> <span class="n">kw</span> <span class="ow">in</span> <span class="n">deprecated_column_kws</span><span class="p">:</span>
<span class="k">if</span> <span class="n">kw</span> <span class="ow">not</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">meta</span> <span class="ow">and</span> <span class="n">kw</span> <span class="ow">not</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">:</span>
<span class="bp">self</span><span class="o">.</span><span class="n">add_column</span><span class="p">(</span><span class="n">deprecated_column_kws</span><span class="p">[</span><span class="n">kw</span><span class="p">],</span> <span class="n">name</span><span class="o">=</span><span class="n">kw</span><span class="p">)</span>
<span class="c1"># Ensure flux and flux_err have the same units</span>
<span class="k">if</span> <span class="bp">self</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span><span class="o">.</span><span class="n">unit</span> <span class="o">!=</span> <span class="bp">self</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span><span class="o">.</span><span class="n">unit</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span><span class="s2">"flux and flux_err must have the same units"</span><span class="p">)</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_new_attributes_relax</span> <span class="o">=</span> <span class="kc">False</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_required_columns_relax</span> <span class="o">=</span> <span class="kc">False</span>
<span class="bp">self</span><span class="o">.</span><span class="n">_check_required_columns</span><span class="p">()</span></div>
<span class="k">def</span> <span class="fm">__getattr__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">name</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""Expose all columns and meta keywords as attributes."""</span>
<span class="k">if</span> <span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">:</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="n">name</span><span class="p">]</span>
<span class="k">elif</span> <span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__class__</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">:</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__class__</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="n">name</span><span class="p">]</span><span class="o">.</span><span class="fm">__get__</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span>
<span class="k">elif</span> <span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">:</span>
<span class="k">return</span> <span class="bp">self</span><span class="p">[</span><span class="n">name</span><span class="p">]</span>
<span class="k">elif</span> <span class="s2">"_meta"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">:</span>
<span class="k">if</span> <span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">]:</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">][</span><span class="n">name</span><span class="p">]</span>
<span class="k">elif</span> <span class="n">name</span><span class="o">.</span><span class="n">upper</span><span class="p">()</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">]:</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">][</span><span class="n">name</span><span class="o">.</span><span class="n">upper</span><span class="p">()]</span>
<span class="k">raise</span> <span class="ne">AttributeError</span><span class="p">(</span><span class="sa">f</span><span class="s2">"object has no attribute </span><span class="si">{</span><span class="n">name</span><span class="si">}</span><span class="s2">"</span><span class="p">)</span>
<span class="k">def</span> <span class="fm">__setattr__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">name</span><span class="p">,</span> <span class="n">value</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""To get copied, attributes have to be stored in the meta dictionary!"""</span>
<span class="n">to_set_as_attr</span> <span class="o">=</span> <span class="kc">False</span>
<span class="k">if</span> <span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">:</span>
<span class="n">to_set_as_attr</span> <span class="o">=</span> <span class="kc">True</span>
<span class="k">elif</span> <span class="n">name</span> <span class="o">==</span> <span class="s2">"time"</span><span class="p">:</span>
<span class="bp">self</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span> <span class="o">=</span> <span class="n">value</span> <span class="c1"># astropy will convert value to Time if needed</span>
<span class="k">elif</span> <span class="p">(</span><span class="s2">"columns"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">)</span> <span class="ow">and</span> <span class="p">(</span><span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"columns"</span><span class="p">]):</span>
<span class="bp">self</span><span class="o">.</span><span class="n">replace_column</span><span class="p">(</span><span class="n">name</span><span class="p">,</span> <span class="n">value</span><span class="p">)</span>
<span class="k">elif</span> <span class="s2">"_meta"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">:</span>
<span class="k">if</span> <span class="n">name</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">]:</span>
<span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">][</span><span class="n">name</span><span class="p">]</span> <span class="o">=</span> <span class="n">value</span>
<span class="k">elif</span> <span class="n">name</span><span class="o">.</span><span class="n">upper</span><span class="p">()</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">]:</span>
<span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span><span class="p">[</span><span class="s2">"_meta"</span><span class="p">][</span><span class="n">name</span><span class="o">.</span><span class="n">upper</span><span class="p">()]</span> <span class="o">=</span> <span class="n">value</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">to_set_as_attr</span> <span class="o">=</span> <span class="kc">True</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">to_set_as_attr</span> <span class="o">=</span> <span class="kc">True</span>
<span class="k">if</span> <span class="n">to_set_as_attr</span><span class="p">:</span>
<span class="k">if</span> <span class="p">(</span>
<span class="n">name</span> <span class="ow">not</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="vm">__dict__</span>
<span class="ow">and</span> <span class="ow">not</span> <span class="n">name</span><span class="o">.</span><span class="n">startswith</span><span class="p">(</span><span class="s2">"_"</span><span class="p">)</span>
<span class="ow">and</span> <span class="ow">not</span> <span class="bp">self</span><span class="o">.</span><span class="n">_new_attributes_relax</span>
<span class="ow">and</span> <span class="n">name</span> <span class="o">!=</span> <span class="s1">'meta'</span>
<span class="p">):</span>
<span class="n">warnings</span><span class="o">.</span><span class="n">warn</span><span class="p">(</span>
<span class="p">(</span>
<span class="s2">"Lightkurve doesn't allow columns or meta values to be created via a new attribute name."</span>
<span class="s2">"A new attribute is created. It will not be carried over when the object is copied."</span>
<span class="s2">" - see https://docs.lightkurve.org/reference/api/lightkurve.LightCurve.html"</span>
<span class="p">),</span>
<span class="ne">UserWarning</span><span class="p">,</span>
<span class="n">stacklevel</span><span class="o">=</span><span class="mi">2</span><span class="p">,</span>
<span class="p">)</span>
<span class="nb">super</span><span class="p">()</span><span class="o">.</span><span class="fm">__setattr__</span><span class="p">(</span><span class="n">name</span><span class="p">,</span> <span class="n">value</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">_repr_simple_</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">-></span> <span class="nb">str</span><span class="p">:</span>
<span class="w"> </span><span class="sd">"""Returns a simple __repr__.</span>
<span class="sd"> Used by `LightCurveCollection`.</span>
<span class="sd"> """</span>
<span class="n">result</span> <span class="o">=</span> <span class="sa">f</span><span class="s2">"<</span><span class="si">{</span><span class="bp">self</span><span class="o">.</span><span class="vm">__class__</span><span class="o">.</span><span class="vm">__name__</span><span class="si">}</span><span class="s2">"</span>
<span class="k">if</span> <span class="s2">"LABEL"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="p">:</span>
<span class="n">result</span> <span class="o">+=</span> <span class="sa">f</span><span class="s2">" LABEL=</span><span class="se">\"</span><span class="si">{</span><span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s1">'LABEL'</span><span class="p">)</span><span class="si">}</span><span class="se">\"</span><span class="s2">"</span>
<span class="k">for</span> <span class="n">kw</span> <span class="ow">in</span> <span class="p">[</span><span class="s2">"QUARTER"</span><span class="p">,</span> <span class="s2">"CAMPAIGN"</span><span class="p">,</span> <span class="s2">"SECTOR"</span><span class="p">,</span> <span class="s2">"AUTHOR"</span><span class="p">,</span> <span class="s2">"FLUX_ORIGIN"</span><span class="p">]:</span>
<span class="k">if</span> <span class="n">kw</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="p">:</span>
<span class="n">result</span> <span class="o">+=</span> <span class="sa">f</span><span class="s2">" </span><span class="si">{</span><span class="n">kw</span><span class="si">}</span><span class="s2">=</span><span class="si">{</span><span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="n">kw</span><span class="p">)</span><span class="si">}</span><span class="s2">"</span>
<span class="n">result</span> <span class="o">+=</span> <span class="s2">">"</span>
<span class="k">return</span> <span class="n">result</span>
<span class="k">def</span> <span class="nf">_base_repr_</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">html</span><span class="o">=</span><span class="kc">False</span><span class="p">,</span> <span class="n">descr_vals</span><span class="o">=</span><span class="kc">None</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""Defines the description shown by `__repr__` and `_html_repr_`."""</span>
<span class="k">if</span> <span class="n">descr_vals</span> <span class="ow">is</span> <span class="kc">None</span><span class="p">:</span>
<span class="n">descr_vals</span> <span class="o">=</span> <span class="p">[</span><span class="bp">self</span><span class="o">.</span><span class="vm">__class__</span><span class="o">.</span><span class="vm">__name__</span><span class="p">]</span>
<span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">masked</span><span class="p">:</span>
<span class="n">descr_vals</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="s2">"masked=True"</span><span class="p">)</span>
<span class="n">descr_vals</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="s2">"length=</span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">)))</span>
<span class="k">if</span> <span class="s2">"LABEL"</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="p">:</span>
<span class="n">descr_vals</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="sa">f</span><span class="s2">"LABEL=</span><span class="se">\"</span><span class="si">{</span><span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="s1">'LABEL'</span><span class="p">)</span><span class="si">}</span><span class="se">\"</span><span class="s2">"</span><span class="p">)</span>
<span class="k">for</span> <span class="n">kw</span> <span class="ow">in</span> <span class="p">[</span><span class="s2">"QUARTER"</span><span class="p">,</span> <span class="s2">"CAMPAIGN"</span><span class="p">,</span> <span class="s2">"SECTOR"</span><span class="p">,</span> <span class="s2">"AUTHOR"</span><span class="p">,</span> <span class="s2">"FLUX_ORIGIN"</span><span class="p">]:</span>
<span class="k">if</span> <span class="n">kw</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="p">:</span>
<span class="n">descr_vals</span><span class="o">.</span><span class="n">append</span><span class="p">(</span><span class="sa">f</span><span class="s2">"</span><span class="si">{</span><span class="n">kw</span><span class="si">}</span><span class="s2">=</span><span class="si">{</span><span class="bp">self</span><span class="o">.</span><span class="n">meta</span><span class="o">.</span><span class="n">get</span><span class="p">(</span><span class="n">kw</span><span class="p">)</span><span class="si">}</span><span class="s2">"</span><span class="p">)</span>
<span class="k">return</span> <span class="nb">super</span><span class="p">()</span><span class="o">.</span><span class="n">_base_repr_</span><span class="p">(</span><span class="n">html</span><span class="o">=</span><span class="n">html</span><span class="p">,</span> <span class="n">descr_vals</span><span class="o">=</span><span class="n">descr_vals</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">)</span>
<span class="c1"># Define `time`, `flux`, `flux_err` as class attributes to enable IDE</span>
<span class="c1"># of these required columns auto-completion.</span>
<span class="nd">@property</span>
<span class="k">def</span> <span class="nf">time</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">-></span> <span class="n">Time</span><span class="p">:</span>
<span class="w"> </span><span class="sd">"""Time values stored as an AstroPy `~astropy.time.Time` object."""</span>
<span class="k">return</span> <span class="bp">self</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span>
<span class="nd">@time</span><span class="o">.</span><span class="n">setter</span>
<span class="k">def</span> <span class="nf">time</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">time</span><span class="p">):</span>
<span class="bp">self</span><span class="p">[</span><span class="s2">"time"</span><span class="p">]</span> <span class="o">=</span> <span class="n">time</span>
<span class="nd">@property</span>
<span class="k">def</span> <span class="nf">flux</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">-></span> <span class="n">Quantity</span><span class="p">:</span>
<span class="w"> </span><span class="sd">"""Brightness values stored as an AstroPy `~astropy.units.Quantity` object."""</span>
<span class="k">return</span> <span class="bp">self</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span>
<span class="nd">@flux</span><span class="o">.</span><span class="n">setter</span>
<span class="k">def</span> <span class="nf">flux</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">flux</span><span class="p">):</span>
<span class="bp">self</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span> <span class="o">=</span> <span class="n">flux</span>
<span class="nd">@property</span>
<span class="k">def</span> <span class="nf">flux_err</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">-></span> <span class="n">Quantity</span><span class="p">:</span>
<span class="w"> </span><span class="sd">"""Brightness uncertainties stored as an AstroPy `~astropy.units.Quantity` object."""</span>
<span class="k">return</span> <span class="bp">self</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span>
<span class="nd">@flux_err</span><span class="o">.</span><span class="n">setter</span>
<span class="k">def</span> <span class="nf">flux_err</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">flux_err</span><span class="p">):</span>
<span class="bp">self</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span> <span class="o">=</span> <span class="n">flux_err</span>
<div class="viewcode-block" id="LightCurve.select_flux"><a class="viewcode-back" href="../../reference/api/lightkurve.LightCurve.select_flux.html#lightkurve.LightCurve.select_flux">[docs]</a> <span class="k">def</span> <span class="nf">select_flux</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">flux_column</span><span class="p">,</span> <span class="n">flux_err_column</span><span class="o">=</span><span class="kc">None</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""Assign a different column to be the flux column.</span>
<span class="sd"> This method returns a copy of the LightCurve in which the ``flux``</span>
<span class="sd"> and ``flux_err`` columns have been replaced by the values contained</span>
<span class="sd"> in a different column.</span>
<span class="sd"> Parameters</span>
<span class="sd"> ----------</span>
<span class="sd"> flux_column : str</span>
<span class="sd"> Name of the column that should become the 'flux' column.</span>
<span class="sd"> flux_err_column : str or `None`</span>
<span class="sd"> Name of the column that should become the 'flux_err' column.</span>
<span class="sd"> By default, the column will be used that is obtained by adding the</span>
<span class="sd"> suffix "_err" to the value of ``flux_column``. If such a</span>
<span class="sd"> column does not exist, ``flux_err`` will be populated with NaN values.</span>
<span class="sd"> Returns</span>
<span class="sd"> -------</span>
<span class="sd"> lc : LightCurve</span>
<span class="sd"> Copy of the ``LightCurve`` object with the new flux values assigned.</span>
<span class="sd"> Examples</span>
<span class="sd"> --------</span>
<span class="sd"> You can use this function to change the flux data on which most Lightkurve</span>
<span class="sd"> features operate. For example, to view a periodogram based on the "sap_flux"</span>
<span class="sd"> column in a TESS light curve, use::</span>
<span class="sd"> >>> lc.select_flux("sap_flux").to_periodogram("lombscargle").plot() # doctest: +SKIP</span>
<span class="sd"> """</span>
<span class="c1"># Input validation</span>
<span class="k">if</span> <span class="n">flux_column</span> <span class="ow">not</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span><span class="sa">f</span><span class="s2">"'</span><span class="si">{</span><span class="n">flux_column</span><span class="si">}</span><span class="s2">' is not a column"</span><span class="p">)</span>
<span class="k">if</span> <span class="n">flux_err_column</span> <span class="ow">and</span> <span class="n">flux_err_column</span> <span class="ow">not</span> <span class="ow">in</span> <span class="bp">self</span><span class="o">.</span><span class="n">columns</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span><span class="sa">f</span><span class="s2">"'</span><span class="si">{</span><span class="n">flux_err_column</span><span class="si">}</span><span class="s2">' is not a column"</span><span class="p">)</span>
<span class="n">lc</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="n">flux_column</span><span class="p">]</span>
<span class="k">if</span> <span class="n">flux_err_column</span><span class="p">:</span> <span class="c1"># not None</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="n">flux_err_column</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="c1"># if `flux_err_column` is unspecified, we attempt to use</span>
<span class="c1"># f"{flux_column}_err" if it exists</span>
<span class="n">flux_err_column</span> <span class="o">=</span> <span class="sa">f</span><span class="s2">"</span><span class="si">{</span><span class="n">flux_column</span><span class="si">}</span><span class="s2">_err"</span>
<span class="k">if</span> <span class="n">flux_err_column</span> <span class="ow">in</span> <span class="n">lc</span><span class="o">.</span><span class="n">columns</span><span class="p">:</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="n">flux_err_column</span><span class="p">]</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">][:]</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">nan</span>
<span class="n">lc</span><span class="o">.</span><span class="n">meta</span><span class="p">[</span><span class="s1">'FLUX_ORIGIN'</span><span class="p">]</span> <span class="o">=</span> <span class="n">flux_column</span>
<span class="n">normalized_new_flux</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span><span class="o">.</span><span class="n">unit</span> <span class="ow">is</span> <span class="kc">None</span> <span class="ow">or</span> <span class="n">lc</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span><span class="o">.</span><span class="n">unit</span> <span class="ow">is</span> <span class="n">u</span><span class="o">.</span><span class="n">dimensionless_unscaled</span>
<span class="c1"># Note: here we assume unitless flux means it's normalized</span>
<span class="c1"># it's not exactly true in many constructed lightcurves in unit test</span>
<span class="c1"># but the assumption should hold for any real world use cases, e.g. TESS QLP</span>
<span class="k">if</span> <span class="n">normalized_new_flux</span><span class="p">:</span>
<span class="n">lc</span><span class="o">.</span><span class="n">meta</span><span class="p">[</span><span class="s2">"NORMALIZED"</span><span class="p">]</span> <span class="o">=</span> <span class="n">normalized_new_flux</span>
<span class="k">else</span><span class="p">:</span>
<span class="c1"># remove it altogether.</span>
<span class="c1"># Setting to False would suffice;</span>
<span class="c1"># but in typical non-normalized LC, the header will not be there at all.</span>
<span class="n">lc</span><span class="o">.</span><span class="n">meta</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s2">"NORMALIZED"</span><span class="p">,</span> <span class="kc">None</span><span class="p">)</span>
<span class="k">return</span> <span class="n">lc</span></div>
<span class="c1"># Define deprecated attributes for compatibility with Lightkurve v1.x:</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span>
<span class="s2">"2.0"</span><span class="p">,</span> <span class="n">alternative</span><span class="o">=</span><span class="s2">"time.format"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span>
<span class="p">)</span>
<span class="k">def</span> <span class="nf">time_format</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">time</span><span class="o">.</span><span class="n">format</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span>
<span class="s2">"2.0"</span><span class="p">,</span> <span class="n">alternative</span><span class="o">=</span><span class="s2">"time.scale"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span>
<span class="p">)</span>
<span class="k">def</span> <span class="nf">time_scale</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">time</span><span class="o">.</span><span class="n">scale</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span><span class="s2">"2.0"</span><span class="p">,</span> <span class="n">alternative</span><span class="o">=</span><span class="s2">"time"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">astropy_time</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">time</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span>
<span class="s2">"2.0"</span><span class="p">,</span> <span class="n">alternative</span><span class="o">=</span><span class="s2">"flux.unit"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span>
<span class="p">)</span>
<span class="k">def</span> <span class="nf">flux_unit</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span><span class="o">.</span><span class="n">unit</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span><span class="s2">"2.0"</span><span class="p">,</span> <span class="n">alternative</span><span class="o">=</span><span class="s2">"flux"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">flux_quantity</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span>
<span class="s2">"2.0"</span><span class="p">,</span>
<span class="n">alternative</span><span class="o">=</span><span class="s2">"fits.open(lc.filename)"</span><span class="p">,</span>
<span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span><span class="p">,</span>
<span class="p">)</span>
<span class="k">def</span> <span class="nf">hdu</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">with</span> <span class="n">fits</span><span class="o">.</span><span class="n">open</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">filename</span><span class="p">)</span> <span class="k">as</span> <span class="n">hdulist</span><span class="p">:</span>
<span class="n">hdulist</span> <span class="o">=</span> <span class="n">hdulist</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="k">return</span> <span class="n">hdulist</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span><span class="s2">"2.0"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">SAP_FLUX</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""A copy of the light curve in which `lc.flux = lc.sap_flux`</span>
<span class="sd"> and `lc.flux_err = lc.sap_flux_err`. It is provided for backwards-</span>
<span class="sd"> compatibility with Lightkurve v1.x and will be removed soon."""</span>
<span class="n">lc</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="s2">"sap_flux"</span><span class="p">]</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="s2">"sap_flux_err"</span><span class="p">]</span>
<span class="k">return</span> <span class="n">lc</span>
<span class="nd">@property</span>
<span class="nd">@deprecated</span><span class="p">(</span><span class="s2">"2.0"</span><span class="p">,</span> <span class="n">warning_type</span><span class="o">=</span><span class="n">LightkurveDeprecationWarning</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">PDCSAP_FLUX</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""A copy of the light curve in which `lc.flux = lc.pdcsap_flux`</span>
<span class="sd"> and `lc.flux_err = lc.pdcsap_flux_err`. It is provided for backwards-</span>
<span class="sd"> compatibility with Lightkurve v1.x and will be removed soon."""</span>
<span class="n">lc</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="s2">"pdcsap_flux"</span><span class="p">]</span>
<span class="n">lc</span><span class="p">[</span><span class="s2">"flux_err"</span><span class="p">]</span> <span class="o">=</span> <span class="n">lc</span><span class="p">[</span><span class="s2">"pdcsap_flux_err"</span><span class="p">]</span>
<span class="k">return</span> <span class="n">lc</span>
<span class="k">def</span> <span class="fm">__add__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="n">newlc</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">other</span><span class="p">,</span> <span class="n">LightCurve</span><span class="p">):</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">!=</span> <span class="nb">len</span><span class="p">(</span><span class="n">other</span><span class="p">):</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span>
<span class="s2">"Cannot add LightCurve objects because "</span>
<span class="s2">"they do not have equal length (</span><span class="si">{}</span><span class="s2"> vs </span><span class="si">{}</span><span class="s2">)."</span>
<span class="s2">""</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">),</span> <span class="nb">len</span><span class="p">(</span><span class="n">other</span><span class="p">))</span>
<span class="p">)</span>
<span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">time</span> <span class="o">!=</span> <span class="n">other</span><span class="o">.</span><span class="n">time</span><span class="p">):</span>
<span class="n">warnings</span><span class="o">.</span><span class="n">warn</span><span class="p">(</span>
<span class="s2">"Two LightCurve objects with inconsistent time "</span>
<span class="s2">"values are being added."</span><span class="p">,</span>
<span class="n">LightkurveWarning</span><span class="p">,</span>
<span class="p">)</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span> <span class="o">+</span> <span class="n">other</span><span class="o">.</span><span class="n">flux</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">hypot</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">flux_err</span><span class="p">,</span> <span class="n">other</span><span class="o">.</span><span class="n">flux_err</span><span class="p">)</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span> <span class="o">+</span> <span class="n">other</span>
<span class="k">return</span> <span class="n">newlc</span>
<span class="k">def</span> <span class="fm">__radd__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="fm">__add__</span><span class="p">(</span><span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="fm">__sub__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="fm">__add__</span><span class="p">(</span><span class="o">-</span><span class="mi">1</span> <span class="o">*</span> <span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="fm">__rsub__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="p">(</span><span class="o">-</span><span class="mi">1</span> <span class="o">*</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="fm">__add__</span><span class="p">(</span><span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="fm">__mul__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="n">newlc</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">other</span><span class="p">,</span> <span class="n">LightCurve</span><span class="p">):</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">!=</span> <span class="nb">len</span><span class="p">(</span><span class="n">other</span><span class="p">):</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span>
<span class="s2">"Cannot multiply LightCurve objects because "</span>
<span class="s2">"they do not have equal length (</span><span class="si">{}</span><span class="s2"> vs </span><span class="si">{}</span><span class="s2">)."</span>
<span class="s2">""</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">),</span> <span class="nb">len</span><span class="p">(</span><span class="n">other</span><span class="p">))</span>
<span class="p">)</span>
<span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">time</span> <span class="o">!=</span> <span class="n">other</span><span class="o">.</span><span class="n">time</span><span class="p">):</span>
<span class="n">warnings</span><span class="o">.</span><span class="n">warn</span><span class="p">(</span>
<span class="s2">"Two LightCurve objects with inconsistent time "</span>
<span class="s2">"values are being multiplied."</span><span class="p">,</span>
<span class="n">LightkurveWarning</span><span class="p">,</span>
<span class="p">)</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span> <span class="o">*</span> <span class="n">other</span><span class="o">.</span><span class="n">flux</span>
<span class="c1"># Applying standard uncertainty propagation, cf.</span>
<span class="c1"># https://en.wikipedia.org/wiki/Propagation_of_uncertainty#Example_formulae</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">=</span> <span class="nb">abs</span><span class="p">(</span><span class="n">newlc</span><span class="o">.</span><span class="n">flux</span><span class="p">)</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">hypot</span><span class="p">(</span>
<span class="bp">self</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">/</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span><span class="p">,</span> <span class="n">other</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">/</span> <span class="n">other</span><span class="o">.</span><span class="n">flux</span>
<span class="p">)</span>
<span class="k">elif</span> <span class="nb">isinstance</span><span class="p">(</span>
<span class="n">other</span><span class="p">,</span> <span class="p">(</span><span class="n">u</span><span class="o">.</span><span class="n">UnitBase</span><span class="p">,</span> <span class="n">u</span><span class="o">.</span><span class="n">FunctionUnitBase</span><span class="p">)</span>
<span class="p">):</span> <span class="c1"># cf. astropy/issues/6517</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="n">other</span> <span class="o">*</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">=</span> <span class="n">other</span> <span class="o">*</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux_err</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="n">other</span> <span class="o">*</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">=</span> <span class="nb">abs</span><span class="p">(</span><span class="n">other</span><span class="p">)</span> <span class="o">*</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux_err</span>
<span class="k">return</span> <span class="n">newlc</span>
<span class="k">def</span> <span class="fm">__rmul__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="fm">__mul__</span><span class="p">(</span><span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="fm">__truediv__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="fm">__mul__</span><span class="p">(</span><span class="mf">1.0</span> <span class="o">/</span> <span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="fm">__rtruediv__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="n">newlc</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span>
<span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">other</span><span class="p">,</span> <span class="n">LightCurve</span><span class="p">):</span>
<span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">)</span> <span class="o">!=</span> <span class="nb">len</span><span class="p">(</span><span class="n">other</span><span class="p">):</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span>
<span class="s2">"Cannot divide LightCurve objects because "</span>
<span class="s2">"they do not have equal length (</span><span class="si">{}</span><span class="s2"> vs </span><span class="si">{}</span><span class="s2">)."</span>
<span class="s2">""</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="bp">self</span><span class="p">),</span> <span class="nb">len</span><span class="p">(</span><span class="n">other</span><span class="p">))</span>
<span class="p">)</span>
<span class="k">if</span> <span class="n">np</span><span class="o">.</span><span class="n">any</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">time</span> <span class="o">!=</span> <span class="n">other</span><span class="o">.</span><span class="n">time</span><span class="p">):</span>
<span class="n">warnings</span><span class="o">.</span><span class="n">warn</span><span class="p">(</span>
<span class="s2">"Two LightCurve objects with inconsistent time "</span>
<span class="s2">"values are being divided."</span><span class="p">,</span>
<span class="n">LightkurveWarning</span><span class="p">,</span>
<span class="p">)</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="n">other</span><span class="o">.</span><span class="n">flux</span> <span class="o">/</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">=</span> <span class="nb">abs</span><span class="p">(</span><span class="n">newlc</span><span class="o">.</span><span class="n">flux</span><span class="p">)</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">hypot</span><span class="p">(</span>
<span class="bp">self</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">/</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span><span class="p">,</span> <span class="n">other</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">/</span> <span class="n">other</span><span class="o">.</span><span class="n">flux</span>
<span class="p">)</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux</span> <span class="o">=</span> <span class="n">other</span> <span class="o">/</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux</span>
<span class="n">newlc</span><span class="o">.</span><span class="n">flux_err</span> <span class="o">=</span> <span class="nb">abs</span><span class="p">((</span><span class="n">other</span> <span class="o">*</span> <span class="bp">self</span><span class="o">.</span><span class="n">flux_err</span><span class="p">)</span> <span class="o">/</span> <span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">flux</span> <span class="o">**</span> <span class="mi">2</span><span class="p">))</span>
<span class="k">return</span> <span class="n">newlc</span>
<span class="k">def</span> <span class="nf">__div__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="fm">__truediv__</span><span class="p">(</span><span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">__rdiv__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">other</span><span class="p">):</span>
<span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="fm">__rtruediv__</span><span class="p">(</span><span class="n">other</span><span class="p">)</span>
<span class="k">def</span> <span class="nf">show_properties</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""Prints a description of all non-callable attributes.</span>
<span class="sd"> Prints in order of type (ints, strings, lists, arrays, others).</span>
<span class="sd"> """</span>
<span class="n">attrs</span> <span class="o">=</span> <span class="p">{}</span>
<span class="n">deprecated_properties</span> <span class="o">=</span> <span class="nb">list</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_deprecated_keywords</span><span class="p">)</span>
<span class="n">deprecated_properties</span> <span class="o">+=</span> <span class="p">[</span>
<span class="s2">"flux_quantity"</span><span class="p">,</span>
<span class="s2">"SAP_FLUX"</span><span class="p">,</span>
<span class="s2">"PDCSAP_FLUX"</span><span class="p">,</span>
<span class="s2">"astropy_time"</span><span class="p">,</span>
<span class="s2">"hdu"</span><span class="p">,</span>
<span class="p">]</span>
<span class="k">for</span> <span class="n">attr</span> <span class="ow">in</span> <span class="nb">dir</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span>
<span class="k">if</span> <span class="ow">not</span> <span class="n">attr</span><span class="o">.</span><span class="n">startswith</span><span class="p">(</span><span class="s2">"_"</span><span class="p">)</span> <span class="ow">and</span> <span class="n">attr</span> <span class="ow">not</span> <span class="ow">in</span> <span class="n">deprecated_properties</span><span class="p">:</span>
<span class="k">try</span><span class="p">:</span>
<span class="n">res</span> <span class="o">=</span> <span class="nb">getattr</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">attr</span><span class="p">)</span>
<span class="k">except</span> <span class="ne">Exception</span><span class="p">:</span>
<span class="k">continue</span>
<span class="k">if</span> <span class="nb">callable</span><span class="p">(</span><span class="n">res</span><span class="p">):</span>
<span class="k">continue</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">]</span> <span class="o">=</span> <span class="p">{</span><span class="s2">"res"</span><span class="p">:</span> <span class="n">res</span><span class="p">}</span>
<span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">res</span><span class="p">,</span> <span class="nb">int</span><span class="p">):</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"</span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">res</span><span class="p">)</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"int"</span>
<span class="k">elif</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">res</span><span class="p">,</span> <span class="n">np</span><span class="o">.</span><span class="n">ndarray</span><span class="p">):</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"array </span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">res</span><span class="o">.</span><span class="n">shape</span><span class="p">)</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"array"</span>
<span class="k">elif</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">res</span><span class="p">,</span> <span class="nb">list</span><span class="p">):</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"list length </span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="nb">len</span><span class="p">(</span><span class="n">res</span><span class="p">))</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"list"</span>
<span class="k">elif</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">res</span><span class="p">,</span> <span class="nb">str</span><span class="p">):</span>
<span class="k">if</span> <span class="n">res</span> <span class="o">==</span> <span class="s2">""</span><span class="p">:</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"</span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="s2">"None"</span><span class="p">)</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"</span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">res</span><span class="p">)</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"str"</span>
<span class="k">elif</span> <span class="n">attr</span> <span class="o">==</span> <span class="s2">"wcs"</span><span class="p">:</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"astropy.wcs.wcs.WCS"</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"other"</span>
<span class="k">else</span><span class="p">:</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"print"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"</span><span class="si">{}</span><span class="s2">"</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="nb">type</span><span class="p">(</span><span class="n">res</span><span class="p">))</span>
<span class="n">attrs</span><span class="p">[</span><span class="n">attr</span><span class="p">][</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">=</span> <span class="s2">"other"</span>
<span class="n">output</span> <span class="o">=</span> <span class="n">Table</span><span class="p">(</span><span class="n">names</span><span class="o">=</span><span class="p">[</span><span class="s2">"Attribute"</span><span class="p">,</span> <span class="s2">"Description"</span><span class="p">],</span> <span class="n">dtype</span><span class="o">=</span><span class="p">[</span><span class="nb">object</span><span class="p">,</span> <span class="nb">object</span><span class="p">])</span>
<span class="n">idx</span> <span class="o">=</span> <span class="mi">0</span>
<span class="n">types</span> <span class="o">=</span> <span class="p">[</span><span class="s2">"int"</span><span class="p">,</span> <span class="s2">"str"</span><span class="p">,</span> <span class="s2">"list"</span><span class="p">,</span> <span class="s2">"array"</span><span class="p">,</span> <span class="s2">"other"</span><span class="p">]</span>
<span class="k">for</span> <span class="n">typ</span> <span class="ow">in</span> <span class="n">types</span><span class="p">:</span>
<span class="k">for</span> <span class="n">attr</span><span class="p">,</span> <span class="n">dic</span> <span class="ow">in</span> <span class="n">attrs</span><span class="o">.</span><span class="n">items</span><span class="p">():</span>
<span class="k">if</span> <span class="n">dic</span><span class="p">[</span><span class="s2">"type"</span><span class="p">]</span> <span class="o">==</span> <span class="n">typ</span><span class="p">:</span>
<span class="n">output</span><span class="o">.</span><span class="n">add_row</span><span class="p">([</span><span class="n">attr</span><span class="p">,</span> <span class="n">dic</span><span class="p">[</span><span class="s2">"print"</span><span class="p">]])</span>
<span class="n">idx</span> <span class="o">+=</span> <span class="mi">1</span>
<span class="n">output</span><span class="o">.</span><span class="n">pprint</span><span class="p">(</span><span class="n">max_lines</span><span class="o">=-</span><span class="mi">1</span><span class="p">,</span> <span class="n">max_width</span><span class="o">=-</span><span class="mi">1</span><span class="p">)</span>
<div class="viewcode-block" id="LightCurve.append"><a class="viewcode-back" href="../../reference/api/lightkurve.LightCurve.append.html#lightkurve.LightCurve.append">[docs]</a> <span class="k">def</span> <span class="nf">append</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">others</span><span class="p">,</span> <span class="n">inplace</span><span class="o">=</span><span class="kc">False</span><span class="p">):</span>
<span class="w"> </span><span class="sd">"""Append one or more other `LightCurve` object(s) to this one.</span>
<span class="sd"> Parameters</span>
<span class="sd"> ----------</span>
<span class="sd"> others : `LightCurve`, or list of `LightCurve`</span>
<span class="sd"> Light curve(s) to be appended to the current one.</span>
<span class="sd"> inplace : bool</span>
<span class="sd"> If True, change the current `LightCurve` instance in place instead</span>
<span class="sd"> of creating and returning a new one. Defaults to False.</span>
<span class="sd"> Returns</span>
<span class="sd"> -------</span>
<span class="sd"> new_lc : `LightCurve`</span>
<span class="sd"> Light curve which has the other light curves appened to it.</span>
<span class="sd"> """</span>
<span class="k">if</span> <span class="n">inplace</span><span class="p">:</span>
<span class="k">raise</span> <span class="ne">ValueError</span><span class="p">(</span>
<span class="s2">"the `inplace` parameter is no longer supported "</span>
<span class="s2">"as of Lightkurve v2.0"</span>
<span class="p">)</span>
<span class="k">if</span> <span class="ow">not</span> <span class="nb">hasattr</span><span class="p">(</span><span class="n">others</span><span class="p">,</span> <span class="s2">"__iter__"</span><span class="p">):</span>
<span class="n">others</span> <span class="o">=</span> <span class="p">(</span><span class="n">others</span><span class="p">,)</span>
<span class="c1"># Re-use LightCurveCollection.stitch() to avoid code duplication</span>
<span class="kn">from</span> <span class="nn">.collections</span> <span class="kn">import</span> <span class="n">LightCurveCollection</span> <span class="c1"># avoid circular import</span>
<span class="k">return</span> <span class="n">LightCurveCollection</span><span class="p">((</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">others</span><span class="p">))</span><span class="o">.</span><span class="n">stitch</span><span class="p">(</span><span class="n">corrector_func</span><span class="o">=</span><span class="kc">None</span><span class="p">)</span></div>
<div class="viewcode-block" id="LightCurve.flatten"><a class="viewcode-back" href="../../reference/api/lightkurve.LightCurve.flatten.html#lightkurve.LightCurve.flatten">[docs]</a> <span class="k">def</span> <span class="nf">flatten</span><span class="p">(</span>
<span class="bp">self</span><span class="p">,</span>
<span class="n">window_length</span><span class="o">=</span><span class="mi">101</span><span class="p">,</span>
<span class="n">polyorder</span><span class="o">=</span><span class="mi">2</span><span class="p">,</span>
<span class="n">return_trend</span><span class="o">=</span><span class="kc">False</span><span class="p">,</span>