Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

from __future__ import division 

import numpy as np 

import palpy 

from lsst.sims.utils.CodeUtilities import _validate_inputs 

from lsst.sims.utils import arcsecFromRadians, cartesianFromSpherical, sphericalFromCartesian 

from lsst.sims.utils import radiansFromArcsec 

from lsst.sims.utils import haversine 

 

__all__ = ["_solarRaDec", "solarRaDec", 

"_distanceToSun", "distanceToSun", 

"applyRefraction", "refractionCoefficients", 

"_applyPrecession", "applyPrecession", 

"_applyProperMotion", "applyProperMotion", 

"_appGeoFromICRS", "appGeoFromICRS", 

"_icrsFromAppGeo", "icrsFromAppGeo", 

"_observedFromAppGeo", "observedFromAppGeo", 

"_appGeoFromObserved", "appGeoFromObserved", 

"_observedFromICRS", "observedFromICRS", 

"_icrsFromObserved", "icrsFromObserved"] 

 

 

def _solarRaDec(mjd, epoch=2000.0): 

""" 

Return the RA and Dec of the Sun in radians 

 

@param [in] mjd is the date represented as a 

ModifiedJulianDate object. 

 

@param [in] epoch is the mean epoch of the coordinate system 

(default is 2000.0) 

 

@param [out] RA of Sun in radians 

 

@param [out] Dec of Sun in radians 

""" 

 

params = palpy.mappa(epoch, mjd.TDB) 

# params[4:7] is a unit vector pointing from the Sun 

# to the Earth (see the docstring for palpy.mappa) 

 

return palpy.dcc2s(-1.0 * params[4:7]) 

 

 

def solarRaDec(mjd, epoch=2000.0): 

""" 

Return the RA and Dec of the Sun in degrees 

 

@param [in] mjd is the date represented as a 

ModifiedJulianDate object. 

 

@param [in] epoch is the mean epoch of the coordinate system 

(default is 2000.0) 

 

@param [out] RA of Sun in degrees 

 

@param [out] Dec of Sun in degress 

""" 

 

solarRA, solarDec = _solarRaDec(mjd, epoch=epoch) 

return np.degrees(solarRA), np.degrees(solarDec) 

 

 

def _distanceToSun(ra, dec, mjd, epoch=2000.0): 

""" 

Calculate the distance from an (ra, dec) point to the Sun (in radians). 

 

@param [in] ra in radians 

 

@param [in] dec in radians 

 

@param [in] mjd is the date represented as a 

ModifiedJulianDate object. 

 

@param [in] epoch is the epoch of the coordinate system 

(default is 2000.0) 

 

@param [out] distance on the sky to the Sun in radians 

""" 

 

sunRa, sunDec = _solarRaDec(mjd, epoch=epoch) 

 

return haversine(ra, dec, sunRa, sunDec) 

 

 

def distanceToSun(ra, dec, mjd, epoch=2000.0): 

""" 

Calculate the distance from an (ra, dec) point to the Sun (in degrees). 

 

@param [in] ra in degrees 

 

@param [in] dec in degrees 

 

@param [in] mjd is the date represented as a 

ModifiedJulianDate object. 

 

@param [in] epoch is the epoch of the coordinate system 

(default is 2000.0) 

 

@param [out] distance on the sky to the Sun in degrees 

""" 

 

return np.degrees(_distanceToSun(np.radians(ra), np.radians(dec), mjd, epoch=epoch)) 

 

 

def refractionCoefficients(wavelength=0.5, site=None): 

""" Calculate the refraction using PAL's refco routine 

 

This calculates the refraction at 2 angles and derives a tanz and tan^3z 

coefficient for subsequent quick calculations. Good for zenith distances < 76 degrees 

 

@param [in] wavelength is effective wavelength in microns (default 0.5) 

 

@param [in] site is an instantiation of the Site class defined in 

sims_utils/../Site.py 

 

One should call PAL refz to apply the coefficients calculated here 

 

""" 

precision = 1.e-10 

 

if site is None: 

raise RuntimeError("Cannot call refractionCoefficients; no site information") 

 

# TODO the latitude in refco needs to be astronomical latitude, 

# not geodetic latitude 

_refcoOutput = palpy.refco(site.height, 

site.temperature_kelvin, 

site.pressure, 

site.humidity, 

wavelength, 

site.latitude_rad, 

site.lapseRate, 

precision) 

 

return _refcoOutput[0], _refcoOutput[1] 

 

 

def applyRefraction(zenithDistance, tanzCoeff, tan3zCoeff): 

""" Calculted refracted Zenith Distance 

 

uses the quick PAL refco routine which approximates the refractin calculation 

 

@param [in] zenithDistance is unrefracted zenith distance of the source in radians. 

Can either be a number or a numpy array (not a list). 

 

@param [in] tanzCoeff is the first output from refractionCoefficients (above) 

 

@param [in] tan3zCoeff is the second output from refractionCoefficients (above) 

 

@param [out] refractedZenith is the refracted zenith distance in radians 

 

""" 

 

if isinstance(zenithDistance, list): 

raise RuntimeError("You passed a list of zenithDistances to " + 

"applyRefraction. The method won't know how to " + 

"handle that. Pass a numpy array.") 

 

