Skip to content

Commit 819b545

Browse files
authored
Merge pull request #118 from karllark/update_calc_av
update A(V) calc to use latest - Decleir et al. (2022)
2 parents de1f372 + c310a85 commit 819b545

File tree

2 files changed

+23
-12
lines changed

2 files changed

+23
-12
lines changed

measure_extinction/extdata.py

Lines changed: 19 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -502,17 +502,18 @@ def calc_EBV(self):
502502
else:
503503
self.columns["EBV"] = (self.exts["BAND"][bindx], self.uncs["BAND"][bindx])
504504

505-
def calc_AV(self, akav=0.112):
505+
def calc_AV(self, akav=0.105):
506506
"""
507507
Calculate A(V) from the observed extinction curve:
508508
- fit a powerlaw to the SpeX extinction curve, if available
509509
- otherwise: extrapolate the K-band extinction
510510
511511
Parameters
512512
----------
513-
akav : float [default = 0.112]
513+
akav : float
514514
Value of A(K)/A(V)
515-
default is from Rieke & Lebofsky (1985)
515+
default is from Decleir et al. (2022)
516+
for Rieke & Lebofsky (1985) use akav=0.112
516517
van de Hulst No. 15 curve has A(K)/A(V) = 0.0885
517518
518519
Returns
@@ -539,23 +540,33 @@ def calc_AV(self, akav=0.112):
539540
avunc = np.absolute(av * (self.uncs["BAND"][kindx] / ekv))
540541
self.columns["AV"] = (av, avunc)
541542

542-
def calc_AV_JHK(self):
543+
def calc_AV_JHK(
544+
self,
545+
ref_waves=np.array([1.25, 1.6, 2.2]) * u.micron,
546+
ref_alav=[0.269, 0.163, 0.105],
547+
):
543548
"""
544549
Calculate A(V) from the observed extinction curve:
545550
- extrapolate from J, H, & K photometry
546-
- assumes functional form from Rieke, Rieke, & Paul (1989)
551+
552+
Parameters
553+
----------
554+
ref_waves : floats
555+
wavelengths for reference values (default = JHK)
556+
ref_alav : floats
557+
A(lambda)/A(V) values for reference
558+
default is for JHK from Decleir et al. (2022)
559+
for Rieke, Rieke, & Paul (1989) use ref_alav=[0.2815534, 0.17475728, 0.11197411],
547560
548561
Returns
549562
-------
550563
Updates self.columns["AV"]
551564
"""
552565
# J, H, K
553-
rrp89_waves = np.array([1.25, 1.6, 2.2]) * u.micron
554-
rrp89_alav = [0.2815534, 0.17475728, 0.11197411]
555566

556567
avs = []
557568
avs_unc = []
558-
for cwave, calav in zip(rrp89_waves, rrp89_alav):
569+
for cwave, calav in zip(ref_waves, ref_alav):
559570
dwaves = np.absolute(self.waves["BAND"] - cwave)
560571
kindx = dwaves.argmin()
561572
if dwaves[kindx] < 0.1 * u.micron:

measure_extinction/tests/test_EbvAvRv.py

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -25,8 +25,8 @@ def test_calc_av():
2525
text.calc_AV()
2626
av = text.columns["AV"]
2727
assert isinstance(av, tuple)
28-
assert np.isclose(av[0], 0.56306, atol=1e-3)
29-
assert np.isclose(av[1], 0.02252, atol=1e-3)
28+
assert np.isclose(av[0], 0.55866, atol=1e-3)
29+
assert np.isclose(av[1], 0.02235, atol=1e-3)
3030

3131

3232
def test_calc_rv_fromelv():
@@ -38,8 +38,8 @@ def test_calc_rv_fromelv():
3838
text.calc_RV()
3939
rv = text.columns["RV"]
4040
assert isinstance(rv, tuple)
41-
assert np.isclose(rv[0], 3.37837, atol=1e-3)
42-
assert np.isclose(rv[1], 0.207647, atol=1e-3)
41+
assert np.isclose(rv[0], 3.351955, atol=1e-3)
42+
assert np.isclose(rv[1], 0.206023, atol=1e-3)
4343

4444

4545
def test_calc_rv_fromebvav():

0 commit comments

Comments
 (0)