API
AsteroidThermoPhysicalModels.CERES
— ConstantCeres
AsteroidThermoPhysicalModels.DIDYMOS
— ConstantAsteroid 65803 Didymos (1996 GT)
Physical parameters
:GM
: GM:M
: Mass:μ
: Standard gravitational parameter about Sun and Ryugu
Orbital elements
:a
: Semi-mojor axis [AU]:e
: Eccentricity [-]:I
: Inclination [deg]:Ω
: Longitude of the ascending node [deg]:ω
: Argument of periapsis [deg]:Φ
: # Mean anomaly [deg]
Spin parameters
:α
: Right ascension (RA) in equatorial coordinate system [deg]:δ
: Declination (Dec) in equatorial coordinate system [deg]:P
: Rotation period [h]
References
- Small Body Database Lookup: https://ssd.jpl.nasa.gov/tools/sbdb_lookup.html#/?sstr=didymos
- Naidu et al. (2020)
AsteroidThermoPhysicalModels.DIMORPHOS
— ConstantDimorphos
:a
: Semi-major axis of the mutual orbit with Dimorphos [m]:P
: Rotation periof of the mutual orbit with Dimorphos [h]
AsteroidThermoPhysicalModels.EARTH
— ConstantEarth
AsteroidThermoPhysicalModels.ERIS
— ConstantEris
AsteroidThermoPhysicalModels.JUPITER
— ConstantJupiter
AsteroidThermoPhysicalModels.MARS
— ConstantMars
AsteroidThermoPhysicalModels.MERCURY
— ConstantMercury
AsteroidThermoPhysicalModels.MOON
— ConstantMoon
AsteroidThermoPhysicalModels.NEPTUNE
— ConstantNeptune
AsteroidThermoPhysicalModels.PLUTO
— ConstantPluto
AsteroidThermoPhysicalModels.RYUGU
— ConstantAsteroid 162173 Ryugu
Physical parameters
:GM
: GM:M
: Mass:μ
: Standard gravitational parameter about Sun and Ryugu
Orbital elements
:a
: Semi-mojor axis [AU]:e
: Eccentricity [-]:I
: Inclination [deg]:Ω
: Longitude of the ascending node [deg]:ω
: Argument of periapsis [deg]:Φ
: # Mean anomaly [deg]
Spin parameters
:α
: Right ascension (RA) in equatorial coordinate system [deg]:δ
: Declination (Dec) in equatorial coordinate system [deg]:P
: Rotation period [h]
AsteroidThermoPhysicalModels.SATURN
— ConstantSaturn
AsteroidThermoPhysicalModels.URANUS
— ConstantUranus
AsteroidThermoPhysicalModels.VENUS
— ConstantVenus
AsteroidThermoPhysicalModels.BackwardEulerSolver
— TypeType of the backward Euler method:
- Implicit in time (Unconditionally stable in the heat conduction equation)
- First order in time
The BackwardEulerSolver
type has vectors for the tridiagonal matrix algorithm.
AsteroidThermoPhysicalModels.BinaryTPM
— Typestruct BinaryTPM <: ThermoPhysicalModel
Fields
pri
: TPM for the primarysec
: TPM for the secondaryMUTUAL_SHADOWING
: Flag to consider mutual shadowingMUTUAL_HEATING
: Flag to consider mutual heating
AsteroidThermoPhysicalModels.BinaryTPM
— MethodBinaryTPM(pri, sec; MUTUAL_SHADOWING=true, MUTUAL_HEATING=true) -> btpm
Construct a thermophysical model for a binary asteroid (BinaryTPM
).
AsteroidThermoPhysicalModels.BinaryTPMResult
— Typestruct BinaryTPMResult
Output data format for BinaryTPM
Fields
pri
: TPM result for the primarysec
: TPM result for the secondary
AsteroidThermoPhysicalModels.BinaryTPMResult
— MethodOuter constructor of BinaryTPMResult
Arguments
btpm
: Thermophysical model for a binary asteroidephem
: Ephemeridestimes_to_save
: Timesteps to save temperature (Common to both the primary and the secondary)face_ID_pri
: Face indices to save subsurface temperature of the primaryface_ID_sec
: Face indices to save subsurface temperature of the secondary
AsteroidThermoPhysicalModels.BoundaryCondition
— TypeAbstract type of a boundary condition for a heat conduction equation
AsteroidThermoPhysicalModels.CrankNicolsonSolver
— TypeType of the Crank-Nicolson method:
- Implicit in time (Unconditionally stable in the heat conduction equation)
- Second order in time
The CrankNicolsonSolver
type has vectors for the tridiagonal matrix algorithm.
References
- https://en.wikipedia.org/wiki/Crank–Nicolson_method
AsteroidThermoPhysicalModels.ForwardEulerSolver
— TypeType of the forward Euler method:
- Explicit in time
- First order in time
The ForwardEulerSolver
type includes a vector for the temperature at the next time step.
AsteroidThermoPhysicalModels.HeatConductionSolver
— TypeAbstract type of a solver for a heat conduction equation
AsteroidThermoPhysicalModels.InsulationBoundaryCondition
— TypeSingleton type of insulation boundary condition
AsteroidThermoPhysicalModels.IsothermalBoundaryCondition
— TypeType of isothermal boundary condition
AsteroidThermoPhysicalModels.NonUniformThermoParams
— Typestruct NonUniformThermoParams
Fields
P
: Cycle of thermal cycle (rotation period) [sec]l
: Thermal skin depth [m]Γ
: Thermal inertia [J ⋅ m⁻² ⋅ K⁻¹ ⋅ s⁻⁰⁵ (tiu)]A_B
: Bond albedoA_TH
: Albedo at thermal radiation wavelengthε
: Emissivityz_max
: Depth of the bottom of a heat conduction equation [m]Δz
: Depth step width [m]Nz
: Number of depth steps
AsteroidThermoPhysicalModels.RadiationBoundaryCondition
— TypeSingleton type of radiation boundary condition
AsteroidThermoPhysicalModels.ShapeModel
— TypeShapeModel
A polyhedral shape model of an asteroid.
Fields
nodes
: Vector of node positionsfaces
: Vector of vertex indices of facesface_centers
: Center position of each faceface_normals
: Normal vector of each faceface_areas
: Area of of each facevisiblefacets
: Vector of vector ofVisibleFacet
AsteroidThermoPhysicalModels.SingleTPM
— Typestruct SingleTPM <: ThermoPhysicalModel
Fields
shape
: Shape modelthermo_params
: Thermophysical parametersflux
: Flux on each face. Matrix of size (Number of faces, 3). Three components are:flux[:, 1]
: F_sun, flux of direct sunlightflux[:, 2]
: F_scat, flux of scattered lightflux[:, 3]
: F_rad, flux of thermal emission from surrounding surface
temperature
: Temperature matrix(Nz, Ns)
according to the number of depth cellsNz
and the number of facesNs
.face_forces
: Thermal force on each faceforce
: Thermal recoil force at body-fixed frame (Yarkovsky effect)torque
: Thermal recoil torque at body-fixed frame (YORP effect)SELF_SHADOWING
: Flag to consider self-shadowingSELF_HEATING
: Flag to consider self-heatingSOLVER
: Solver of heat conduction equationBC_UPPER
: Boundary condition at the upper boundaryBC_LOWER
: Boundary condition at the lower boundary
TO DO:
- roughness_maps ::ShapeModel[]
AsteroidThermoPhysicalModels.SingleTPM
— MethodSingleTPM(shape, thermo_params; SELF_SHADOWING=true, SELF_HEATING=true) -> stpm
Construct a thermophysical model for a single asteroid (SingleTPM
).
Arguments
shape
: Shape modelthermo_params
: Thermophysical parameters
Keyword arguments
SELF_SHADOWING
: Flag to consider self-shadowingSELF_HEATING
: Flag to consider self-heatingSOLVER
: Solver of heat conduction equationBC_UPPER
: Boundary condition at the upper boundaryBC_LOWER
: Boundary condition at the lower boundary
AsteroidThermoPhysicalModels.SingleTPMResult
— Typestruct SingleTPMResult
Output data format for SingleTPM
Fields
Saved at all time steps
times
: Timesteps, given the same vector asephem.time
[s]E_in
: Input energy per second on the whole surface [W]E_out
: Output enegey per second from the whole surface [W]E_cons
: Energy conservation ratio [-], ratio of total energy going out to total energy coming in in the last rotation cycleforce
: Thermal force on the asteroid [N]torque
: Thermal torque on the asteroid [N ⋅ m]
Saved only at the time steps desired by the user
times_to_save
: Timesteps to save temperature [s]depth_nodes
: Depths of the calculation nodes for 1-D heat conduction [m], a vector of sizeNz
surface_temperature
: Surface temperature [K], a matrix in size of(Ns, Nt)
.Ns
: Number of facesNt
: Number of time steps to save surface temperature
subsurface_temperature
: Temperature [K] as a function of depth [m] and time [s],Dict
with face ID as key and a matrix(Nz, Nt)
as an entry.Nz
: The number of the depth nodesNt
: The number of time steps to save temperature
AsteroidThermoPhysicalModels.SingleTPMResult
— MethodOuter constructor of SingleTPMResult
Arguments
stpm
: Thermophysical model for a single asteroidephem
: Ephemeridestimes_to_save
: Timesteps to save temperatureface_ID
: Face indices to save subsurface temperature
AsteroidThermoPhysicalModels.ThermoPhysicalModel
— TypeAbstract type of a thermophysical model
AsteroidThermoPhysicalModels.UniformThermoParams
— Typestruct UniformThermoParams
Fields
P
: Thermal cycle (rotation period) [sec]l
: Thermal skin depth [m]Γ
: Thermal inertia [J ⋅ m⁻² ⋅ K⁻¹ ⋅ s⁻⁰⁵ (tiu)]A_B
: Bond albedoA_TH
: Albedo at thermal radiation wavelengthε
: Emissivityz_max
: Depth of the bottom of a heat conduction equation [m]Δz
: Depth step width [m]Nz
: Number of depth steps
AsteroidThermoPhysicalModels.VisibleFacet
— Typestruct VisibleFacet
Index of an interfacing facet and its view factor
Fields
id
: Index of the interfacing facetf
: View factor from facet i to jd
: Distance from facet i to jd̂
: Normal vector from facet i to j
AsteroidThermoPhysicalModels.backward_euler!
— Methodbackward_euler!(stpm::SingleTPM, Δt)
Predict the temperature at the next time step by the backward Euler method.
- Implicit in time (Unconditionally stable in the heat conduction equation)
- First order in time
- Second order in space
In this function, the heat conduction equation is non-dimensionalized in time and length.
AsteroidThermoPhysicalModels.concave_spherical_segment
— Methodconcave_spherical_segment(r, h, xc, yc, x, y) -> z
Return the z-coordinate of a concave spherical segment.
Arguments
r
: Crater radiush
: Crater depthxc
: x-coordinate of crater centeryc
: y-coordinate of crater centerx
: x-coordinate where to calculate zy
: y-coordinate where to calculate z
AsteroidThermoPhysicalModels.concave_spherical_segment
— Methodconcave_spherical_segment(r, h; Nx=2^5, Ny=2^5, xc=0.5, yc=0.5) -> xs, ys, zs
Return (x, y, z) grid of a concave spherical segment.
Arguments
r
: Crater radiush
: Crater depthNx
: Number of nodes in the x-directionNy
: Number of nodes in the y-directionxc
: x-coordinate of crater centeryc
: y-coordinate of crater center
AsteroidThermoPhysicalModels.crank_nicolson!
— Methodcrank_nicolson!(stpm::SingleTPM, Δt)
Predict the temperature at the next time step by the Crank-Nicolson method.
- Implicit in time (Unconditionally stable in the heat conduction equation)
- Second order in time
- Second order in space
In this function, the heat conduction equation is non-dimensionalized in time and length.
AsteroidThermoPhysicalModels.crater_curvature_radius
— Methodcrater_curvature_radius(r, h) -> R
Return the curvature radius of a concave spherical segment.
Arguments
r
: Crater radiush
: Crater depth
AsteroidThermoPhysicalModels.energy_in
— Methodenergy_in(stpm::SingleTPM) -> E_in
Input energy per second on the whole surface [W]
AsteroidThermoPhysicalModels.energy_out
— Methodenergy_out(stpm::SingleTPM) -> E_out
Output enegey per second from the whole surface [W]
Arguments
stpm
: Thermophysical model for a single asteroid
AsteroidThermoPhysicalModels.export_TPM_results
— Methodexport_TPM_results(dirpath, result::BinaryTPMResult)
Export the result of BinaryTPM
to CSV files. The output files are saved in the following directory structure:
dirpath
├── pri
│ ├── physical_quantities.csv
│ ├── subsurface_temperature.csv
│ └── surface_temperature.csv
└── sec
├── physical_quantities.csv
├── subsurface_temperature.csv
└── surface_temperature.csv
Arguments
dirpath
: Path to the directory to save CSV files.result
: Output data format forBinaryTPM
AsteroidThermoPhysicalModels.export_TPM_results
— Methodexport_TPM_results(dirpath, result::SingleTPMResult)
Export the result of SingleTPM
to CSV files. The output files are saved in the following directory structure:
dirpath
├── physical_quantities.csv
├── subsurface_temperature.csv
└── surface_temperature.csv
Arguments
dirpath
: Path to the directory to save CSV files.result
: Output data format forSingleTPM
AsteroidThermoPhysicalModels.find_visiblefacets!
— Methodfind_visiblefacets!(obs::Facet, facets)
Find facets that is visible from the facet where the observer is located.
Parameters
obs
: Facet where the observer standsfacets
: Array ofFacet
AsteroidThermoPhysicalModels.flux_total
— Methodflux_total(A_B, A_TH, F_sun, F_scat, F_rad) -> F_total
Total energy absorbed by the surface [W/m²]
Arguments
A_B
: Bond albedo [-]A_TH
: Albedo at thermal infrared wavelength [-]F_sun
: Flux of direct sunlight [W/m²]F_scat
: Flux of scattered light [W/m²]F_rad
: Flux of thermal radiation from surrounding surface [W/m²]
AsteroidThermoPhysicalModels.forward_euler!
— Methodforward_euler!(stpm::SingleTPM, Δt)
Predict the temperature at the next time step by the forward Euler method.
- Explicit in time
- First order in time
In this function, the heat conduction equation is non-dimensionalized in time and length.
Arguments
stpm
: Thermophysical model for a single asteroidΔt
: Time step [sec]
AsteroidThermoPhysicalModels.grid_to_faces
— Methodgrid_to_faces(xs::AbstractVector, ys::AbstractVector, zs::AbstractMatrix) -> nodes, faces
Convert a regular grid (x, y) and corresponding z-coordinates to triangular facets
| ⧹| ⧹| ⧹|
j+1 ・–C–D–・ |⧹ |⧹ |⧹ | | ⧹| ⧹| ⧹| j ・–A–B–・ |⧹ |⧹ |⧹ | i i+1
Arguments
xs::AbstractVector
: x-coordinates of grid points (should be sorted)ys::AbstractVector
: y-coordinates of grid points (should be sorted)zs::AbstractMatrix
: z-coordinates of grid points
AsteroidThermoPhysicalModels.init_temperature!
— Methodinit_temperature!(btpm::BinaryTPM, T₀::Real)
Initialize all temperature cells at the given temperature T₀
Arguments
btpm
: Thermophysical model for a binary asteroidT₀
: Initial temperature of all cells [K]
AsteroidThermoPhysicalModels.init_temperature!
— Methodinit_temperature!(stpm::SingleTPM, T₀::Real)
Initialize all temperature cells at the given temperature T₀
Arguments
stpm
: Thermophysical model for a single asteroidT₀
: Initial temperature of all cells [K]
AsteroidThermoPhysicalModels.isilluminated
— Methodisilluminated(shape::ShapeModel, r☉::StaticVector{3}, i::Integer) -> Bool
Return if the i
-th face of the shape
model is illuminated by the direct sunlight or not
Arguments
shape
: Shape model of an asteroidr☉
: Sun's position in the asteroid-fixed frame, which doesn't have to be normalized.i
: Index of the face to be checked
AsteroidThermoPhysicalModels.load_shape_grid
— Methodload_shape_grid(xs, ys, zs; scale=1.0, find_visible_facets=false) -> shape
Convert a regular grid (x, y) to a shape model
Arguments
xs::AbstractVector
: x-coordinates of grid pointsys::AbstractVector
: y-coordinates of grid pointszs::AbstractMatrix
: z-coordinates of grid points
AsteroidThermoPhysicalModels.loadobj
— Methodloadobj(shapepath::String; scale=1, message=true) -> nodes, faces
AsteroidThermoPhysicalModels.mutual_heating!
— Methodmutual_heating!(btpm::BinaryTPM, rₛ, R₂₁)
Calculate the mutual heating between the primary and secondary asteroids.
Arguments
btpm
: Thermophysical model for a binary asteroidrₛ
: Position of the secondary relative to the primary (NOT normalized)R₂₁
: Rotation matrix from secondary to primary
TO DO
- Need to consider local horizon?
AsteroidThermoPhysicalModels.mutual_shadowing!
— Methodmutual_shadowing!(btpm::BinaryTPM, r☉, rₛ, R₂₁)
Detect eclipse events between the primary and secondary, and update the solar fluxes of the faces.
Arguments
btpm
: Thermophysical model for a binary asteroidr☉
: Position of the sun relative to the primary (NOT normalized)rₛ
: Position of the secondary relative to the primary (NOT normalized)R₂₁
: Rotation matrix from secondary to primary
AsteroidThermoPhysicalModels.polyhedron_volume
— Methodpolyhedron_volume(nodes, faces) -> vol
polyhedron_volume(shape::ShapeModel) -> vol
Calculate volume of a polyhedral
AsteroidThermoPhysicalModels.raycast
— Methodraycast(A, B, C, R) -> Bool
Intersection detection between ray R and triangle ABC. Note that the starting point of the ray is the origin (0, 0, 0).
AsteroidThermoPhysicalModels.raycast
— Methodraycast(A, B, C, R, O) -> Bool
Intersection detection between ray R and triangle ABC. Use when the starting point of the ray is an arbitrary point O
.
AsteroidThermoPhysicalModels.run_TPM!
— Methodrun_TPM!(btpm::BinaryTPM, ephem, savepath)
Run TPM for a binary asteroid.
Arguments
btpm
: Thermophysical model for a binary asteroidephem
: Ephemeridestime
: Ephemeris timessun1
: Sun's position in the primary's framesun2
: Sun's position in the secondary's framesec
: Secondary's position in the primary's frameP2S
: Rotation matrix from primary to secondary framesS2P
: Rotation matrix from secondary to primary frames
times_to_save
: Timesteps to save temperatureface_ID_pri
: Face indices where to save subsurface termperature for the primaryface_ID_sec
: Face indices where to save subsurface termperature for the secondary
Keyword arguments
show_progress
: Flag to show the progress meter
AsteroidThermoPhysicalModels.run_TPM!
— Methodrun_TPM!(stpm::SingleTPM, ephem, savepath)
Run TPM for a single asteroid.
Arguments
stpm
: Thermophysical model for a single asteroidephem
: Ephemeridesephem.time
: Ephemeris timesephem.sun
: Sun's position in the asteroid-fixed frame (Not normalized)
times_to_save
: Timesteps to save temperatureface_ID
: Face indices where to save subsurface termperature
Keyword arguments
show_progress
: Flag to show the progress meter
AsteroidThermoPhysicalModels.subsolar_temperature
— Methodsubsolar_temperature(r☉) -> Tₛₛ
Subsolar temperature [K] on an asteroid at a heliocentric distance r☉
[m], assuming radiative equilibrium with zero conductivity.
AsteroidThermoPhysicalModels.surface_temperature
— MethodReturn surface temperature of a single asteroid corrsponding to each face.
AsteroidThermoPhysicalModels.thermal_inertia
— Methodthermal_inertia(k, ρ, Cp) -> Γ
Arguments
k
: Thermal conductivity [W/m/K]ρ
: Material density [kg/m³]Cₚ
: Heat capacity [J/kg/K]
Return
Γ
: Thermal inertia [J ⋅ m⁻² ⋅ K⁻¹ ⋅ s⁻⁰⁵ (tiu)]
AsteroidThermoPhysicalModels.thermal_skin_depth
— Methodthermal_skin_depth(P, k, ρ, Cp) -> l_2π
Arguments
P
: Cycle of thermal cycle [sec]k
: Thermal conductivity [W/m/K]ρ
: Material density [kg/m³]Cₚ
: Heat capacity [J/kg/K]
Return
l_2π
: Thermal skin depth [m], as defined in Rozitis & Green (2011).
AsteroidThermoPhysicalModels.thermoparams
— Methodthermoparams(; A_B, A_TH, k, ρ, Cp, ε, t_begin, t_end, Nt, z_max, Nz, P)
AsteroidThermoPhysicalModels.tridiagonal_matrix_algorithm!
— Methodtridiagonal_matrix_algorithm!(a, b, c, d, x)
tridiagonal_matrix_algorithm!(stpm::SingleTPM)
Tridiagonal matrix algorithm to solve the heat conduction equation by the backward Euler and Crank-Nicolson methods.
| b₁ c₁ 0 ⋯ 0 | | x₁ | | d₁ |
| a₂ b₂ c₂ ⋯ 0 | | x₂ | | d₂ |
| 0 a₃ b₃ ⋯ 0 | | x₃ | = | d₃ |
| ⋮ ⋮ ⋮ ⋱ cₙ₋₁| | ⋮ | | ⋮ |
| 0 0 0 aₙ bₙ | | xₙ | | dₙ |
References
- https://en.wikipedia.org/wiki/Tridiagonalmatrixalgorithm
AsteroidThermoPhysicalModels.update_TPM_result!
— Methodupdate_TPM_result!(result::BinaryTPMResult, btpm::BinaryTPM, ephem, nₜ::Integer)
Save the results of TPM at the time step nₜ
to result
.
Arguments
result
: Output data format forBinaryTPM
btpm
: Thermophysical model for a binary asteroidephem
: Ephemeridesnₜ
: Time step
AsteroidThermoPhysicalModels.update_TPM_result!
— Methodupdate_TPM_result!(result::SingleTPMResult, stpm::SingleTPM, nₜ::Integer)
Save the results of TPM at the time step nₜ
to result
.
Arguments
result
: Output data format forSingleTPM
stpm
: Thermophysical model for a single asteroidnₜ
: Time step to save data
AsteroidThermoPhysicalModels.update_flux_rad_single!
— Methodupdate_flux_rad_single!(btpm::BinaryTPM)
Update flux of absorption of thermal radiation from surrounding surface. Single radiation-absorption is only considered, assuming albedo is close to zero at thermal infrared wavelength.
Arguments
btpm
: Thermophysical model for a binary asteroid
AsteroidThermoPhysicalModels.update_flux_rad_single!
— Methodupdate_flux_rad_single!(stpm::SingleTPM)
Update flux of absorption of thermal radiation from surrounding surface. Single radiation-absorption is only considered, assuming albedo is close to zero at thermal infrared wavelength.
Arguments
stpm
: Thermophysical model for a single asteroid
AsteroidThermoPhysicalModels.update_flux_scat_single!
— Methodupdate_flux_scat_single!(btpm::BinaryTPM)
Update flux of scattered sunlight, only considering single scattering.
Arguments
btpm
: Thermophysical model for a binary asteroid
AsteroidThermoPhysicalModels.update_flux_scat_single!
— Methodupdate_flux_scat_single!(stpm::SingleTPM)
Update flux of scattered sunlight, only considering single scattering.
Arguments
stpm
: Thermophysical model for a single asteroid
AsteroidThermoPhysicalModels.update_flux_sun!
— Methodupdate_flux_sun!(btpm::BinaryTPM, r☉₁::StaticVector{3}, r☉₂::StaticVector{3})
Arguments
btpm
: Thermophysical model for a binary asteroidr☉₁
: Sun's position in the body-fixed frame of the primary, which is not normalized.r☉₂
: Sun's position in the body-fixed frame of the secondary, which is not normalized.
AsteroidThermoPhysicalModels.update_flux_sun!
— Methodupdate_flux_sun!(stpm::SingleTPM, r̂☉::StaticVector{3}, F☉::Real)
Update solar irradiation flux on every face of a shape model.
shape
: Shape modelr̂☉
: Normalized vector indicating the direction of the sun in the body-fixed frameF☉
: Solar radiation flux [W/m²]
AsteroidThermoPhysicalModels.update_flux_sun!
— Methodupdate_flux_sun!(stpm::SingleTPM, r☉::StaticVector{3})
Update solar irradiation flux on every face of a shape model.
Arguments
stpm
: Thermophysical model for a single asteroidr☉
: Position of the sun in the body-fixed frame (NOT normalized)
AsteroidThermoPhysicalModels.update_lower_temperature!
— Methodupdate_bottom_temperature!(shape::ShapeModel)
Update the temperature of the bottom surface based on the boundary condition stpm.BC_LOWER
.
Arguments
stpm
: Thermophysical model for a single asteroid
AsteroidThermoPhysicalModels.update_surface_temperature!
— Methodupdate_surface_temperature!(T::AbstractVector, F_total::Real, k::Real, l::Real, Δz::Real, ε::Real)
Newton's method to update the surface temperature under radiation boundary condition.
Arguments
T
: 1-D array of temperaturesF_total
: Total energy absorbed by the facetΓ
: Thermal inertia [tiu]P
: Period of thermal cycle [sec]Δz̄
: Non-dimensional step in depth, normalized by thermal skin depthl
ε
: Emissivity
AsteroidThermoPhysicalModels.update_temperature!
— Methodupdate_temperature!(btpm::BinaryTPM, Δt)
Calculate the temperature for the next time step based on 1D heat conductivity equation.
Arguments
btpm
: Thermophysical model for a binary asteroidΔt
: Time step [sec]
AsteroidThermoPhysicalModels.update_temperature!
— Methodupdate_temperature!(stpm::SingleTPM, Δt)
Calculate the temperature for the next time step based on 1D heat conduction equation.
Arguments
stpm
: Thermophysical model for a single asteroidΔt
: Time step [sec]
AsteroidThermoPhysicalModels.update_thermal_force!
— Methodupdate_thermal_force!(btpm::BinaryTPM)
Calculate the thermal force and torque on every face and integrate them over all faces.
Arguments
btpm
: Thermophysical model for a binary asteroid
AsteroidThermoPhysicalModels.update_thermal_force!
— Methodupdate_thermal_force!(stpm::SingleTPM)
Calculate the thermal force and torque on every face and integrate them over all faces.
Arguments
stpm
: Thermophysical model for a single asteroid
AsteroidThermoPhysicalModels.update_upper_temperature!
— Methodupdate_upper_temperature!(stpm::SingleTPM, nₛ::Integer)
Update the temperature of the upper surface based on the boundary condition stpm.BC_UPPER
.
Arguments
stpm
: Thermophysical model for a single asteroidnₛ
: Index of the face of the shape model
AsteroidThermoPhysicalModels.view_factor
— Methodview_factor(cᵢ, cⱼ, n̂ᵢ, n̂ⱼ, aⱼ) -> fᵢⱼ, dᵢⱼ, d̂ᵢⱼ
View factor from facet i to j, assuming Lambertian emission.
- (i) fᵢⱼ (j)
- △ –> △
- cᵢ cⱼ : Center of each face
- n̂ᵢ n̂ⱼ : Normal vector of each face
aⱼ : Area of j-th face