if isinstance(zenithDistance, np.ndarray): 

refractedZenith = palpy.refzVector( 

zenithDistance, tanzCoeff, tan3zCoeff) 

else: 

refractedZenith = palpy.refz(zenithDistance, tanzCoeff, tan3zCoeff) 

 

return refractedZenith 

 

 

def applyPrecession(ra, dec, epoch=2000.0, mjd=None): 

""" 

applyPrecession() applies precesion and nutation to coordinates between two epochs. 

Accepts RA and dec as inputs. Returns corrected RA and dec (in degrees). 

 

Assumes FK5 as the coordinate system 

units: ra_in (degrees), dec_in (degrees) 

 

The precession-nutation matrix is calculated by the palpy.prenut method 

which uses the IAU 2006/2000A model 

 

@param [in] ra in degrees 

 

@param [in] dec in degrees 

 

@param [in] epoch is the epoch of the mean equinox (in years; default 2000) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the RA 

corrected for precession and nutation and the second row is the 

Dec corrected for precession and nutation (both in degrees) 

 

""" 

 

output = _applyPrecession(np.radians(ra), np.radians(dec), 

epoch=epoch, mjd=mjd) 

 

return np.degrees(output) 

 

 

def _applyPrecession(ra, dec, epoch=2000.0, mjd=None): 

""" 

_applyPrecession() applies precesion and nutation to coordinates between two epochs. 

Accepts RA and dec as inputs. Returns corrected RA and dec (in radians). 

 

Assumes FK5 as the coordinate system 

units: ra_in (radians), dec_in (radians) 

 

The precession-nutation matrix is calculated by the palpy.prenut method 

which uses the IAU 2006/2000A model 

 

@param [in] ra in radians 

 

@param [in] dec in radians 

 

@param [in] epoch is the epoch of the mean equinox (in years; default 2000) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the RA 

corrected for precession and nutation and the second row is the 

Dec corrected for precession and nutation (both in radians) 

""" 

 

if hasattr(ra, '__len__'): 

if len(ra) != len(dec): 

raise RuntimeError("You supplied %d RAs but %d Decs to applyPrecession" % 

(len(ra), len(dec))) 

 

if mjd is None: 

raise RuntimeError("You need to supply applyPrecession with an mjd") 

 

# Determine the precession and nutation 

# palpy.prenut takes the julian epoch for the mean coordinates 

# and the MJD for the the true coordinates 

# 

# TODO it is not specified what this MJD should be (i.e. in which 

# time system it should be reckoned) 

rmat = palpy.prenut(epoch, mjd.TT) 

 

# Apply rotation matrix 

xyz = cartesianFromSpherical(ra, dec) 

xyz = np.dot(rmat, xyz.transpose()).transpose() 

 

raOut, decOut = sphericalFromCartesian(xyz) 

return np.array([raOut, decOut]) 

 

 

def applyProperMotion(ra, dec, pm_ra, pm_dec, parallax, v_rad, 

epoch=2000.0, mjd=None): 

"""Applies proper motion between two epochs. 

 

units: ra (degrees), dec (degrees), pm_ra (arcsec/year), pm_dec 

(arcsec/year), parallax (arcsec), v_rad (km/sec, positive if receding), 

epoch (Julian years) 

 

Returns corrected ra and dec (in radians) 

 

The function palpy.pm does not work properly if the parallax is below 

0.00045 arcseconds 

 

@param [in] ra in degrees. Can be a number or a numpy array (not a list). 

 

@param [in] dec in degrees. Can be a number or a numpy array (not a list). 

 

@param [in] pm_ra is ra proper motion multiplied by cos(Dec) in arcsec/year. 

Can be a number or a numpy array (not a list). 

 

@param [in] pm_dec is dec proper motion in arcsec/year. 

Can be a number or a numpy array (not a list). 

 

@param [in] parallax in arcsec. Can be a number or a numpy array (not a list). 

 

@param [in] v_rad is radial velocity in km/sec (positive if the object is receding). 

Can be a number or a numpy array (not a list). 

 

@param [in] epoch is epoch in Julian years (default: 2000.0) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the RA corrected 

for proper motion and the second row is the Dec corrected for proper motion 

(both in degrees) 

""" 

 

output = _applyProperMotion(np.radians(ra), np.radians(dec), 

radiansFromArcsec(pm_ra), 

radiansFromArcsec(pm_dec), 

radiansFromArcsec(parallax), 

v_rad, epoch=epoch, mjd=mjd) 

 

return np.degrees(output) 

 

 

def _applyProperMotion(ra, dec, pm_ra, pm_dec, parallax, v_rad, 

epoch=2000.0, mjd=None): 

