Geophysical Models

This example demonstrates gravity and magnetic field models essential for high-precision navigation, geodesy, and aerospace applications.

Overview

Geophysical models are critical for:

  • Inertial navigation: Gravity compensation in INS

  • Geodesy: Height reference and surveying

  • Aerospace: Satellite orbit determination

  • Geophysics: Subsurface exploration

Gravity Models

Normal Gravity (Somigliana)
  • Gravity variation with latitude

  • ~0.5% increase from equator to poles

  • Due to Earth’s rotation and flattening

WGS84 Gravity Model
  • Full gravity vector computation

  • Includes deflection of vertical

  • Essential for inertial navigation

J2 Gravity Model
  • Simplified model using J2 oblateness term

  • Adequate for many applications

  • Faster computation than full models

Geoid Height
  • Separation between geoid and ellipsoid

  • J2 approximation captures main flattening

  • Full models (EGM96/2008) for precision

Gravity Anomalies
  • Free-air anomaly

  • Gravity disturbance

  • Used for geophysical exploration

Tidal Effects
  • Solid Earth tide displacement (~30 cm)

  • Tidal gravity variations (~300 uGal)

  • Essential for precision gravimetry

Magnetic Field Models

World Magnetic Model (WMM2020)
  • Standard model for navigation

  • Updated every 5 years

  • Valid 2020-2025

IGRF-13
  • International Geomagnetic Reference Field

  • Historical and predictive coefficients

  • Used for scientific applications

Key Parameters
  • Declination: angle between true and magnetic north

  • Inclination: dip angle of field lines

  • Total intensity: field strength in nT

Notable Features
  • South Atlantic Anomaly: weak field region

  • Magnetic poles: ~11° offset from geographic

  • Secular variation: field changes over time

Code Highlights

The example demonstrates:

  • Normal gravity with normal_gravity_somigliana()

  • WGS84 gravity with gravity_wgs84()

  • Geoid height with geoid_height_j2()

  • Free-air anomaly with free_air_anomaly()

  • Solid Earth tides with solid_earth_tide_displacement()

  • Magnetic field with wmm() and igrf()

  • Magnetic declination with magnetic_declination()

Source Code

  1"""
  2Geophysical Models Example
  3==========================
  4
  5This example demonstrates the geophysical models in PyTCL:
  6
  7Gravity Models:
  8- Normal gravity (Somigliana formula)
  9- WGS84 and J2 gravity models
 10- Spherical harmonic expansions
 11- Geoid height computation
 12- Gravity anomalies and disturbances
 13- Tidal effects (solid Earth, ocean loading)
 14
 15Magnetic Field Models:
 16- World Magnetic Model (WMM2020)
 17- International Geomagnetic Reference Field (IGRF-13)
 18- Enhanced Magnetic Model (EMM)
 19- Magnetic declination, inclination, and intensity
 20
 21These models are essential for high-precision navigation, geodesy,
 22and aerospace applications.
 23"""
 24
 25from pathlib import Path
 26
 27import numpy as np
 28import plotly.graph_objects as go
 29from plotly.subplots import make_subplots
 30
 31# Output directory for generated plots
 32OUTPUT_DIR = Path(__file__).parent.parent / "docs" / "_static" / "images" / "examples"
 33OUTPUT_DIR.mkdir(parents=True, exist_ok=True)
 34
 35# Global flag to control plotting
 36SHOW_PLOTS = True
 37
 38
 39from pytcl.gravity import (  # Normal gravity; Gravity models; Geoid; Anomalies; Tidal effects; Constants
 40    GRS80,
 41    WGS84,
 42    free_air_anomaly,
 43    geoid_height,
 44    geoid_height_j2,
 45    gravity_anomaly,
 46    gravity_disturbance,
 47    gravity_j2,
 48    gravity_wgs84,
 49    normal_gravity,
 50    normal_gravity_somigliana,
 51    solid_earth_tide_displacement,
 52    solid_earth_tide_gravity,
 53    tidal_gravity_correction,
 54)
 55from pytcl.magnetism import (  # World Magnetic Model; IGRF
 56    dipole_moment,
 57    igrf,
 58    igrf_declination,
 59    igrf_inclination,
 60    magnetic_declination,
 61    magnetic_field_intensity,
 62    magnetic_inclination,
 63    magnetic_north_pole,
 64    wmm,
 65)
 66
 67
 68def demo_normal_gravity():
 69    """Demonstrate normal gravity computation."""
 70    print("=" * 70)
 71    print("Normal Gravity Demo")
 72    print("=" * 70)
 73
 74    # Test locations at different latitudes
 75    locations = [
 76        ("Equator", 0.0, 0.0, 0.0),
 77        ("Washington DC", 38.9, -77.0, 0.0),
 78        ("North Pole", 90.0, 0.0, 0.0),
 79        ("Mount Everest", 27.99, 86.93, 8848.0),
 80        ("Dead Sea", 31.5, 35.5, -430.0),
 81    ]
 82
 83    print("\nNormal gravity at various locations:")
 84    print("-" * 60)
 85    print(f"{'Location':<20} {'Lat':>8} {'Alt (m)':>10} {'g (m/s²)':>12}")
 86    print("-" * 60)
 87
 88    for name, lat, lon, alt in locations:
 89        # Normal gravity at sea level
 90        g_somigliana = normal_gravity_somigliana(np.radians(lat))
 91
 92        # Gravity at altitude (free-air correction)
 93        g_at_alt = normal_gravity(np.radians(lat), alt)
 94
 95        print(f"{name:<20} {lat:>8.2f} {alt:>10.1f} {g_at_alt:>12.6f}")
 96
 97    # Show latitude variation
 98    print("\n--- Latitude Variation ---")
 99    lats = np.array([0, 15, 30, 45, 60, 75, 90])
