• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

dsavransky / EXOSIMS / 14536111075

18 Apr 2025 01:44PM UTC coverage: 65.708% (+0.2%) from 65.466%
14536111075

push

github

web-flow
Merge pull request #402 from dsavransky/exozodi_fix

Exozodi fix. Extensive re-plumbing of all things related to exozodi.

210 of 269 new or added lines in 26 files covered. (78.07%)

14 existing lines in 7 files now uncovered.

9665 of 14709 relevant lines covered (65.71%)

0.66 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

85.42
/EXOSIMS/Prototypes/SimulatedUniverse.py
1
import astropy.constants as const
1✔
2
import astropy.units as u
1✔
3
import numpy as np
1✔
4

5
from EXOSIMS.util.deltaMag import deltaMag
1✔
6
from EXOSIMS.util.eccanom import eccanom
1✔
7
from EXOSIMS.util.get_dirs import get_cache_dir
1✔
8
from EXOSIMS.util.get_module import get_module
1✔
9
from EXOSIMS.util.keplerSTM import planSys
1✔
10
from EXOSIMS.util.vprint import vprint
1✔
11

12

13
class SimulatedUniverse(object):
1✔
14
    r""":ref:`SimulatedUniverse` Prototype
15

16
    Args:
17
        fixedPlanPerStar (int, optional):
18
            If set, every system will have the same number of planets.
19
            Defaults to None
20
        Min (float, optional):
21
            Initial mean anomaly for all planets.  If set, every planet
22
            has the same mean anomaly at mission start. Defaults to None
23
        cachedir (str, optional):
24
            Full path to cachedir.
25
            If None (default) use default (see :ref:`EXOSIMSCACHE`)
26
        lucky_planets (bool):
27
            Used downstream in survey simulation. If True, planets are
28
            observed at optimal times. Defaults to False
29
        commonSystemPlane (bool):
30
            Planet inclinations are sampled as normally distributed about a
31
            common system plane. Defaults to False
32
        commonSystemPlaneParams (list(float)):
33
            [inclination mean, inclination standard deviation, Omega mean, Omega
34
            standard deviation] defining the normal distribution of
35
            inclinations and longitudes of the ascending node about a common
36
            system plane in units of degrees.  Ignored if commonSystemPlane is
37
            False. Defaults to [0 2.25, 0, 2.25], where the standard deviation
38
            is approximately the standard deviation of solar system planet
39
            inclinations.
40
        commonSystemnEZ (bool):
41
            Assume same zodi for planets in the same system. Defaults to False.
42
        **specs:
43
            :ref:`sec:inputspec`
44

45
    Attributes:
46
        _outspec (dict):
47
            :ref:`sec:outspec`
48
        a (astropy.units.quantity.Quantity):
49
             Planet semi-major axis (length units)
50
        BackgroundSources (:ref:`BackgroundSources`):
51
            BackgroundSources object
52
        cachedir (str):
53
            Path to the EXOSIMS cache directory (see :ref:`EXOSIMSCACHE`)
54
        commonSystemPlaneParams (list):
55
            2 element list of [mean, standard deviation] in units of degrees,
56
            describing the distribution of inclinations relative to a common orbital
57
            plane.  Ignored if commonSystemPlane is False.
58
        commonSystemPlane (bool):
59
            If False, planet inclinations are independently drawn for all planets,
60
            including those in the same target system.  If True, inclinations will be
61
            drawn from a normal distribution defined by
62
            commonSystemPlaneParams and added to a single inclination value drawn
63
            for each system.
64
        Completeness (:ref:`Completeness`):
65
            Completeness object
66
        d (astropy.units.quantity.Quantity):
67
            Current orbital radius magnitude (length units)
68
        dMag (numpy.ndarray):
69
            Current planet :math:`\Delta\mathrm{mag}`
70
        e (numpy.ndarray):
71
            Planet eccentricity
72
        fixedPlanPerStar (int or None):
73
            If set, every system has the same number of planets, given by
74
            this attribute
75
        I (astropy.units.quantity.Quantity):
76
            Planet inclinations (angle units)
77
        lucky_planets (bool):
78
            If True, planets are observed at optimal times.
79
        M0 (astropy.units.quantity.Quantity):
80
            Initial planet mean anomaly (at mission start time).
81
        Min (float or None):
82
            Input constant initial mean anomaly.  If none, initial
83
            mean anomaly is randomly distributed from a uniform distribution in
84
            [0, 360] degrees.
85
        Mp (astropy.units.quantity.Quantity):
86
            Planet mass.
87
        nPlans (int):
88
            Number of planets in all target systems.
89
        O (astropy.units.quantity.Quantity):
90
            Planet longitude of the ascending node (angle units)
91
        OpticalSystem (:ref:`OpticalSystem`):
92
            Optical System object
93
        p (numpy.ndarray):
94
            Planet geometric albedo
95
        phi (numpy.ndarray):
96
            Current value of planet phase function.
97
        phiIndex (numpy.ndarray):
98
            Intended for use with input
99
            'whichPlanetPhaseFunction'='realSolarSystemPhaseFunc'
100
            When None, the default is the phi_lambert function, otherwise it is Solar
101
            System Phase Functions
102
        plan2star (numpy.ndarray):
103
            Index of host star or each planet.  Indexes attributes of TargetsList.
104
        planet_atts (list):
105
            List of planet attributes
106
        PlanetPhysicalModel (:ref:`PlanetPhysicalModel`):
107
            Planet physical model object.
108
        PlanetPopulation (:ref:`PlanetPopulation`):
109
            Planet population object.
110
        PostProcessing (:ref:`PostProcessing`):
111
            Postprocessing object.
112
        r (astropy.units.quantity.Quantity):
113
            Current planet orbital radius (3xnPlans). Length units.
114
        Rp (astropy.units.quantity.Quantity):
115
            Planet radius (length units).
116
        s (astropy.units.quantity.Quantity):
117
            Current planet projected separation. Length units.
118
        sInds (numpy.ndarray):
119
            Indices of stars with planets.  Equivalent to unique entries of
120
            ``plan2star``.
121
        nEZ (numpy.ndarray):
122
            Number of exozodi for each planet.
123
        TargetList (:ref:`TargetList`):
124
            Target list object.
125
        v (astropy.units.quantity.Quantity):
126
            Current orbital velocity vector (3xnPlans). Velocity units.
127
        w (astropy.units.quantity.Quantity):
128
            Planet argument of periapsis.
129
        WA (astropy.units.quantity.Quantity):
130
            Current planet angular separation (angle units)
131
        ZodiacalLight (:ref:`ZodiacalLight`):
132
            Zodiacal light object.
133

134

135
    .. note::
136

137
        When generating planets, :ref:`PlanetPopulation` attribute ``eta`` is
138
        treated as the rate parameter of a Poisson distribution.
139
        Each target's number of planets is a Poisson random variable
140
        sampled with :math:`\lambda\equiv\eta`.
141

142
    .. warning::
143

144
        All attributes described as 'current' are updated only when planets are
145
        observed.  As such, during mission simulations, these values for different
146
        planets correspond to different times (bookkept in the survey simulation
147
        object).
148

149
    """