"""Applies proper motion between two epochs. 

 

units: ra (radians), dec (radians), pm_ra (radians/year), pm_dec 

(radians/year), parallax (radians), v_rad (km/sec, positive if receding), 

epoch (Julian years) 

 

Returns corrected ra and dec (in radians) 

 

The function palpy.pm does not work properly if the parallax is below 

0.00045 arcseconds 

 

@param [in] ra in radians. Can be a number or a numpy array (not a list). 

 

@param [in] dec in radians. Can be a number or a numpy array (not a list). 

 

@param [in] pm_ra is ra proper motion multiplied by cos(Dec) in radians/year. 

Can be a number or a numpy array (not a list). 

 

@param [in] pm_dec is dec proper motion in radians/year. 

Can be a number or a numpy array (not a list). 

 

@param [in] parallax in radians. Can be a number or a numpy array (not a list). 

 

@param [in] v_rad is radial velocity in km/sec (positive if the object is receding). 

Can be a number or a numpy array (not a list). 

 

@param [in] epoch is epoch in Julian years (default: 2000.0) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the RA corrected 

for proper motion and the second row is the Dec corrected for proper motion 

(both in radians) 

 

""" 

 

if (isinstance(ra, list) or isinstance(dec, list) or 

isinstance(pm_ra, list) or isinstance(pm_dec, list) or 

isinstance(parallax, list) or isinstance(v_rad, list)): 

 

raise RuntimeError("You tried to pass lists to applyPm. " + 

"The method does not know how to handle lists. " + 

"Use numpy arrays.") 

 

if mjd is None: 

raise RuntimeError("cannot call applyProperMotion; mjd is None") 

 

parallaxArcsec = arcsecFromRadians(parallax) 

# convert to Arcsec because that is what PALPY expects 

 

# Generate Julian epoch from MJD 

# 

# 19 November 2015 

# I am assuming here that the time scale should be 

# Terrestrial Dynamical Time (TT), since that is used 

# as the independent variable for apparent geocentric 

# ephemerides 

julianEpoch = palpy.epj(mjd.TT) 

 

# because PAL and ERFA expect proper motion in terms of "coordinate 

# angle; not true angle" (as stated in erfa/starpm.c documentation) 

pm_ra_corrected = pm_ra / np.cos(dec) 

 

if isinstance(ra, np.ndarray): 

if ((len(ra) != len(dec) or 

len(ra) != len(pm_ra) or 

len(ra) != len(pm_dec) or 

len(ra) != len(parallaxArcsec)) or 

len(ra) != len(v_rad)): 

 

raise RuntimeError("You passed: " + 

"%d RAs, " % len(ra) + 

"%d Dec, " % len(dec) + 

"%d pm_ras, " % len(pm_ra) + 

"%d pm_decs, " % len(pm_dec) + 

"%d parallaxes, " % len(parallaxArcsec) + 

"%d v_rads " % len(v_rad) + 

"to applyPm; those numbers need to be identical.") 

 

raOut, decOut = palpy.pmVector(ra, dec, pm_ra_corrected, pm_dec, 

parallaxArcsec, v_rad, epoch, julianEpoch) 

else: 

raOut, decOut = palpy.pm(ra, dec, pm_ra_corrected, pm_dec, parallaxArcsec, v_rad, epoch, julianEpoch) 

 

return np.array([raOut, decOut]) 

 

 

def appGeoFromICRS(ra, dec, pm_ra=None, pm_dec=None, parallax=None, 

v_rad=None, epoch=2000.0, mjd=None): 

""" 

Convert the mean position (RA, Dec) in the International Celestial Reference 

System (ICRS) to the mean apparent geocentric position 

 

units: ra (degrees), dec (degrees), pm_ra (arcsec/year), pm_dec 

(arcsec/year), parallax (arcsec), v_rad (km/sec; positive if receding), 

epoch (Julian years) 

 

@param [in] ra in degrees (ICRS). Can be a numpy array or a number. 

 

@param [in] dec in degrees (ICRS). Can be a numpy array or a number. 

 

@param [in] pm_ra is ra proper motion multiplied by cos(Dec) in arcsec/year 

 

@param [in] pm_dec is dec proper motion in arcsec/year 

 

@param [in] parallax in arcsec 

 

@param [in] v_rad is radial velocity in km/sec (positive if the object is receding) 

 

@param [in] epoch is the julian epoch (in years) of the equinox against which to 

measure RA (default: 2000.0) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the apparent 

geocentric RA and the second row is the apparent geocentric Dec (both in degrees) 

""" 

 

if pm_ra is not None: 

pm_ra_in = radiansFromArcsec(pm_ra) 

else: 

pm_ra_in = None 

 

if pm_dec is not None: 

pm_dec_in = radiansFromArcsec(pm_dec) 

else: 

pm_dec_in = None 

 

if parallax is not None: 

px_in = radiansFromArcsec(parallax) 

else: 

px_in = None 

 

output = _appGeoFromICRS(np.radians(ra), np.radians(dec), 

pm_ra=pm_ra_in, pm_dec=pm_dec_in, 

parallax=px_in, v_rad=v_rad, epoch=epoch, mjd=mjd) 

 

return np.degrees(output) 

 

 

def _appGeoFromICRS(ra, dec, pm_ra=None, pm_dec=None, parallax=None, 

v_rad=None, epoch=2000.0, mjd=None): 

