DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.Domain
DiffusionGarnet.D_update!
DiffusionGarnet.D_update!
DiffusionGarnet.IC1DMajor
DiffusionGarnet.IC1DTrace
DiffusionGarnet.IC2DMajor
DiffusionGarnet.IC2DTrace
DiffusionGarnet.IC3DMajor
DiffusionGarnet.IC3DTrace
DiffusionGarnet.ICSphMajor
DiffusionGarnet.ICSphTrace
DiffusionGarnet.save_data
DiffusionGarnet.save_data_paraview
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.simulate
DiffusionGarnet.update_diffusion_coef
DiffusionGarnet.Domain Type
Domain(IC::InitialConditionsSphericalMajor, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to spherical initial conditions, define corresponding spherical domain. Assume that the center of the grain is on the left side.
sourceDiffusionGarnet.Domain Type
Domain(IC::IC3DMajor, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to 3D initial conditions, define corresponding 3D domain.
sourceDiffusionGarnet.Domain Type
Domain(IC::InitialConditions1D, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr"; bc_neumann::Tuple=(false, false))
When applied to 1D initial conditions, define corresponding 1D domain. bc_neumann
can be used to define Neumann boundary conditions on the left or right side of the domain if set to true.
DiffusionGarnet.Domain Type
Domain(IC::InitialConditionsSphericalTrace, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to spherical initial conditions, define corresponding spherical domain for trace element diffusion.
sourceDiffusionGarnet.Domain Type
Domain(IC::InitialConditions3DTrace, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to 3D initial conditions, define corresponding 3D domain for trace element diffusion.
sourceDiffusionGarnet.Domain Type
Domain(IC::IC2DMajor, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to 2D initial conditions, define corresponding 2D domain.
sourceDiffusionGarnet.Domain Type
Domain(IC::InitialConditions2DTrace, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to 2D initial conditions, define corresponding 2D domain for trace element diffusion.
sourceDiffusionGarnet.Domain Type
Domain(IC, T, P, time_update=0u"Myr")
Return a struct containing the struct IC
, the PT conditions T
and P
and the nondimensionalised parameters of the problem. time_update
can be used to update the PT conditions during the simulation. If so, T
and P
have to be of the same size as time_update
and correspond to the values of PT to update to.
DiffusionGarnet.Domain Type
Domain(IC::InitialConditions1DTrace, T::Union{Unitful.Temperature,Array{<:Unitful.Temperature{<:Real}, 1}}, P::Union{Unitful.Pressure,Array{<:Unitful.Pressure{<:Real}, 1}}, time_update::Union{Unitful.Time,Array{<:Unitful.Time{<:Real}, 1}}=0u"Myr")
When applied to 1D initial conditions, define corresponding 1D domain for trace element diffusion. bc_neumann
can be used to define Neumann boundary conditions on the left or right side of the domain if set to true.
DiffusionGarnet.D_update! Method
D_update!(Domain, IC::T, T_K, P_kbar, fO2) where T <: InitialConditionsMajor
Update the diffusion coefficients D0
based on the temperature T_K
, pressure P_kbar
and oxygen fugacity fO2
DiffusionGarnet.D_update! Method
D_update!(Domain, IC::T, T_K, P_kbar, fO2) where T <: InitialConditionsTrace
Update the diffusion coefficients D0
based on the temperature T_K
, pressure P_kbar
and oxygen fugacity fO2
DiffusionGarnet.IC1DMajor Method
IC1DMajor(CMg0::Array{<:Real}, CFe0::Array{<:Real}, CMn0::Array{<:Real}, Lx::Unitful.Length, x::AbstractArray{<:Unitful.Length}=range(0u"µm", length=size(CMg0, 1), stop=Lx), tfinal::Unitful.Time)
Return a structure containing the initial conditions for a 1D diffusion problem. CMg0, CFe0 and CMn0 need to be in molar fraction. Convert Lx
, x
, and tfinal
to µm, µm and Myr, respectively.
DiffusionGarnet.IC1DTrace Method
IC1DTrace(;C0::Array{<:Real, 1}, D, Lx::Unitful.Length, tfinal::Unitful.Time)
Return a structure containing the initial conditions for a 1D diffusion problem. C needs to be in µg/g. Convert the Lx
and tfinal
to cm and Myr respectively.
DiffusionGarnet.IC2DMajor Method
IC2DMajor(;CMg0::AbstractArray{<:Real, 2}, CFe0::AbstractArray{<:Real, 2}, CMn0::AbstractArray{<:Real, 2}, Lx::Unitful.Length, Ly::Unitful.Length, tfinal::Unitful.Time, grt_boundary::AbstractArray{<:Real, 2}=zeros(eltype(CMg0), size(CMg0)...))
Return a structure containing the initial conditions for a 2D diffusion problem. CMg0, CFe0 and CMn0 need to be in molar fraction. Convert Lx
, Ly
and tfinal
to µm, µm and Myr respectively. grt_boundary is a matrix of the same size as CMg0, CFe0 and CMn0. It is a binary matrix with 1 where the contour of the garnet is present and 0 otherwise. It is used to apply the Dirichlet boundary condition in the model. If not provided, it is equal to zero everywhere.
DiffusionGarnet.IC2DTrace Method
IC2DTrace(;C0::Array{<:Real, 2}, D, Lx::Unitful.Length, Ly::Unitful.Length, tfinal::Unitful.Time)
Return a structure containing the initial conditions for a 2D diffusion problem. C needs to be in µg/g. Convert Lx
, Ly
and tfinal
to µm, µm and Myr respectively.
DiffusionGarnet.IC3DMajor Method
IC3DMajor(;CMg0::AbstractArray{<:Real, 3}, CFe0::AbstractArray{<:Real, 3}, CMn0::AbstractArray{<:Real, 3}, Lx::Unitful.Length, Ly::Unitful.Length, Lz::Unitful.Length, tfinal::Unitful.Time, grt_boundary::Union{AbstractArray{<:Real, 3}, AbstractArray{<:Bool, 3}}=zeros(Bool, size(CMg0)...))
Return a structure containing the initial conditions for a 3D diffusion problem. CMg0, CFe0 and CMn0 need to be in molar fraction. Convert Lx
, Ly
, Lz
and tfinal
to µm, µm, µm and Myr respectively. grt_boundary is a matrix of the same size as CMg0, CFe0 and CMn0. It is a binary matrix with 1 where the contour of the garnet is present and 0 otherwise. It is used to apply the Dirichlet boundary condition in the model. If not provided, it is equal to zero everywhere.
DiffusionGarnet.IC3DTrace Method
IC3DTrace(;C0::Array{<:Real, 3}, D, Lx::Unitful.Length, Ly::Unitful.Length, Lz::Unitful.Length, tfinal::Unitful.Time)
Return a structure containing the initial conditions for a 3D diffusion problem. C needs to be in µg/g. Convert Lx
, Ly
, Lz
and tfinal
to µm, µm, µm and Myr respectively.
DiffusionGarnet.ICSphMajor Method
ICSphMajor(CMg0::AbstractArray{<:Real, 1}, CFe0::AbstractArray{<:Real, 1}, CMn0::AbstractArray{<:Real, 1}, Lr::Unitful.Length, tfinal::Unitful.Time)
Return a structure containing the initial conditions for a spherical diffusion problem. CMg0, CFe0 and CMn0 need to be in molar fraction. Convert Lr
and tfinal
to µm and Myr respectively.
DiffusionGarnet.ICSphTrace Method
ICSphTrace(;C0::Array{<:Real, 1}, D, Lr::Unitful.Length, tfinal::Unitful.Time)
Return a structure containing the initial conditions for a spherical diffusion problem. C needs to be in µg/g. Convert Lr
and tfinal
to µm and Myr respectively.
DiffusionGarnet.save_data Method
save_data(integrator)
Callback function used to save major element compositions to an HDF5 file at a specific timestep.
sourceDiffusionGarnet.save_data_paraview Method
save_data_paraview(integrator)
Callback function used to save data to an HDF5 file and produce an XMDF file. XMDF files can be used to read data on software like Paraview.
Note that data are saved in row major format. For column major, use save_data()
instead.
DiffusionGarnet.simulate Function
simulate(domain; path_save=nothing, solver=ROCK2(), p_extra = nothing, kwargs...)
Solve the coupled major element diffusion equations or the linear diffusion equation for a given domain using finite differences for the discretisation in space and return a solution type variable.
The default time discretisation is based on the ROCK2 method, a stabilized explicit method (Adbdulle and Medovikov, 2001) using OrdinaryDiffEq.
The solution type variable is following the format of OrdinaryDiffEq.jl (see https://docs.sciml.ai/DiffEqDocs/stable/basics/solution/), and can be used to plot the solution, and to extract the solution at a given time. The time of the solution is non-dimensional but can be converted back using the characteristic time (t_charact
contained in the Domain
structure).
path_save
is an optional argument, which can be used to define the path of the HDF5 output file. Default is to nothing.
solver
is an optional argument, which can be used to define the solver to use for the time discretisation. Default is the ROCK2 method. All other ODE solvers accepted as the ones from DifferentialEquations (see https://docs.sciml.ai/DiffEqDocs/stable/solvers/ode_solve/).
p_extra
is an optional argument that can be used to pass additional parameters to the ODEProblem. It is not used by default, but can be useful to write custom callbacks.
All other accepted arguments such as callback
or progress
are the same as those of the solve
function from DifferentialEquations (see https://docs.sciml.ai/DiffEqDocs/stable/basics/common_solver_opts/).
DiffusionGarnet.simulate Method
simulate(domain::Domain1DMajor; path_save=nothing, solver=ROCK2(), abstol=1e-8, reltol=1e-4, p_extra = nothing, kwargs...)
Solve the coupled major element diffusion equations in 1D. Save all timesteps in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::Domain1DTrace; path_save=nothing, solver=ROCK2(), abstol=1e-8, reltol=1e-6, kwargs...)
Solve the linear diffusion equation in 1D. Save all timesteps in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::Domain2DMajor; path_save=nothing, solver=ROCK2(), p_extra = nothing, kwargs...)
Solve the coupled major element diffusion equations in 2D. Save only the first and last timestep in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::Domain2DTrace; path_save=nothing, solver=ROCK2(), p_extra = nothing, kwargs...)
Solve the linear diffusion equations for trace elements in 2D. Save only the first and last timestep in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::Domain3DMajor; path_save=nothing, solver=ROCK2(), p_extra = nothing, kwargs...)
Solve the coupled major element diffusion equations in 3D. Save only the first and last timestep in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::Domain3DTrace; path_save=nothing, solver=ROCK2(), p_extra = nothing, kwargs...)
Solve the linear diffusion equations in 3D. Save only the first and last timestep in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::DomainSphericalMajor; path_save=nothing, solver=ROCK2(), abstol=1e-8, reltol=1e-4, p_extra = nothing, kwargs...)
Solve the coupled major element diffusion equations in spherical coordinates. Save all timesteps in the output solution type variable by default.
sourceDiffusionGarnet.simulate Method
simulate(domain::DomainSphericalTrace; path_save=nothing, solver=ROCK2(), abstol=1e-8, reltol=1e-6, p_extra=nothing, kwargs...)
Solve the linear diffusion equation for trace elements in spherical coordinates. Save all timesteps in the output solution type variable by default.
sourceDiffusionGarnet.update_diffusion_coef Method
update_diffusion_coef(integrator)
Callback function to update the diffusion coefficients during the simulation based on the time steps defined in time_update_ad
. Works for both major and trace elements.