100    for lat in lats:
101        g = normal_gravity_somigliana(np.radians(lat))
102        print(f"  {lat:>2}°: {g:.6f} m/s²")
103
104    print("\nNote: Gravity increases from equator to poles due to")
105    print("Earth's rotation (centrifugal) and flattening effects.")
106
107    # Plot gravity variation with latitude
108    if SHOW_PLOTS:
109        # High-resolution latitude array
110        lats_fine = np.linspace(0, 90, 181)
111        g_values = np.array(
112            [normal_gravity_somigliana(np.radians(lat)) for lat in lats_fine]
113        )
114
115        fig = make_subplots(
116            rows=1,
117            cols=2,
118            subplot_titles=(
119                "Normal Gravity vs Latitude (Somigliana)",
120                "Gravity Increase from Equator",
121            ),
122        )
123
124        # Left plot: Gravity vs latitude
125        fig.add_trace(
126            go.Scatter(
127                x=lats_fine,
128                y=g_values,
129                mode="lines",
130                name="Gravity",
131                line=dict(color="blue", width=2),
132            ),
133            row=1,
134            col=1,
135        )
136        # Mark equator and pole values
137        fig.add_hline(
138            y=g_values[0],
139            line_dash="dash",
140            line_color="red",
141            annotation_text=f"Equator: {g_values[0]:.4f}",
142            row=1,
143            col=1,
144        )
145        fig.add_hline(
146            y=g_values[-1],
147            line_dash="dash",
148            line_color="green",
149            annotation_text=f"Pole: {g_values[-1]:.4f}",
150            row=1,
151            col=1,
152        )
153
154        # Right plot: Gravity difference from equator
155        g_diff = (g_values - g_values[0]) * 1000  # mGal
156        fig.add_trace(
157            go.Scatter(
158                x=lats_fine,
159                y=g_diff,
160                mode="lines",
161                name="Δg",
162                line=dict(color="blue", width=2),
163            ),
164            row=1,
165            col=2,
166        )
167
168        fig.update_xaxes(title_text="Latitude (°)", range=[0, 90], row=1, col=1)
169        fig.update_yaxes(title_text="Normal Gravity (m/s²)", row=1, col=1)
170        fig.update_xaxes(title_text="Latitude (°)", range=[0, 90], row=1, col=2)
171        fig.update_yaxes(title_text="Δg from Equator (mGal)", row=1, col=2)
172
173        fig.update_layout(height=500, width=1200, showlegend=False)
174        fig.write_html(str(OUTPUT_DIR / "geophysical_gravity_latitude.html"))
175        print("\n  [Plot saved to geophysical_gravity_latitude.html]")
176
177
178def demo_gravity_models():
179    """Demonstrate different gravity models."""
180    print("\n" + "=" * 70)
181    print("Gravity Models Comparison Demo")
182    print("=" * 70)
183
184    # Test point: Washington DC
185    lat = np.radians(38.9)
186    lon = np.radians(-77.0)
187    alt = 0.0  # Sea level
188
189    print(f"\nLocation: Washington DC (38.9°N, 77.0°W)")
190    print("-" * 50)
191
192    # WGS84 gravity model
193    g_wgs84 = gravity_wgs84(lat, lon, alt)
194    print(f"\nWGS84 gravity model:")
195    print(f"  Total gravity: {g_wgs84.magnitude:.6f} m/s²")
196    print(f"  Down component: {g_wgs84.g_down:.6f} m/s²")
197    print(f"  North component: {g_wgs84.g_north:.9f} m/s²")
198
199    # J2 gravity model (simpler, uses only J2 term)
200    g_j2 = gravity_j2(lat, lon, alt)
201    print(f"\nJ2 gravity model:")
202    print(f"  Total gravity: {g_j2.magnitude:.6f} m/s²")
203    print(
204        f"  Difference from WGS84: {(g_j2.magnitude - g_wgs84.magnitude)*1e6:.2f} µGal"
205    )
206
207    # Compare at different altitudes
208    print("\n--- Gravity vs Altitude ---")
209    altitudes = [0, 1000, 10000, 100000, 400000]  # meters
210    for alt in altitudes:
211        g = gravity_wgs84(lat, lon, alt)
212        print(f"  {alt:>7} m: {g.magnitude:.6f} m/s²")
213
214    print("\nNote: Gravity decreases with altitude approximately as")
215    print("g(h) ≈ g₀(1 - 2h/R), about 3.1 mGal per meter at low altitudes.")
216
217    # Plot gravity vs altitude
218    if SHOW_PLOTS:
219        fig = make_subplots(
220            rows=1,
221            cols=2,
222            subplot_titles=(
223                "Gravity vs Altitude (WGS84)",
224                "Gravity Reduction with Altitude",
225            ),
226        )
227
228        # Altitude range from surface to ISS altitude
229        alts = np.linspace(0, 500000, 500)  # 0 to 500 km
230        g_values = np.array([gravity_wgs84(lat, lon, a).magnitude for a in alts])
231
232        # Left plot: Gravity vs altitude
233        fig.add_trace(
234            go.Scatter(
235                x=alts / 1000,
236                y=g_values,
237                mode="lines",
238                name="Gravity",
239                line=dict(color="blue", width=2),
240            ),
241            row=1,
242            col=1,
243        )
244
245        # Right plot: Gravity reduction rate
246        g_reduction = (g_values[0] - g_values) / g_values[0] * 100  # percent
247        fig.add_trace(
248            go.Scatter(
249                x=alts / 1000,
250                y=g_reduction,
251                mode="lines",
252                name="Reduction",
253                line=dict(color="red", width=2),
254            ),
255            row=1,
256            col=2,
257        )
258
259        # Mark ISS at ~400 km
260        fig.add_hline(
261            y=g_reduction[400],
262            line_dash="dash",
263            line_color="green",
264            annotation_text=f"ISS (~400 km): {g_reduction[400]:.1f}% reduction",
265            row=1,
266            col=2,
267        )
268
269        fig.update_xaxes(title_text="Altitude (km)", row=1, col=1)
270        fig.update_yaxes(title_text="Gravity (m/s²)", row=1, col=1)
271        fig.update_xaxes(title_text="Altitude (km)", row=1, col=2)
272        fig.update_yaxes(title_text="Gravity Reduction (%)", row=1, col=2)
273
274        fig.update_layout(height=500, width=1200, showlegend=False)
275        fig.write_html(str(OUTPUT_DIR / "geophysical_gravity_altitude.html"))
276        print("\n  [Plot saved to geophysical_gravity_altitude.html]")
277
278
279def demo_geoid():
280    """Demonstrate geoid height computation."""
281    print("\n" + "=" * 70)
282    print("Geoid Height Demo")
283    print("=" * 70)
284
285    # Various latitudes (J2 geoid approximation is zonal - only latitude dependent)
286    locations = [
287        ("Equator", 0.0),
288        ("Mid-latitude (30N)", 30.0),
289        ("Mid-latitude (45N)", 45.0),
290        ("High latitude (60N)", 60.0),
291        ("Pole", 90.0),
292    ]
293
294    print("\nJ2 Geoid heights (zonal approximation - latitude dependent only):")
295    print("-" * 40)
296    print(f"{'Location':<25} {'Lat':>7} {'N (m)':>10}")
297    print("-" * 40)
298
299    for name, lat in locations:
300        # J2 geoid approximation - only depends on latitude
301        N = geoid_height_j2(np.radians(lat))
302        print(f"{name:<25} {lat:>7.1f} {N:>10.2f}")
303
304    print("\nNote: The J2 geoid model captures the main flattening effect.")
305    print("Real geoid variations are ±100m from the ellipsoid (need EGM96/2008).")
306
307
308def demo_gravity_anomalies():
309    """Demonstrate gravity anomaly computation."""
310    print("\n" + "=" * 70)
311    print("Gravity Anomalies Demo")
312    print("=" * 70)
313
314    # Simulated gravity survey
315    np.random.seed(42)
316
317    # Locations along a survey line
318    n_points = 5
319    lats = np.linspace(38.0, 39.0, n_points)
320    lons = np.full(n_points, -77.0)
321    alts = np.array([100, 150, 200, 180, 120])  # meters
322
323    # Simulated observed gravity (with anomaly)
324    g_normal = np.array(
325        [normal_gravity(np.radians(lat), alt) for lat, alt in zip(lats, alts)]
326    )
327    # Add a gravity anomaly (e.g., from subsurface density variation)
328    anomaly_true = np.array([0.0, 10.0, 25.0, 15.0, 5.0]) * 1e-5  # m/s²
329    g_observed = g_normal + anomaly_true
330
331    print("\nGravity survey along profile:")
332    print("-" * 65)
333    print(f"{'Point':>5} {'Lat':>7} {'Alt(m)':>7} {'g_obs':>12} {'FAA':>12}")
334    print("-" * 65)
335
336    for i in range(n_points):
337        lat_rad = np.radians(lats[i])
338        faa = free_air_anomaly(g_observed[i], lat_rad, alts[i])
339        print(
340            f"{i+1:>5} {lats[i]:>7.2f} {alts[i]:>7.0f} "
341            f"{g_observed[i]:>12.6f} {faa*1e5:>12.2f} mGal"
342        )
343
344    print("\nNote: Free-air anomaly removes the effect of elevation")
345    print("to reveal subsurface density variations.")
346
347
348def demo_tidal_effects():
349    """Demonstrate tidal effects on gravity and position."""
350    print("\n" + "=" * 70)
351    print("Tidal Effects Demo")
352    print("=" * 70)
353
354    # Observer location
355    lat = np.radians(38.9)  # Washington DC
356    lon = np.radians(-77.0)
357
358    # Time (Julian date) - approximate
359    # Let's use a few different times
360    times = [
361        ("2025-01-01 00:00 UTC", 2460676.5),
362        ("2025-01-01 06:00 UTC", 2460676.75),
363        ("2025-01-01 12:00 UTC", 2460677.0),
364        ("2025-01-01 18:00 UTC", 2460677.25),
365    ]
366
367    print(f"\nLocation: Washington DC (38.9°N, 77.0°W)")
368    print("\nSolid Earth tide effects (displacement and gravity):")
369    print("-" * 60)
370
371    for time_str, jd in times:
372        # Solid Earth tide displacement
373        disp = solid_earth_tide_displacement(lat, lon, jd)
374
375        # Tidal gravity correction
376        dg = tidal_gravity_correction(lat, lon, jd)
377
378        print(f"\n{time_str}:")
379        print(
380            f"  Displacement: ({disp[0]*1000:.1f}, {disp[1]*1000:.1f}, "
381            f"{disp[2]*1000:.1f}) mm (N, E, Up)"
382        )
383        print(f"  Gravity change: {dg*1e8:.2f} µGal")
384
385    print("\nNote: Solid Earth tides cause displacements of ~30 cm")
386    print("and gravity changes of ~300 µGal peak-to-peak.")
387
388    # Plot tidal effects over 24 hours
389    if SHOW_PLOTS:
390        # Time array: 24 hours at 15-minute intervals
391        hours = np.linspace(0, 24, 97)
392        jd_start = 2460676.5  # 2025-01-01 00:00 UTC
393        jds = jd_start + hours / 24
394
395        # Compute displacements and gravity changes
396        disp_n = np.zeros_like(hours)
397        disp_e = np.zeros_like(hours)
398        disp_u = np.zeros_like(hours)
399        dg = np.zeros_like(hours)
400
401        for i, jd in enumerate(jds):
402            disp = solid_earth_tide_displacement(lat, lon, jd)
403            disp_n[i] = disp[0] * 1000  # mm
404            disp_e[i] = disp[1] * 1000
405            disp_u[i] = disp[2] * 1000
406            dg[i] = tidal_gravity_correction(lat, lon, jd) * 1e8  # µGal
407
408        fig = make_subplots(
409            rows=2,
410            cols=1,
411            subplot_titles=(
412                "Solid Earth Tide - Washington DC (2025-01-01)",
413                "Tidal Gravity Variation",
414            ),
415            shared_xaxes=True,
416        )
417
418        # Top plot: Displacement components
419        fig.add_trace(
420            go.Scatter(
421                x=hours,
422                y=disp_n,
423                mode="lines",
424                name="North",
425                line=dict(color="blue", width=1.5),
426            ),
427            row=1,
428            col=1,
429        )
430        fig.add_trace(
431            go.Scatter(
432                x=hours,
433                y=disp_e,
434                mode="lines",
435                name="East",
436                line=dict(color="green", width=1.5),
437            ),
438            row=1,
439            col=1,
440        )
441        fig.add_trace(
442            go.Scatter(
443                x=hours,
444                y=disp_u,
445                mode="lines",
446                name="Up",
447                line=dict(color="red", width=2),
448            ),
449            row=1,
450            col=1,
451        )
452        fig.add_hline(y=0, line_color="black", line_width=0.5, row=1, col=1)
453
454        # Bottom plot: Gravity change
455        fig.add_trace(
456            go.Scatter(
457                x=hours,
458                y=dg,
459                mode="lines",
460                name="Gravity",
461                line=dict(color="purple", width=2),
462                fill="tozeroy",
463                fillcolor="rgba(128,0,128,0.3)",
464            ),
465            row=2,
466            col=1,
467        )
468        fig.add_hline(y=0, line_color="black", line_width=0.5, row=2, col=1)
469
470        fig.update_xaxes(title_text="Hour (UTC)", range=[0, 24], row=2, col=1)
471        fig.update_yaxes(title_text="Displacement (mm)", row=1, col=1)
472        fig.update_yaxes(title_text="Gravity Change (µGal)", row=2, col=1)
473
474        fig.update_layout(height=600, width=1000)
475        fig.write_html(str(OUTPUT_DIR / "geophysical_tides.html"))
476        print("\n  [Plot saved to geophysical_tides.html]")
477
478
479def demo_magnetic_field():
480    """Demonstrate magnetic field computation."""
481    print("\n" + "=" * 70)
482    print("Magnetic Field Models Demo")
483    print("=" * 70)
484
485    # Test locations
486    locations = [
487        ("Washington DC", 38.9, -77.0, 0.0),
488        ("London", 51.5, -0.1, 0.0),
489        ("Sydney", -33.9, 151.2, 0.0),
490        ("Magnetic North", 86.5, -175.3, 0.0),
491        ("Magnetic Equator", 0.0, 0.0, 0.0),
492        ("South Atlantic Anomaly", -25.0, -50.0, 0.0),
493    ]
494
495    # Use 2025 epoch
496    decimal_year = 2025.0
497
498    print(f"\nMagnetic field at various locations (WMM2020, {decimal_year}):")
499    print("-" * 75)
500    print(f"{'Location':<25} {'Dec':>8} {'Inc':>8} {'F (nT)':>10}")
501    print("-" * 75)
502
503    for name, lat, lon, alt in locations:
504        # Compute magnetic field using WMM
505        result = wmm(np.radians(lat), np.radians(lon), alt, decimal_year)
506
507        dec = np.degrees(result.D)  # Declination
508        inc = np.degrees(result.I)  # Inclination
509        F = result.F  # Total intensity
510
511        print(f"{name:<25} {dec:>8.2f}° {inc:>8.2f}° {F:>10.0f}")
512
513    # Show magnetic declination map concept
514    print("\n--- Magnetic Declination Grid ---")
515    lats_grid = np.array([-60, -30, 0, 30, 60])
516    lons_grid = np.array([-120, -60, 0, 60, 120])
517
518    print(f"{'Lat\\Lon':<8}", end="")
519    for lon in lons_grid:
520        print(f"{lon:>8}°", end="")
521    print()
522
523    for lat in lats_grid:
524        print(f"{lat:>6}°  ", end="")
525        for lon in lons_grid:
526            dec = magnetic_declination(
527                np.radians(lat), np.radians(lon), 0.0, decimal_year
528            )
529            print(f"{np.degrees(dec):>8.1f}", end="")
530        print()
531
532    # Plot magnetic declination and field intensity maps
533    if SHOW_PLOTS:
534        # Create grid for plotting
535        lat_grid = np.linspace(-80, 80, 33)
536        lon_grid = np.linspace(-180, 180, 73)
537        LAT, LON = np.meshgrid(lat_grid, lon_grid)
538
539        # Compute declination and intensity on grid
540        DEC = np.zeros_like(LAT)
541        F = np.zeros_like(LAT)
542        for i in range(LAT.shape[0]):
543            for j in range(LAT.shape[1]):
544                result = wmm(
545                    np.radians(LAT[i, j]), np.radians(LON[i, j]), 0.0, decimal_year
546                )
547                DEC[i, j] = np.degrees(result.D)
548                F[i, j] = result.F
549
550        fig = make_subplots(
551            rows=1,
552            cols=2,
553            subplot_titles=(
554                f"Magnetic Declination (WMM {decimal_year:.0f})",
555                f"Magnetic Field Intensity (WMM {decimal_year:.0f})",
556            ),
557        )
558
559        # Left plot: Magnetic declination
560        fig.add_trace(
561            go.Contour(
562                x=lon_grid,
563                y=lat_grid,
564                z=DEC.T,
565                colorscale="RdBu_r",
566                contours=dict(showlines=True),
567                colorbar=dict(title="Declination (°)", x=0.45),
568            ),
569            row=1,
570            col=1,
571        )
572
573        # Right plot: Total field intensity
574        fig.add_trace(
575            go.Contour(
576                x=lon_grid,
577                y=lat_grid,
578                z=F.T,
579                colorscale="Viridis",
580                colorbar=dict(title="Intensity (nT)", x=1.0),
581            ),
582            row=1,
583            col=2,
584        )
585
586        # Mark South Atlantic Anomaly region
587        fig.add_trace(
588            go.Scatter(
589                x=[-50],
590                y=[-25],
591                mode="markers",
592                marker=dict(symbol="star", size=15, color="red"),
593                name="South Atlantic Anomaly",
594                showlegend=True,
595            ),
596            row=1,
597            col=2,
598        )
599
600        fig.update_xaxes(title_text="Longitude (°)", range=[-180, 180], row=1, col=1)
601        fig.update_yaxes(title_text="Latitude (°)", range=[-80, 80], row=1, col=1)
602        fig.update_xaxes(title_text="Longitude (°)", range=[-180, 180], row=1, col=2)
603        fig.update_yaxes(title_text="Latitude (°)", range=[-80, 80], row=1, col=2)
604
605        fig.update_layout(height=500, width=1400)
606        fig.write_html(str(OUTPUT_DIR / "geophysical_magnetic_field.html"))
607        print("\n  [Plot saved to geophysical_magnetic_field.html]")
608
609
610def demo_magnetic_properties():
611    """Demonstrate magnetic field properties and variations."""
612    print("\n" + "=" * 70)
613    print("Magnetic Field Properties Demo")
614    print("=" * 70)
615
616    decimal_year = 2025.0
617
618    # Magnetic pole location
619    pole_lat, pole_lon = magnetic_north_pole(decimal_year)
620    print(f"\nMagnetic North Pole location ({decimal_year}):")
621    print(f"  Latitude: {np.degrees(pole_lat):.2f}°N")
622    print(f"  Longitude: {np.degrees(pole_lon):.2f}°E")
623
624    # Earth's dipole moment (using WMM2020 coefficients)
625    moment = dipole_moment()  # Uses default WMM2020 coefficients
626    print(f"\nEarth's dipole moment (WMM2020): {moment:.4e} A·m²")
627
628    # Compare WMM and IGRF at a test location
629    lat, lon, alt = np.radians(40.0), np.radians(-100.0), 0.0
630
631    wmm_result = wmm(lat, lon, alt, decimal_year)
632    igrf_result = igrf(lat, lon, alt, decimal_year)
633
634    print(f"\nModel comparison at (40°N, 100°W):")
635    print("-" * 40)
636    print(f"{'Component':<12} {'WMM':>12} {'IGRF':>12}")
637    print("-" * 40)
638    print(f"{'X (nT)':<12} {wmm_result.X:>12.1f} {igrf_result.X:>12.1f}")
639    print(f"{'Y (nT)':<12} {wmm_result.Y:>12.1f} {igrf_result.Y:>12.1f}")
640    print(f"{'Z (nT)':<12} {wmm_result.Z:>12.1f} {igrf_result.Z:>12.1f}")
641    print(f"{'F (nT)':<12} {wmm_result.F:>12.1f} {igrf_result.F:>12.1f}")
642
643    # Secular variation
644    print("\n--- Secular Variation ---")
645    print("Magnetic field changes over time due to core dynamics.")
646    years = [2020, 2022, 2024, 2025]
647    for year in years:
648        dec = magnetic_declination(lat, lon, alt, float(year))
649        print(f"  {year}: declination = {np.degrees(dec):.2f}°")
650
651
652def demo_navigation_application():
653    """Demonstrate application to navigation."""
654    print("\n" + "=" * 70)
655    print("Navigation Application Demo")
656    print("=" * 70)
657
658    # Aircraft navigation scenario
659    lat = np.radians(40.0)
660    lon = np.radians(-74.0)
661    alt = 10000.0  # meters (cruise altitude)
662    decimal_year = 2025.0
663
664    print("\nAircraft navigation correction scenario:")
665    print(f"  Position: 40°N, 74°W")
666    print(f"  Altitude: {alt:.0f} m")
667
668    # Magnetic declination for compass correction
669    dec = magnetic_declination(lat, lon, alt, decimal_year)
670    print(f"\n  Magnetic declination: {np.degrees(dec):.2f}°")
671    print(f"  → Add {np.degrees(dec):.2f}° to magnetic heading for true heading")
672
673    # Gravity for inertial navigation
674    g_result = gravity_wgs84(lat, lon, alt)
675    print(f"\n  Local gravity: {g_result.magnitude:.6f} m/s²")
676    print(f"  → Used for vertical channel in INS")
677
678    # Deflection of vertical (simplified)
679    print(f"\n  Gravity deflection (N): {g_result.g_north*1e6:.2f} µrad")
680    print(f"  → Correction for inertial alignment")
681
682    # Geoid for altitude reference (J2 approximation - latitude only)
683    N = geoid_height_j2(lat)
684    print(f"\n  Geoid undulation (J2 approx): {N:.1f} m")
685    print(f"  → Ellipsoid alt = GPS alt, Orthometric alt = GPS alt - N")
686
687
688def main():
689    """Run all demonstrations."""
690    print("\n" + "#" * 70)
691    print("# PyTCL Geophysical Models Example")
692    print("#" * 70)
693
694    # Gravity models
695    demo_normal_gravity()
696    demo_gravity_models()
697    demo_geoid()
698    demo_gravity_anomalies()
699    demo_tidal_effects()
700
701    # Magnetic field models
702    demo_magnetic_field()
703    demo_magnetic_properties()
704
705    # Application
706    demo_navigation_application()
707
708    print("\n" + "=" * 70)
709    print("Example complete!")
710    if SHOW_PLOTS:
711        print("Plots saved: geophysical_gravity_latitude.html,")
712        print("             geophysical_gravity_altitude.html,")
713        print("             geophysical_tides.html,")
714        print("             geophysical_magnetic_field.html")
715    print("=" * 70)
716
717
718if __name__ == "__main__":
719    main()

Running the Example

python examples/geophysical_models.py

See Also

  • ins_gnss_navigation - INS/GNSS integration

  • navigation_geodesy - Geodetic calculations

  • Magnetism Models - Magnetic field details