""" 

Convert the mean position (RA, Dec) in the International Celestial Reference 

System (ICRS) to the mean apparent geocentric position 

 

units: ra (radians), dec (radians), pm_ra (radians/year), pm_dec 

(radians/year), parallax (radians), v_rad (km/sec; positive if receding), 

epoch (Julian years) 

 

@param [in] ra in radians (ICRS). Can be a numpy array or a number. 

 

@param [in] dec in radians (ICRS). Can be a numpy array or a number. 

 

@param [in] pm_ra is ra proper motion multiplied by cos(Dec) in radians/year. 

Can be a numpy array or a number or None. 

 

@param [in] pm_dec is dec proper motion in radians/year. 

Can be a numpy array or a number or None. 

 

@param [in] parallax in radians. Can be a numpy array or a number or None. 

 

@param [in] v_rad is radial velocity in km/sec (positive if the object is receding). 

Can be a numpy array or a number or None. 

 

@param [in] epoch is the julian epoch (in years) of the equinox against which to 

measure RA (default: 2000.0) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the apparent 

geocentric RAand the second row is the apparent geocentric Dec (both in radians) 

""" 

 

if mjd is None: 

raise RuntimeError("cannot call appGeoFromICRS; mjd is None") 

 

include_px = False 

 

if (pm_ra is not None or pm_dec is not None or 

v_rad is not None or parallax is not None): 

 

include_px = True 

 

if isinstance(ra, np.ndarray): 

fill_value = np.zeros(len(ra), dtype=float) 

else: 

fill_value = 0.0 

 

if pm_ra is None: 

pm_ra = fill_value 

 

if pm_dec is None: 

pm_dec = fill_value 

 

if v_rad is None: 

v_rad = fill_value 

 

if parallax is None: 

parallax = fill_value 

 

are_arrays = _validate_inputs([ra, dec, pm_ra, pm_dec, v_rad, parallax], 

['ra', 'dec', 'pm_ra', 'pm_dec', 'v_rad', 

'parallax'], 

"appGeoFromICRS") 

else: 

are_arrays = _validate_inputs([ra, dec], ['ra', 'dec'], "appGeoFromICRS") 

 

 

# Define star independent mean to apparent place parameters 

# palpy.mappa calculates the star-independent parameters 

# needed to correct RA and Dec 

# e.g the Earth barycentric and heliocentric position and velocity, 

# the precession-nutation matrix, etc. 

# 

# arguments of palpy.mappa are: 

# epoch of mean equinox to be used (Julian) 

# 

# date (MJD) 

prms = palpy.mappa(epoch, mjd.TDB) 

 

# palpy.mapqk does a quick mean to apparent place calculation using 

# the output of palpy.mappa 

# 

# Taken from the palpy source code (palMap.c which calls both palMappa and palMapqk): 

# The accuracy is sub-milliarcsecond, limited by the 

# precession-nutation model (see palPrenut for details). 

 

if include_px: 

# because PAL and ERFA expect proper motion in terms of "coordinate 

# angle; not true angle" (as stated in erfa/starpm.c documentation) 

pm_ra_corrected = pm_ra / np.cos(dec) 

 

if are_arrays: 

if include_px: 

raOut, decOut = palpy.mapqkVector(ra, dec, pm_ra_corrected, pm_dec, 

arcsecFromRadians(parallax), v_rad, prms) 

else: 

raOut, decOut = palpy.mapqkzVector(ra, dec, prms) 

else: 

if include_px: 

raOut, decOut = palpy.mapqk(ra, dec, pm_ra_corrected, pm_dec, 

arcsecFromRadians(parallax), v_rad, prms) 

else: 

raOut, decOut = palpy.mapqkz(ra, dec, prms) 

 

return np.array([raOut, decOut]) 

 

 

def _icrsFromAppGeo(ra, dec, epoch=2000.0, mjd=None): 

""" 

Convert the apparent geocentric position in (RA, Dec) to 

the mean position in the International Celestial Reference 

System (ICRS) 

 

This method undoes the effects of precession, annual aberration, 

and nutation. It is meant for mapping pointing RA and Dec (which 

presumably include the above effects) back to mean ICRS RA and Dec 

so that the user knows how to query a database of mean RA and Decs 

for objects observed at a given telescope pointing. 

 

WARNING: This method does not account for apparent motion due to parallax. 

This means it should not be used to invert the ICRS-to-apparent geocentric 

transformation for actual celestial objects. This method is only useful 

for mapping positions on a theoretical celestial sphere. 

 

This method works in radians. 

 

@param [in] ra in radians (apparent geocentric). Can be a numpy array or a number. 

 

@param [in] dec in radians (apparent geocentric). Can be a numpy array or a number. 

 

@param [in] epoch is the julian epoch (in years) of the equinox against which to 

measure RA (default: 2000.0) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the mean ICRS RA and 

the second row is the mean ICRS Dec (both in radians) 

""" 

 

are_arrays = _validate_inputs([ra, dec], ['ra', 'dec'], "icrsFromAppGeo") 

 

# Define star independent mean to apparent place parameters 

# palpy.mappa calculates the star-independent parameters 

# needed to correct RA and Dec 

# e.g the Earth barycentric and heliocentric position and velocity, 

# the precession-nutation matrix, etc. 

# 

# arguments of palpy.mappa are: 

# epoch of mean equinox to be used (Julian) 

# 

# date (MJD) 

params = palpy.mappa(epoch, mjd.TDB) 

 