150

151
    _modtype = "SimulatedUniverse"
1✔
152

153
    def __init__(
1✔
154
        self,
155
        fixedPlanPerStar=None,
156
        Min=None,
157
        cachedir=None,
158
        lucky_planets=False,
159
        commonSystemPlane=False,
160
        commonSystemPlaneParams=[0, 2.25, 0, 2.25],
161
        commonSystemnEZ=False,
162
        **specs,
163
    ):
164
        # start the outspec
165
        self._outspec = {}
1✔
166

167
        # load the vprint function (same line in all prototype module constructors)
168
        self.vprint = vprint(specs.get("verbose", True))
1✔
169
        self.lucky_planets = lucky_planets
1✔
170
        self._outspec["lucky_planets"] = lucky_planets
1✔
171
        self.commonSystemPlane = bool(commonSystemPlane)
1✔
172
        self._outspec["commonSystemPlane"] = commonSystemPlane
1✔
173
        assert (
1✔
174
            len(commonSystemPlaneParams) == 4
175
        ), "commonSystemPlaneParams must be a four-element list"
176
        self.commonSystemPlaneParams = commonSystemPlaneParams
1✔
177
        self._outspec["commonSystemPlaneParams"] = commonSystemPlaneParams
1✔
178

179
        # Set the number of exozodi
180
        self.commonSystemnEZ = commonSystemnEZ
1✔
181

182
        # save fixed number of planets to generate
183
        self.fixedPlanPerStar = fixedPlanPerStar
1✔
184
        self._outspec["fixedPlanPerStar"] = fixedPlanPerStar
1✔
185

186
        # get cache directory
187
        self.cachedir = get_cache_dir(cachedir)
1✔
188
        self._outspec["cachedir"] = self.cachedir
1✔
189
        specs["cachedir"] = self.cachedir
1✔
190

191
        # check if KnownRVPlanetsUniverse has correct input modules
192
        if specs["modules"]["SimulatedUniverse"] == "KnownRVPlanetsUniverse":
1✔
193
            val = (
1✔
194
                specs["modules"]["TargetList"] == "KnownRVPlanetsTargetList"
195
                and specs["modules"]["PlanetPopulation"] == "KnownRVPlanets"
196
            )
197
            assert val, (
1✔
198
                "KnownRVPlanetsUniverse must use KnownRVPlanetsTargetList "
199
                "and KnownRVPlanets"
200
            )
201
        else:
