Skip to content

Commit

Permalink
Zero-conductivity case (#128)
Browse files Browse the repository at this point in the history
* Update heat_conduction.jl

* Update .gitignore

* Update runtests.jl

* Create TPM_zero-conductivity.jl

* Update TPM_zero-conductivity.jl

* Update TPM_zero-conductivity.jl

* Update TPM_zero-conductivity.jl

* Update heat_conduction.jl
  • Loading branch information
MasanoriKanamaru authored Mar 29, 2024
1 parent c5a0d17 commit 4327ce3
Show file tree
Hide file tree
Showing 4 changed files with 107 additions and 12 deletions.
1 change: 1 addition & 0 deletions .gitignore
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,7 @@

/test/TPM_Ryugu/output
/test/TPM_non-uniform_thermoparams/output
/test/TPM_zero-conductivity/output
/test/TPM_Didymos/output

# Keep these shape models
Expand Down
42 changes: 30 additions & 12 deletions src/heat_conduction.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,8 @@
update_temperature!(stpm::SingleTPM, Δt)
Calculate the temperature for the next time step based on 1D heat conduction equation.
If the thermal inertia (conductivity) is zero, omit to solve the heat conduction equation.
The surface termperature is determined only by radiative equilibrium.
# Arguments
- `stpm` : Thermophysical model for a single asteroid
Expand Down Expand Up @@ -62,23 +64,39 @@ function forward_euler!(stpm::SingleTPM, Δt)
n_depth = size(T, 1)
n_face = size(T, 2)

for i_face in 1:n_face
P = stpm.thermo_params.P
Δz = stpm.thermo_params.Δz
l = (stpm.thermo_params.l isa Real ? stpm.thermo_params.l : stpm.thermo_params.l[i_face])
## Zero-conductivity (thermal inertia) case
if iszero(stpm.thermo_params.Γ)
for i_face in 1:n_face
A_B = (stpm.thermo_params.A_B isa Real ? stpm.thermo_params.A_B : stpm.thermo_params.A_B[i_face] )
A_TH = (stpm.thermo_params.A_TH isa Real ? stpm.thermo_params.A_TH : stpm.thermo_params.A_TH[i_face])
ε = (stpm.thermo_params.ε isa Real ? stpm.thermo_params.ε : stpm.thermo_params.ε[i_face] )
εσ = ε * σ_SB

λ = (Δt/P) / (Δz/l)^2 / 4π
λ 0.5 && error("The forward Euler method is unstable because λ = . This should be less than 0.5.")
F_sun, F_scat, F_rad = stpm.flux[i_face, :]
F_total = flux_total(A_B, A_TH, F_sun, F_scat, F_rad)

for i_depth in 2:(n_depth-1)
stpm.SOLVER.T[i_depth] = (1-2λ)*T[i_depth, i_face] + λ*(T[i_depth+1, i_face] + T[i_depth-1, i_face]) # Predict temperature at next time step
stpm.temperature[begin, i_face] = (F_total / εσ)^(1/4)
end
## Non-zero-conductivity (thermal inertia) case
else
for i_face in 1:n_face
P = stpm.thermo_params.P
Δz = stpm.thermo_params.Δz
l = (stpm.thermo_params.l isa Real ? stpm.thermo_params.l : stpm.thermo_params.l[i_face])

λ = (Δt/P) / (Δz/l)^2 / 4π
λ 0.5 && error("The forward Euler method is unstable because λ = . This should be less than 0.5.")

## Apply boundary conditions
update_upper_temperature!(stpm, i_face)
update_lower_temperature!(stpm)
for i_depth in 2:(n_depth-1)
stpm.SOLVER.T[i_depth] = (1-2λ)*T[i_depth, i_face] + λ*(T[i_depth+1, i_face] + T[i_depth-1, i_face]) # Predict temperature at next time step
end

T[:, i_face] .= stpm.SOLVER.T # Copy temperature at next time step
## Apply boundary conditions
update_upper_temperature!(stpm, i_face)
update_lower_temperature!(stpm)

T[:, i_face] .= stpm.SOLVER.T # Copy temperature at next time step
end
end
end

Expand Down
75 changes: 75 additions & 0 deletions test/TPM_zero-conductivity/TPM_zero-conductivity.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,75 @@
# See https://github.com/Astroshaper/Astroshaper-examples/tree/main/TPM_Ryugu for more information.
@testset "TPM_zero-conductivity" begin
DIR_OUTPUT = joinpath(@__DIR__, "output")
rm(DIR_OUTPUT; recursive=true, force=true)
mkpath(DIR_OUTPUT)

msg = """\n
┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
| Test: TPM_zero-conductivity |
┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
"""
println(msg)

##= Ephemerides =##
P = SPICE.convrt(8, "hours", "seconds") # Rotation period of the asteroid [s]
ncycles = 2 # Number of cycles to perform TPM
nsteps_in_cycle = 72 # Number of time steps in one rotation period

et_begin = 0.0 # Start time of TPM
et_end = et_begin + P * ncycles # End time of TPM
et_range = range(et_begin, et_end; length=nsteps_in_cycle*ncycles+1)

"""
- `time` : Ephemeris times
- `sun` : Sun's position in the asteroid-fixed frame
"""
ephem = (
time = collect(et_range),
sun = [inv(RotZ(2π * et / P)) * [SPICE.convrt(1, "au", "m"), 0, 0] for et in et_range],
)

##= Load obj file =##
path_obj = joinpath("shape", "ryugu_test.obj")
shape = AsteroidThermoPhysicalModels.load_shape_obj(path_obj; scale=1000, find_visible_facets=true)

##= Thermal properties: zero-conductivity case =##
thermo_params = AsteroidThermoPhysicalModels.thermoparams(
P = P,
l = 0.0,
Γ = 0.0,
A_B = 0.1,
A_TH = 0.0,
ε = 1.0,
z_max = 0.6,
n_depth = 11,
)

println(thermo_params)

##= Setting of TPM =##
stpm = AsteroidThermoPhysicalModels.SingleTPM(shape, thermo_params;
SELF_SHADOWING = true,
SELF_HEATING = false,
SOLVER = AsteroidThermoPhysicalModels.ForwardEulerSolver(thermo_params),
BC_UPPER = AsteroidThermoPhysicalModels.RadiationBoundaryCondition(),
BC_LOWER = AsteroidThermoPhysicalModels.InsulationBoundaryCondition(),
)
AsteroidThermoPhysicalModels.init_temperature!(stpm, 0)

##= Run TPM =##
times_to_save = ephem.time[end-nsteps_in_cycle:end] # Save temperature during the final rotation
face_ID = [1, 2] # Face indices to save subsurface temperature

result = AsteroidThermoPhysicalModels.run_TPM!(stpm, ephem, times_to_save, face_ID)

##= Save TPM result =##
@testset "Save TPM result" begin
AsteroidThermoPhysicalModels.export_TPM_results(DIR_OUTPUT, result)

@test isfile(joinpath(DIR_OUTPUT, "physical_quantities.csv"))
@test isfile(joinpath(DIR_OUTPUT, "subsurface_temperature.csv"))
@test isfile(joinpath(DIR_OUTPUT, "surface_temperature.csv"))
@test isfile(joinpath(DIR_OUTPUT, "thermal_force.csv"))
end
end
1 change: 1 addition & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -15,5 +15,6 @@ Aqua.test_all(AsteroidThermoPhysicalModels, ambiguities=false)
include("find_visiblefacets.jl")
include("TPM_Ryugu/TPM_Ryugu.jl")
include("TPM_non-uniform_thermoparams/TPM_non-uniform_thermoparams.jl")
include("TPM_zero-conductivity/TPM_zero-conductivity.jl")
include("TPM_Didymos/TPM_Didymos.jl")
include("heat_conduction_1D.jl")

0 comments on commit 4327ce3

Please sign in to comment.