if are_arrays: 

raOut, decOut = palpy.ampqkVector(ra, dec, params) 

else: 

raOut, decOut = palpy.ampqk(ra, dec, params) 

 

return np.array([raOut, decOut]) 

 

 

def icrsFromAppGeo(ra, dec, epoch=2000.0, mjd=None): 

""" 

Convert the apparent geocentric position in (RA, Dec) to 

the mean position in the International Celestial Reference 

System (ICRS) 

 

This method undoes the effects of precession, annual aberration, 

and nutation. It is meant for mapping pointing RA and Dec (which 

presumably include the above effects) back to mean ICRS RA and Dec 

so that the user knows how to query a database of mean RA and Decs 

for objects observed at a given telescope pointing. 

 

WARNING: This method does not account for apparent motion due to parallax. 

This means it should not be used to invert the ICRS-to-apparent geocentric 

transformation for actual celestial objects. This method is only useful 

for mapping positions on a theoretical celestial sphere. 

 

This method works in degrees. 

 

@param [in] ra in degrees (apparent geocentric). Can be a numpy array or a number. 

 

@param [in] dec in degrees (apparent geocentric). Can be a numpy array or a number. 

 

@param [in] epoch is the julian epoch (in years) of the equinox against which to 

measure RA (default: 2000.0) 

 

@param [in] mjd is an instantiation of the ModifiedJulianDate class 

representing the date of the observation 

 

@param [out] a 2-D numpy array in which the first row is the mean ICRS RA and 

the second row is the mean ICRS Dec (both in degrees) 

""" 

 

raOut, decOut = _icrsFromAppGeo(np.radians(ra), np.radians(dec), 

epoch=epoch, mjd=mjd) 

 

return np.array([np.degrees(raOut), np.degrees(decOut)]) 

 

 

def observedFromAppGeo(ra, dec, includeRefraction=True, 

altAzHr=False, wavelength=0.5, obs_metadata=None): 

""" 

Convert apparent geocentric (RA, Dec) to observed (RA, Dec). More 

specifically: apply refraction and diurnal aberration. 

 

This method works in degrees. 

 

@param [in] ra is geocentric apparent RA (degrees). Can be a numpy array or a number. 

 

@param [in] dec is geocentric apparent Dec (degrees). Can be a numpy array or a number. 

 

@param [in] includeRefraction is a boolean to turn refraction on and off 

 

@param [in] altAzHr is a boolean indicating whether or not to return altitude 

and azimuth 

 

@param [in] wavelength is effective wavelength in microns (default: 0.5) 

 

@param [in] obs_metadata is an ObservationMetaData characterizing the 

observation. 

 

@param [out] a 2-D numpy array in which the first row is the observed RA 

and the second row is the observed Dec (both in degrees) 

 

@param [out] a 2-D numpy array in which the first row is the altitude 

and the second row is the azimuth (both in degrees). Only returned 

if altAzHr == True. 

""" 

 

if altAzHr: 

raDec, altAz = _observedFromAppGeo(np.radians(ra), np.radians(dec), 

includeRefraction=includeRefraction, 

altAzHr=altAzHr, wavelength=wavelength, 

obs_metadata=obs_metadata) 

 

return np.degrees(raDec), np.degrees(altAz) 

 

else: 

output = _observedFromAppGeo(np.radians(ra), np.radians(dec), 

includeRefraction=includeRefraction, 

altAzHr=altAzHr, wavelength=wavelength, 

obs_metadata=obs_metadata) 

 

return np.degrees(output) 

 

 

def _calculateObservatoryParameters(obs_metadata, wavelength, includeRefraction): 

""" 

Computer observatory-based parameters using palpy.aoppa 

 

@param [in] obs_metadata is an ObservationMetaData characterizing 

the specific telescope site and pointing 

 

@param [in] wavelength is the effective wavelength in microns 

 

@param [in] includeRefraction is a boolean indicating whether or not 

to include the effects of refraction 

 

@param [out] the numpy array of observatory paramters calculated by 

palpy.aoppa 

""" 

 

# Correct site longitude for polar motion slaPolmo 

# 

# 5 January 2016 

# palAop.c (which calls Aoppa and Aopqk, as we do here) says 

# * - The azimuths etc produced by the present routine are with 

# * respect to the celestial pole. Corrections to the terrestrial 

# * pole can be computed using palPolmo. 

# 

# As a future issue, we should figure out how to incorporate polar motion 

# into these calculations. For now, we will set polar motion to zero. 

xPolar = 0.0 

yPolar = 0.0 

 

# 

# palpy.aoppa computes star-independent parameters necessary for 

# converting apparent place into observed place 

# i.e. it calculates geodetic latitude, magnitude of diurnal aberration, 

# refraction coefficients and the like based on data about the observation site 

if includeRefraction: 

obsPrms = palpy.aoppa(obs_metadata.mjd.UTC, obs_metadata.mjd.dut1, 

obs_metadata.site.longitude_rad, 

obs_metadata.site.latitude_rad, 

obs_metadata.site.height, 

xPolar, 

yPolar, 

obs_metadata.site.temperature_kelvin, 

obs_metadata.site.pressure, 

obs_metadata.site.humidity, 

wavelength, 

obs_metadata.site.lapseRate) 