202
            val = (
1✔
203
                specs["modules"]["TargetList"] == "KnownRVPlanetsTargetList"
204
                or specs["modules"]["PlanetPopulation"] == "KnownRVPlanets"
205
            )
206
            assert not (val), (
1✔
207
                "KnownRVPlanetsTargetList or KnownRVPlanets should not be used "
208
                "with this SimulatedUniverse"
209
            )
210

211
        # import TargetList class
212
        self.TargetList = get_module(specs["modules"]["TargetList"], "TargetList")(
1✔
213
            **specs
214
        )
215

216
        # bring inherited class objects to top level of Simulated Universe
217
        TL = self.TargetList
1✔
218
        self.StarCatalog = TL.StarCatalog
1✔
219
        self.PlanetPopulation = TL.PlanetPopulation
1✔
220
        self.PlanetPhysicalModel = TL.PlanetPhysicalModel
1✔
221
        self.OpticalSystem = TL.OpticalSystem
1✔
222
        self.ZodiacalLight = TL.ZodiacalLight
1✔
223
        self.BackgroundSources = TL.BackgroundSources
1✔
224
        self.PostProcessing = TL.PostProcessing
1✔
225
        self.Completeness = TL.Completeness
1✔
226

227
        # initial constant mean anomaly
228
        assert isinstance(Min, (int, float)) or (
1✔
229
            Min is None
230
        ), "Min may be int, float, or None"
231
        if Min is not None:
1✔
232
            self.Min = float(Min) * u.deg
1✔
233
        else:
234
            self.Min = Min
1✔
235
        self._outspec["Min"] = Min
1✔
236

237
        # list of possible planet attributes
238
        self.planet_atts = [
1✔
239
            "plan2star",
240
            "a",
241
            "e",
242
            "I",
243
            "O",
244
            "w",
245
            "M0",
246
            "Min",
247
            "Rp",
248
            "Mp",
249
            "p",
250
            "r",
251
            "v",
252
            "d",
253
            "s",
254
            "phi",
255
            "nEZ",
256
            "dMag",
257
            "WA",
258
        ]
259

260
        self.phiIndex = (
1✔
261
            None  # Used to switch select specific phase function for each planet
262
        )
263

264
        # generate orbital elements, albedos, radii, and masses
265
        self.gen_physical_properties(**specs)
1✔
266

267
        # find initial position-related parameters: position, velocity, planet-star
268
        # distance, apparent separation, number of zodis
269
        self.init_systems()
1✔
270

271
    def __str__(self):
1✔
272
        """String representation of Simulated Universe object
273

274
        When the command 'print' is used on the Simulated Universe object,
275
        this method will return the values contained in the object
276

277
        """
278

279
        for att in self.__dict__:
1✔
280
            print("%s: %r" % (att, getattr(self, att)))
1✔
281

282
        return "Simulated Universe class object attributes"
1✔
283

284
    def gen_physical_properties(self, **specs):
1✔
285
        """Generates the planetary systems' physical properties.
286

287
        Populates arrays of the orbital elements, albedos, masses and radii
288
        of all planets, and generates indices that map from planet to parent star.
289

290
        Args:
291
            **specs:
292
                :ref:`sec:inputspec`
293

294
        Returns:
295
            None
296

297
        """
298

299
        PPop = self.PlanetPopulation
1✔
300
        ZL = self.ZodiacalLight
1✔
301
        TL = self.TargetList
1✔
302

303
        if self.fixedPlanPerStar is not None:  # Must be an int for fixedPlanPerStar
1✔
304
            # Create array of length TL.nStars each w/ value ppStar
305
            targetSystems = np.ones(TL.nStars).astype(int) * self.fixedPlanPerStar
1✔
306
        else:
307
            # treat eta as the rate parameter of a Poisson distribution
308
            targetSystems = np.random.poisson(lam=PPop.eta, size=TL.nStars)
1✔
309

310
        plan2star = []
1✔
311
        for j, n in enumerate(targetSystems):
1✔
312
            plan2star = np.hstack((plan2star, [j] * n))
1✔
313
        self.plan2star = plan2star.astype(int)
1✔
314
        self.sInds = np.unique(self.plan2star)
1✔
315
        self.nPlans = len(self.plan2star)
1✔
316

317
        # The prototype StarCatalog module is made of one single G star at 1pc.
318
        # In that case, the SimulatedUniverse prototype generates one Jupiter
319
        # at 5 AU to allow for characterization testing.
320
        # Also generates at least one Jupiter if no planet was generated.
321
        if (TL.Name[0].startswith("Prototype") and (TL.nStars == 1)) or (
1✔
322
            self.nPlans == 0
323
        ):
324
            if self.nPlans == 0:
1✔
325
                self.vprint("No planets were generated. Creating single fake planet.")
1✔
326
            else:
