aboutsummaryrefslogtreecommitdiffstats
path: root/Helpers/kepler_helpers.py
blob: 1058870f4a1e88530afea1aa4e1e74b09784c949 (plain)
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
"""
Various Helper Function for Kepler Orbits.

They usually break down for inc=0 and/or ecc=0.
Might need to upgrade to Regular Equations some day
cf. http://server.faia.upm.es/moda/curso1112/kepler.pdf
"""

import numpy as np
import vector_helpers as vh

def helio2bary(r_h, v_h, m, return_sun=False, central_mass=1.0):
    """
    Convert Helio- to Barycentric Coordinates.
    """

    total_mass = np.sum(m) + central_mass

    # Barycenter Offsets
    r_b = np.sum(r_h * m) / total_mass
    v_b = np.sum(v_h * m) / total_mass

    # New Coordinates
    r = r_h - r_b
    v = v_h - v_b

    if return_sun:
        return r, v, -r_b, -v_b
    else:
        return r, v

def bary2helio(r_b, v_b, m, return_barycenter=False, central_mass=1.0):
    """
    Convert Bary- to Heliocentric Coordinates.
    """

    # Heliocenter Offsets
    r_h = - np.sum(r_b * m) / central_mass
    v_h = - np.sum(v_b * m) / central_mass

    # New Coordinates
    r = r_b - r_h
    v = v_b - v_h

    # Return
    if return_barycenter:
        return r, v, -r_h, -v_h
    else:
        return r, v

def cart2kepX(x, y, z, vx, vy, vz, mass, central_mass=1.0):
    """
    Vectorized version of cart2kep. Much love for many particles.

    @params
    r - (x,y,z) Cartesian Positions
    v - (vx,vy,vz) Cartesian Velocities
    mass - Particle Mass
    central_mass - Mass of Central Object

    @returns
    a - Semi-Major Axis
    ecc - Eccentricity
    inc - Inclination
    Omega - Longitude of the Ascending Node
    omega - Argument of Periapsis
    M - Mean Anomaly at Epoch

    Cf. http://www.bruce-shapiro.com/pair/ElementConversionRecipes.pdf
    """

    # Gravitational Parameter
    G = 1.0
    mu = G * ( central_mass + mass )

    # Angular Momentum Vector
    hx, hy, hz = vh.cross(x, y, z, vx, vy, vz)

    # Laplace-Runge-Lenz Vector
    # Scalar Eccentricity
    tmp_x, tmp_y, tmp_z = vh.cross(vx, vy, vz, hx, hy, hz)
    r_norm = vh.norm(x, y, z)
    lrl_x = tmp_x / mu - x / r_norm
    lrl_y = tmp_y / mu - y / r_norm
    lrl_z = tmp_z / mu - z / r_norm
    ecc = vh.norm(lrl_x, lrl_y, lrl_z)

    # Semi-Major Axis
    a = vh.dot(hx, hy, hz, hx, hy, hz) / ( mu * ( 1.0 - ecc**2.0 ) )

    # Inclination
    h_norm = vh.norm(hx, hy, hz)
    inc = np.arccos(vh.dot(0.0, 0.0, 1.0, hx, hy, hz) / h_norm)

    # Longitude of the Ascending Node
    nx, ny, nz = vh.cross(0.0, 0.0, 1.0, hx, hy, hz)
    n_norm = vh.norm(nx, ny, nz)
    tmp0 = 0.0
    tmpX = np.arccos(vh.dot(1.0, 0.0, 0.0, nx, ny, nz) / n_norm)
    Omega = np.where(inc==0.0, tmp0, tmpX)
    Omega[vh.dot(nx, ny, nz, 0.0, 1.0, 0.0) < 0.0] = \
        2.0 * np.pi - Omega[vh.dot(nx, ny, nz, 0.0, 1.0, 0.0) < 0.0]

    # Argument of Perigee
    # For Zero Inclination, Fall Back to 2D Case
    # http://en.wikipedia.org/wiki/Argument_of_periapsis
    tmp0 = np.arctan2(lrl_y/ecc, lrl_x/ecc)
    tmpX = np.arccos(vh.dot(nx, ny, nz, lrl_x, lrl_y, lrl_z) / (n_norm * ecc))
    omega = np.where(inc==0.0, tmp0, tmpX)
    omega[vh.dot(lrl_x, lrl_y, lrl_z, 0.0, 0.0, 1.0) < 0.0] = \
        2.0 * np.pi - omega[vh.dot(lrl_x, lrl_y, lrl_z, 0.0, 0.0, 1.0) < 0.0]

    # True Anomaly
    theta = np.arccos(vh.dot(lrl_x, lrl_y, lrl_z, x, y, z) / (ecc * r_norm))
    theta[vh.dot(x, y, z, vx, vy, vz) < 0.0] = \
        2.0 * np.pi - theta[vh.dot(x, y, z, vx, vy, vz) < 0.0]

    # Eccentric Anomaly
    E = np.arccos((ecc + np.cos(theta)) / (1 + ecc * np.cos(theta)))
    E[np.logical_and(np.pi < theta, theta < 2.0 * np.pi)] = \
        2.0 * np.pi - E[np.logical_and(np.pi < theta, theta < 2.0 * np.pi)]

    # Mean Anomaly
    M = E - ecc * np.sin(E)

    # Return Set
    return a, ecc, inc, Omega, omega, M