else: 

# we can discard refraction by setting pressure and humidity to zero 

obsPrms = palpy.aoppa(obs_metadata.mjd.UTC, obs_metadata.mjd.dut1, 

obs_metadata.site.longitude_rad, 

obs_metadata.site.latitude_rad, 

obs_metadata.site.height, 

xPolar, 

yPolar, 

obs_metadata.site.temperature, 

0.0, 

0.0, 

wavelength, 

obs_metadata.site.lapseRate) 

 

return obsPrms 

 

 

def _observedFromAppGeo(ra, dec, includeRefraction=True, 

altAzHr=False, wavelength=0.5, obs_metadata=None): 

""" 

Convert apparent geocentric (RA, Dec) to observed (RA, Dec). More specifically: 

apply refraction and diurnal aberration. 

 

This method works in radians. 

 

@param [in] ra is geocentric apparent RA (radians). Can be a numpy array or a number. 

 

@param [in] dec is geocentric apparent Dec (radians). Can be a numpy array or a number. 

 

@param [in] includeRefraction is a boolean to turn refraction on and off 

 

@param [in] altAzHr is a boolean indicating whether or not to return altitude 

and azimuth 

 

@param [in] wavelength is effective wavelength in microns (default: 0.5) 

 

@param [in] obs_metadata is an ObservationMetaData characterizing the 

observation. 

 

@param [out] a 2-D numpy array in which the first row is the observed RA 

and the second row is the observed Dec (both in radians) 

 

@param [out] a 2-D numpy array in which the first row is the altitude 

and the second row is the azimuth (both in radians). Only returned 

if altAzHr == True. 

 

""" 

 

are_arrays = _validate_inputs( 

[ra, dec], ['ra', 'dec'], "observedFromAppGeo") 

 

if obs_metadata is None: 

raise RuntimeError( 

"Cannot call observedFromAppGeo without an obs_metadata") 

 

if obs_metadata.site is None: 

raise RuntimeError( 

"Cannot call observedFromAppGeo: obs_metadata has no site info") 

 

if obs_metadata.mjd is None: 

raise RuntimeError( 

"Cannot call observedFromAppGeo: obs_metadata has no mjd") 

 

obsPrms = _calculateObservatoryParameters( 

obs_metadata, wavelength, includeRefraction) 

 

# palpy.aopqk does an apparent to observed place 

# correction 

# 

# it corrects for diurnal aberration and refraction 

# (using a fast algorithm for refraction in the case of 

# a small zenith distance and a more rigorous algorithm 

# for a large zenith distance) 

# 

 

if are_arrays: 

azimuth, zenith, hourAngle, decOut, raOut = palpy.aopqkVector( 

ra, dec, obsPrms) 

else: 

azimuth, zenith, hourAngle, decOut, raOut = palpy.aopqk( 

ra, dec, obsPrms) 

 

# 

# Note: this is a choke point. Even the vectorized version of aopqk 

# is expensive (it takes about 0.006 seconds per call) 

# 

# Actually, this is only a choke point if you are dealing with zenith 

# distances of greater than about 70 degrees 

 

if altAzHr: 

# 

# palpy.de2h converts equatorial to horizon coordinates 

# 

if are_arrays: 

az, alt = palpy.de2hVector( 

hourAngle, decOut, obs_metadata.site.latitude_rad) 

else: 

az, alt = palpy.de2h( 

hourAngle, decOut, obs_metadata.site.latitude_rad) 

 

return np.array([raOut, decOut]), np.array([alt, az]) 

return np.array([raOut, decOut]) 

 

 

def appGeoFromObserved(ra, dec, includeRefraction=True, 

wavelength=0.5, obs_metadata=None): 

""" 

Convert observed (RA, Dec) to apparent geocentric (RA, Dec). More 

specifically: undo the effects of refraction and diurnal aberration. 

 

Note: This method is only accurate at zenith distances less than ~75 degrees. 

 

This method works in degrees. 

 

@param [in] ra is observed RA (degrees). Can be a numpy array or a number. 

 

@param [in] dec is observed Dec (degrees). Can be a numpy array or a number. 

 

@param [in] includeRefraction is a boolean to turn refraction on and off 

 

@param [in] wavelength is effective wavelength in microns (default: 0.5) 

 

@param [in] obs_metadata is an ObservationMetaData characterizing the 

observation. 

 

@param [out] a 2-D numpy array in which the first row is the apparent 

geocentric RA and the second row is the apparentGeocentric Dec (both 

in degrees) 

""" 

 

raOut, decOut = _appGeoFromObserved(np.radians(ra), np.radians(dec), 

includeRefraction=includeRefraction, 

wavelength=wavelength, 

obs_metadata=obs_metadata) 

 

return np.array([np.degrees(raOut), np.degrees(decOut)]) 

 

 

def _appGeoFromObserved(ra, dec, includeRefraction=True, 

wavelength=0.5, obs_metadata=None): 