327
                self.vprint(
1✔
328
                    (
329
                        "Prototype StarCatalog with 1 target. "
330
                        "Creating single fake planet."
331
                    )
332
                )
333
            self.plan2star = np.array([0], dtype=int)
1✔
334
            self.sInds = np.unique(self.plan2star)
1✔
335
            self.nPlans = len(self.plan2star)
1✔
336
            self.a = np.array([5.0]) * u.AU
1✔
337
            self.e = np.array([0.0])
1✔
338
            self.I = np.array([0.0]) * u.deg
1✔
339
            self.O = np.array([0.0]) * u.deg
1✔
340
            self.w = np.array([0.0]) * u.deg
1✔
341
            self.gen_M0()
1✔
342
            self.Rp = np.array([10.0]) * u.earthRad
1✔
343
            self.Mp = np.array([300.0]) * u.earthMass
1✔
344
            self.p = np.array([0.6])
1✔
345
        else:
346
            # sample all of the orbital and physical parameters
347
            self.I, self.O, self.w = PPop.gen_angles(
1✔
348
                self.nPlans,
349
                commonSystemPlane=self.commonSystemPlane,
350
                commonSystemPlaneParams=self.commonSystemPlaneParams,
351
            )
352
            self.setup_system_planes()
1✔
353

354
            self.a, self.e, self.p, self.Rp = PPop.gen_plan_params(self.nPlans)
1✔
355
            if PPop.scaleOrbits:
1✔
356
                self.a *= np.sqrt(TL.L[self.plan2star])
1✔
357
            self.gen_M0()  # initial mean anomaly
1✔
358
            self.Mp = PPop.gen_mass(self.nPlans)  # mass
1✔
359

360
        if self.commonSystemnEZ:
1✔
361
            # Assign the same nEZ to all planets in the system
NEW
362
            self.nEZ = ZL.gen_systemnEZ(TL.nStars)[self.plan2star]
×
363
        else:
364
            # Assign a unique nEZ to each planet
365
            self.nEZ = ZL.gen_systemnEZ(self.nPlans)
1✔
366

367
        self.phiIndex = np.asarray(
1✔
368
            []
369
        )  # Used to switch select specific phase function for each planet
370

371
    def gen_M0(self):
1✔
372
        """Set initial mean anomaly for each planet"""
373
        if self.Min is not None:
1✔
374
            self.M0 = np.ones((self.nPlans,)) * self.Min
1✔
375
        else:
376
            self.M0 = np.random.uniform(360, size=int(self.nPlans)) * u.deg
1✔
377

378
    def init_systems(self):
1✔
379
        """Finds initial time-dependant parameters. Assigns each planet an
380
        initial position, velocity, planet-star distance, apparent separation,
381
        phase function, surface brightness of exo-zodiacal light, delta
382
        magnitude, and working angle.
383

384
        This method makes use of the systems' physical properties (masses,
385
        distances) and their orbital elements (a, e, I, O, w, M0).
386
        """
387

388
        PPMod = self.PlanetPhysicalModel
1✔
389
        TL = self.TargetList
1✔
390

391
        a = self.a.to("AU").value  # semi-major axis
1✔
392
        e = self.e  # eccentricity
1✔
393
        I = self.I.to("rad").value  # inclinations #noqa: E741
1✔
394
        O = self.O.to("rad").value  # right ascension of the ascending node #noqa: E741
1✔
395
        w = self.w.to("rad").value  # argument of perigee
1✔
396
        M0 = self.M0.to("rad").value  # initial mean anomany
1✔
397
        E = eccanom(M0, e)  # eccentric anomaly
1✔
398
        Mp = self.Mp  # planet masses
1✔
399

400
        # This is the a1 a2 a3 and b1 b2 b3 are the euler angle transformation from
401
        # the I,J,K refernece frame to an x,y,z frame
402
        a1 = np.cos(O) * np.cos(w) - np.sin(O) * np.cos(I) * np.sin(w)
1✔
403
        a2 = np.sin(O) * np.cos(w) + np.cos(O) * np.cos(I) * np.sin(w)
1✔
404
        a3 = np.sin(I) * np.sin(w)
1✔
405
        A = a * np.vstack((a1, a2, a3)) * u.AU
1✔
406
        b1 = -np.sqrt(1 - e**2) * (
1✔
407
            np.cos(O) * np.sin(w) + np.sin(O) * np.cos(I) * np.cos(w)
408
        )
409
        b2 = np.sqrt(1 - e**2) * (
1✔
410
            -np.sin(O) * np.sin(w) + np.cos(O) * np.cos(I) * np.cos(w)
411
        )
412
        b3 = np.sqrt(1 - e**2) * np.sin(I) * np.cos(w)
1✔
413
        B = a * np.vstack((b1, b2, b3)) * u.AU
1✔
414
        r1 = np.cos(E) - e
1✔
415
        r2 = np.sin(E)