def cart2kep(r, v, mass, central_mass=1.0):
    """
    @params
    r - (x,y,z) Cartesian Positions
    v - (vx,vy,vz) Cartesian Velocities
    mass - Particle Mass
    central_mass - Mass of Central Object

    @returns
    a - Semi-Major Axis
    ecc - Eccentricity
    inc - Inclination
    Omega - Longitude of the Ascending Node
    omega - Argument of Periapsis
    M - Mean Anomaly at Epoch

    Cf. http://www.bruce-shapiro.com/pair/ElementConversionRecipes.pdf
    """

    # Gravitational Parameter
    G = 1.0
    mu = G * ( central_mass + mass )

    # Cartesian Unit Vectors
    iHat = np.array([1., 0., 0.])
    jHat = np.array([0., 1., 0.])
    kHat = np.array([0., 0., 1.])

    # Eccentricity
    h = np.cross(r, v)
    evec = 1. / mu * np.cross(v, h) - r / np.linalg.norm(r)
    ecc = np.linalg.norm(evec)

    # Semi Major Axis
    a = np.dot(h,h) / ( mu * ( 1. - ecc**2. ))

    # Inclination
    inc = np.arccos(np.dot(kHat, h) / np.linalg.norm(h))

    # Longitude of the Ascending Node
    n = np.cross(kHat, h)
    if inc == 0.0:
        Omega = 0.0
    else:
        Omega = np.arccos(np.dot(iHat, n) / np.linalg.norm(n))
        if np.dot(n, jHat) < 0:
            Omega = 2. * np.pi - Omega

    # Argument of Perigee
    # For Zero Inclination, Fall Back to 2D Case
    # http://en.wikipedia.org/wiki/Argument_of_periapsis
    if inc == 0.0:
        omega = np.arctan2(evec[1]/ecc,evec[0]/ecc)
    else:
        omega = np.arccos(np.dot(n, evec) / (np.linalg.norm(n) * ecc))
        if np.dot(evec, kHat) < 0:
            omega = 2. * np.pi - omega

    # True Anomaly
    theta = np.arccos(np.dot(evec, r) / (ecc * np.linalg.norm(r)))
    if np.dot(r,v) < 0:
        theta = 2. * np.pi - theta

    # Eccentric, Circular
    if ecc < 1.0:
        
        # Eccentric Anomaly
        E = np.arccos((ecc + np.cos(theta)) / (1 + ecc * np.cos(theta)))
        if np.pi < theta and theta < 2. * np.pi:
            E = 2. * np.pi - E

        # Mean Anomaly
        M = E - ecc * np.sin(E)
    
    # Hyperbolic
    elif ecc > 1.0:
        
        # Hyperbolic Anomaly
        H = np.arccosh((ecc + np.cos(theta)) / (1.0 + ecc * np.cos(theta)))
        
        # Mean Anomlay
        M = ecc * np.sinh(H) - H

        # Flip Around Focus?
        if np.pi < theta and theta < 2. * np.pi:
            M *= -1.0
    
    # Parabolic
    elif ecc == 1.0:
        raise Exception("Parabolic Orbits Unsupported.")

    # Return Set
    return a, ecc, inc, Omega, omega, M

