/
example_build_planet.py
181 lines (146 loc) · 8 KB
/
example_build_planet.py
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
# This file is part of BurnMan - a thermoelastic and thermodynamic toolkit for the Earth and Planetary Sciences
# Copyright (C) 2012 - 2015 by the BurnMan team, released under the GNU
# GPL v2 or later.
"""
example_build_planet
--------------------
For Earth we have well-constrained one-dimensional density models. This allows us to
calculate pressure as a function of depth. Furthermore, petrologic data and assumptions
regarding the convective state of the planet allow us to estimate the temperature.
For planets other than Earth we have much less information, and in particular we
know almost nothing about the pressure and temperature in the interior. Instead, we tend
to have measurements of things like mass, radius, and moment-of-inertia. We would like
to be able to make a model of the planet's interior that is consistent with those
measurements.
However, there is a difficulty with this. In order to know the density of the planetary
material, we need to know the pressure and temperature. In order to know the pressure,
we need to know the gravity profile. And in order to the the gravity profile, we need
to know the density. This is a nonlinear problem which requires us to iterate to find
a self-consistent solution.
This example allows the user to define layers of planets of known outer radius and self-
consistently solve for the density, pressure and gravity profiles. The calculation will
iterate until the difference between central pressure calculations are less than 1e-5.
The planet class in BurnMan (../burnman/planet.py) allows users to call multiple
properties of the model planet after calculations, such as the mass of an individual layer,
the total mass of the planet and the moment if inertia. See planets.py for information
on each of the parameters which can be called.
*Uses:*
* :doc:`mineral_database`
* :class:`burnman.Planet`
* :class:`burnman.Layer`
*Demonstrates:*
* setting up a planet
* computing its self-consistent state
* computing various parameters for the planet
* seismic comparison
"""
from __future__ import absolute_import
from __future__ import print_function
import matplotlib.pyplot as plt
import numpy as np
import burnman_path # adds the local burnman directory to the path
import burnman
import warnings
assert burnman_path # silence pyflakes warning
if __name__ == '__main__':
# FIRST: we must define the composition of the planet as individual layers.
# A layer is defined by 4 parameters: Name, min_depth, max_depth,and number of slices within the layer.
# Separately the composition and the temperature_mode need to set.
radius_planet = 6371.e3
# inner_core
inner_core = burnman.Layer('inner core', radii = np.linspace(0,1220.e3,10))
inner_core.set_material(burnman.minerals.other.Fe_Dewaele())
# The minerals that make up our core do not currently implement the thermal equation of state, so we will set the temperature at 300 K.
inner_core.set_temperature_mode('user-defined',
300.*np.ones_like(inner_core.radii))
# outer_core
outer_core = burnman.Layer('outer core', radii=np.linspace(1220.e3,3480.e3,10))
outer_core.set_material(burnman.minerals.other.Liquid_Fe_Anderson())
# The minerals that make up our core do not currently implement the thermal equation of state, so we will define the temperature at 300 K.
outer_core.set_temperature_mode('user-defined', 300.*np.ones_like(outer_core.radii))
# Next the Mantle.
lower_mantle_radii = np.linspace(3480.e3, 5711.e3, 101)
lower_mantle = burnman.Layer('lower mantle', radii=lower_mantle_radii)
lower_mantle.set_material(burnman.minerals.SLB_2011.mg_bridgmanite())
# Here we add a thermal boundary layer perturbation, assuming that the
# lower mantle has a Rayleigh number of 1.e7, and that there
# is an 840 K jump across the basal thermal boundary layer and a
# 0 K jump at the top of the lower mantle.
tbl_perturbation = burnman.BoundaryLayerPerturbation(radius_bottom=3480.e3,
radius_top=5711.e3,
rayleigh_number=1.e7,
temperature_change=840.,
boundary_layer_ratio=0.)
lower_mantle_tbl = tbl_perturbation.temperature(lower_mantle_radii)
lower_mantle.set_temperature_mode('perturbed-adiabatic',
temperatures=lower_mantle_tbl)
upper_mantle = burnman.Layer('upper mantle', radii = np.linspace(5711.e3, 6371e3, 10))
upper_mantle.set_material(burnman.minerals.SLB_2011.forsterite())
upper_mantle.set_temperature_mode('adiabatic', temperature_top = 1200.)
# Now we calculate the planet.
planet_zog = burnman.Planet('Planet Zog',
[inner_core, outer_core, lower_mantle, upper_mantle],
verbose=True)
print(planet_zog)
# Here we compute its state. Go BurnMan Go!
# (If we were to change composition of one of the layers, we would have to
# recompute the state)
planet_zog.make()
# Now we output the mass of the planet and moment of inertia
print('\nmass/Earth= {0:.3f}, moment of inertia factor= {1:.4f}'.format(planet_zog.mass / 5.97e24,
planet_zog.moment_of_inertia_factor))
# And here's the mass of the individual layers:
for layer in planet_zog.layers:
print('{0} mass fraction of planet {1:.3f}'.format(layer.name, layer.mass / planet_zog.mass))
print('')
# Let's get PREM to compare everything to as we are trying
# to imitate Earth
prem = burnman.seismic.PREM()
premradii = 6371.e3 - prem.internal_depth_list()
with warnings.catch_warnings(record=True) as w:
eval = prem.evaluate(['density', 'pressure', 'gravity', 'v_s', 'v_p'])
premdensity, prempressure, premgravity, premvs, premvp = eval
print(w[-1].message)
# Now let's plot everything up
# Optional prettier plotting
# plt.style.use('ggplot')
figure = plt.figure(figsize=(10, 12))
figure.suptitle(
'{0} has a mass {1:.3f} times that of Earth,\n'
'has an average density of {2:.1f} kg/m$^3$,\n'
'and a moment of inertia factor of {3:.4f}'.format(
planet_zog.name,
planet_zog.mass/5.97e24,
planet_zog.average_density,
planet_zog.moment_of_inertia_factor),
fontsize=20)
ax = [figure.add_subplot(4, 1, i) for i in range(1, 5)]
ax[0].plot(planet_zog.radii / 1.e3, planet_zog.density / 1.e3, 'k', linewidth=2.,
label=planet_zog.name)
ax[0].plot( premradii / 1.e3, premdensity / 1.e3, '--k', linewidth=1.,
label='PREM')
ax[0].set_ylim(0., (max(planet_zog.density) / 1.e3) + 1.)
ax[0].set_ylabel('Density ($\cdot 10^3$ kg/m$^3$)')
ax[0].legend()
# Make a subplot showing the calculated pressure profile
ax[1].plot(planet_zog.radii / 1.e3, planet_zog.pressure / 1.e9, 'b', linewidth=2.)
ax[1].plot(premradii / 1.e3, prempressure / 1.e9, '--b', linewidth=1.)
ax[1].set_ylim(0., (max(planet_zog.pressure) / 1e9) + 10.)
ax[1].set_ylabel('Pressure (GPa)')
# Make a subplot showing the calculated gravity profile
ax[2].plot(planet_zog.radii / 1.e3, planet_zog.gravity, 'g', linewidth=2.)
ax[2].plot(premradii / 1.e3, premgravity, '--g', linewidth=1.)
ax[2].set_ylabel('Gravity (m/s$^2)$')
ax[2].set_ylim(0., max(planet_zog.gravity) + 0.5)
# Make a subplot showing the calculated temperature profile
mask = planet_zog.temperature > 301.
ax[3].plot(planet_zog.radii[mask] / 1.e3,
planet_zog.temperature[mask], 'r', linewidth=2.)
ax[3].set_ylabel('Temperature ($K$)')
ax[3].set_xlabel('Radius (km)')
ax[3].set_ylim(0., max(planet_zog.temperature) + 100)
for i in range(3):
ax[i].set_xticklabels([])
for i in range(4):
ax[i].set_xlim(0., max(planet_zog.radii) / 1.e3)
plt.show()