1✔
416

417
        mu = const.G * (Mp + TL.MsTrue[self.plan2star])
1✔
418
        v1 = np.sqrt(mu / self.a**3) / (1 - e * np.cos(E))
1✔
419
        v2 = np.cos(E)
1✔
420

421
        self.r = (A * r1 + B * r2).T.to("AU")  # position
1✔
422
        self.v = (v1 * (-A * r2 + B * v2)).T.to("AU/day")  # velocity
1✔
423
        self.s = np.linalg.norm(self.r[:, 0:2], axis=1)  # apparent separation
1✔
424
        self.d = np.linalg.norm(self.r, axis=1)  # planet-star distance
1✔
425
        try:
1✔
426
            self.phi = PPMod.calc_Phi(
1✔
427
                np.arccos(self.r[:, 2] / self.d), phiIndex=self.phiIndex
428
            )  # planet phase
429
        except u.UnitTypeError:
×
430
            self.d = self.d * self.r.unit  # planet-star distance
×
431
            self.phi = PPMod.calc_Phi(
×
432
                np.arccos(self.r[:, 2] / self.d), phiIndex=self.phiIndex
433
            )  # planet phase
434

435
        self.dMag = deltaMag(self.p, self.Rp, self.d, self.phi)  # delta magnitude
1✔
436
        try:
1✔
437
            self.WA = np.arctan(self.s / TL.dist[self.plan2star]).to(
1✔
438
                "arcsec"
439
            )  # working angle
440
        except u.UnitTypeError:
×
441
            self.s = self.s * self.r.unit
×
442
            self.WA = np.arctan(self.s / TL.dist[self.plan2star]).to(
×
443
                "arcsec"
444
            )  # working angle
445

446
    def propag_system(self, sInd, dt):
1✔
447
        """Propagates planet time-dependant parameters: position, velocity,
448
        planet-star distance, apparent separation, phase function, surface brightness
449
        of exo-zodiacal light, delta magnitude, and working angle.
450

451
        This method uses the Kepler state transition matrix to propagate a
452
        planet's state (position and velocity vectors) forward in time using
453
        the Kepler state transition matrix.
454

455
        Args:
456
            sInd (int):
457
                Index of the target system of interest
458
            dt (~astropy.units.Quantity(float)):
459
                Time increment in units of day, for planet position propagation
460

461
        Returns:
462
            None
463
        """
464

465
        PPMod = self.PlanetPhysicalModel
1✔
466
        TL = self.TargetList
1✔
467

468
        assert np.isscalar(
1✔
469
            sInd
470
        ), "Can only propagate one system at a time, sInd must be scalar."
471
        # check for planets around this target
472
        pInds = np.where(self.plan2star == sInd)[0]
1✔
473
        if len(pInds) == 0:
1✔
474
            return
×
475
        # check for positive time increment
476
        assert dt >= 0, "Time increment (dt) to propagate a planet must be positive."
1✔
477
        if dt == 0:
1✔
478
            return
×
479

480
        # Calculate initial positions in AU and velocities in AU/day
481
        r0 = self.r[pInds].to("AU").value
1✔
482
        v0 = self.v[pInds].to("AU/day").value
1✔
483
        # stack dimensionless positions and velocities
484
        nPlans = pInds.size
1✔
485
        x0 = np.reshape(np.concatenate((r0, v0), axis=1), nPlans * 6)
1✔
486

487
        # Calculate vector of gravitational parameter in AU3/day2
488
        Ms = TL.MsTrue[[sInd]]
1✔
489
        Mp = self.Mp[pInds]
1✔
490
        mu = (const.G * (Mp + Ms)).to("AU3/day2").value
1✔
491

492
        # use keplerSTM.py to propagate the system
493
        prop = planSys(x0, mu, epsmult=10.0)
1✔
494
        try:
1✔
495
            prop.takeStep(dt.to("day").value)
1✔
496
        except ValueError:
×
497
            # try again with larger epsmult and two steps to force convergence
498
            prop = planSys(x0, mu, epsmult=100.0)
×
499
            try:
×
500
                prop.takeStep(dt.to("day").value / 2.0)
×
501
                prop.takeStep(dt.to("day").value / 2.0)
×
502
            except ValueError:
×
503
                raise ValueError("planSys error")
×
504

505
        # split off position and velocity vectors
506
        x1 = np.array(np.hsplit(prop.x0, 2 * nPlans))
1✔
507
        rind = np.array(range(0, len(x1), 2))  # even indices
1✔
508
        vind = np.array(range(1, len(x1), 2))  # odd indices
1✔
509

510
        # update planets' position, velocity, planet-star distance, apparent,
511
        # separation, phase function, exozodi surface brightness, delta magnitude and
512
        # working angle
513
        self.r[pInds] = x1[rind] * u.AU
