DiffusionGarnet.Domain โ€” Type
Domain(IC::InitialConditions2D, 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.

DiffusionGarnet.Domain โ€” Type
Domain(IC::InitialConditions3D, 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.

DiffusionGarnet.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::InitialConditionsSpherical, 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.

DiffusionGarnet.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.InitialConditions1D โ€” Method
InitialConditions1D(CMg0::Array{<:Real, 1}, CFe0::Array{<:Real, 1}, CMn0::Array{<:Real, 1}, Lx::Unitful.Length, 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 the Lxandtfinal`` to ยตm and Myr respectively.

DiffusionGarnet.InitialConditions2D โ€” Method
InitialConditions2D(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.InitialConditions3D โ€” Method
InitialConditions3D(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.InitialConditionsSpherical โ€” Method
InitialConditionsSpherical(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.save_data โ€” Method
save_data(integrator)

Callback function used to save major element compositions to an HDF5 file at a specific timestep.

DiffusionGarnet.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(), kwargs...)

Solve the coupled major element diffusion equations for a given domain using finite differences for the discretisation in space and return a solution type variable.

The time discretisation is based on the ROCK2 method, a stabilized explicit method (Adbdulle and Medovikov, 2001 ; https://doi.org/10.1007/s002110100292) using OrdinaryDiffEq.jl.

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 in Myr.

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/).

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/commonsolveropts/).

DiffusionGarnet.simulate โ€” Method
simulate(domain::Domain1D; path_save=nothing, solver=ROCK2(), abstol=1e-8, reltol=1e-6, kwargs...)

Solve the coupled major element diffusion equations in 1D. Save all timesteps in the output solution type variable by default.

DiffusionGarnet.simulate โ€” Method
simulate(domain::Domain2D; path_save=nothing, solver=ROCK2(), kwargs...)

Solve the coupled major element diffusion equations in 2D. By default, save only the first and last timestep in the output solution type variable by default.

DiffusionGarnet.simulate โ€” Method
simulate(domain::Domain3D; path_save=nothing, solver=ROCK2(), 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.

DiffusionGarnet.simulate โ€” Method
simulate(domain::DomainSpherical; path_save=nothing, solver=ROCK2(), abstol=1e-8, reltol=1e-6, kwargs...)

Solve the coupled major element diffusion equations in spherical coordinates. Save all timesteps in the output solution type variable by default.

DiffusionGarnet.update_diffusion_coef โ€” Method
update_diffusion_coef(integrator)

Callback function to update the diffusion coefficients at a given time from a new pressure and temperature. To use with the callback PresetTimeCallback (https://docs.sciml.ai/stable/basics/callbacks/#PresetTimeCallback-1).

Follows the syntax of callback functions defined by DiffEqCallbacks.jl (https://docs.sciml.ai/DiffEqCallbacks/stable/).