def kep2cart(a, ecc, inc, Omega, omega, M, mass, central_mass=1.0):
    """
    @params
    a - Semi-Major Axis
    ecc - Eccentricity
    inc - Inclination
    Omega - Longitude of the Ascending Node
    omega - Argument of Periapsis
    M - Mean Anomaly at Epoch

    @returns
    r - (x,y,z) Cartesian Positions
    v - (vx,vy,vz) Cartesian Velocities

    Cf. http://www.bruce-shapiro.com/pair/ElementConversionRecipes.pdf
    """

    # Zero Inclination = No Line of Nodes
    if inc == 0.0:
        Omega = 0.0

    # Ciruclar Orbit = No Argument of Perigee
    if ecc == 0.0:
        omega = 0.0

    # Mean Anomaly -> Eccentric/Hyperbolic Anomaly
    E = nr(M, ecc)
        
    # Eccentric, Circular
    if ecc < 1.0:

        # Compute XY in Orbit Frame
        X = a * (np.cos(E) - ecc)
        Y = a * np.sqrt(1.0 - ecc**2.) * np.sin(E)

        # Compute VX,VY in Orbit Frame
        G = 1.0; mu = G * ( central_mass + mass )
        n = np.sqrt(mu / a**3.0)
        Edot = n / ( 1.0 - ecc * np.cos(E) )
        Vx = - a * np.sin(E) * Edot
        Vy =   a * np.sqrt(1.0 - ecc**2.0) * Edot * np.cos(E)

    # Hyperbolic
    elif ecc > 1.0:

        # Convenience
        H = E
        
        # Compute XY in Orbit Frame
        X = np.abs(a) * ( ecc - np.cosh(H) )
        Y = np.abs(a) * np.sqrt(ecc**2.0 - 1.0) * np.sinh(H)

        # Compute VX,VY in Orbit Frame
        G = 1.0; mu = G * ( central_mass + mass )
        n = np.sqrt(mu / np.abs(a)**3.0)
        Hdot = n / ( ecc * np.cosh(H) - 1.0 )
        Vx = - np.abs(a) * Hdot * np.sinh(H)
        Vy =   np.abs(a) * np.sqrt(ecc**2.0 - 1.0) * Hdot * np.cosh(H)

    # Parabolic
    elif ecc == 1.0:
        raise Exception("Parabolic Orbits Unsupported.")

    # Rotation Matrix Components
    Px, Py, Pz, Qx, Qy, Qz = PQW(Omega, omega, inc)

    # Rotate Positions
    x = X * Px + Y * Qx
    y = X * Py + Y * Qy
    z = X * Pz + Y * Qz

    # Rotate Velocities
    vx = Vx * Px + Vy * Qx
    vy = Vx * Py + Vy * Qy
    vz = Vx * Pz + Vy * Qz

    # Build Arrays
    r = np.array([x,y,z])
    v = np.array([vx,vy,vz])

    # Return
    return r, v

def PQW(Omega, omega, inc):
    """
    Rotation Matrix Components (Orbit Frame => Inertial Frame)
    http://biomathman.com/pair/KeplerElements.pdf (End)
    http://astro.geo.tu-dresden.de/~klioner/celmech.pdf (Eqn. 2.30)
    NB: Shapiro Notation Tricky; Multiplies x = R.T * X, Dresden x = R * X
    """

    Px = np.cos(omega) * np.cos(Omega) - \
         np.sin(omega) * np.cos(inc) * np.sin(Omega)
    Py = np.cos(omega) * np.sin(Omega) + \
         np.sin(omega) * np.cos(inc) * np.cos(Omega)
    Pz = np.sin(omega) * np.sin(inc)

    Qx = - np.sin(omega) * np.cos(Omega) - \
           np.cos(omega) * np.cos(inc) * np.sin(Omega)
    Qy = - np.sin(omega) * np.sin(Omega) + \
           np.cos(omega) * np.cos(inc) * np.cos(Omega)
    Qz =   np.sin(inc) * np.cos(omega)

    return Px, Py, Pz, Qx, Qy, Qz

def nr(M, ecc, epsilon_target=1.0e-5):
    """
    Newton-Raphson Iteration.
    Computes Eccentric/Hyperbolic Anomaly from Mean Anomaly.
    Cf. Slide 13/26
    http://mmae.iit.edu/~mpeet/Classes/MMAE441/Spacecraft/441Lecture17.pdf
    """

    Ei = M; ii = 1
    # print "Running Newton-Raphson for M=%.2e, e=%.2f" % ( M, ecc )
    while True:
        # Eccentric Anomaly
        if ecc < 1.0:
            Ei1 = Ei - \
                ( Ei - ecc * np.sin(Ei) - M ) / (  1.0 - ecc * np.cos(Ei)  )
        # Hyperbolic Anomaly
        elif ecc > 1.0:
            Ei1 = Ei + \
                ( M - ecc * np.sinh(Ei) + Ei ) / ( ecc * np.cosh(Ei) - 1.0 )
        epsilon = np.abs(Ei1 - Ei)
        Ei = Ei1
        # print "Iteration %i, Residual %.2e" % ( ii, epsilon )
        if epsilon < epsilon_target:
            break
        if ii > 100:
            raise Exception("NR Iteration Failed To Converge.")
        ii += 1
    # print "Found E=%.2f" % Ei1
    return Ei1