1✔
514
        self.v[pInds] = x1[vind] * u.AU / u.day
1✔
515

516
        try:
1✔
517
            self.d[pInds] = np.linalg.norm(self.r[pInds], axis=1)
1✔
518
            if len(self.phiIndex) == 0:
1✔
519
                self.phi[pInds] = PPMod.calc_Phi(
1✔
520
                    np.arccos(self.r[pInds, 2] / self.d[pInds]), phiIndex=self.phiIndex
521
                )
522
            else:
523
                self.phi[pInds] = PPMod.calc_Phi(
×
524
                    np.arccos(self.r[pInds, 2] / self.d[pInds]),
525
                    phiIndex=self.phiIndex[pInds],
526
                )
527
        except u.UnitTypeError:
×
528
            self.d[pInds] = np.linalg.norm(self.r[pInds], axis=1) * self.r.unit
×
529
            if len(self.phiIndex) == 0:
×
530
                self.phi[pInds] = PPMod.calc_Phi(
×
531
                    np.arccos(self.r[pInds, 2] / self.d[pInds]), phiIndex=self.phiIndex
532
                )
533
            else:
534
                self.phi[pInds] = PPMod.calc_Phi(
×
535
                    np.arccos(self.r[pInds, 2] / self.d[pInds]),
536
                    phiIndex=self.phiIndex[pInds],
537
                )
538

539
        self.dMag[pInds] = deltaMag(
1✔
540
            self.p[pInds], self.Rp[pInds], self.d[pInds], self.phi[pInds]
541
        )
542
        try:
1✔
543
            self.s[pInds] = np.linalg.norm(self.r[pInds, 0:2], axis=1)
1✔
544
            self.WA[pInds] = np.arctan(self.s[pInds] / TL.dist[sInd]).to("arcsec")
1✔
545
        except u.UnitTypeError:
×
546
            self.s[pInds] = np.linalg.norm(self.r[pInds, 0:2], axis=1) * self.r.unit
×
547
            self.WA[pInds] = np.arctan(self.s[pInds] / TL.dist[sInd]).to("arcsec")
×
548

549
    def scale_JEZ(self, sInd, mode, pInds=None):
1✔
550
        """Scales the exozodi intensity by the inverse square of the planet-star distance.
551

552
        Args:
553
            sInd (int):
554
                Index of the target system of interest
555
            mode (dict):
556
                Selected observing mode
557
            pInds (numpy.ndarray):
558
                Planet indices. Default value (None) will return all planets in the system
559

560
        Returns:
561
            numpy.ndarray:
562
                Scaled exozodi intensity for each planet in the system in units of
563
                photon/s/m2/arcsec2
564
        """
565
        # Get the 1 AU value of JEZ for the system
566
        JEZ0 = self.TargetList.JEZ0[mode["hex"]][sInd]
1✔
567

568
        # Scale JEZ by nEZ and the inverse square of the planet-star distance
569
        all_pinds = np.where(self.plan2star == sInd)[0]
1✔
570
        if pInds is None:
1✔
571
            pinds = all_pinds
1✔
572
        else:
573
            pinds = np.intersect1d(all_pinds, pInds)
1✔
574
        JEZ = JEZ0 * self.nEZ[pinds] * (1 / self.d[pinds].to("AU").value) ** 2
1✔
575
        return JEZ
1✔
576

577
    def set_planet_phase(self, beta=np.pi / 2):
1✔
578
        """Positions all planets at input star-planet-observer phase angle
579
        where possible. For systems where the input phase angle is not achieved,
580
        planets are positioned at quadrature (phase angle of 90 deg).
581

582
        The position found here is not unique. The desired phase angle will be
583
        achieved at two points on the planet's orbit (for non-face on orbits).
584

585
        Args:
586
            beta (float):
587
                star-planet-observer phase angle in radians.
588

589
        """
590

591
        PPMod = self.PlanetPhysicalModel
1✔
592
        TL = self.TargetList
1✔
593

594
        a = self.a.to("AU").value  # semi-major axis
1✔
595
        e = self.e  # eccentricity
1✔
596
        I = self.I.to("rad").value  # noqa: E741 # inclinations
1✔
597
        O = self.O.to("rad").value  # noqa: E741 # long. of the ascending node
1✔
598
        w = self.w.to("rad").value  # argument of perigee
1✔
599
        Mp = self.Mp  # planet masses
1✔
600

601
        # make list of betas
602
        betas = beta * np.ones(w.shape)
1✔
603
        mask = np.cos(betas) / np.sin(I) > 1.0
1✔
604
        num = len(np.where(mask)[0])
1✔
605
        betas[mask] = np.pi / 2.0
1✔
606
        mask = np.cos(betas) / np.sin(I) < -1.0
1✔
607
        num += len(np.where(mask)[0])
1✔
608
        betas[mask] = np.pi / 2.0