""" 

Convert observed (RA, Dec) to apparent geocentric (RA, Dec). 

More specifically: undo the effects of refraction and diurnal aberration. 

 

Note: This method is only accurate at zenith distances less than ~ 75 degrees. 

 

This method works in radians. 

 

@param [in] ra is observed RA (radians). Can be a numpy array or a number. 

 

@param [in] dec is observed Dec (radians). Can be a numpy array or a number. 

 

@param [in] includeRefraction is a boolean to turn refraction on and off 

 

@param [in] wavelength is effective wavelength in microns (default: 0.5) 

 

@param [in] obs_metadata is an ObservationMetaData characterizing the 

observation. 

 

@param [out] a 2-D numpy array in which the first row is the apparent 

geocentric RA and the second row is the apparentGeocentric Dec (both 

in radians) 

""" 

 

are_arrays = _validate_inputs([ra, dec], ['ra', 'dec'], "appGeoFromObserved") 

 

if obs_metadata is None: 

raise RuntimeError("Cannot call appGeoFromObserved without an obs_metadata") 

 

if obs_metadata.site is None: 

raise RuntimeError("Cannot call appGeoFromObserved: obs_metadata has no site info") 

 

if obs_metadata.mjd is None: 

raise RuntimeError("Cannot call appGeoFromObserved: obs_metadata has no mjd") 

 

obsPrms = _calculateObservatoryParameters(obs_metadata, wavelength, includeRefraction) 

 

if are_arrays: 

raOut, decOut = palpy.oapqkVector('r', ra, dec, obsPrms) 

else: 

raOut, decOut = palpy.oapqk('r', ra, dec, obsPrms) 

 

return np.array([raOut, decOut]) 

 

 

def observedFromICRS(ra, dec, pm_ra=None, pm_dec=None, parallax=None, v_rad=None, 

obs_metadata=None, epoch=None, includeRefraction=True): 

""" 

Convert mean position (RA, Dec) in the International Celestial Reference Frame 

to observed (RA, Dec). 

 

included are precession-nutation, aberration, proper motion, parallax, refraction, 

radial velocity, diurnal aberration. 

 

This method works in degrees. 

 

@param [in] ra is the unrefracted RA in degrees (ICRS). Can be a numpy array or a number. 

 

@param [in] dec is the unrefracted Dec in degrees (ICRS). Can be a numpy array or a number. 

 

@param [in] pm_ra is proper motion in RA multiplied by cos(Dec) (arcsec/yr) 

Can be a numpy array or a number or None (default=None). 

 

@param [in] pm_dec is proper motion in dec (arcsec/yr) 

Can be a numpy array or a number or None (default=None). 

 

@param [in] parallax is parallax in arcsec 

Can be a numpy array or a number or None (default=None). 

 

@param [in] v_rad is radial velocity (km/s) 

Can be a numpy array or a number or None (default=None). 

 

@param [in] obs_metadata is an ObservationMetaData object describing the 

telescope pointing. 

 

@param [in] epoch is the julian epoch (in years) against which the mean 

equinoxes are measured. 

 

@param [in] includeRefraction toggles whether or not to correct for refraction 

 

@param [out] a 2-D numpy array in which the first row is the observed 

RA and the second row is the observed Dec (both in degrees) 

""" 

 

if pm_ra is not None: 

pm_ra_in = radiansFromArcsec(pm_ra) 

else: 

pm_ra_in = None 

 

if pm_dec is not None: 

pm_dec_in = radiansFromArcsec(pm_dec) 

else: 

pm_dec_in = None 

 

if parallax is not None: 

parallax_in = radiansFromArcsec(parallax) 

else: 

parallax_in = None 

 

output = _observedFromICRS(np.radians(ra), np.radians(dec), 

pm_ra=pm_ra_in, pm_dec=pm_dec_in, parallax=parallax_in, 

v_rad=v_rad, obs_metadata=obs_metadata, epoch=epoch, 

includeRefraction=includeRefraction) 

 

return np.degrees(output) 

 

 

def _observedFromICRS(ra, dec, pm_ra=None, pm_dec=None, parallax=None, v_rad=None, 

obs_metadata=None, epoch=None, includeRefraction=True): 

""" 

Convert mean position (RA, Dec) in the International Celestial Reference Frame 

to observed (RA, Dec)-like coordinates. 

 

included are precession-nutation, aberration, proper motion, parallax, refraction, 

radial velocity, diurnal aberration. 

 

This method works in radians. 

 

@param [in] ra is the unrefracted RA in radians (ICRS). Can be a numpy array or a number. 

 

@param [in] dec is the unrefracted Dec in radians (ICRS). Can be a numpy array or a number. 

 

@param [in] pm_ra is proper motion in RA multiplied by cos(Dec) (radians/yr) 

Can be a numpy array or a number or None (default=None). 

 

@param [in] pm_dec is proper motion in dec (radians/yr) 

Can be a numpy array or a number or None (default=None). 

 

@param [in] parallax is parallax in radians 

Can be a numpy array or a number or None (default=None). 

 

@param [in] v_rad is radial velocity (km/s) 

Can be a numpy array or a number or None (default=None). 

 

@param [in] obs_metadata is an ObservationMetaData object describing the 

telescope pointing. 

 

@param [in] epoch is the julian epoch (in years) against which the mean 

equinoxes are measured. 

 

@param [in] includeRefraction toggles whether or not to correct for refraction 

 

@param [out] a 2-D numpy array in which the first row is the observed 

RA and the second row is the observed Dec (both in radians) 

 

""" 

 