def compute_ellipseX(a, ecc, inc, Omega, omega):
    """
    Compute XYZ Sequence for a Kepler Ellipse.
    Vectorized Version. Expects 1D Arrays Passed.
    """

    # Some Reshaping
    a = a[:,np.newaxis]
    ecc = ecc[:,np.newaxis]
    inc = inc[:,np.newaxis]
    Omega = Omega[:,np.newaxis]
    omega = omega[:,np.newaxis]

    # Eccentric Anomaly
    E = np.linspace(0.0, 2.*np.pi, 128)
    E = E[np.newaxis,:]

    # Rotation Matrix Components
    Px, Py, Pz, Qx, Qy, Qz = PQW(Omega, omega, inc)

    # Compute Ellipse
    x = a * (np.cos(E) - ecc) * Px + \
        a * np.sqrt(1.0 - ecc**2.) * np.sin(E) * Qx
    y = a * (np.cos(E) - ecc) * Py + \
        a * np.sqrt(1.0 - ecc**2.) * np.sin(E) * Qy
    z = a * (np.cos(E) - ecc) * Pz + \
        a * np.sqrt(1.0 - ecc**2.) * np.sin(E) * Qz

    # Return
    return x, y, z

def compute_ellipse(a, ecc, inc, Omega, omega):
    """
    Compute XYZ Sequence for a Kepler Ellipses/Hyperbolas.
    """
  
    # Eccentric, Circular
    if ecc < 1.0:
        
        # Eccentric Anomaly
        E = np.linspace(0.0, 2.*np.pi, 128)
        
        # Compute XY in Orbit Frame
        X = a * (np.cos(E) - ecc)
        Y = a * np.sqrt(1.0 - ecc**2.) * np.sin(E)
    
    # Hyperbolic
    elif ecc > 1.0:
        
        # Hyperbolic Anomaly
        H = np.linspace(-np.pi/2.0, np.pi/2.0, 128)
        
        # Compute XY in Orbit Frame
        X = np.abs(a) * ( ecc - np.cosh(H) )
        Y = np.abs(a) * np.sqrt(ecc**2.0 - 1.0) * np.sinh(H)
    
    # Parabolic
    elif ecc == 1.0:
        raise Exception("Parabolic Orbits Unsupported.")

    # Rotation Matrix Components
    Px, Py, Pz, Qx, Qy, Qz = PQW(Omega, omega, inc)

    # Rotate Positions
    x = X * Px + Y * Qx
    y = X * Py + Y * Qy
    z = X * Pz + Y * Qz
    
    # Return
    return x, y, z

def kep2del(a, e, i, Omega, omega, M):
    """
    Keper to Delaunay Elements.
    Cf. (a) Joachim's Unpublished Paper
        (b) http://www.bourbaphy.fr/chenciner.pdf
    """

    Msolar = 1.99e30 # kg
    m = 5.0 * 5.97e24 / 2048.0
    mu = m**2.0 * Msolar
    Lambda = np.sqrt(mu * a)
    L = np.sqrt(mu * a * ( 1.0 - e**2.0 ))
    Lz = L * np.cos(i)
    return Lambda, L, Lz, Omega, omega, M

def del2poi(Lambda, L, Lz, Omega, omega, M):
    """
    Delaunay to Poincare Elements.
    Cf. (a) Joachim's Unpublished Paper
        (b) http://www.bourbaphy.fr/chenciner.pdf
    """

    lambda_small = M + omega + Omega
    lambda_small = lambda_small % ( 2.0 * np.pi )
    xi_real = np.sqrt(2.0 * ( Lambda - L )) * np.cos(omega + Omega)
    xi_imag = np.sqrt(2.0 * ( Lambda - L )) * np.sin(omega + Omega)
    eta_real = np.sqrt(2.0 * ( L - Lz )) * np.cos(Omega)
    eta_imag = np.sqrt(2.0 * ( L - Lz )) * np.sin(Omega)
    return Lambda, lambda_small, xi_real, xi_imag, eta_real, eta_imag