1✔
609
        if num > 0:
1✔
610
            self.vprint("***Warning***")
1✔
611
            self.vprint(
1✔
612
                (
613
                    "{} planets out of {} could not be set to phase angle {} radians."
614
                ).format(num, self.nPlans, beta)
615
            )
616
            self.vprint("These planets are set to quadrature (phase angle pi/2)")
1✔
617

618
        # solve for true anomaly
619
        nu = np.arcsin(np.cos(betas) / np.sin(I)) - w
1✔
620

621
        # setup for position and velocity
622
        a1 = np.cos(O) * np.cos(w) - np.sin(O) * np.cos(I) * np.sin(w)
1✔
623
        a2 = np.sin(O) * np.cos(w) + np.cos(O) * np.cos(I) * np.sin(w)
1✔
624
        a3 = np.sin(I) * np.sin(w)
1✔
625
        A = np.vstack((a1, a2, a3))
1✔
626

627
        b1 = -(np.cos(O) * np.sin(w) + np.sin(O) * np.cos(I) * np.cos(w))
1✔
628
        b2 = -np.sin(O) * np.sin(w) + np.cos(O) * np.cos(I) * np.cos(w)
1✔
629
        b3 = np.sin(I) * np.cos(w)
1✔
630
        B = np.vstack((b1, b2, b3))
1✔
631

632
        r = a * (1.0 - e**2) / (1.0 - e * np.cos(nu))
1✔
633
        mu = const.G * (Mp + TL.MsTrue[self.plan2star])
1✔
634
        v1 = -np.sqrt(mu / (self.a * (1.0 - self.e**2))) * np.sin(nu)
1✔
635
        v2 = np.sqrt(mu / (self.a * (1.0 - self.e**2))) * (self.e + np.cos(nu))
1✔
636

637
        self.r = (A * r * np.cos(nu) + B * r * np.sin(nu)).T * u.AU  # position
1✔
638
        self.v = (A * v1 + B * v2).T.to("AU/day")  # velocity
1✔
639

640
        try:
1✔
641
            self.d = np.linalg.norm(self.r, axis=1)  # planet-star distance
1✔
642
            self.phi = PPMod.calc_Phi(
1✔
643
                np.arccos(self.r[:, 2].to("AU").value / self.d.to("AU").value) * u.rad,
644
                phiIndex=self.phiIndex,
645
            )  # planet phase
646
        except u.UnitTypeError:
×
647
            self.d = (
×
648
                np.linalg.norm(self.r, axis=1) * self.r.unit
649
            )  # planet-star distance
650
            self.phi = PPMod.calc_Phi(
×
651
                np.arccos(self.r[:, 2].to("AU").value / self.d.to("AU").value) * u.rad,
652
                phiIndex=self.phiIndex,
653
            )  # planet phase
654

655
        self.dMag = deltaMag(self.p, self.Rp, self.d, self.phi)  # delta magnitude
1✔
656

657
        try:
1✔
658
            self.s = np.linalg.norm(self.r[:, 0:2], axis=1)  # apparent separation
1✔
659
            self.WA = np.arctan(self.s / TL.dist[self.plan2star]).to(
1✔
660
                "arcsec"
661
            )  # working angle
662
        except u.UnitTypeError:
×
663
            self.s = (
×
664
                np.linalg.norm(self.r[:, 0:2], axis=1) * self.r.unit
665
            )  # apparent separation
666
            self.WA = np.arctan(self.s / TL.dist[self.plan2star]).to(
×
667
                "arcsec"
668
            )  # working angle
669

670
    def dump_systems(self):
1✔
671
        """Create a dictionary of planetary properties for archiving use.
672

673
        Args:
674
            None
675

676
        Returns:
677
            dict:
678
                Dictionary of planetary properties
679

680
        """
681

682
        systems = {
1✔
683
            "a": self.a,
684
            "e": self.e,
685
            "I": self.I,
686
            "O": self.O,
687
            "w": self.w,
688
            "M0": self.M0,
689
            "Mp": self.Mp,
690
            "mu": (
691
                const.G * (self.Mp + self.TargetList.MsTrue[self.plan2star])
692
            ).decompose(),
693
            "Rp": self.Rp,
694
            "p": self.p,
695
            "nEZ": self.nEZ,
696
            "plan2star": self.plan2star,
697
            "star": self.TargetList.Name[self.plan2star],
698
        }
699
        if self.commonSystemPlane:
1✔
700
            systems["systemInclination"] = self.TargetList.systemInclination
×
701
            systems["systemOmega"] = self.TargetList.systemOmega
×
702

703
        return systems
1✔
704

705
    def load_systems(self, systems):