if obs_metadata is None: 

raise RuntimeError("Cannot call observedFromICRS; obs_metadata is none") 

 

if obs_metadata.mjd is None: 

raise RuntimeError("Cannot call observedFromICRS; obs_metadata.mjd is none") 

 

if epoch is None: 

raise RuntimeError("Cannot call observedFromICRS; you have not specified an epoch") 

 

ra_apparent, dec_apparent = _appGeoFromICRS(ra, dec, pm_ra=pm_ra, 

pm_dec=pm_dec, parallax=parallax, 

v_rad=v_rad, epoch=epoch, mjd=obs_metadata.mjd) 

 

ra_out, dec_out = _observedFromAppGeo(ra_apparent, dec_apparent, obs_metadata=obs_metadata, 

includeRefraction=includeRefraction) 

 

return np.array([ra_out, dec_out]) 

 

 

def icrsFromObserved(ra, dec, obs_metadata=None, epoch=None, includeRefraction=True): 

""" 

Convert observed RA, Dec into mean International Celestial Reference Frame (ICRS) 

RA, Dec. This method undoes the effects of precession, nutation, aberration (annual 

and diurnal), and refraction. It is meant so that users can take pointing RA and Decs, 

which will be in the observed reference system, and transform them into ICRS for 

purposes of querying database tables (likely to contain mean ICRS RA, Dec) for objects 

visible from a given pointing. 

 

Note: This method is only accurate at angular distances from the sun of greater 

than 45 degrees and zenith distances of less than 75 degrees. 

 

WARNING: This method does not account for apparent motion due to parallax. 

This means it should not be used to invert the ICRS-to-observed coordinates 

transformation for actual celestial objects. This method is only useful 

for mapping positions on a theoretical celestial sphere. 

 

This method works in degrees. 

 

@param [in] ra is the observed RA in degrees. Can be a numpy array or a number. 

 

@param [in] dec is the observed Dec in degrees. Can be a numpy array or a number. 

 

@param [in] obs_metadata is an ObservationMetaData object describing the 

telescope pointing. 

 

@param [in] epoch is the julian epoch (in years) against which the mean 

equinoxes are measured. 

 

@param [in] includeRefraction toggles whether or not to correct for refraction 

 

@param [out] a 2-D numpy array in which the first row is the mean ICRS 

RA and the second row is the mean ICRS Dec (both in degrees) 

""" 

 

ra_out, dec_out = _icrsFromObserved(np.radians(ra), np.radians(dec), 

obs_metadata=obs_metadata, 

epoch=epoch, includeRefraction=includeRefraction) 

 

return np.array([np.degrees(ra_out), np.degrees(dec_out)]) 

 

 

def _icrsFromObserved(ra, dec, obs_metadata=None, epoch=None, includeRefraction=True): 

""" 

Convert observed RA, Dec into mean International Celestial Reference Frame (ICRS) 

RA, Dec. This method undoes the effects of precession, nutation, aberration (annual 

and diurnal), and refraction. It is meant so that users can take pointing RA and Decs, 

which will be in the observed reference system, and transform them into ICRS for 

purposes of querying database tables (likely to contain mean ICRS RA, Dec) for objects 

visible from a given pointing. 

 

Note: This method is only accurate at angular distances from the sun of greater 

than 45 degrees and zenith distances of less than 75 degrees. 

 

WARNING: This method does not account for apparent motion due to parallax. 

This means it should not be used to invert the ICRS-to-observed coordinates 

transformation for actual celestial objects. This method is only useful 

for mapping positions on a theoretical celestial sphere. 

 

This method works in radians. 

 

@param [in] ra is the observed RA in radians. Can be a numpy array or a number. 

 

@param [in] dec is the observed Dec in radians. Can be a numpy array or a number. 

 

@param [in] obs_metadata is an ObservationMetaData object describing the 

telescope pointing. 

 

@param [in] epoch is the julian epoch (in years) against which the mean 

equinoxes are measured. 

 

@param [in] includeRefraction toggles whether or not to correct for refraction 

 

@param [out] a 2-D numpy array in which the first row is the mean ICRS 

RA and the second row is the mean ICRS Dec (both in radians) 

""" 

 

_validate_inputs([ra, dec], ['ra', 'dec'], "icrsFromObserved") 

 

if obs_metadata is None: 

raise RuntimeError("Cannot call icrsFromObserved; obs_metadata is None") 

 

if obs_metadata.mjd is None: 

raise RuntimeError("Cannot call icrsFromObserved; obs_metadata.mjd is None") 

 

if epoch is None: 

raise RuntimeError("Cannot call icrsFromObserved; you have not specified an epoch") 

 

ra_app, dec_app = _appGeoFromObserved(ra, dec, obs_metadata=obs_metadata, 

includeRefraction=includeRefraction) 

 

ra_icrs, dec_icrs = _icrsFromAppGeo(ra_app, dec_app, epoch=epoch, 

mjd=obs_metadata.mjd) 

 

return np.array([ra_icrs, dec_icrs])