def kep2metric(a1, a2, e1, e2, i1, i2, \
               Omega1, Omega2, omega1, omega2):
    """
    Computes Metric Distance from Keplerian Elements.
    Cf. http://adsabs.harvard.edu/abs/2008CeMDA.100..169K
    Eq. 21

    This breaks down for small inclinations, eccentricities.
    Here, small changes in coordinates cause omega, Omega to vary hugely.
    Thus, we cannot compute small differences anymore.
     
    Calculate directly from state vector (x,v) in this case.
    See below in cart2metic, cart2metricX.
    """

    # Some sort of normalization
    L = 1.0
    # Help Me 01
    cos_xi = np.cos(i1) * np.cos(i2) + np.sin(i1) * np.sin(i2) * np.cos(Omega1 - Omega2)
    # Help Me 02
    cos_eta = ( np.cos(omega1) * np.cos(omega2) + np.cos(i1) * np.cos(i2) * np.sin(omega1) * np.sin(omega2) ) * np.cos(Omega1 - Omega2) \
            + ( np.cos(i2) * np.cos(omega1) * np.sin(omega2) - np.cos(i1) * np.sin(omega1) * np.cos(omega2) ) * np.sin(Omega1 - Omega2) \
            + np.sin(i1) * np.sin(i2) * np.sin(omega1) * np.sin(omega2)
    # One More Conversion
    p1 = a1 * ( 1.0 - e1**2.0 )
    p2 = a2 * ( 1.0 - e2**2.0 )
    # Hit Me
    rho2 = 1.0 / L * ( p1 + p2 - 2.0 * np.sqrt(p1*p2) * cos_xi ) + ( e1**2.0 + e2**2.0 - 2.0 * e1 * e2 * cos_eta )
    return rho2

def cart2metric(r1, v1, r2, v2):
    """
    Computes Metric Distance.
    Cf. Kholshevnikov 2007, Sec. 3, Eq. (4)
    """
    # mu = G * ( central_mass + mass )
    mu = 1.0
    L = 1.0
    # Specific Angular Momentum Vectors
    h1 = np.cross(r1,v1)
    h2 = np.cross(r2,v2)
    # Laplace-Runge-Lenz Vectors / Eccentricity Vectors
    e1 = np.cross(v1,h1)/mu - r1/np.linalg.norm(r1)
    e2 = np.cross(v2,h2)/mu - r2/np.linalg.norm(r2)
    # Distance (Squared)
    dh = h1 - h2
    de = e1 - e2
    rho2 = 1.0/mu/L * np.dot(dh,dh) + np.dot(de,de)
    return rho2

def cart2metricX(x1, y1, z1, vx1, vy1, vz1, \
                 x2, y2, z2, vx2, vy2, vz2, \
                 vanilla=True):
    """
    Computes Metric Distance.
    Cf. Kholshevnikov 2007, Sec. 3, Eq. (4)
    """
    # Gravitational Constant
    mu = 1.0
    # Scale Factors
    L = 1.0
    L1 = 1.0
    # Specific Angular Momentum Vectors
    hx1, hy1, hz1 = vh.cross(x1, y1, z1, vx1, vy1, vz1)
    hx2, hy2, hz2 = vh.cross(x2, y2, z2, vx2, vy2, vz2)
    # Laplace-Runge-Lenz Vectors / Eccentricity Vectors
    vx1, vy1, vz1 = vh.cross(vx1, vy1, vz1, hx1, hy1, hz1)
    vx2, vy2, vz2 = vh.cross(vx2, vy2, vz2, hx2, hy2, hz2)
    ex1 = vx1 / mu - x1 / vh.norm(x1, y1, z1)
    ey1 = vy1 / mu - y1 / vh.norm(x1, y1, z1)
    ez1 = vz1 / mu - z1 / vh.norm(x1, y1, z1)
    ex2 = vx2 / mu - x2 / vh.norm(x2, y2, z2)
    ey2 = vy2 / mu - y2 / vh.norm(x2, y2, z2)
    ez2 = vz2 / mu - z2 / vh.norm(x2, y2, z2)
    # Energy Constants
    E1 = vh.dot(vx1, vy1, vz1, vx1, vy1, vz1) / 2.0 - mu / vh.norm(x1, y1, z1)
    E2 = vh.dot(vx2, vy2, vz2, vx2, vy2, vz2) / 2.0 - mu / vh.norm(x2, y2, z2)
    # Coordinate Distances (Squared)
    dh2 = (hx1 - hx2)**2.0 + (hy1 - hy2)**2.0 + (hz1 - hz2)**2.0
    de2 = (ex1 - ex2)**2.0 + (ey1 - ey2)**2.0 + (ez1 - ez2)**2.0
    dE2 = (E1 - E2)**2.0
    # Total Distances (Squared)
    rho2 = 1.0/mu/L * dh2 + de2
    rho2e = 1.0/mu/L * dh2 + de2 + L1**2.0 / mu**2.0 * dE2
    # Return Values
    if vanilla:
        return rho2
    else:
        return rho2, rho2e, dh2, de2, dE2