1✔
706
        """Load a dictionary of planetary properties (nominally created by dump_systems)
707

708
        Args:
709
            systems (dict):
710
                Dictionary of planetary properties corresponding to the output of
711
                dump_systems.
712

713
            Returns:
714
                None
715

716
        .. note::
717

718
            If keyword ``systemInclination`` is present in the dictionary, it
719
            is assumed that it was generated with ``commonSystemPlane`` set to
720
            True.
721

722
        .. warning::
723

724
            This method assumes that the exact same targetlist is being used as in the
725
            object that generated the systems dictionary.  If this assumption is
726
            violated unexpected results may occur.
727
        """
728

729
        self.a = systems["a"]
1✔
730
        self.e = systems["e"]
1✔
731
        self.I = systems["I"]  # noqa: E741
1✔
732
        self.O = systems["O"]  # noqa: E741
1✔
733
        self.w = systems["w"]
1✔
734
        self.M0 = systems["M0"]
1✔
735
        self.Mp = systems["Mp"]
1✔
736
        self.Rp = systems["Rp"]
1✔
737
        self.p = systems["p"]
1✔
738
        self.nEZ = systems["nEZ"]
1✔
739
        self.plan2star = systems["plan2star"]
1✔
740

741
        if "systemInclination" in systems:
1✔
742
            self.TargetList.systemInclination = systems["systemInclination"]
×
743
            self.commonSystemPlane = True
×
744
            if "systemOmega" in systems:
×
745
                # leaving as if for backwards compatibility with old dumped
746
                # params for now
747
                self.TargetList.systemOmega = systems["systemOmega"]
×
748

749
        self.init_systems()
1✔
750

751
    def dump_system_params(self, sInd=None):
1✔
752
        """Create a dictionary of time-dependant planet properties for a specific target
753

754
        Args:
755
            sInd (int):
756
                Index of the target system of interest. Default value (None) will
757
                return an empty dictionary with the selected parameters and their units.
758

759
        Returns:
760
            dict:
761
                Dictionary of time-dependant planet properties
762

763
        """
764

765
        # get planet indices
766
        if sInd is None:
1✔
767
            pInds = np.array([], dtype=int)
×
768
        else:
769
            pInds = np.where(self.plan2star == sInd)[0]
1✔
770

771
        # build dictionary
772
        system_params = {
1✔
773
            "d": self.d[pInds],
774
            "phi": self.phi[pInds],
775
            "dMag": self.dMag[pInds],
776
            "WA": self.WA[pInds],
777
        }
778

779
        return system_params
1✔
780

781
    def revise_planets_list(self, pInds):
1✔
782
        """Replaces Simulated Universe planet attributes with filtered values,
783
        and updates the number of planets.
784

785
        Args:
786
            pInds (~numpy.ndarray(int)):
787
                Planet indices to keep
788

789
        Returns:
790
            None
791

792
        .. warning::
793

794
            Throws AssertionError if all planets are removed
795

796
        """
797

798
        # planet attributes which are floats and should not be filtered
799
        bad_atts = ["Min"]
1✔
800

801
        if len(pInds) == 0:
1✔
802
            raise IndexError("Planets list filtered to empty.")
×
803

804
        for att in self.planet_atts:
1✔
805
            if att not in bad_atts:
1✔
806
                if getattr(self, att).size != 0:
1✔
807
                    setattr(self, att, getattr(self, att)[pInds])
1✔
808
        self.nPlans = len(pInds)
1✔
809
        assert self.nPlans, "Planets list is empty: nPlans = %r" % self.nPlans
1✔
810

811
    def revise_stars_list(self, sInds):
1✔
812
        """Revises the TargetList with filtered values, and updates the
813
        planets list accordingly.
814

815
        Args:
816
            sInds (~numpy.ndarray(int)):
817
                Star indices to keep
818

819
        Returns:
820
            None
821

822
        """
823
        self.TargetList.revise_lists(sInds)
1✔
824
        pInds = np.sort(
1✔
825
            np.concatenate([np.where(self.plan2star == x)[0] for x in sInds])
826
        )
827
        self.revise_planets_list(pInds)
1✔
828
        for i, ind in enumerate(sInds):
1✔
829
            self.plan2star[np.where(self.plan2star == ind)[0]] = i
1✔
830

831
    def setup_system_planes(self):
1✔
832
        """
833
        Helper function that augments the system planes if
834
        commonSystemPlane is true
835

836
        Args:
837
            None
838

839
        Returns:
840
            None
841
        """
842
        if self.commonSystemPlane:
1✔
843
            self.I += self.TargetList.systemInclination[self.plan2star]
×
844
            # Ensure all inclinations are in [0, pi]
845
            self.I = (self.I.to(u.deg).value % 180) * u.deg
×
846

847
            self.O += self.TargetList.systemOmega[self.plan2star]
×
848
            # Cut longitude of the ascending nodes to [0, 2pi]
849
            self.O = (self.O.to(u.deg).value % 360) * u.deg
×
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2026 Coveralls, Inc